Nektar++
ExpList.h
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: ExpList.h
4//
5// For more information, please see: http://www.nektar.info
6//
7// The MIT License
8//
9// Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10// Department of Aeronautics, Imperial College London (UK), and Scientific
11// Computing and Imaging Institute, University of Utah (USA).
12//
13// Permission is hereby granted, free of charge, to any person obtaining a
14// copy of this software and associated documentation files (the "Software"),
15// to deal in the Software without restriction, including without limitation
16// the rights to use, copy, modify, merge, publish, distribute, sublicense,
17// and/or sell copies of the Software, and to permit persons to whom the
18// Software is furnished to do so, subject to the following conditions:
19//
20// The above copyright notice and this permission notice shall be included
21// in all copies or substantial portions of the Software.
22//
23// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29// DEALINGS IN THE SOFTWARE.
30//
31// Description: Expansion list top class definition
32//
33///////////////////////////////////////////////////////////////////////////////
34
35#ifndef NEKTAR_LIBS_MULTIREGIONS_EXPLIST_H
36#define NEKTAR_LIBS_MULTIREGIONS_EXPLIST_H
37
52#include <tinyxml.h>
53
55{
56
57// Forward declarations
58class ExpList;
59class GlobalLinSys;
60class AssemblyMapDG;
61class AssemblyMapCG;
62class GlobalLinSysKey;
63class GlobalMatrix;
64
66{
71 eN
72};
73
75{
84};
85
87
88/// A map between global matrix keys and their associated block
89/// matrices.
90typedef std::map<GlobalMatrixKey, DNekScalBlkMatSharedPtr> BlockMatrixMap;
91/// A shared pointer to a BlockMatrixMap.
92typedef std::shared_ptr<BlockMatrixMap> BlockMatrixMapShPtr;
93/// Shared pointer to an ExpList object.
94typedef std::shared_ptr<ExpList> ExpListSharedPtr;
95
96/// Base class for all multi-elemental spectral/hp expansions.
97class ExpList : public std::enable_shared_from_this<ExpList>
98{
99public:
100 /// The default constructor using a type
102
103 /// The copy constructor.
105 const bool DeclareCoeffPhysArrays = true);
106
107 /// Constructor copying only elements defined in eIds.
109 const std::vector<unsigned int> &eIDs,
110 const bool DeclareCoeffPhysArrays = true,
111 const Collections::ImplementationType ImpType =
113
114 /// Generate an ExpList from a meshgraph \a graph and session file
118 const bool DeclareCoeffPhysArrays = true,
119 const std::string &var = "DefaultVar",
120 const Collections::ImplementationType ImpType =
122
123 /// Sets up a list of local expansions based on an expansion Map
126 const SpatialDomains::ExpansionInfoMap &expansions,
127 const bool DeclareCoeffPhysArrays = true,
128 const Collections::ImplementationType ImpType =
130
131 //---------------------------------------------------------
132 // Specialised constructors in ExpListConstructor.cpp
133 //---------------------------------------------------------
134 /// Specialised constructors for 0D Expansions
135 /// Wrapper around LocalRegion::PointExp - used in PrePacing.cpp
138
139 /// Generate expansions for the trace space expansions used in
140 /// DisContField.
143 const Array<OneD, const ExpListSharedPtr> &bndConstraint,
145 &bndCond,
146 const LocalRegions::ExpansionVector &locexp,
148 const LibUtilities::CommSharedPtr &comm,
149 const bool DeclareCoeffPhysArrays = true,
150 const std::string variable = "DefaultVar",
151 const Collections::ImplementationType ImpType =
153
154 /// Generate an trace ExpList from a meshgraph \a graph and session file
157 const LocalRegions::ExpansionVector &locexp,
159 const bool DeclareCoeffPhysArrays, const std::string variable,
160 const Collections::ImplementationType ImpType =
162
163 /// Constructor based on domain information only for 1D &
164 /// 2D boundary conditions
167 const SpatialDomains::CompositeMap &domain,
169 const bool DeclareCoeffPhysArrays = true,
170 const std::string variable = "DefaultVar",
171 bool SetToOneSpaceDimension = false,
173 const Collections::ImplementationType ImpType =
175
176 /// The default destructor.
178
179 /// Returns the total number of local degrees of freedom
180 /// \f$N_{\mathrm{eof}}=\sum_{e=1}^{{N_{\mathrm{el}}}}N^{e}_m\f$.
181 inline int GetNcoeffs(void) const;
182
183 /// Returns the total number of local degrees of freedom
184 /// for element eid
185 MULTI_REGIONS_EXPORT int GetNcoeffs(const int eid) const;
186
187 /// Returns the type of the expansion
189
190 /// Returns the type of the expansion
192
193 /// Evaulates the maximum number of modes in the elemental basis
194 /// order over all elements
195 inline int EvalBasisNumModesMax(void) const;
196
197 /// Returns the vector of the number of modes in the elemental
198 /// basis order over all elements.
200 void) const;
201
202 /// Returns the total number of quadrature points #m_npoints
203 /// \f$=Q_{\mathrm{tot}}\f$.
204 inline int GetTotPoints(void) const;
205
206 /// Returns the total number of quadrature points for eid's element
207 /// \f$=Q_{\mathrm{tot}}\f$.
208 inline int GetTotPoints(const int eid) const;
209
210 /// Returns the total number of quadrature points #m_npoints
211 /// \f$=Q_{\mathrm{tot}}\f$.
212 inline int GetNpoints(void) const;
213
214 /// Returns the total number of qudature points scaled by
215 /// the factor scale on each 1D direction
216 inline int Get1DScaledTotPoints(const NekDouble scale) const;
217
218 /// Sets the wave space to the one of the possible configuration
219 /// true or false
220 inline void SetWaveSpace(const bool wavespace);
221
222 /// Set Modified Basis for the stability analysis
223 inline void SetModifiedBasis(const bool modbasis);
224
225 /// This function returns the third direction expansion condition,
226 /// which can be in wave space (coefficient) or not
227 /// It is stored in the variable m_WaveSpace.
228 inline bool GetWaveSpace(void) const;
229
230 /// Set the \a i th value of \a m_phys to value \a val
231 inline void SetPhys(int i, NekDouble val);
232 /// Fills the array #m_phys
233 inline void SetPhys(const Array<OneD, const NekDouble> &inarray);
234 /// Sets the array #m_phys
235 inline void SetPhysArray(Array<OneD, NekDouble> &inarray);
236 /// This function manually sets whether the array of physical
237 /// values \f$\boldsymbol{u}_l\f$ (implemented as #m_phys) is
238 /// filled or not.
239 inline void SetPhysState(const bool physState);
240 /// This function indicates whether the array of physical values
241 /// \f$\boldsymbol{u}_l\f$ (implemented as #m_phys) is filled or
242 /// not.
243 inline bool GetPhysState(void) const;
244
245 /// multiply the metric jacobi and quadrature weights
247 const Array<OneD, const NekDouble> &inarray,
248 Array<OneD, NekDouble> &outarray);
249 /// Divided by the metric jacobi and quadrature weights
251 const Array<OneD, const NekDouble> &inarray,
252 Array<OneD, NekDouble> &outarray);
253 /// This function calculates the inner product of a function
254 /// \f$f(\boldsymbol{x})\f$ with respect to all \em local
255 /// expansion modes \f$\phi_n^e(\boldsymbol{x})\f$.
256 inline void IProductWRTBase(const Array<OneD, const NekDouble> &inarray,
257 Array<OneD, NekDouble> &outarray);
258 /// This function calculates the inner product of a function
259 /// \f$f(\boldsymbol{x})\f$ with respect to the derivative (in
260 /// direction \param dir) of all \em local expansion modes
261 /// \f$\phi_n^e(\boldsymbol{x})\f$.
263 const int dir, const Array<OneD, const NekDouble> &inarray,
264 Array<OneD, NekDouble> &outarray);
265
267 const Array<OneD, const NekDouble> &direction,
268 const Array<OneD, const NekDouble> &inarray,
269 Array<OneD, NekDouble> &outarray);
270
271 /// This function calculates the inner product of a
272 /// function \f$f(\boldsymbol{x})\f$ with respect to the
273 /// derivative of all \em local expansion modes
274 /// \f$\phi_n^e(\boldsymbol{x})\f$.
276 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
277 Array<OneD, NekDouble> &outarray);
278 /// This function elementally evaluates the forward transformation
279 /// of a function \f$u(\boldsymbol{x})\f$ onto the global
280 /// spectral/hp expansion.
281 inline void FwdTransLocalElmt(const Array<OneD, const NekDouble> &inarray,
282 Array<OneD, NekDouble> &outarray);
283 ///
284 inline void FwdTrans(const Array<OneD, const NekDouble> &inarray,
285 Array<OneD, NekDouble> &outarray);
287 const NekDouble alpha,
288 const NekDouble exponent,
289 const NekDouble cutoff);
290 /// This function elementally mulplies the coefficient space of
291 /// Sin my the elemental inverse of the mass matrix.
293 const Array<OneD, const NekDouble> &inarray,
294 Array<OneD, NekDouble> &outarray);
295 inline void MultiplyByInvMassMatrix(
296 const Array<OneD, const NekDouble> &inarray,
297 Array<OneD, NekDouble> &outarray);
299 const Array<OneD, const NekDouble> &inarray,
300 Array<OneD, NekDouble> &outarray)
301 {
303 BwdTrans(inarray, tmp);
304 IProductWRTBase(tmp, outarray);
305 }
306 /// Smooth a field across elements
307 inline void SmoothField(Array<OneD, NekDouble> &field);
308
309 /// Solve helmholtz problem
311 const Array<OneD, const NekDouble> &inarray,
312 Array<OneD, NekDouble> &outarray,
315 const MultiRegions::VarFactorsMap &varfactors =
318 const bool PhysSpaceForcing = true);
319
320 /// Solve Advection Diffusion Reaction
322 const Array<OneD, const NekDouble> &inarray,
323 Array<OneD, NekDouble> &outarray,
326 const MultiRegions::VarFactorsMap &varfactors =
329 const bool PhysSpaceForcing = true);
330
331 /// Solve Advection Diffusion Reaction
334 const Array<OneD, const NekDouble> &inarray,
335 Array<OneD, NekDouble> &outarray, const NekDouble lambda,
337 ///
339 const Array<OneD, const NekDouble> &inarray,
340 Array<OneD, NekDouble> &outarray);
341 /// This function elementally evaluates the backward transformation
342 /// of the global spectral/hp element expansion.
343 inline void BwdTrans(const Array<OneD, const NekDouble> &inarray,
344 Array<OneD, NekDouble> &outarray);
345 /// This function calculates the coordinates of all the elemental
346 /// quadrature points \f$\boldsymbol{x}_i\f$.
347 inline void GetCoords(
348 Array<OneD, NekDouble> &coord_0,
351
352 inline void GetCoords(
353 const int eid, Array<OneD, NekDouble> &coord_0,
356
357 // Homogeneous transforms
358 inline void HomogeneousFwdTrans(const int npts,
359 const Array<OneD, const NekDouble> &inarray,
360 Array<OneD, NekDouble> &outarray,
361 bool Shuff = true, bool UnShuff = true);
362 inline void HomogeneousBwdTrans(const int npts,
363 const Array<OneD, const NekDouble> &inarray,
364 Array<OneD, NekDouble> &outarray,
365 bool Shuff = true, bool UnShuff = true);
366 inline void DealiasedProd(const int num_dofs,
367 const Array<OneD, NekDouble> &inarray1,
368 const Array<OneD, NekDouble> &inarray2,
369 Array<OneD, NekDouble> &outarray);
370 inline void DealiasedDotProd(
371 const int num_dofs, const Array<OneD, Array<OneD, NekDouble>> &inarray1,
372 const Array<OneD, Array<OneD, NekDouble>> &inarray2,
373 Array<OneD, Array<OneD, NekDouble>> &outarray);
374 inline void GetBCValues(Array<OneD, NekDouble> &BndVals,
375 const Array<OneD, NekDouble> &TotField, int BndID);
378 Array<OneD, NekDouble> &outarray,
379 int BndID);
380 inline void NormVectorIProductWRTBase(
382 Array<OneD, NekDouble> &outarray);
383 /// Apply geometry information to each expansion.
385 /// Reset geometry information and reset matrices
387 {
388 v_Reset();
389 }
390 // not sure we really need these in ExpList
391 void WriteTecplotHeader(std::ostream &outfile, std::string var = "")
392 {
393 v_WriteTecplotHeader(outfile, var);
394 }
395 void WriteTecplotZone(std::ostream &outfile, int expansion = -1)
396 {
397 v_WriteTecplotZone(outfile, expansion);
398 }
399 void WriteTecplotField(std::ostream &outfile, int expansion = -1)
400 {
401 v_WriteTecplotField(outfile, expansion);
402 }
403 void WriteTecplotConnectivity(std::ostream &outfile, int expansion = -1)
404 {
405 v_WriteTecplotConnectivity(outfile, expansion);
406 }
407 MULTI_REGIONS_EXPORT void WriteVtkHeader(std::ostream &outfile);
408 MULTI_REGIONS_EXPORT void WriteVtkFooter(std::ostream &outfile);
409 void WriteVtkPieceHeader(std::ostream &outfile, int expansion,
410 int istrip = 0)
411 {
412 v_WriteVtkPieceHeader(outfile, expansion, istrip);
413 }
414 MULTI_REGIONS_EXPORT void WriteVtkPieceFooter(std::ostream &outfile,
415 int expansion);
416 void WriteVtkPieceData(std::ostream &outfile, int expansion,
417 std::string var = "v")
418 {
419 v_WriteVtkPieceData(outfile, expansion, var);
420 }
421
422 /// This function returns the dimension of the coordinates of the
423 /// element \a eid.
424 // inline
425 MULTI_REGIONS_EXPORT int GetCoordim(int eid);
426
427 /// Set the \a i th coefficiient in \a m_coeffs to value \a val
428 inline void SetCoeff(int i, NekDouble val);
429 /// Set the \a i th coefficiient in #m_coeffs to value \a val
430 inline void SetCoeffs(int i, NekDouble val);
431 /// Set the #m_coeffs array to inarray
432 inline void SetCoeffsArray(Array<OneD, NekDouble> &inarray);
433 /// This function returns the dimension of the shape of the
434 /// element \a eid.
435 // inline
437 /// This function returns (a reference to) the array
438 /// \f$\boldsymbol{\hat{u}}_l\f$ (implemented as #m_coeffs)
439 /// containing all local expansion coefficients.
440 inline const Array<OneD, const NekDouble> &GetCoeffs() const;
441 /// Impose Dirichlet Boundary Conditions onto Array
443 /// Fill Bnd Condition expansion from the values stored in expansion
444 inline void FillBndCondFromField(const Array<OneD, NekDouble> coeffs);
445 /// Fill Bnd Condition expansion in nreg from the values
446 /// stored in expansion
447 inline void FillBndCondFromField(const int nreg,
448 const Array<OneD, NekDouble> coeffs);
449 /// Gathers the global coefficients \f$\boldsymbol{\hat{u}}_g\f$
450 /// from the local coefficients \f$\boldsymbol{\hat{u}}_l\f$.
451 // inline
452 MULTI_REGIONS_EXPORT inline void LocalToGlobal(bool useComm = true);
454 const Array<OneD, const NekDouble> &inarray,
455 Array<OneD, NekDouble> &outarray, bool useComm = true);
456 /// Scatters from the global coefficients
457 /// \f$\boldsymbol{\hat{u}}_g\f$ to the local coefficients
458 /// \f$\boldsymbol{\hat{u}}_l\f$.
459 // inline
460 MULTI_REGIONS_EXPORT inline void GlobalToLocal(void);
461 /**
462 * This operation is evaluated as:
463 * \f{tabbing}
464 * \hspace{1cm} \= Do \= $e=$ $1, N_{\mathrm{el}}$ \\
465 * > > Do \= $i=$ $0,N_m^e-1$ \\
466 * > > > $\boldsymbol{\hat{u}}^{e}[i] = \mbox{sign}[e][i] \cdot
467 * \boldsymbol{\hat{u}}_g[\mbox{map}[e][i]]$ \\
468 * > > continue \\
469 * > continue
470 * \f}
471 * where \a map\f$[e][i]\f$ is the mapping array and \a
472 * sign\f$[e][i]\f$ is an array of similar dimensions ensuring the
473 * correct modal connectivity between the different elements (both
474 * these arrays are contained in the data member #m_locToGloMap). This
475 * operation is equivalent to the scatter operation
476 * \f$\boldsymbol{\hat{u}}_l=\mathcal{A}\boldsymbol{\hat{u}}_g\f$,
477 * where \f$\mathcal{A}\f$ is the
478 * \f$N_{\mathrm{eof}}\times N_{\mathrm{dof}}\f$ permutation matrix.
479 *
480 * @param inarray An array of size \f$N_\mathrm{dof}\f$
481 * containing the global degrees of freedom
482 * \f$\boldsymbol{x}_g\f$.
483 * @param outarray The resulting local degrees of freedom
484 * \f$\boldsymbol{x}_l\f$ will be stored in this
485 * array of size \f$N_\mathrm{eof}\f$.
486 */
488 const Array<OneD, const NekDouble> &inarray,
489 Array<OneD, NekDouble> &outarray);
490 /// Get the \a i th value (coefficient) of #m_coeffs
491 inline NekDouble GetCoeff(int i);
492 /// Get the \a i th value (coefficient) of #m_coeffs
493 inline NekDouble GetCoeffs(int i);
494 /// This function returns (a reference to) the array
495 /// \f$\boldsymbol{u}_l\f$ (implemented as #m_phys) containing the
496 /// function \f$u^{\delta}(\boldsymbol{x})\f$ evaluated at the
497 /// quadrature points.
498 // inline
500 /// This function calculates the \f$L_\infty\f$ error of the global
501 /// spectral/hp element approximation.
503 Linf(const Array<OneD, const NekDouble> &inarray,
505 /// This function calculates the \f$L_\infty\f$ error of the global
506 /// This function calculates the \f$L_2\f$ error with
507 /// respect to soln of the global
508 /// spectral/hp element approximation.
510 const Array<OneD, const NekDouble> &inarray,
512 {
513 return v_L2(inarray, soln);
514 }
515 /// Calculates the \f$H^1\f$ error of the global spectral/hp
516 /// element approximation.
518 H1(const Array<OneD, const NekDouble> &inarray,
520 /// Calculates the \f$H^1\f$ error of the global spectral/hp
521 /// element approximation.
522 /**
523 * The integration is evaluated locally, that is
524 * \f[\int
525 * f(\boldsymbol{x})d\boldsymbol{x}=\sum_{e=1}^{{N_{\mathrm{el}}}}
526 * \left\{\int_{\Omega_e}f(\boldsymbol{x})d\boldsymbol{x}\right\}, \f]
527 * where the integration over the separate elements is done by the
528 * function StdRegions#StdExpansion#Integral, which discretely
529 * evaluates the integral using Gaussian quadrature.
530 *
531 * Note that the array #m_phys should be filled with the values of the
532 * function \f$f(\boldsymbol{x})\f$ at the quadrature points
533 * \f$\boldsymbol{x}_i\f$.
534 *
535 * @return The value of the discretely evaluated integral
536 * \f$\int f(\boldsymbol{x})d\boldsymbol{x}\f$.
537 */
539 {
540 ASSERTL1(m_physState == true, "local physical space is not true ");
541 return Integral(m_phys);
542 }
543 /**
544 * The integration is evaluated locally, that is
545 * \f[\int
546 * f(\boldsymbol{x})d\boldsymbol{x}=\sum_{e=1}^{{N_{\mathrm{el}}}}
547 * \left\{\int_{\Omega_e}f(\boldsymbol{x})d\boldsymbol{x}\right\}, \f]
548 * where the integration over the separate elements is done by the
549 * function StdRegions#StdExpansion#Integral, which discretely
550 * evaluates the integral using Gaussian quadrature.
551 *
552 * @param inarray An array of size \f$Q_{\mathrm{tot}}\f$
553 * containing the values of the function
554 * \f$f(\boldsymbol{x})\f$ at the quadrature
555 * points \f$\boldsymbol{x}_i\f$.
556 * @return The value of the discretely evaluated integral
557 * \f$\int f(\boldsymbol{x})d\boldsymbol{x}\f$.
558 */
560 {
561 return v_Integral(inarray);
562 }
564 {
565 return v_VectorFlux(inarray);
566 }
567 /// This function calculates the energy associated with
568 /// each one of the modesof a 3D homogeneous nD expansion
570 {
571 return v_HomogeneousEnergy();
572 }
573
574 /// This function sets the Spectral Vanishing Viscosity
575 /// in homogeneous1D expansion.
577 {
579 }
580
581 /// This function returns a vector containing the wave
582 /// numbers in z-direction associated
583 /// with the 3D homogenous expansion. Required if a
584 /// parellelisation is applied in the Fourier direction
586 {
587 return v_GetZIDs();
588 }
589
590 /// This function returns the transposition class
591 /// associated with the homogeneous expansion.
593 {
594 return v_GetTransposition();
595 }
596
597 /// This function returns the Width of homogeneous direction
598 /// associated with the homogeneous expansion.
600 {
601 return v_GetHomoLen();
602 }
603
604 /// This function sets the Width of homogeneous direction
605 /// associated with the homogeneous expansion.
606 void SetHomoLen(const NekDouble lhom)
607 {
608 return v_SetHomoLen(lhom);
609 }
610
611 /// This function returns a vector containing the wave
612 /// numbers in y-direction associated
613 /// with the 3D homogenous expansion. Required if a
614 /// parellelisation is applied in the Fourier direction
616 {
617 return v_GetYIDs();
618 }
619
620 /// This function interpolates the physical space points in
621 /// \a inarray to \a outarray using the same points defined in the
622 /// expansion but where the number of points are rescaled
623 /// by \a 1DScale
625 const Array<OneD, NekDouble> &inarray,
626 Array<OneD, NekDouble> &outarray)
627 {
628 v_PhysInterp1DScaled(scale, inarray, outarray);
629 }
630
631 /// This function Galerkin projects the physical space points in
632 /// \a inarray to \a outarray where inarray is assumed to
633 /// be defined in the expansion but where the number of
634 /// points are rescaled by \a 1DScale
636 const Array<OneD, NekDouble> &inarray,
637 Array<OneD, NekDouble> &outarray)
638 {
639 v_PhysGalerkinProjection1DScaled(scale, inarray, outarray);
640 }
641
642 /// This function returns the number of elements in the expansion.
643 inline int GetExpSize(void);
644
645 /// This function returns the number of elements in the
646 /// expansion which may be different for a homogeoenous extended
647 /// expansionp.
648 inline size_t GetNumElmts(void)
649 {
650 return v_GetNumElmts();
651 }
652
653 /// This function returns the vector of elements in the expansion.
654 inline const std::shared_ptr<LocalRegions::ExpansionVector> GetExp() const;
655
656 /// This function returns (a shared pointer to) the local elemental
657 /// expansion of the \f$n^{\mathrm{th}}\f$ element.
658 inline LocalRegions::ExpansionSharedPtr &GetExp(int n) const;
659
660 /// This function returns (a shared pointer to) the local elemental
661 /// expansion of the \f$n^{\mathrm{th}}\f$ element given a global
662 /// geometry ID.
664
665 /// This function returns (a shared pointer to) the local elemental
666 /// expansion containing the arbitrary point given by \a gloCoord.
668 const Array<OneD, const NekDouble> &gloCoord);
669
670 /**
671 * @brief This function returns the index of the local elemental
672 * expansion containing the arbitrary point given by \a gloCoord,
673 * within a distance tolerance of tol.
674 *
675 * If returnNearestElmt is true and no element contains the point,
676 * this function returns the nearest element whose bounding box contains
677 * the point. The bounding box has a 10% margin in each direction.
678 *
679 * @param gloCoord (input) coordinate in physics space
680 * @param locCoords (output) local coordinate xi in the returned element
681 * @param tol distance tolerance to judge if a point lies in an
682 * element
683 * @param returnNearestElmt if true and no element contains this point, the
684 * nearest element whose bounding box contains this
685 * point is returned
686 * @param cachedId an initial guess of the most possible element index
687 * @param maxDistance if returnNearestElmt is set as true, the nearest
688 * element will be returned. But the distance of the
689 * nearest element and this point should be <=
690 * maxDistance.
691 *
692 * @return element index; if no element is found, -1 is returned.
693 */
695 const Array<OneD, const NekDouble> &gloCoord, NekDouble tol = 0.0,
696 bool returnNearestElmt = false, int cachedId = -1,
697 NekDouble maxDistance = 1e6);
698
699 /** This function returns the index and the Local
700 * Cartesian Coordinates \a locCoords of the local
701 * elemental expansion containing the arbitrary point
702 * given by \a gloCoords.
703 **/
705 const Array<OneD, const NekDouble> &gloCoords,
706 Array<OneD, NekDouble> &locCoords, NekDouble tol = 0.0,
707 bool returnNearestElmt = false, int cachedId = -1,
708 NekDouble maxDistance = 1e6);
709
710 /** This function return the expansion field value
711 * at the coordinates given as input.
712 **/
715 const Array<OneD, const NekDouble> &phys);
716
717 /// Get the start offset position for a local contiguous list of
718 /// coeffs correspoinding to element n.
719 inline int GetCoeff_Offset(int n) const;
720
721 /// Get the start offset position for a local contiguous list of
722 /// quadrature points in a full array correspoinding to element n.
723 inline int GetPhys_Offset(int n) const;
724
725 /// This function returns (a reference to) the array
726 /// \f$\boldsymbol{\hat{u}}_l\f$ (implemented as #m_coeffs)
727 /// containing all local expansion coefficients.
729 /// This function returns (a reference to) the array
730 /// \f$\boldsymbol{u}_l\f$ (implemented as #m_phys) containing the
731 /// function \f$u^{\delta}(\boldsymbol{x})\f$ evaluated at the
732 /// quadrature points.
734 inline void PhysDeriv(Direction edir,
735 const Array<OneD, const NekDouble> &inarray,
737 /// This function discretely evaluates the derivative of a function
738 /// \f$f(\boldsymbol{x})\f$ on the domain consisting of all
739 /// elements of the expansion.
740 inline void PhysDeriv(
741 const Array<OneD, const NekDouble> &inarray,
745 inline void PhysDeriv(const int dir,
746 const Array<OneD, const NekDouble> &inarray,
748
749 inline void Curl(Array<OneD, Array<OneD, NekDouble>> &Vel,
751
752 inline void CurlCurl(Array<OneD, Array<OneD, NekDouble>> &Vel,
754 inline void PhysDirectionalDeriv(
755 const Array<OneD, const NekDouble> &direction,
756 const Array<OneD, const NekDouble> &inarray,
757 Array<OneD, NekDouble> &outarray);
758 inline void GetMovingFrames(const SpatialDomains::GeomMMF MMFdir,
759 const Array<OneD, const NekDouble> &CircCentre,
760 Array<OneD, Array<OneD, NekDouble>> &outarray);
761 // functions associated with DisContField
764 /// Get the weight value for boundary conditions
766 /// Set the weight value for boundary conditions
767 inline void SetBndCondBwdWeight(const int index, const NekDouble value);
768 inline std::shared_ptr<ExpList> &UpdateBndCondExpansion(int i);
769 inline void Upwind(const Array<OneD, const NekDouble> &Vn,
773 inline void Upwind(const Array<OneD, const Array<OneD, NekDouble>> &Vec,
777 /**
778 * Return a reference to the trace space associated with this
779 * expansion list.
780 */
781 inline std::shared_ptr<ExpList> &GetTrace();
782 inline std::shared_ptr<AssemblyMapDG> &GetTraceMap(void);
783 inline const Array<OneD, const int> &GetTraceBndMap(void);
784 inline void GetNormals(Array<OneD, Array<OneD, NekDouble>> &normals);
785 /// Get the length of elements in boundary normal direction
787 Array<OneD, NekDouble> &lengthsFwd, Array<OneD, NekDouble> &lengthsBwd);
788 /// Get the weight value for boundary conditions
789 /// for boundary average and jump calculations
791 Array<OneD, NekDouble> &weightJump);
793 Array<OneD, NekDouble> &outarray);
796 Array<OneD, NekDouble> &outarray);
799 inline void GetFwdBwdTracePhys(const Array<OneD, const NekDouble> &field,
802 bool FillBnd = true,
803 bool PutFwdInBwdOnBCs = false,
804 bool DoExchange = true);
805 inline void FillBwdWithBoundCond(const Array<OneD, NekDouble> &Fwd,
807 bool PutFwdInBwdOnBCs = false);
808 /// Add Fwd and Bwd value to field,
809 /// a reverse procedure of GetFwdBwdTracePhys
816 {
817 v_AddTraceQuadPhysToOffDiag(Fwd, Bwd, field);
818 }
821 Array<OneD, NekDouble> &locTraceFwd,
822 Array<OneD, NekDouble> &locTraceBwd);
823 /// Fill Bwd with boundary conditions
824 inline void FillBwdWithBwdWeight(Array<OneD, NekDouble> &weightave,
825 Array<OneD, NekDouble> &weightjmp);
826 /// Copy and fill the Periodic boundaries
827 inline void PeriodicBwdCopy(const Array<OneD, const NekDouble> &Fwd,
829 inline const std::vector<bool> &GetLeftAdjacentFaces(void) const;
830 inline void ExtractTracePhys(Array<OneD, NekDouble> &outarray);
831 inline void ExtractTracePhys(const Array<OneD, const NekDouble> &inarray,
832 Array<OneD, NekDouble> &outarray);
837 inline void EvaluateBoundaryConditions(
838 const NekDouble time = 0.0, const std::string varName = "",
841 // Routines for continous matrix solution
842 /// This function calculates the result of the multiplication of a
843 /// matrix of type specified by \a mkey with a vector given by \a
844 /// inarray.
846 const GlobalMatrixKey &gkey,
847 const Array<OneD, const NekDouble> &inarray,
848 Array<OneD, NekDouble> &outarray);
849 inline void SetUpPhysNormals();
850 inline void GetBoundaryToElmtMap(Array<OneD, int> &ElmtID,
851 Array<OneD, int> &EdgeID);
852 virtual void GetBndElmtExpansion(int i, std::shared_ptr<ExpList> &result,
853 const bool DeclareCoeffPhysArrays = true);
854
855 inline void ExtractElmtToBndPhys(int i, const Array<OneD, NekDouble> &elmt,
856 Array<OneD, NekDouble> &boundary);
857
858 inline void ExtractPhysToBndElmt(int i,
860 Array<OneD, NekDouble> &bndElmt);
861
862 inline void ExtractPhysToBnd(int i,
865
866 inline void GetBoundaryNormals(
867 int i, Array<OneD, Array<OneD, NekDouble>> &normals);
868
870 std::vector<LibUtilities::FieldDefinitionsSharedPtr> &fielddef,
871 int NumHomoDir = 0,
874 std::vector<NekDouble> &HomoLen = LibUtilities::NullNekDoubleVector,
875 bool homoStrips = false,
876 std::vector<unsigned int> &HomoSIDs =
878 std::vector<unsigned int> &HomoZIDs =
880 std::vector<unsigned int> &HomoYIDs =
882
883 std::map<int, RobinBCInfoSharedPtr> GetRobinBCInfo()
884 {
885 return v_GetRobinBCInfo();
886 }
887
888 void GetPeriodicEntities(PeriodicMap &periodicVerts,
889 PeriodicMap &periodicEdges,
890 PeriodicMap &periodicFaces = NullPeriodicMap)
891 {
892 v_GetPeriodicEntities(periodicVerts, periodicEdges, periodicFaces);
893 }
894
895 std::vector<LibUtilities::FieldDefinitionsSharedPtr> GetFieldDefinitions()
896 {
897 return v_GetFieldDefinitions();
898 }
899
901 std::vector<LibUtilities::FieldDefinitionsSharedPtr> &fielddef)
902 {
903 v_GetFieldDefinitions(fielddef);
904 }
905
906 /// Append the element data listed in elements
907 /// fielddef->m_ElementIDs onto fielddata
909 std::vector<NekDouble> &fielddata)
910 {
911 v_AppendFieldData(fielddef, fielddata);
912 }
913 /// Append the data in coeffs listed in elements
914 /// fielddef->m_ElementIDs onto fielddata
916 std::vector<NekDouble> &fielddata,
918 {
919 v_AppendFieldData(fielddef, fielddata, coeffs);
920 }
921 /** \brief Extract the data in fielddata into the coeffs
922 * using the basic ExpList Elemental expansions rather
923 * than planes in homogeneous case
924 */
927 std::vector<NekDouble> &fielddata, std::string &field,
928 Array<OneD, NekDouble> &coeffs);
929 /** \brief Extract the data from fromField using
930 * fromExpList the coeffs using the basic ExpList
931 * Elemental expansions rather than planes in homogeneous
932 * case
933 */
935 const std::shared_ptr<ExpList> &fromExpList,
936 const Array<OneD, const NekDouble> &fromCoeffs,
937 Array<OneD, NekDouble> &toCoeffs);
938 // Extract data in fielddata into the m_coeffs_list for the 3D stability
939 // analysis (base flow is 2D)
942 std::vector<NekDouble> &fielddata, std::string &field,
944 std::unordered_map<int, int> zIdToPlane =
945 std::unordered_map<int, int>());
946 // Extract data from file fileName and put coefficents into array coefffs
948 const std::string &fileName, LibUtilities::CommSharedPtr comm,
949 const std::string &varName, Array<OneD, NekDouble> &coeffs);
951 const int ElementID, const NekDouble scalar1, const NekDouble scalar2,
952 Array<OneD, NekDouble> &outarray);
953 /// Returns a shared pointer to the current object.
954 std::shared_ptr<ExpList> GetSharedThisPtr()
955 {
956 return shared_from_this();
957 }
958 /// Returns the session object
959 std::shared_ptr<LibUtilities::SessionReader> GetSession() const
960 {
961 return m_session;
962 }
963 /// Returns the comm object
964 std::shared_ptr<LibUtilities::Comm> GetComm() const
965 {
966 return m_comm;
967 }
969 {
970 return m_graph;
971 }
972 // Wrapper functions for Homogeneous Expansions
974 {
975 return v_GetHomogeneousBasis();
976 }
977 std::shared_ptr<ExpList> &GetPlane(int n)
978 {
979 return v_GetPlane(n);
980 }
984
985 MULTI_REGIONS_EXPORT int GetPoolCount(std::string);
986
988
991
992 /// Get m_coeffs to elemental value map
994 &GetCoeffsToElmt() const;
1000 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
1001 const int nDirctn, Array<OneD, DNekMatSharedPtr> &mtxPerVar);
1003 const TensorOfArray3D<NekDouble> &inarray,
1006 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
1009 const Array<OneD, const NekDouble> &FwdFlux,
1010 const Array<OneD, const NekDouble> &BwdFlux,
1011 Array<OneD, NekDouble> &outarray)
1012 {
1013 v_AddTraceIntegralToOffDiag(FwdFlux, BwdFlux, outarray);
1014 }
1016 const int dir, const Array<OneD, const DNekMatSharedPtr> ElmtJacQuad,
1017 Array<OneD, DNekMatSharedPtr> ElmtJacCoef);
1019 const Array<OneD, const DNekMatSharedPtr> ElmtJacQuad,
1020 Array<OneD, DNekMatSharedPtr> ElmtJacCoef);
1022 GetLocTraceToTraceMap() const;
1023 // Return the internal vector which identifieds if trace
1024 // is left adjacent definiing which trace the normal
1025 // points otwards from
1026 MULTI_REGIONS_EXPORT std::vector<bool> &GetLeftAdjacentTraces(void);
1027
1028 /// This function returns the map of index inside m_exp to geom id
1029 MULTI_REGIONS_EXPORT inline const std::unordered_map<int, int>
1031 {
1032 return m_elmtToExpId;
1033 }
1034
1035 /// This function returns the index inside m_exp for a given geom id
1037 {
1038 auto it = m_elmtToExpId.find(elmtId);
1039 ASSERTL0(it != m_elmtToExpId.end(), "Global geometry ID " +
1040 std::to_string(elmtId) +
1041 " not found in element ID to "
1042 "expansion ID map.")
1043 return it->second;
1044 }
1045
1046protected:
1047 /// Exapnsion type
1049 std::shared_ptr<DNekMat> GenGlobalMatrixFull(
1050 const GlobalLinSysKey &mkey,
1051 const std::shared_ptr<AssemblyMapCG> &locToGloMap);
1052 /// Communicator
1054 /// Session
1056 /// Mesh associated with this expansion list.
1058 /// The total number of local degrees of freedom. #m_ncoeffs
1059 /// \f$=N_{\mathrm{eof}}=\sum_{e=1}^{{N_{\mathrm{el}}}}N^{e}_l\f$
1061 /** The total number of quadrature points. #m_npoints
1062 *\f$=Q_{\mathrm{tot}}=\sum_{e=1}^{{N_{\mathrm{el}}}}N^{e}_Q\f$
1063 **/
1065 /**
1066 * \brief Concatenation of all local expansion coefficients.
1067 *
1068 * The array of length #m_ncoeffs\f$=N_{\mathrm{eof}}\f$ which is
1069 * the concatenation of the local expansion coefficients
1070 * \f$\hat{u}_n^e\f$ over all \f$N_{\mathrm{el}}\f$ elements
1071 * \f[\mathrm{\texttt{m\_coeffs}}=\boldsymbol{\hat{u}}_{l} =
1072 * \underline{\boldsymbol{\hat{u}}}^e = \left [ \begin{array}{c}
1073 * \boldsymbol{\hat{u}}^{1} \ \
1074 * \boldsymbol{\hat{u}}^{2} \ \
1075 * \vdots \ \
1076 * \boldsymbol{\hat{u}}^{{{N_{\mathrm{el}}}}} \end{array} \right ],
1077 * \quad
1078 * \mathrm{where}\quad \boldsymbol{\hat{u}}^{e}[n]=\hat{u}_n^{e}\f]
1079 */
1081 /**
1082 * \brief The global expansion evaluated at the quadrature points
1083 *
1084 * The array of length #m_npoints\f$=Q_{\mathrm{tot}}\f$ containing
1085 * the evaluation of \f$u^{\delta}(\boldsymbol{x})\f$ at the
1086 * quadrature points \f$\boldsymbol{x}_i\f$.
1087 * \f[\mathrm{\texttt{m\_phys}}=\boldsymbol{u}_{l} =
1088 * \underline{\boldsymbol{u}}^e = \left [ \begin{array}{c}
1089 * \boldsymbol{u}^{1} \ \
1090 * \boldsymbol{u}^{2} \ \
1091 * \vdots \ \
1092 * \boldsymbol{u}^{{{N_{\mathrm{el}}}}} \end{array} \right ],\quad
1093 * \mathrm{where}\quad
1094 * \boldsymbol{u}^{e}[i]=u^{\delta}(\boldsymbol{x}_i)\f]
1095 */
1097 /**
1098 * \brief The state of the array #m_phys.
1099 *
1100 * Indicates whether the array #m_phys, created to contain the
1101 * evaluation of \f$u^{\delta}(\boldsymbol{x})\f$ at the quadrature
1102 * points \f$\boldsymbol{x}_i\f$, is filled with these values.
1103 */
1105 /**
1106 * \brief The list of local expansions.
1107 *
1108 * The (shared pointer to the) vector containing (shared pointers
1109 * to) all local expansions. The fact that the local expansions are
1110 * all stored as a (pointer to a) #StdExpansion, the abstract base
1111 * class for all local expansions, allows a general implementation
1112 * where most of the routines for the derived classes are defined
1113 * in the #ExpList base class.
1114 */
1115 std::shared_ptr<LocalRegions::ExpansionVector> m_exp;
1117 /// Vector of bools to act as an initialise on first call flag
1118 std::vector<bool> m_collectionsDoInit;
1119 /// Offset of elemental data into the array #m_coeffs
1120 std::vector<int> m_coll_coeff_offset;
1121 /// Offset of elemental data into the array #m_phys
1122 std::vector<int> m_coll_phys_offset;
1123 /// Offset of elemental data into the array #m_coeffs
1125 /// Offset of elemental data into the array #m_phys
1127 /// m_coeffs to elemental value map
1130 //@todo should this be in ExpList or
1131 // ExpListHomogeneous1D.cpp it's a bool which determine if
1132 // the expansion is in the wave space (coefficient space)
1133 // or not
1135 /// Mapping from geometry ID of element to index inside #m_exp
1136 std::unordered_map<int, int> m_elmtToExpId;
1137 /// This function assembles the block diagonal matrix of local
1138 /// matrices of the type \a mtype.
1141 void MultiplyByBlockMatrix(const GlobalMatrixKey &gkey,
1142 const Array<OneD, const NekDouble> &inarray,
1143 Array<OneD, NekDouble> &outarray);
1144 /// Generates a global matrix from the given key and map.
1145 std::shared_ptr<GlobalMatrix> GenGlobalMatrix(
1146 const GlobalMatrixKey &mkey,
1147 const std::shared_ptr<AssemblyMapCG> &locToGloMap);
1149 const std::shared_ptr<DNekMat> &Gmat,
1150 Array<OneD, NekDouble> &EigValsReal,
1151 Array<OneD, NekDouble> &EigValsImag,
1153 /// This operation constructs the global linear system of type \a
1154 /// mkey.
1155 std::shared_ptr<GlobalLinSys> GenGlobalLinSys(
1156 const GlobalLinSysKey &mkey,
1157 const std::shared_ptr<AssemblyMapCG> &locToGloMap);
1158 /// Generate a GlobalLinSys from information provided by the key
1159 /// "mkey" and the mapping provided in LocToGloBaseMap.
1160 std::shared_ptr<GlobalLinSys> GenGlobalBndLinSys(
1161 const GlobalLinSysKey &mkey, const AssemblyMapSharedPtr &locToGloMap);
1162 // Virtual prototypes
1163 virtual size_t v_GetNumElmts(void)
1164 {
1165 return (*m_exp).size();
1166 }
1170 virtual void v_SetBndCondBwdWeight(const int index, const NekDouble value);
1171 virtual std::shared_ptr<ExpList> &v_UpdateBndCondExpansion(int i);
1172 virtual void v_Upwind(const Array<OneD, const Array<OneD, NekDouble>> &Vec,
1176 virtual void v_Upwind(const Array<OneD, const NekDouble> &Vn,
1180 virtual std::shared_ptr<ExpList> &v_GetTrace();
1181 virtual std::shared_ptr<AssemblyMapDG> &v_GetTraceMap();
1183 virtual const std::shared_ptr<LocTraceToTraceMap> &v_GetLocTraceToTraceMap(
1184 void) const;
1185 virtual std::vector<bool> &v_GetLeftAdjacentTraces(void);
1186 /// Populate \a normals with the normals of all expansions.
1187 virtual void v_GetNormals(Array<OneD, Array<OneD, NekDouble>> &normals);
1188 virtual void v_AddTraceIntegral(const Array<OneD, const NekDouble> &Fn,
1189 Array<OneD, NekDouble> &outarray);
1190 virtual void v_AddFwdBwdTraceIntegral(
1193 Array<OneD, NekDouble> &outarray);
1196 virtual void v_GetFwdBwdTracePhys(const Array<OneD, const NekDouble> &field,
1199 bool FillBnd = true,
1200 bool PutFwdInBwdOnBCs = false,
1201 bool DoExchange = true);
1202 virtual void v_FillBwdWithBoundCond(const Array<OneD, NekDouble> &Fwd,
1204 bool PutFwdInBwdOnBCs);
1205 virtual void v_AddTraceQuadPhysToField(
1208 virtual void v_AddTraceQuadPhysToOffDiag(
1211 virtual void v_GetLocTraceFromTracePts(
1214 Array<OneD, NekDouble> &locTraceFwd,
1215 Array<OneD, NekDouble> &locTraceBwd);
1216 virtual void v_FillBwdWithBwdWeight(Array<OneD, NekDouble> &weightave,
1217 Array<OneD, NekDouble> &weightjmp);
1218 virtual void v_PeriodicBwdCopy(const Array<OneD, const NekDouble> &Fwd,
1220 virtual const std::vector<bool> &v_GetLeftAdjacentFaces(void) const;
1221 virtual void v_ExtractTracePhys(Array<OneD, NekDouble> &outarray);
1222 virtual void v_ExtractTracePhys(const Array<OneD, const NekDouble> &inarray,
1223 Array<OneD, NekDouble> &outarray);
1224 virtual void v_MultiplyByInvMassMatrix(
1225 const Array<OneD, const NekDouble> &inarray,
1226 Array<OneD, NekDouble> &outarray);
1227
1229 const Array<OneD, const NekDouble> &inarray,
1230 Array<OneD, NekDouble> &outarray,
1232 const StdRegions::VarCoeffMap &varcoeff,
1233 const MultiRegions::VarFactorsMap &varfactors,
1234 const Array<OneD, const NekDouble> &dirForcing,
1235 const bool PhysSpaceForcing);
1236
1238 const Array<OneD, const NekDouble> &inarray,
1239 Array<OneD, NekDouble> &outarray,
1241 const StdRegions::VarCoeffMap &varcoeff,
1242 const MultiRegions::VarFactorsMap &varfactors,
1243 const Array<OneD, const NekDouble> &dirForcing,
1244 const bool PhysSpaceForcing);
1245
1246 virtual void v_LinearAdvectionReactionSolve(
1248 const Array<OneD, const NekDouble> &inarray,
1249 Array<OneD, NekDouble> &outarray, const NekDouble lambda,
1251 // wrapper functions about virtual functions
1253 virtual void v_FillBndCondFromField(const Array<OneD, NekDouble> coeffs);
1254 virtual void v_FillBndCondFromField(const int nreg,
1255 const Array<OneD, NekDouble> coeffs);
1256 virtual void v_Reset();
1257 virtual void v_LocalToGlobal(bool UseComm);
1258 virtual void v_LocalToGlobal(const Array<OneD, const NekDouble> &inarray,
1259 Array<OneD, NekDouble> &outarray,
1260 bool UseComm);
1261 virtual void v_GlobalToLocal(void);
1262 virtual void v_GlobalToLocal(const Array<OneD, const NekDouble> &inarray,
1263 Array<OneD, NekDouble> &outarray);
1264 virtual void v_BwdTrans(const Array<OneD, const NekDouble> &inarray,
1265 Array<OneD, NekDouble> &outarray);
1266 virtual void v_FwdTrans(const Array<OneD, const NekDouble> &inarray,
1267 Array<OneD, NekDouble> &outarray);
1268 virtual void v_FwdTransLocalElmt(
1269 const Array<OneD, const NekDouble> &inarray,
1270 Array<OneD, NekDouble> &outarray);
1271 virtual void v_FwdTransBndConstrained(
1272 const Array<OneD, const NekDouble> &inarray,
1273 Array<OneD, NekDouble> &outarray);
1274 virtual void v_SmoothField(Array<OneD, NekDouble> &field);
1275 virtual void v_IProductWRTBase(const Array<OneD, const NekDouble> &inarray,
1276 Array<OneD, NekDouble> &outarray);
1277
1278 // Define ExpList::IProductWRTDerivBase as virtual function
1279 virtual void v_IProductWRTDerivBase(
1280 const int dir, const Array<OneD, const NekDouble> &inarray,
1281 Array<OneD, NekDouble> &outarray);
1282
1283 virtual void v_IProductWRTDerivBase(
1284 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
1285 Array<OneD, NekDouble> &outarray);
1286
1287 virtual void v_GetCoords(
1288 Array<OneD, NekDouble> &coord_0,
1291
1292 virtual void v_GetCoords(const int eid, Array<OneD, NekDouble> &xc0,
1295
1296 virtual void v_PhysDeriv(const Array<OneD, const NekDouble> &inarray,
1297 Array<OneD, NekDouble> &out_d0,
1298 Array<OneD, NekDouble> &out_d1,
1299 Array<OneD, NekDouble> &out_d2);
1300 virtual void v_PhysDeriv(const int dir,
1301 const Array<OneD, const NekDouble> &inarray,
1302 Array<OneD, NekDouble> &out_d);
1303 virtual void v_PhysDeriv(Direction edir,
1304 const Array<OneD, const NekDouble> &inarray,
1305 Array<OneD, NekDouble> &out_d);
1306
1307 virtual void v_Curl(Array<OneD, Array<OneD, NekDouble>> &Vel,
1309
1310 virtual void v_CurlCurl(Array<OneD, Array<OneD, NekDouble>> &Vel,
1312
1313 virtual void v_PhysDirectionalDeriv(
1314 const Array<OneD, const NekDouble> &direction,
1315 const Array<OneD, const NekDouble> &inarray,
1316 Array<OneD, NekDouble> &outarray);
1317 virtual void v_GetMovingFrames(
1318 const SpatialDomains::GeomMMF MMFdir,
1319 const Array<OneD, const NekDouble> &CircCentre,
1320 Array<OneD, Array<OneD, NekDouble>> &outarray);
1321 virtual void v_HomogeneousFwdTrans(
1322 const int npts, const Array<OneD, const NekDouble> &inarray,
1323 Array<OneD, NekDouble> &outarray, bool Shuff = true,
1324 bool UnShuff = true);
1325 virtual void v_HomogeneousBwdTrans(
1326 const int npts, const Array<OneD, const NekDouble> &inarray,
1327 Array<OneD, NekDouble> &outarray, bool Shuff = true,
1328 bool UnShuff = true);
1329 virtual void v_DealiasedProd(const int num_dofs,
1330 const Array<OneD, NekDouble> &inarray1,
1331 const Array<OneD, NekDouble> &inarray2,
1332 Array<OneD, NekDouble> &outarray);
1333 virtual void v_DealiasedDotProd(
1334 const int num_dofs, const Array<OneD, Array<OneD, NekDouble>> &inarray1,
1335 const Array<OneD, Array<OneD, NekDouble>> &inarray2,
1336 Array<OneD, Array<OneD, NekDouble>> &outarray);
1337 virtual void v_GetBCValues(Array<OneD, NekDouble> &BndVals,
1338 const Array<OneD, NekDouble> &TotField,
1339 int BndID);
1342 Array<OneD, NekDouble> &outarray,
1343 int BndID);
1344 virtual void v_NormVectorIProductWRTBase(
1346 Array<OneD, NekDouble> &outarray);
1347 virtual void v_SetUpPhysNormals();
1348 virtual void v_GetBoundaryToElmtMap(Array<OneD, int> &ElmtID,
1349 Array<OneD, int> &EdgeID);
1350 virtual void v_GetBndElmtExpansion(int i, std::shared_ptr<ExpList> &result,
1351 const bool DeclareCoeffPhysArrays);
1352
1353 virtual void v_ExtractElmtToBndPhys(const int i,
1354 const Array<OneD, NekDouble> &elmt,
1355 Array<OneD, NekDouble> &boundary);
1356
1357 virtual void v_ExtractPhysToBndElmt(
1358 const int i, const Array<OneD, const NekDouble> &phys,
1359 Array<OneD, NekDouble> &bndElmt);
1360
1361 virtual void v_ExtractPhysToBnd(const int i,
1362 const Array<OneD, const NekDouble> &phys,
1364
1365 virtual void v_GetBoundaryNormals(
1366 int i, Array<OneD, Array<OneD, NekDouble>> &normals);
1367
1368 virtual std::vector<LibUtilities::FieldDefinitionsSharedPtr>
1370
1371 virtual void v_GetFieldDefinitions(
1372 std::vector<LibUtilities::FieldDefinitionsSharedPtr> &fielddef);
1373
1374 virtual void v_AppendFieldData(
1376 std::vector<NekDouble> &fielddata);
1377 virtual void v_AppendFieldData(
1379 std::vector<NekDouble> &fielddata, Array<OneD, NekDouble> &coeffs);
1380 virtual void v_ExtractDataToCoeffs(
1382 std::vector<NekDouble> &fielddata, std::string &field,
1383 Array<OneD, NekDouble> &coeffs,
1384 std::unordered_map<int, int> zIdToPlane);
1385 virtual void v_ExtractCoeffsToCoeffs(
1386 const std::shared_ptr<ExpList> &fromExpList,
1387 const Array<OneD, const NekDouble> &fromCoeffs,
1388 Array<OneD, NekDouble> &toCoeffs);
1389 virtual void v_WriteTecplotHeader(std::ostream &outfile,
1390 std::string var = "");
1391 virtual void v_WriteTecplotZone(std::ostream &outfile, int expansion);
1392 virtual void v_WriteTecplotField(std::ostream &outfile, int expansion);
1393 virtual void v_WriteTecplotConnectivity(std::ostream &outfile,
1394 int expansion);
1395 virtual void v_WriteVtkPieceData(std::ostream &outfile, int expansion,
1396 std::string var);
1397 virtual void v_WriteVtkPieceHeader(std::ostream &outfile, int expansion,
1398 int istrip);
1399
1400 virtual NekDouble v_L2(
1401 const Array<OneD, const NekDouble> &phys,
1403
1404 virtual NekDouble v_Integral(const Array<OneD, const NekDouble> &inarray);
1405 virtual NekDouble v_VectorFlux(
1406 const Array<OneD, Array<OneD, NekDouble>> &inarray);
1407
1409
1411 virtual NekDouble v_GetHomoLen(void);
1412 virtual void v_SetHomoLen(const NekDouble lhom);
1415
1416 // 1D Scaling and projection
1417 virtual void v_PhysInterp1DScaled(const NekDouble scale,
1418 const Array<OneD, NekDouble> &inarray,
1419 Array<OneD, NekDouble> &outarray);
1420
1422 const NekDouble scale, const Array<OneD, NekDouble> &inarray,
1423 Array<OneD, NekDouble> &outarray);
1424
1425 virtual void v_ClearGlobalLinSysManager(void);
1426
1427 virtual int v_GetPoolCount(std::string);
1428
1429 virtual void v_UnsetGlobalLinSys(GlobalLinSysKey, bool);
1430
1433
1434 void ExtractFileBCs(const std::string &fileName,
1436 const std::string &varName,
1437 const std::shared_ptr<ExpList> locExpList);
1438
1439 // Utility function for a common case of retrieving a
1440 // BoundaryCondition from a boundary condition collection.
1444 unsigned int index, const std::string &variable);
1445
1448
1451
1452 virtual void v_EvaluateBoundaryConditions(
1453 const NekDouble time = 0.0, const std::string varName = "",
1456
1457 virtual std::map<int, RobinBCInfoSharedPtr> v_GetRobinBCInfo(void);
1458
1459 virtual void v_GetPeriodicEntities(PeriodicMap &periodicVerts,
1460 PeriodicMap &periodicEdges,
1461 PeriodicMap &periodicFaces);
1462
1463 // Homogeneous direction wrapper functions.
1465 {
1467 "This method is not defined or valid for this class type");
1469 }
1470
1471 // wrapper function to set viscosity for Homo1D expansion
1473 [[maybe_unused]] Array<OneD, NekDouble> visc)
1474 {
1476 "This method is not defined or valid for this class type");
1477 }
1478
1479 virtual std::shared_ptr<ExpList> &v_GetPlane(int n);
1480
1481 virtual void v_AddTraceIntegralToOffDiag(
1482 const Array<OneD, const NekDouble> &FwdFlux,
1483 const Array<OneD, const NekDouble> &BwdFlux,
1484 Array<OneD, NekDouble> &outarray);
1485
1486private:
1487 /// Definition of the total number of degrees of freedom and
1488 /// quadrature points and offsets to access data
1489 void SetupCoeffPhys(bool DeclareCoeffPhysArrays = true,
1490 bool SetupOffsets = true);
1491
1492 /// Define a list of elements using the geometry and basis
1493 /// key information in expmap;
1495};
1496
1497/// An empty ExpList object.
1500
1501// An empty GlobaLinSysManager and GlobalLinSysKey object
1505
1506// Inline routines follow.
1507
1508/**
1509 * Returns the total number of local degrees of freedom
1510 * \f$N_{\mathrm{eof}}=\sum_{e=1}^{{N_{\mathrm{el}}}}N^{e}_m\f$.
1511 */
1512inline int ExpList::GetNcoeffs() const
1513{
1514 return m_ncoeffs;
1515}
1516
1517inline int ExpList::GetNcoeffs(const int eid) const
1518{
1519 return (*m_exp)[eid]->GetNcoeffs();
1520}
1521
1522/**
1523 * Evaulates the maximum number of modes in the elemental basis
1524 * order over all elements
1525 */
1527{
1528 int returnval = 0;
1529
1530 for (size_t i = 0; i < (*m_exp).size(); ++i)
1531 {
1532 returnval = (std::max)(returnval, (*m_exp)[i]->EvalBasisNumModesMax());
1533 }
1534
1535 return returnval;
1536}
1537
1538/**
1539 *
1540 */
1542{
1543 Array<OneD, int> returnval((*m_exp).size(), 0);
1544
1545 for (size_t i = 0; i < (*m_exp).size(); ++i)
1546 {
1547 returnval[i] =
1548 (std::max)(returnval[i], (*m_exp)[i]->EvalBasisNumModesMax());
1549 }
1550
1551 return returnval;
1552}
1553
1554/**
1555 *
1556 */
1557inline int ExpList::GetTotPoints() const
1558{
1559 return m_npoints;
1560}
1561
1562inline int ExpList::GetTotPoints(const int eid) const
1563{
1564 return (*m_exp)[eid]->GetTotPoints();
1565}
1566
1567inline int ExpList::Get1DScaledTotPoints(const NekDouble scale) const
1568{
1569 int returnval = 0;
1570 size_t cnt;
1571 size_t nbase = (*m_exp)[0]->GetNumBases();
1572
1573 for (size_t i = 0; i < (*m_exp).size(); ++i)
1574 {
1575 cnt = 1;
1576 for (size_t j = 0; j < nbase; ++j)
1577 {
1578 cnt *= scale * ((*m_exp)[i]->GetNumPoints(j));
1579 }
1580 returnval += cnt;
1581 }
1582 return returnval;
1583}
1584
1585/**
1586 *
1587 */
1588inline int ExpList::GetNpoints() const
1589{
1590 return m_npoints;
1591}
1592
1593/**
1594 *
1595 */
1596inline void ExpList::SetWaveSpace(const bool wavespace)
1597{
1598 m_WaveSpace = wavespace;
1599}
1600
1601/**
1602 *
1603 */
1604inline bool ExpList::GetWaveSpace() const
1605{
1606 return m_WaveSpace;
1607}
1608
1609/// Set the \a i th value of\a m_phys to value \a val
1610inline void ExpList::SetPhys(int i, NekDouble val)
1611{
1612 m_phys[i] = val;
1613}
1614/**
1615 * This function fills the array \f$\boldsymbol{u}_l\f$, the evaluation
1616 * of the expansion at the quadrature points (implemented as #m_phys),
1617 * with the values of the array \a inarray.
1618 *
1619 * @param inarray The array containing the values where
1620 * #m_phys should be filled with.
1621 */
1623{
1624 ASSERTL0((int)inarray.size() == m_npoints,
1625 "Input array does not have correct number of elements.");
1626 Vmath::Vcopy(m_npoints, &inarray[0], 1, &m_phys[0], 1);
1627 m_physState = true;
1628}
1630{
1631 m_phys = inarray;
1632}
1633/**
1634 * @param physState \a true (=filled) or \a false (=not filled).
1635 */
1636inline void ExpList::SetPhysState(const bool physState)
1637{
1638 m_physState = physState;
1639}
1640/**
1641 * @return physState \a true (=filled) or \a false (=not filled).
1642 */
1643inline bool ExpList::GetPhysState() const
1644{
1645 return m_physState;
1646}
1647/**
1648 *
1649 */
1651 const Array<OneD, const NekDouble> &inarray,
1652 Array<OneD, NekDouble> &outarray)
1653{
1654 v_IProductWRTBase(inarray, outarray);
1655}
1656/**
1657 *
1658 */
1660 const int dir, const Array<OneD, const NekDouble> &inarray,
1661 Array<OneD, NekDouble> &outarray)
1662{
1663 v_IProductWRTDerivBase(dir, inarray, outarray);
1664}
1665
1666/**
1667 *
1668 */
1670 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
1671 Array<OneD, NekDouble> &outarray)
1672{
1673 v_IProductWRTDerivBase(inarray, outarray);
1674}
1675
1676/**
1677 *
1678 */
1680 Array<OneD, NekDouble> &outarray)
1681{
1682 v_FwdTrans(inarray, outarray);
1683}
1684/**
1685 *
1686 */
1688 const Array<OneD, const NekDouble> &inarray,
1689 Array<OneD, NekDouble> &outarray)
1690{
1691 v_FwdTransLocalElmt(inarray, outarray);
1692}
1693/**
1694 *
1695 */
1697 const Array<OneD, const NekDouble> &inarray,
1698 Array<OneD, NekDouble> &outarray)
1699{
1700 v_FwdTransBndConstrained(inarray, outarray);
1701}
1702/**
1703 *
1704 */
1706{
1707 v_SmoothField(field);
1708}
1709/**
1710 *
1711 */
1712/**
1713 * Given the coefficients of an expansion, this function evaluates the
1714 * spectral/hp expansion \f$u^{\delta}(\boldsymbol{x})\f$ at the
1715 * quadrature points \f$\boldsymbol{x}_i\f$.
1716 */
1718 Array<OneD, NekDouble> &outarray)
1719{
1720 v_BwdTrans(inarray, outarray);
1721}
1722/**
1723 *
1724 */
1726 const Array<OneD, const NekDouble> &inarray,
1727 Array<OneD, NekDouble> &outarray)
1728{
1729 v_MultiplyByInvMassMatrix(inarray, outarray);
1730}
1731/**
1732 *
1733 */
1735 const Array<OneD, const NekDouble> &inarray,
1737 const StdRegions::VarCoeffMap &varcoeff,
1738 const MultiRegions::VarFactorsMap &varfactors,
1739 const Array<OneD, const NekDouble> &dirForcing, const bool PhysSpaceForcing)
1740{
1741 return v_HelmSolve(inarray, outarray, factors, varcoeff, varfactors,
1742 dirForcing, PhysSpaceForcing);
1743}
1744/**
1745 *
1746 */
1748 const Array<OneD, const NekDouble> &inarray,
1750 const StdRegions::VarCoeffMap &varcoeff,
1751 const MultiRegions::VarFactorsMap &varfactors,
1752 const Array<OneD, const NekDouble> &dirForcing, const bool PhysSpaceForcing)
1753{
1755 inarray, outarray, factors, varcoeff, varfactors, dirForcing,
1756 PhysSpaceForcing);
1757}
1760 const Array<OneD, const NekDouble> &inarray,
1761 Array<OneD, NekDouble> &outarray, const NekDouble lambda,
1762 const Array<OneD, const NekDouble> &dirForcing)
1763{
1764 v_LinearAdvectionReactionSolve(velocity, inarray, outarray, lambda,
1765 dirForcing);
1766}
1767/**
1768 *
1769 */
1771 Array<OneD, NekDouble> &coord_1,
1772 Array<OneD, NekDouble> &coord_2)
1773{
1774 v_GetCoords(coord_0, coord_1, coord_2);
1775}
1776
1777inline void ExpList::GetCoords(const int eid, Array<OneD, NekDouble> &xc0,
1780{
1781 v_GetCoords(eid, xc0, xc1, xc2);
1782}
1783
1784/**
1785 *
1786 */
1788 const SpatialDomains::GeomMMF MMFdir,
1789 const Array<OneD, const NekDouble> &CircCentre,
1790 Array<OneD, Array<OneD, NekDouble>> &outarray)
1791{
1792 v_GetMovingFrames(MMFdir, CircCentre, outarray);
1793}
1794/**
1795 *
1796 */
1798 Array<OneD, NekDouble> &out_d0,
1799 Array<OneD, NekDouble> &out_d1,
1800 Array<OneD, NekDouble> &out_d2)
1801{
1802 v_PhysDeriv(inarray, out_d0, out_d1, out_d2);
1803}
1804/**
1805 *
1806 */
1807inline void ExpList::PhysDeriv(const int dir,
1808 const Array<OneD, const NekDouble> &inarray,
1810{
1811 v_PhysDeriv(dir, inarray, out_d);
1812}
1814 const Array<OneD, const NekDouble> &inarray,
1816{
1817 v_PhysDeriv(edir, inarray, out_d);
1818}
1819/**
1820 *
1821 */
1823 const Array<OneD, const NekDouble> &direction,
1824 const Array<OneD, const NekDouble> &inarray,
1825 Array<OneD, NekDouble> &outarray)
1826{
1827 v_PhysDirectionalDeriv(direction, inarray, outarray);
1828}
1829/**
1830 *
1831 */
1834{
1835 v_Curl(Vel, Q);
1836}
1837/**
1838 *
1839 */
1842{
1843 v_CurlCurl(Vel, Q);
1844}
1845/**
1846 *
1847 */
1849 const int npts, const Array<OneD, const NekDouble> &inarray,
1850 Array<OneD, NekDouble> &outarray, bool Shuff, bool UnShuff)
1851{
1852 v_HomogeneousFwdTrans(npts, inarray, outarray, Shuff, UnShuff);
1853}
1854/**
1855 *
1856 */
1858 const int npts, const Array<OneD, const NekDouble> &inarray,
1859 Array<OneD, NekDouble> &outarray, bool Shuff, bool UnShuff)
1860{
1861 v_HomogeneousBwdTrans(npts, inarray, outarray, Shuff, UnShuff);
1862}
1863/**
1864 *
1865 */
1866inline void ExpList::DealiasedProd(const int num_dofs,
1867 const Array<OneD, NekDouble> &inarray1,
1868 const Array<OneD, NekDouble> &inarray2,
1869 Array<OneD, NekDouble> &outarray)
1870{
1871 v_DealiasedProd(num_dofs, inarray1, inarray2, outarray);
1872}
1873/**
1874 *
1875 */
1877 const int num_dofs, const Array<OneD, Array<OneD, NekDouble>> &inarray1,
1878 const Array<OneD, Array<OneD, NekDouble>> &inarray2,
1879 Array<OneD, Array<OneD, NekDouble>> &outarray)
1880{
1881 v_DealiasedDotProd(num_dofs, inarray1, inarray2, outarray);
1882}
1883/**
1884 *
1885 */
1887 const Array<OneD, NekDouble> &TotField,
1888 int BndID)
1889{
1890 v_GetBCValues(BndVals, TotField, BndID);
1891}
1892/**
1893 *
1894 */
1897 Array<OneD, NekDouble> &outarray,
1898 int BndID)
1899{
1900 v_NormVectorIProductWRTBase(V1, V2, outarray, BndID);
1901}
1904{
1905 v_NormVectorIProductWRTBase(V, outarray);
1906}
1907/**
1908 * @param eid The index of the element to be checked.
1909 * @return The dimension of the coordinates of the specific element.
1910 */
1911inline int ExpList::GetCoordim(int eid)
1912{
1913 ASSERTL2(eid <= (*m_exp).size(), "eid is larger than number of elements");
1914 return (*m_exp)[eid]->GetCoordim();
1915}
1916/**
1917 * @param eid The index of the element to be checked.
1918 * @return The dimension of the shape of the specific element.
1919 */
1921{
1922 return (*m_exp)[0]->GetShapeDimension();
1923}
1924/**
1925 * @param i The index of m_coeffs to be set
1926 * @param val The value which m_coeffs[i] is to be set to.
1927 */
1928inline void ExpList::SetCoeff(int i, NekDouble val)
1929{
1930 m_coeffs[i] = val;
1931}
1932/**
1933 * @param i The index of #m_coeffs to be set.
1934 * @param val The value which #m_coeffs[i] is to be set to.
1935 */
1936inline void ExpList::SetCoeffs(int i, NekDouble val)
1937{
1938 m_coeffs[i] = val;
1939}
1941{
1942 m_coeffs = inarray;
1943}
1944/**
1945 * As the function returns a constant reference to a
1946 * <em>const Array</em>, it is not possible to modify the
1947 * underlying data of the array #m_coeffs. In order to do
1948 * so, use the function #UpdateCoeffs instead.
1949 *
1950 * @return (A constant reference to) the array #m_coeffs.
1951 */
1953{
1954 return m_coeffs;
1955}
1957{
1959}
1961{
1962 v_FillBndCondFromField(coeffs);
1963}
1964inline void ExpList::FillBndCondFromField(const int nreg,
1965 const Array<OneD, NekDouble> coeffs)
1966{
1967 v_FillBndCondFromField(nreg, coeffs);
1968}
1969inline void ExpList::LocalToGlobal(bool useComm)
1970{
1971 v_LocalToGlobal(useComm);
1972}
1974 Array<OneD, NekDouble> &outarray,
1975 bool useComm)
1976{
1977 v_LocalToGlobal(inarray, outarray, useComm);
1978}
1979inline void ExpList::GlobalToLocal(void)
1980{
1982}
1983/**
1984 * This operation is evaluated as:
1985 * \f{tabbing}
1986 * \hspace{1cm} \= Do \= $e=$ $1, N_{\mathrm{el}}$ \\
1987 * > > Do \= $i=$ $0,N_m^e-1$ \\
1988 * > > > $\boldsymbol{\hat{u}}^{e}[i] = \mbox{sign}[e][i] \cdot
1989 * \boldsymbol{\hat{u}}_g[\mbox{map}[e][i]]$ \\
1990 * > > continue \\
1991 * > continue
1992 * \f}
1993 * where \a map\f$[e][i]\f$ is the mapping array and \a
1994 * sign\f$[e][i]\f$ is an array of similar dimensions ensuring the
1995 * correct modal connectivity between the different elements (both
1996 * these arrays are contained in the data member #m_locToGloMap). This
1997 * operation is equivalent to the scatter operation
1998 * \f$\boldsymbol{\hat{u}}_l=\mathcal{A}\boldsymbol{\hat{u}}_g\f$, where
1999 * \f$\mathcal{A}\f$ is the
2000 * \f$N_{\mathrm{eof}}\times N_{\mathrm{dof}}\f$ permutation matrix.
2001 *
2002 * @param inarray An array of size \f$N_\mathrm{dof}\f$
2003 * containing the global degrees of freedom
2004 * \f$\boldsymbol{x}_g\f$.
2005 * @param outarray The resulting local degrees of freedom
2006 * \f$\boldsymbol{x}_l\f$ will be stored in this
2007 * array of size \f$N_\mathrm{eof}\f$.
2008 */
2010 Array<OneD, NekDouble> &outarray)
2011{
2012 v_GlobalToLocal(inarray, outarray);
2013}
2014/**
2015 * @param i The index of #m_coeffs to be returned
2016 * @return The NekDouble held in #m_coeffs[i].
2017 */
2019{
2020 return m_coeffs[i];
2021}
2022/**
2023 * @param i The index of #m_coeffs to be returned
2024 * @return The NekDouble held in #m_coeffs[i].
2025 */
2027{
2028 return m_coeffs[i];
2029}
2030/**
2031 * As the function returns a constant reference to a
2032 * <em>const Array</em> it is not possible to modify the
2033 * underlying data of the array #m_phys. In order to do
2034 * so, use the function #UpdatePhys instead.
2035 *
2036 * @return (A constant reference to) the array #m_phys.
2037 */
2039{
2040 return m_phys;
2041}
2042/**
2043 * @return \f$N_{\mathrm{el}}\f$, the number of elements in the
2044 * expansion.
2045 */
2046inline int ExpList::GetExpSize(void)
2047{
2048 return (*m_exp).size();
2049}
2050/**
2051 * @param n The index of the element concerned.
2052 *
2053 * @return (A shared pointer to) the local expansion of the
2054 * \f$n^{\mathrm{th}}\f$ element.
2055 */
2057{
2058 return (*m_exp)[n];
2059}
2060/**
2061 * @param n The global id of the element concerned.
2062 *
2063 * @return (A shared pointer to) the local expansion of the
2064 * \f$n^{\mathrm{th}}\f$ element.
2065 */
2067{
2068 auto it = m_elmtToExpId.find(n);
2069 ASSERTL0(it != m_elmtToExpId.end(), "Global geometry ID " +
2070 std::to_string(n) +
2071 " not found in element ID to "
2072 "expansion ID map.")
2073 return (*m_exp)[it->second];
2074}
2075
2076/**
2077 * @return (A const shared pointer to) the local expansion vector #m_exp
2078 */
2079inline const std::shared_ptr<LocalRegions::ExpansionVector> ExpList::GetExp(
2080 void) const
2081{
2082 return m_exp;
2083}
2084/**
2085 *
2086 */
2087inline int ExpList::GetCoeff_Offset(int n) const
2088{
2089 return m_coeff_offset[n];
2090}
2091/**
2092 *
2093 */
2094inline int ExpList::GetPhys_Offset(int n) const
2095{
2096 return m_phys_offset[n];
2097}
2098/**
2099 * If one wants to get hold of the underlying data without modifying
2100 * them, rather use the function #GetCoeffs instead.
2101 *
2102 * @return (A reference to) the array #m_coeffs.
2103 */
2105{
2106 return m_coeffs;
2107}
2108/**
2109 * If one wants to get hold of the underlying data without modifying
2110 * them, rather use the function #GetPhys instead.
2111 *
2112 * @return (A reference to) the array #m_phys.
2113 */
2115{
2116 m_physState = true;
2117 return m_phys;
2118}
2119// functions associated with DisContField
2122{
2123 return v_GetBndCondExpansions();
2124}
2125/// Get m_coeffs to elemental value map
2128{
2129 return m_coeffsToElmt;
2130}
2133{
2134 return v_GetLocTraceToTraceMap();
2135}
2137{
2138 return v_GetBndCondBwdWeight();
2139}
2140inline void ExpList::SetBndCondBwdWeight(const int index, const NekDouble value)
2141{
2142 v_SetBndCondBwdWeight(index, value);
2143}
2144inline std::shared_ptr<ExpList> &ExpList::UpdateBndCondExpansion(int i)
2145{
2146 return v_UpdateBndCondExpansion(i);
2147}
2149 const Array<OneD, const Array<OneD, NekDouble>> &Vec,
2152{
2153 v_Upwind(Vec, Fwd, Bwd, Upwind);
2154}
2158 Array<OneD, NekDouble> &Upwind)
2159{
2160 v_Upwind(Vn, Fwd, Bwd, Upwind);
2161}
2162inline std::shared_ptr<ExpList> &ExpList::GetTrace()
2163{
2164 return v_GetTrace();
2165}
2166inline std::shared_ptr<AssemblyMapDG> &ExpList::GetTraceMap()
2167{
2168 return v_GetTraceMap();
2169}
2171{
2172 return v_GetTraceBndMap();
2173}
2175{
2176 v_GetNormals(normals);
2177}
2179 Array<OneD, NekDouble> &outarray)
2180{
2181 v_AddTraceIntegral(Fn, outarray);
2182}
2186{
2187 v_AddFwdBwdTraceIntegral(Fwd, Bwd, outarray);
2188}
2191{
2192 v_GetFwdBwdTracePhys(Fwd, Bwd);
2193}
2196 Array<OneD, NekDouble> &Bwd, bool FillBnd, bool PutFwdInBwdOnBCs,
2197 bool DoExchange)
2198{
2199 v_GetFwdBwdTracePhys(field, Fwd, Bwd, FillBnd, PutFwdInBwdOnBCs,
2200 DoExchange);
2201}
2204 bool PutFwdInBwdOnBCs)
2205{
2206 v_FillBwdWithBoundCond(Fwd, Bwd, PutFwdInBwdOnBCs);
2207}
2211{
2212 v_AddTraceQuadPhysToField(Fwd, Bwd, field);
2213}
2217 Array<OneD, NekDouble> &locTraceFwd, Array<OneD, NekDouble> &locTraceBwd)
2218{
2219 v_GetLocTraceFromTracePts(Fwd, Bwd, locTraceFwd, locTraceBwd);
2220}
2222 Array<OneD, NekDouble> &weightjmp)
2223{
2224 v_FillBwdWithBwdWeight(weightave, weightjmp);
2225}
2228{
2229 v_PeriodicBwdCopy(Fwd, Bwd);
2230}
2231inline const std::vector<bool> &ExpList::GetLeftAdjacentFaces(void) const
2232{
2233 return v_GetLeftAdjacentFaces();
2234}
2236{
2237 v_ExtractTracePhys(outarray);
2238}
2240 const Array<OneD, const NekDouble> &inarray,
2241 Array<OneD, NekDouble> &outarray)
2242{
2243 v_ExtractTracePhys(inarray, outarray);
2244}
2247{
2248 return v_GetBndConditions();
2249}
2252{
2253 return v_UpdateBndConditions();
2254}
2256 const std::string varName,
2257 const NekDouble x2_in,
2258 const NekDouble x3_in)
2259{
2260 v_EvaluateBoundaryConditions(time, varName, x2_in, x3_in);
2261}
2263{
2265}
2267 Array<OneD, int> &EdgeID)
2268{
2269 v_GetBoundaryToElmtMap(ElmtID, EdgeID);
2270}
2272 std::shared_ptr<ExpList> &result,
2273 const bool DeclareCoeffPhysArrays)
2274{
2275 v_GetBndElmtExpansion(i, result, DeclareCoeffPhysArrays);
2276}
2277
2279 const Array<OneD, NekDouble> &elmt,
2280 Array<OneD, NekDouble> &boundary)
2281{
2282 v_ExtractElmtToBndPhys(i, elmt, boundary);
2283}
2284
2286 int i, const Array<OneD, const NekDouble> &phys,
2287 Array<OneD, NekDouble> &bndElmt)
2288{
2289 v_ExtractPhysToBndElmt(i, phys, bndElmt);
2290}
2291
2293 const Array<OneD, const NekDouble> &phys,
2295{
2296 v_ExtractPhysToBnd(i, phys, bnd);
2297}
2298
2300 int i, Array<OneD, Array<OneD, NekDouble>> &normals)
2301{
2302 v_GetBoundaryNormals(i, normals);
2303}
2304
2305inline std::vector<bool> &ExpList::GetLeftAdjacentTraces(void)
2306{
2307 return v_GetLeftAdjacentTraces();
2308}
2309
2311
2312} // namespace Nektar::MultiRegions
2313
2314#endif // EXPLIST_H
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:202
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
Definition: ErrorUtil.hpp:265
#define MULTI_REGIONS_EXPORT
Base class for all multi-elemental spectral/hp expansions.
Definition: ExpList.h:98
virtual void v_Reset()
Reset geometry information, metrics, matrix managers and geometry information.
Definition: ExpList.cpp:3043
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Definition: ExpList.h:1080
Array< OneD, NekDouble > & UpdateCoeffs()
This function returns (a reference to) the array (implemented as m_coeffs) containing all local expa...
Definition: ExpList.h:2104
void AddTraceJacToElmtJac(const Array< OneD, const DNekMatSharedPtr > &FwdMat, const Array< OneD, const DNekMatSharedPtr > &BwdMat, Array< OneD, DNekMatSharedPtr > &fieldMat)
inverse process of v_GetFwdBwdTracePhys. Given Trace integration of Fwd and Bwd Jacobian,...
Definition: ExpList.cpp:5922
void LinearAdvectionReactionSolve(const Array< OneD, Array< OneD, NekDouble > > &velocity, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const NekDouble lambda, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
Solve Advection Diffusion Reaction.
Definition: ExpList.h:1758
void IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function calculates the inner product of a function with respect to all local expansion modes .
Definition: ExpList.h:1650
const Array< OneD, const std::shared_ptr< ExpList > > & GetBndCondExpansions()
Definition: ExpList.h:2121
std::shared_ptr< ExpList > & GetTrace()
Definition: ExpList.h:2162
void WriteTecplotConnectivity(std::ostream &outfile, int expansion=-1)
Definition: ExpList.h:403
int GetNcoeffs(void) const
Returns the total number of local degrees of freedom .
Definition: ExpList.h:1512
virtual void v_ExtractCoeffsToCoeffs(const std::shared_ptr< ExpList > &fromExpList, const Array< OneD, const NekDouble > &fromCoeffs, Array< OneD, NekDouble > &toCoeffs)
Definition: ExpList.cpp:4042
const DNekScalBlkMatSharedPtr & GetBlockMatrix(const GlobalMatrixKey &gkey)
Definition: ExpList.cpp:2374
virtual void v_AddFwdBwdTraceIntegral(const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:4748
LibUtilities::NekManager< GlobalLinSysKey, GlobalLinSys > & GetGlobalLinSysManager(void)
Definition: ExpList.cpp:5646
ExpList(const ExpansionType Type=eNoType)
The default constructor using a type.
Definition: ExpList.cpp:90
void SetCoeffs(int i, NekDouble val)
Set the i th coefficiient in m_coeffs to value val.
Definition: ExpList.h:1936
void CreateCollections(Collections::ImplementationType ImpType=Collections::eNoImpType)
Construct collections of elements containing a single element type and polynomial order from the list...
Definition: ExpList.cpp:5461
virtual void v_GetBoundaryToElmtMap(Array< OneD, int > &ElmtID, Array< OneD, int > &EdgeID)
Definition: ExpList.cpp:5354
void ExponentialFilter(Array< OneD, NekDouble > &array, const NekDouble alpha, const NekDouble exponent, const NekDouble cutoff)
Definition: ExpList.cpp:2109
std::shared_ptr< GlobalMatrix > GenGlobalMatrix(const GlobalMatrixKey &mkey, const std::shared_ptr< AssemblyMapCG > &locToGloMap)
Generates a global matrix from the given key and map.
Definition: ExpList.cpp:2497
virtual void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:5069
void MultiplyByMassMatrix(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:298
void GlobalEigenSystem(const std::shared_ptr< DNekMat > &Gmat, Array< OneD, NekDouble > &EigValsReal, Array< OneD, NekDouble > &EigValsImag, Array< OneD, NekDouble > &EigVecs=NullNekDouble1DArray)
void FillBwdWithBoundCond(const Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd, bool PutFwdInBwdOnBCs=false)
Definition: ExpList.h:2202
void AddRightIPTPhysDerivBase(const int dir, const Array< OneD, const DNekMatSharedPtr > ElmtJacQuad, Array< OneD, DNekMatSharedPtr > ElmtJacCoef)
Definition: ExpList.cpp:6158
virtual void v_ExtractTracePhys(Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:4834
void ExtractPhysToBndElmt(int i, const Array< OneD, const NekDouble > &phys, Array< OneD, NekDouble > &bndElmt)
Definition: ExpList.h:2285
void Reset()
Reset geometry information and reset matrices.
Definition: ExpList.h:386
int GetExpIndex(const Array< OneD, const NekDouble > &gloCoord, NekDouble tol=0.0, bool returnNearestElmt=false, int cachedId=-1, NekDouble maxDistance=1e6)
This function returns the index of the local elemental expansion containing the arbitrary point given...
Definition: ExpList.cpp:2892
virtual GlobalLinSysKey v_HelmSolve(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::ConstFactorMap &factors, const StdRegions::VarCoeffMap &varcoeff, const MultiRegions::VarFactorsMap &varfactors, const Array< OneD, const NekDouble > &dirForcing, const bool PhysSpaceForcing)
Definition: ExpList.cpp:4857
void GetFieldDefinitions(std::vector< LibUtilities::FieldDefinitionsSharedPtr > &fielddef)
Definition: ExpList.h:900
virtual void v_WriteTecplotConnectivity(std::ostream &outfile, int expansion)
Definition: ExpList.cpp:3201
BlockMatrixMapShPtr m_blockMat
Definition: ExpList.h:1129
virtual Array< OneD, const unsigned int > v_GetYIDs(void)
Definition: ExpList.cpp:3663
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
multiply the metric jacobi and quadrature weights
Definition: ExpList.cpp:1647
virtual void v_LinearAdvectionReactionSolve(const Array< OneD, Array< OneD, NekDouble > > &velocity, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const NekDouble lambda, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
Definition: ExpList.cpp:4884
void SetPhys(int i, NekDouble val)
Set the i th value of m_phys to value val.
Definition: ExpList.h:1610
void ExtractTracePhys(Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:2235
void GetBoundaryToElmtMap(Array< OneD, int > &ElmtID, Array< OneD, int > &EdgeID)
Definition: ExpList.h:2266
void PeriodicBwdCopy(const Array< OneD, const NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
Copy and fill the Periodic boundaries.
Definition: ExpList.h:2226
void WriteVtkPieceFooter(std::ostream &outfile, int expansion)
Definition: ExpList.cpp:3461
virtual void v_ExtractPhysToBnd(const int i, const Array< OneD, const NekDouble > &phys, Array< OneD, NekDouble > &bnd)
Definition: ExpList.cpp:5273
void GetBwdWeight(Array< OneD, NekDouble > &weightAver, Array< OneD, NekDouble > &weightJump)
Get the weight value for boundary conditions for boundary average and jump calculations.
Definition: ExpList.cpp:4071
virtual void v_AddTraceIntegralToOffDiag(const Array< OneD, const NekDouble > &FwdFlux, const Array< OneD, const NekDouble > &BwdFlux, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:5700
virtual std::vector< bool > & v_GetLeftAdjacentTraces(void)
Definition: ExpList.cpp:4304
void InitialiseExpVector(const SpatialDomains::ExpansionInfoMap &expmap)
Define a list of elements using the geometry and basis key information in expmap;.
Definition: ExpList.cpp:1446
void SetupCoeffPhys(bool DeclareCoeffPhysArrays=true, bool SetupOffsets=true)
Definition of the total number of degrees of freedom and quadrature points and offsets to access data...
Definition: ExpList.cpp:1403
virtual LibUtilities::BasisSharedPtr v_GetHomogeneousBasis(void)
Definition: ExpList.h:1464
static SpatialDomains::BoundaryConditionShPtr GetBoundaryCondition(const SpatialDomains::BoundaryConditionCollection &collection, unsigned int index, const std::string &variable)
Definition: ExpList.cpp:5429
const Array< OneD, const NekDouble > & GetCoeffs() const
This function returns (a reference to) the array (implemented as m_coeffs) containing all local expa...
Definition: ExpList.h:1952
virtual void GetBndElmtExpansion(int i, std::shared_ptr< ExpList > &result, const bool DeclareCoeffPhysArrays=true)
Definition: ExpList.h:2271
void SetPhysArray(Array< OneD, NekDouble > &inarray)
Sets the array m_phys.
Definition: ExpList.h:1629
virtual void v_PhysDirectionalDeriv(const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:2079
virtual void v_GetBCValues(Array< OneD, NekDouble > &BndVals, const Array< OneD, NekDouble > &TotField, int BndID)
Definition: ExpList.cpp:4935
void FillBndCondFromField(const Array< OneD, NekDouble > coeffs)
Fill Bnd Condition expansion from the values stored in expansion.
Definition: ExpList.h:1960
const Array< OneD, const std::pair< int, int > > & GetCoeffsToElmt() const
Get m_coeffs to elemental value map.
Definition: ExpList.h:2127
virtual void v_GetLocTraceFromTracePts(const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &locTraceFwd, Array< OneD, NekDouble > &locTraceBwd)
Definition: ExpList.cpp:4801
const Array< OneD, const NekDouble > & GetPhys() const
This function returns (a reference to) the array (implemented as m_phys) containing the function ev...
Definition: ExpList.h:2038
std::shared_ptr< ExpList > & UpdateBndCondExpansion(int i)
Definition: ExpList.h:2144
virtual void v_WriteTecplotHeader(std::ostream &outfile, std::string var="")
Definition: ExpList.cpp:3070
virtual std::shared_ptr< ExpList > & v_GetPlane(int n)
Definition: ExpList.cpp:5450
void NormVectorIProductWRTBase(Array< OneD, const NekDouble > &V1, Array< OneD, const NekDouble > &V2, Array< OneD, NekDouble > &outarray, int BndID)
Definition: ExpList.h:1895
int GetCoeff_Offset(int n) const
Get the start offset position for a local contiguous list of coeffs correspoinding to element n.
Definition: ExpList.h:2087
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:2046
virtual void v_AddTraceQuadPhysToOffDiag(const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &field)
Definition: ExpList.cpp:4792
std::shared_ptr< LibUtilities::SessionReader > GetSession() const
Returns the session object.
Definition: ExpList.h:959
virtual void v_ClearGlobalLinSysManager(void)
Definition: ExpList.cpp:3671
virtual void v_UnsetGlobalLinSys(GlobalLinSysKey, bool)
Definition: ExpList.cpp:3683
void GeneralGetFieldDefinitions(std::vector< LibUtilities::FieldDefinitionsSharedPtr > &fielddef, int NumHomoDir=0, Array< OneD, LibUtilities::BasisSharedPtr > &HomoBasis=LibUtilities::NullBasisSharedPtr1DArray, std::vector< NekDouble > &HomoLen=LibUtilities::NullNekDoubleVector, bool homoStrips=false, std::vector< unsigned int > &HomoSIDs=LibUtilities::NullUnsignedIntVector, std::vector< unsigned int > &HomoZIDs=LibUtilities::NullUnsignedIntVector, std::vector< unsigned int > &HomoYIDs=LibUtilities::NullUnsignedIntVector)
Definition: ExpList.cpp:3771
virtual void v_GetPeriodicEntities(PeriodicMap &periodicVerts, PeriodicMap &periodicEdges, PeriodicMap &periodicFaces)
Definition: ExpList.cpp:5421
void IProductWRTDerivBase(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function calculates the inner product of a function with respect to the derivative (in directio...
Definition: ExpList.h:1659
void SetCoeffsArray(Array< OneD, NekDouble > &inarray)
Set the m_coeffs array to inarray.
Definition: ExpList.h:1940
virtual void v_CurlCurl(Array< OneD, Array< OneD, NekDouble > > &Vel, Array< OneD, Array< OneD, NekDouble > > &Q)
Definition: ExpList.cpp:2011
LibUtilities::TranspositionSharedPtr GetTransposition(void)
This function returns the transposition class associated with the homogeneous expansion.
Definition: ExpList.h:592
virtual void v_LocalToGlobal(bool UseComm)
Definition: ExpList.cpp:5023
void SetBndCondBwdWeight(const int index, const NekDouble value)
Set the weight value for boundary conditions.
Definition: ExpList.h:2140
virtual const Array< OneD, const SpatialDomains::BoundaryConditionShPtr > & v_GetBndConditions()
Definition: ExpList.cpp:5378
virtual LibUtilities::NekManager< GlobalLinSysKey, GlobalLinSys > & v_GetGlobalLinSysManager(void)
Definition: ExpList.cpp:3691
void SetHomoLen(const NekDouble lhom)
This function sets the Width of homogeneous direction associated with the homogeneous expansion.
Definition: ExpList.h:606
virtual void v_DealiasedDotProd(const int num_dofs, const Array< OneD, Array< OneD, NekDouble > > &inarray1, const Array< OneD, Array< OneD, NekDouble > > &inarray2, Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: ExpList.cpp:4925
int GetNpoints(void) const
Returns the total number of quadrature points m_npoints .
Definition: ExpList.h:1588
SpatialDomains::MeshGraphSharedPtr GetGraph()
Definition: ExpList.h:968
virtual ~ExpList()
The default destructor.
Definition: ExpList.cpp:1615
void GetFwdBwdTracePhys(Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
Definition: ExpList.h:2189
int GetElmtToExpId(int elmtId)
This function returns the index inside m_exp for a given geom id.
Definition: ExpList.h:1036
const Array< OneD, const int > & GetTraceBndMap(void)
Definition: ExpList.h:2170
virtual void v_FwdTransLocalElmt(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:2168
int GetPoolCount(std::string)
Definition: ExpList.cpp:5635
virtual void v_ImposeDirichletConditions(Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:4997
void FwdTransLocalElmt(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function elementally evaluates the forward transformation of a function onto the global spectra...
Definition: ExpList.h:1687
Array< OneD, const unsigned int > GetZIDs(void)
This function returns a vector containing the wave numbers in z-direction associated with the 3D homo...
Definition: ExpList.h:585
virtual void v_EvaluateBoundaryConditions(const NekDouble time=0.0, const std::string varName="", const NekDouble x2_in=NekConstants::kNekUnsetDouble, const NekDouble x3_in=NekConstants::kNekUnsetDouble)
Definition: ExpList.cpp:5399
void FwdTransBndConstrained(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:1696
virtual size_t v_GetNumElmts(void)
Definition: ExpList.h:1163
virtual const std::vector< bool > & v_GetLeftAdjacentFaces(void) const
Definition: ExpList.cpp:4826
Array< OneD, int > m_coeff_offset
Offset of elemental data into the array m_coeffs.
Definition: ExpList.h:1124
int EvalBasisNumModesMax(void) const
Evaulates the maximum number of modes in the elemental basis order over all elements.
Definition: ExpList.h:1526
void AddFwdBwdTraceIntegral(const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:2183
std::shared_ptr< AssemblyMapDG > & GetTraceMap(void)
Definition: ExpList.h:2166
void WriteTecplotField(std::ostream &outfile, int expansion=-1)
Definition: ExpList.h:399
virtual void v_GetMovingFrames(const SpatialDomains::GeomMMF MMFdir, const Array< OneD, const NekDouble > &CircCentre, Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: ExpList.cpp:4084
size_t GetNumElmts(void)
This function returns the number of elements in the expansion which may be different for a homogeoeno...
Definition: ExpList.h:648
virtual std::shared_ptr< AssemblyMapDG > & v_GetTraceMap()
Definition: ExpList.cpp:4291
Array< OneD, std::pair< int, int > > m_coeffsToElmt
m_coeffs to elemental value map
Definition: ExpList.h:1128
virtual void v_GetFwdBwdTracePhys(Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
Definition: ExpList.cpp:4757
NekDouble PhysEvaluate(const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &phys)
Definition: ExpList.cpp:3008
const std::unordered_map< int, int > & GetElmtToExpId(void)
This function returns the map of index inside m_exp to geom id.
Definition: ExpList.h:1030
virtual GlobalLinSysKey v_LinearAdvectionDiffusionReactionSolve(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::ConstFactorMap &factors, const StdRegions::VarCoeffMap &varcoeff, const MultiRegions::VarFactorsMap &varfactors, const Array< OneD, const NekDouble > &dirForcing, const bool PhysSpaceForcing)
Definition: ExpList.cpp:4870
virtual NekDouble v_Integral(const Array< OneD, const NekDouble > &inarray)
Definition: ExpList.cpp:3591
virtual void v_ExtractDataToCoeffs(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, std::string &field, Array< OneD, NekDouble > &coeffs, std::unordered_map< int, int > zIdToPlane)
Extract data from raw field data into expansion list.
Definition: ExpList.cpp:3954
void PhysDirectionalDeriv(const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:1822
virtual std::map< int, RobinBCInfoSharedPtr > v_GetRobinBCInfo(void)
Definition: ExpList.cpp:5411
virtual const Array< OneD, const int > & v_GetTraceBndMap()
Definition: ExpList.cpp:4299
virtual void v_GetBoundaryNormals(int i, Array< OneD, Array< OneD, NekDouble > > &normals)
Definition: ExpList.cpp:5310
void GetBoundaryNormals(int i, Array< OneD, Array< OneD, NekDouble > > &normals)
Definition: ExpList.h:2299
virtual Array< OneD, const NekDouble > v_HomogeneousEnergy(void)
Definition: ExpList.cpp:3625
virtual NekDouble v_GetHomoLen(void)
Definition: ExpList.cpp:3641
void ExtractFileBCs(const std::string &fileName, LibUtilities::CommSharedPtr comm, const std::string &varName, const std::shared_ptr< ExpList > locExpList)
void GetMovingFrames(const SpatialDomains::GeomMMF MMFdir, const Array< OneD, const NekDouble > &CircCentre, Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: ExpList.h:1787
const LocTraceToTraceMapSharedPtr & GetLocTraceToTraceMap() const
Definition: ExpList.h:2132
virtual Array< OneD, const unsigned int > v_GetZIDs(void)
Definition: ExpList.cpp:3655
Array< OneD, const NekDouble > HomogeneousEnergy(void)
This function calculates the energy associated with each one of the modesof a 3D homogeneous nD expan...
Definition: ExpList.h:569
void SetHomo1DSpecVanVisc(Array< OneD, NekDouble > visc)
This function sets the Spectral Vanishing Viscosity in homogeneous1D expansion.
Definition: ExpList.h:576
void ExtractCoeffsToCoeffs(const std::shared_ptr< ExpList > &fromExpList, const Array< OneD, const NekDouble > &fromCoeffs, Array< OneD, NekDouble > &toCoeffs)
Extract the data from fromField using fromExpList the coeffs using the basic ExpList Elemental expans...
Definition: ExpList.cpp:3938
void GetCoords(Array< OneD, NekDouble > &coord_0, Array< OneD, NekDouble > &coord_1=NullNekDouble1DArray, Array< OneD, NekDouble > &coord_2=NullNekDouble1DArray)
This function calculates the coordinates of all the elemental quadrature points .
Definition: ExpList.h:1770
const Array< OneD, int > EvalBasisNumModesMaxPerExp(void) const
Returns the vector of the number of modes in the elemental basis order over all elements.
Definition: ExpList.h:1541
std::shared_ptr< GlobalLinSys > GenGlobalBndLinSys(const GlobalLinSysKey &mkey, const AssemblyMapSharedPtr &locToGloMap)
Generate a GlobalLinSys from information provided by the key "mkey" and the mapping provided in LocTo...
Definition: ExpList.cpp:2800
virtual void v_SetBndCondBwdWeight(const int index, const NekDouble value)
Definition: ExpList.cpp:4819
void UnsetGlobalLinSys(GlobalLinSysKey, bool)
Definition: ExpList.cpp:5640
void DealiasedProd(const int num_dofs, const Array< OneD, NekDouble > &inarray1, const Array< OneD, NekDouble > &inarray2, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:1866
void GetNormals(Array< OneD, Array< OneD, NekDouble > > &normals)
Definition: ExpList.h:2174
void ExtractPhysToBnd(int i, const Array< OneD, const NekDouble > &phys, Array< OneD, NekDouble > &bnd)
Definition: ExpList.h:2292
virtual void v_SmoothField(Array< OneD, NekDouble > &field)
Definition: ExpList.cpp:2200
virtual void v_AppendFieldData(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata)
Definition: ExpList.cpp:3898
std::vector< bool > m_collectionsDoInit
Vector of bools to act as an initialise on first call flag.
Definition: ExpList.h:1118
NekDouble GetCoeff(int i)
Get the i th value (coefficient) of m_coeffs.
Definition: ExpList.h:2018
bool m_physState
The state of the array m_phys.
Definition: ExpList.h:1104
virtual void v_HomogeneousFwdTrans(const int npts, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true)
Definition: ExpList.cpp:4895
void ExtractCoeffsFromFile(const std::string &fileName, LibUtilities::CommSharedPtr comm, const std::string &varName, Array< OneD, NekDouble > &coeffs)
Definition: ExpList.cpp:3698
NekDouble Integral()
Calculates the error of the global spectral/hp element approximation.
Definition: ExpList.h:538
std::shared_ptr< ExpList > & GetPlane(int n)
Definition: ExpList.h:977
virtual void v_DealiasedProd(const int num_dofs, const Array< OneD, NekDouble > &inarray1, const Array< OneD, NekDouble > &inarray2, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:4915
void AddRightIPTBaseMatrix(const Array< OneD, const DNekMatSharedPtr > ElmtJacQuad, Array< OneD, DNekMatSharedPtr > ElmtJacCoef)
Definition: ExpList.cpp:6215
void GetLocTraceFromTracePts(const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &locTraceFwd, Array< OneD, NekDouble > &locTraceBwd)
Definition: ExpList.h:2214
virtual void v_WriteVtkPieceData(std::ostream &outfile, int expansion, std::string var)
Definition: ExpList.cpp:3468
virtual void v_WriteVtkPieceHeader(std::ostream &outfile, int expansion, int istrip)
Definition: ExpList.cpp:3328
std::shared_ptr< ExpList > GetSharedThisPtr()
Returns a shared pointer to the current object.
Definition: ExpList.h:954
int GetShapeDimension()
This function returns the dimension of the shape of the element eid.
Definition: ExpList.h:1920
void ExtractElmtDataToCoeffs(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, std::string &field, Array< OneD, NekDouble > &coeffs)
Extract the data in fielddata into the coeffs using the basic ExpList Elemental expansions rather tha...
NekDouble H1(const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &soln=NullNekDouble1DArray)
Calculates the error of the global spectral/hp element approximation.
Definition: ExpList.cpp:3753
std::vector< int > m_coll_phys_offset
Offset of elemental data into the array m_phys.
Definition: ExpList.h:1122
LibUtilities::CommSharedPtr m_comm
Communicator.
Definition: ExpList.h:1053
std::shared_ptr< GlobalLinSys > GenGlobalLinSys(const GlobalLinSysKey &mkey, const std::shared_ptr< AssemblyMapCG > &locToGloMap)
This operation constructs the global linear system of type mkey.
Definition: ExpList.cpp:2782
void WriteVtkHeader(std::ostream &outfile)
Definition: ExpList.cpp:3314
bool GetPhysState(void) const
This function indicates whether the array of physical values (implemented as m_phys) is filled or no...
Definition: ExpList.h:1643
std::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:1115
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
Definition: ExpList.h:1060
void GenerateElementVector(const int ElementID, const NekDouble scalar1, const NekDouble scalar2, Array< OneD, NekDouble > &outarray)
Generate vector v such that v[i] = scalar1 if i is in the element < ElementID. Otherwise,...
Definition: ExpList.cpp:4128
NekDouble VectorFlux(const Array< OneD, Array< OneD, NekDouble > > &inarray)
Definition: ExpList.h:563
void AddTraceIntegralToOffDiag(const Array< OneD, const NekDouble > &FwdFlux, const Array< OneD, const NekDouble > &BwdFlux, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:1008
std::shared_ptr< DNekMat > GenGlobalMatrixFull(const GlobalLinSysKey &mkey, const std::shared_ptr< AssemblyMapCG > &locToGloMap)
Definition: ExpList.cpp:2630
const Array< OneD, const SpatialDomains::BoundaryConditionShPtr > & GetBndConditions()
Definition: ExpList.h:2246
LocalRegions::ExpansionSharedPtr & GetExpFromGeomId(int n)
This function returns (a shared pointer to) the local elemental expansion of the element given a glo...
Definition: ExpList.h:2066
void SetPhysState(const bool physState)
This function manually sets whether the array of physical values (implemented as m_phys) is filled o...
Definition: ExpList.h:1636
void PhysGalerkinProjection1DScaled(const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function Galerkin projects the physical space points in inarray to outarray where inarray is ass...
Definition: ExpList.h:635
virtual void v_NormVectorIProductWRTBase(Array< OneD, const NekDouble > &V1, Array< OneD, const NekDouble > &V2, Array< OneD, NekDouble > &outarray, int BndID)
Definition: ExpList.cpp:4944
virtual void v_ExtractElmtToBndPhys(const int i, const Array< OneD, NekDouble > &elmt, Array< OneD, NekDouble > &boundary)
Definition: ExpList.cpp:5196
virtual void v_ExtractPhysToBndElmt(const int i, const Array< OneD, const NekDouble > &phys, Array< OneD, NekDouble > &bndElmt)
Definition: ExpList.cpp:5234
void CurlCurl(Array< OneD, Array< OneD, NekDouble > > &Vel, Array< OneD, Array< OneD, NekDouble > > &Q)
Definition: ExpList.h:1840
std::map< int, RobinBCInfoSharedPtr > GetRobinBCInfo()
Definition: ExpList.h:883
GlobalLinSysKey HelmSolve(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::ConstFactorMap &factors, const StdRegions::VarCoeffMap &varcoeff=StdRegions::NullVarCoeffMap, const MultiRegions::VarFactorsMap &varfactors=MultiRegions::NullVarFactorsMap, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray, const bool PhysSpaceForcing=true)
Solve helmholtz problem.
Definition: ExpList.h:1734
virtual void v_Upwind(const Array< OneD, const Array< OneD, NekDouble > > &Vec, const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &Upwind)
Definition: ExpList.cpp:4185
void EvaluateBoundaryConditions(const NekDouble time=0.0, const std::string varName="", const NekDouble=NekConstants::kNekUnsetDouble, const NekDouble=NekConstants::kNekUnsetDouble)
Definition: ExpList.h:2255
std::vector< LibUtilities::FieldDefinitionsSharedPtr > GetFieldDefinitions()
Definition: ExpList.h:895
virtual void v_FillBwdWithBwdWeight(Array< OneD, NekDouble > &weightave, Array< OneD, NekDouble > &weightjmp)
Definition: ExpList.cpp:5361
Array< OneD, NekDouble > & UpdatePhys()
This function returns (a reference to) the array (implemented as m_phys) containing the function ev...
Definition: ExpList.h:2114
virtual void v_Curl(Array< OneD, Array< OneD, NekDouble > > &Vel, Array< OneD, Array< OneD, NekDouble > > &Q)
Definition: ExpList.cpp:1962
virtual void v_PhysGalerkinProjection1DScaled(const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:6272
void GeneralMatrixOp(const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function calculates the result of the multiplication of a matrix of type specified by mkey with ...
Definition: ExpList.cpp:2418
virtual void v_SetUpPhysNormals()
: Set up a normal along the trace elements between two elements at elemental level
Definition: ExpList.cpp:5173
virtual void v_GlobalToLocal(void)
Definition: ExpList.cpp:5038
virtual LibUtilities::TranspositionSharedPtr v_GetTransposition(void)
Definition: ExpList.cpp:3633
void ExtractElmtToBndPhys(int i, const Array< OneD, NekDouble > &elmt, Array< OneD, NekDouble > &boundary)
Definition: ExpList.h:2278
void HomogeneousBwdTrans(const int npts, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true)
Definition: ExpList.h:1857
void SetWaveSpace(const bool wavespace)
Sets the wave space to the one of the possible configuration true or false.
Definition: ExpList.h:1596
const std::shared_ptr< LocalRegions::ExpansionVector > GetExp() const
This function returns the vector of elements in the expansion.
Definition: ExpList.h:2079
SpatialDomains::MeshGraphSharedPtr m_graph
Mesh associated with this expansion list.
Definition: ExpList.h:1057
virtual void v_MultiplyByInvMassMatrix(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:4849
virtual void v_AddTraceQuadPhysToField(const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &field)
Definition: ExpList.cpp:4783
virtual const Array< OneD, const NekDouble > & v_GetBndCondBwdWeight()
Definition: ExpList.cpp:4811
void GlobalToLocal(void)
Scatters from the global coefficients to the local coefficients .
Definition: ExpList.h:1979
virtual NekDouble v_L2(const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &soln=NullNekDouble1DArray)
Definition: ExpList.cpp:3546
const Array< OneD, const NekDouble > & GetBndCondBwdWeight()
Get the weight value for boundary conditions.
Definition: ExpList.h:2136
void WriteVtkFooter(std::ostream &outfile)
Definition: ExpList.cpp:3322
void MultiplyByBlockMatrix(const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:1627
virtual std::shared_ptr< ExpList > & v_GetTrace()
Definition: ExpList.cpp:4283
void GetElmtNormalLength(Array< OneD, NekDouble > &lengthsFwd, Array< OneD, NekDouble > &lengthsBwd)
Get the length of elements in boundary normal direction.
Definition: ExpList.cpp:4593
int Get1DScaledTotPoints(const NekDouble scale) const
Returns the total number of qudature points scaled by the factor scale on each 1D direction.
Definition: ExpList.h:1567
GlobalLinSysKey LinearAdvectionDiffusionReactionSolve(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::ConstFactorMap &factors, const StdRegions::VarCoeffMap &varcoeff=StdRegions::NullVarCoeffMap, const MultiRegions::VarFactorsMap &varfactors=MultiRegions::NullVarFactorsMap, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray, const bool PhysSpaceForcing=true)
Solve Advection Diffusion Reaction.
Definition: ExpList.h:1747
void PhysInterp1DScaled(const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function interpolates the physical space points in inarray to outarray using the same points def...
Definition: ExpList.h:624
virtual NekDouble v_VectorFlux(const Array< OneD, Array< OneD, NekDouble > > &inarray)
Definition: ExpList.cpp:3605
NekDouble L2(const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &soln=NullNekDouble1DArray)
This function calculates the error of the global This function calculates the error with respect to...
Definition: ExpList.h:509
virtual std::vector< LibUtilities::FieldDefinitionsSharedPtr > v_GetFieldDefinitions(void)
Definition: ExpList.cpp:3883
std::vector< int > m_coll_coeff_offset
Offset of elemental data into the array m_coeffs.
Definition: ExpList.h:1120
void IProductWRTDirectionalDerivBase(const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Directional derivative along a given direction.
Definition: ExpList.cpp:1710
void WriteVtkPieceData(std::ostream &outfile, int expansion, std::string var="v")
Definition: ExpList.h:416
void FillBwdWithBwdWeight(Array< OneD, NekDouble > &weightave, Array< OneD, NekDouble > &weightjmp)
Fill Bwd with boundary conditions.
Definition: ExpList.h:2221
virtual void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:2836
virtual Array< OneD, SpatialDomains::BoundaryConditionShPtr > & v_UpdateBndConditions()
Definition: ExpList.cpp:5389
void SetCoeff(int i, NekDouble val)
Set the i th coefficiient in m_coeffs to value val.
Definition: ExpList.h:1928
void GetMatIpwrtDeriveBase(const Array< OneD, const Array< OneD, NekDouble > > &inarray, const int nDirctn, Array< OneD, DNekMatSharedPtr > &mtxPerVar)
Definition: ExpList.cpp:5708
void LocalToGlobal(bool useComm=true)
Gathers the global coefficients from the local coefficients .
Definition: ExpList.h:1969
void GetPeriodicEntities(PeriodicMap &periodicVerts, PeriodicMap &periodicEdges, PeriodicMap &periodicFaces=NullPeriodicMap)
Definition: ExpList.h:888
NekDouble Integral(const Array< OneD, const NekDouble > &inarray)
Definition: ExpList.h:559
virtual void v_GetCoords(Array< OneD, NekDouble > &coord_0, Array< OneD, NekDouble > &coord_1=NullNekDouble1DArray, Array< OneD, NekDouble > &coord_2=NullNekDouble1DArray)
Definition: ExpList.cpp:5113
virtual void v_WriteTecplotField(std::ostream &outfile, int expansion)
Definition: ExpList.cpp:3281
virtual void v_FillBndCondFromField(const Array< OneD, NekDouble > coeffs)
Definition: ExpList.cpp:5006
Collections::CollectionVector m_collections
Definition: ExpList.h:1116
void AddTraceQuadPhysToField(const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &field)
Add Fwd and Bwd value to field, a reverse procedure of GetFwdBwdTracePhys.
Definition: ExpList.h:2208
virtual int v_GetPoolCount(std::string)
Definition: ExpList.cpp:3677
virtual void v_PhysInterp1DScaled(const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:5651
virtual void v_HomogeneousBwdTrans(const int npts, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true)
Definition: ExpList.cpp:4905
void ApplyGeomInfo()
Apply geometry information to each expansion.
Definition: ExpList.cpp:3031
virtual void v_PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2)
Definition: ExpList.cpp:1864
virtual void v_FwdTransBndConstrained(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:2177
void AddTraceIntegral(const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:2178
std::unordered_map< int, int > m_elmtToExpId
Mapping from geometry ID of element to index inside m_exp.
Definition: ExpList.h:1136
void FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:1679
void SmoothField(Array< OneD, NekDouble > &field)
Smooth a field across elements.
Definition: ExpList.h:1705
virtual void v_WriteTecplotZone(std::ostream &outfile, int expansion)
Definition: ExpList.cpp:3108
Array< OneD, const unsigned int > GetYIDs(void)
This function returns a vector containing the wave numbers in y-direction associated with the 3D homo...
Definition: ExpList.h:615
LibUtilities::SessionReaderSharedPtr m_session
Session.
Definition: ExpList.h:1055
virtual void v_PeriodicBwdCopy(const Array< OneD, const NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
Definition: ExpList.cpp:5368
void MultiplyByInvMassMatrix(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:1725
void AppendFieldData(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, Array< OneD, NekDouble > &coeffs)
Append the data in coeffs listed in elements fielddef->m_ElementIDs onto fielddata.
Definition: ExpList.h:915
virtual void v_FillBwdWithBoundCond(const Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd, bool PutFwdInBwdOnBCs)
Definition: ExpList.cpp:4774
void DivideByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Divided by the metric jacobi and quadrature weights.
Definition: ExpList.cpp:1664
virtual const Array< OneD, const std::shared_ptr< ExpList > > & v_GetBndCondExpansions(void)
Definition: ExpList.cpp:4158
void ImposeDirichletConditions(Array< OneD, NekDouble > &outarray)
Impose Dirichlet Boundary Conditions onto Array.
Definition: ExpList.h:1956
Array< OneD, int > m_phys_offset
Offset of elemental data into the array m_phys.
Definition: ExpList.h:1126
virtual void v_SetHomoLen(const NekDouble lhom)
Definition: ExpList.cpp:3649
void Upwind(const Array< OneD, const NekDouble > &Vn, const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &Upwind)
Definition: ExpList.h:2155
virtual void v_GetNormals(Array< OneD, Array< OneD, NekDouble > > &normals)
Populate normals with the normals of all expansions.
Definition: ExpList.cpp:4384
Array< OneD, SpatialDomains::BoundaryConditionShPtr > & UpdateBndConditions()
Definition: ExpList.h:2251
const DNekScalBlkMatSharedPtr GenBlockMatrix(const GlobalMatrixKey &gkey)
This function assembles the block diagonal matrix of local matrices of the type mtype.
Definition: ExpList.cpp:2239
void Curl(Array< OneD, Array< OneD, NekDouble > > &Vel, Array< OneD, Array< OneD, NekDouble > > &Q)
Definition: ExpList.h:1832
void PhysDeriv(Direction edir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d)
Definition: ExpList.h:1813
void MultiplyByElmtInvMass(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function elementally mulplies the coefficient space of Sin my the elemental inverse of the mass ...
Definition: ExpList.cpp:2130
void BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function elementally evaluates the backward transformation of the global spectral/hp element exp...
Definition: ExpList.h:1717
NekDouble GetHomoLen(void)
This function returns the Width of homogeneous direction associated with the homogeneous expansion.
Definition: ExpList.h:599
virtual void v_GetBndElmtExpansion(int i, std::shared_ptr< ExpList > &result, const bool DeclareCoeffPhysArrays)
Definition: ExpList.cpp:5186
virtual void v_FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:5052
bool GetWaveSpace(void) const
This function returns the third direction expansion condition, which can be in wave space (coefficien...
Definition: ExpList.h:1604
virtual std::shared_ptr< ExpList > & v_UpdateBndCondExpansion(int i)
Definition: ExpList.cpp:4166
int GetPhys_Offset(int n) const
Get the start offset position for a local contiguous list of quadrature points in a full array corres...
Definition: ExpList.h:2094
virtual const std::shared_ptr< LocTraceToTraceMap > & v_GetLocTraceToTraceMap(void) const
Definition: ExpList.cpp:6347
virtual void v_IProductWRTDerivBase(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:1691
void ExtractDataToCoeffs(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, std::string &field, Array< OneD, NekDouble > &coeffs, std::unordered_map< int, int > zIdToPlane=std::unordered_map< int, int >())
Extract the data in fielddata into the coeffs.
Definition: ExpList.cpp:3930
void SetModifiedBasis(const bool modbasis)
Set Modified Basis for the stability analysis.
ExpansionType m_expType
Exapnsion type.
Definition: ExpList.h:1048
void WriteTecplotHeader(std::ostream &outfile, std::string var="")
Definition: ExpList.h:391
const std::vector< bool > & GetLeftAdjacentFaces(void) const
Definition: ExpList.h:2231
void DealiasedDotProd(const int num_dofs, const Array< OneD, Array< OneD, NekDouble > > &inarray1, const Array< OneD, Array< OneD, NekDouble > > &inarray2, Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: ExpList.h:1876
void AppendFieldData(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata)
Append the element data listed in elements fielddef->m_ElementIDs onto fielddata.
Definition: ExpList.h:908
std::shared_ptr< LibUtilities::Comm > GetComm() const
Returns the comm object.
Definition: ExpList.h:964
virtual void v_AddTraceIntegral(const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:4740
Array< OneD, NekDouble > m_phys
The global expansion evaluated at the quadrature points.
Definition: ExpList.h:1096
std::vector< bool > & GetLeftAdjacentTraces(void)
Definition: ExpList.h:2305
ExpansionType GetExpType(void)
Returns the type of the expansion.
Definition: ExpList.cpp:1610
void AddTraceQuadPhysToOffDiag(const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &field)
Definition: ExpList.h:813
NekDouble Linf(const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &soln=NullNekDouble1DArray)
This function calculates the error of the global spectral/hp element approximation.
Definition: ExpList.cpp:3508
void HomogeneousFwdTrans(const int npts, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true)
Definition: ExpList.h:1848
void WriteVtkPieceHeader(std::ostream &outfile, int expansion, int istrip=0)
Definition: ExpList.h:409
int GetTotPoints(void) const
Returns the total number of quadrature points m_npoints .
Definition: ExpList.h:1557
void SetExpType(ExpansionType Type)
Returns the type of the expansion.
LibUtilities::BasisSharedPtr GetHomogeneousBasis(void)
Definition: ExpList.h:973
void WriteTecplotZone(std::ostream &outfile, int expansion=-1)
Definition: ExpList.h:395
ExpList(const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, const ExpListSharedPtr > &bndConstraint, const Array< OneD, const SpatialDomains ::BoundaryConditionShPtr > &bndCond, const LocalRegions::ExpansionVector &locexp, const SpatialDomains::MeshGraphSharedPtr &graph, const LibUtilities::CommSharedPtr &comm, const bool DeclareCoeffPhysArrays=true, const std::string variable="DefaultVar", const Collections::ImplementationType ImpType=Collections::eNoImpType)
Generate expansions for the trace space expansions used in DisContField.
virtual void v_SetHomo1DSpecVanVisc(Array< OneD, NekDouble > visc)
Definition: ExpList.h:1472
void GetBCValues(Array< OneD, NekDouble > &BndVals, const Array< OneD, NekDouble > &TotField, int BndID)
Definition: ExpList.h:1886
int GetCoordim(int eid)
This function returns the dimension of the coordinates of the element eid.
Definition: ExpList.h:1911
void GetDiagMatIpwrtBase(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, DNekMatSharedPtr > &mtxPerVar)
Definition: ExpList.cpp:5851
Describes a matrix with ordering defined by a local to global map.
std::vector< Collection > CollectionVector
Definition: Collection.h:118
static BasisSharedPtr NullBasisSharedPtr
Definition: Basis.h:350
static std::vector< unsigned int > NullUnsignedIntVector
std::shared_ptr< Basis > BasisSharedPtr
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< FieldDefinitions > FieldDefinitionsSharedPtr
Definition: FieldIO.h:184
static std::vector< NekDouble > NullNekDoubleVector
std::shared_ptr< Transposition > TranspositionSharedPtr
static Array< OneD, BasisSharedPtr > NullBasisSharedPtr1DArray
Definition: Basis.h:351
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:55
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:66
std::vector< ExpansionSharedPtr > ExpansionVector
Definition: Expansion.h:68
const std::vector< NekDouble > velocity
static ExpList NullExpList
An empty ExpList object.
Definition: ExpList.h:1498
static PeriodicMap NullPeriodicMap
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:86
std::shared_ptr< BlockMatrixMap > BlockMatrixMapShPtr
A shared pointer to a BlockMatrixMap.
Definition: ExpList.h:92
static VarFactorsMap NullVarFactorsMap
static ExpListSharedPtr NullExpListSharedPtr
Definition: ExpList.h:1499
static LibUtilities::NekManager< GlobalLinSysKey, GlobalLinSys > NullGlobalLinSysManager
Definition: ExpList.h:1503
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
static GlobalLinSysKey NullGlobalLinSysKey(StdRegions::eNoMatrixType)
std::shared_ptr< LocTraceToTraceMap > LocTraceToTraceMapSharedPtr
static const Array< OneD, ExpListSharedPtr > NullExpListSharedPtrArray
Definition: ExpList.h:2310
std::map< StdRegions::ConstFactorType, Array< OneD, NekDouble > > VarFactorsMap
std::shared_ptr< AssemblyMap > AssemblyMapSharedPtr
Definition: AssemblyMap.h:50
std::map< int, std::vector< PeriodicEntity > > PeriodicMap
std::map< GlobalMatrixKey, DNekScalBlkMatSharedPtr > BlockMatrixMap
A map between global matrix keys and their associated block matrices.
Definition: ExpList.h:90
static const NekDouble kNekUnsetDouble
GeomMMF
Principle direction for MMF.
std::shared_ptr< BoundaryConditionBase > BoundaryConditionShPtr
Definition: Conditions.h:212
std::map< int, BoundaryConditionMapShPtr > BoundaryConditionCollection
Definition: Conditions.h:219
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174
std::shared_ptr< PointGeom > PointGeomSharedPtr
Definition: Geometry.h:57
std::map< int, ExpansionInfoShPtr > ExpansionInfoMap
Definition: MeshGraph.h:141
std::map< int, CompositeSharedPtr > CompositeMap
Definition: MeshGraph.h:136
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:402
static VarCoeffMap NullVarCoeffMap
Definition: StdRegions.hpp:347
std::map< StdRegions::VarCoeffType, VarCoeffEntry > VarCoeffMap
Definition: StdRegions.hpp:346
StdRegions::ConstFactorMap factors
std::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
Definition: NekTypeDefs.hpp:79
static Array< OneD, NekDouble > NullNekDouble1DArray
double NekDouble
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825