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 
38 #include <boost/core/ignore_unused.hpp>
39 
40 #include <Collections/Collection.h>
45 #include <LocalRegions/Expansion.h>
54 #include <tinyxml.h>
55 
56 namespace Nektar
57 {
58 namespace MultiRegions
59 {
60 // Forward declarations
61 class ExpList;
62 class GlobalLinSys;
63 class AssemblyMapDG;
64 class AssemblyMapCG;
65 class GlobalLinSysKey;
66 class GlobalMatrix;
67 
69 {
70  eX,
71  eY,
72  eZ,
73  eS,
74  eN
75 };
76 
78 {
79  e0D,
80  e1D,
81  e2D,
85  e3D,
86  eNoType
87 };
88 
90 
91 /// A map between global matrix keys and their associated block
92 /// matrices.
93 typedef std::map<GlobalMatrixKey, DNekScalBlkMatSharedPtr> BlockMatrixMap;
94 /// A shared pointer to a BlockMatrixMap.
95 typedef std::shared_ptr<BlockMatrixMap> BlockMatrixMapShPtr;
96 /// Shared pointer to an ExpList object.
97 typedef std::shared_ptr<ExpList> ExpListSharedPtr;
98 
99 /// Base class for all multi-elemental spectral/hp expansions.
100 class ExpList : public std::enable_shared_from_this<ExpList>
101 {
102 public:
103  /// The default constructor using a type
105 
106  /// The copy constructor.
108  const bool DeclareCoeffPhysArrays = true);
109 
110  /// Constructor copying only elements defined in eIds.
112  const std::vector<unsigned int> &eIDs,
113  const bool DeclareCoeffPhysArrays = true,
114  const Collections::ImplementationType ImpType =
116 
117  /// Generate an ExpList from a meshgraph \a graph and session file
119  const LibUtilities::SessionReaderSharedPtr &pSession,
121  const bool DeclareCoeffPhysArrays = true,
122  const std::string &var = "DefaultVar",
123  const Collections::ImplementationType ImpType =
125 
126  /// Sets up a list of local expansions based on an expansion Map
128  const LibUtilities::SessionReaderSharedPtr &pSession,
129  const SpatialDomains::ExpansionInfoMap &expansions,
130  const bool DeclareCoeffPhysArrays = true,
131  const Collections::ImplementationType ImpType =
133 
134  //---------------------------------------------------------
135  // Specialised constructors in ExpListConstructor.cpp
136  //---------------------------------------------------------
137  /// Specialised constructors for 0D Expansions
138  /// Wrapper around LocalRegion::PointExp - used in PrePacing.cpp
141 
142  /// Generate expansions for the trace space expansions used in
143  /// DisContField.
145  const LibUtilities::SessionReaderSharedPtr &pSession,
146  const Array<OneD, const ExpListSharedPtr> &bndConstraint,
148  &bndCond,
149  const LocalRegions::ExpansionVector &locexp,
151  const LibUtilities::CommSharedPtr &comm,
152  const bool DeclareCoeffPhysArrays = true,
153  const std::string variable = "DefaultVar",
154  const Collections::ImplementationType ImpType =
156 
157  /// Generate an trace ExpList from a meshgraph \a graph and session file
159  const LibUtilities::SessionReaderSharedPtr &pSession,
160  const LocalRegions::ExpansionVector &locexp,
162  const bool DeclareCoeffPhysArrays, const std::string variable,
163  const Collections::ImplementationType ImpType =
165 
166  /// Constructor based on domain information only for 1D &
167  /// 2D boundary conditions
169  const LibUtilities::SessionReaderSharedPtr &pSession,
170  const SpatialDomains::CompositeMap &domain,
172  const bool DeclareCoeffPhysArrays = true,
173  const std::string variable = "DefaultVar",
174  bool SetToOneSpaceDimension = false,
176  const Collections::ImplementationType ImpType =
178 
179  /// The default destructor.
180  MULTI_REGIONS_EXPORT virtual ~ExpList();
181 
182  /// Returns the total number of local degrees of freedom
183  /// \f$N_{\mathrm{eof}}=\sum_{e=1}^{{N_{\mathrm{el}}}}N^{e}_m\f$.
184  inline int GetNcoeffs(void) const;
185 
186  /// Returns the total number of local degrees of freedom
187  /// for element eid
188  MULTI_REGIONS_EXPORT int GetNcoeffs(const int eid) const;
189 
190  /// Returns the type of the expansion
192 
193  /// Returns the type of the expansion
195 
196  /// Evaulates the maximum number of modes in the elemental basis
197  /// order over all elements
198  inline int EvalBasisNumModesMax(void) const;
199 
200  /// Returns the vector of the number of modes in the elemental
201  /// basis order over all elements.
203  void) const;
204 
205  /// Returns the total number of quadrature points #m_npoints
206  /// \f$=Q_{\mathrm{tot}}\f$.
207  inline int GetTotPoints(void) const;
208 
209  /// Returns the total number of quadrature points for eid's element
210  /// \f$=Q_{\mathrm{tot}}\f$.
211  inline int GetTotPoints(const int eid) const;
212 
213  /// Returns the total number of quadrature points #m_npoints
214  /// \f$=Q_{\mathrm{tot}}\f$.
215  inline int GetNpoints(void) const;
216 
217  /// Returns the total number of qudature points scaled by
218  /// the factor scale on each 1D direction
219  inline int Get1DScaledTotPoints(const NekDouble scale) const;
220 
221  /// Sets the wave space to the one of the possible configuration
222  /// true or false
223  inline void SetWaveSpace(const bool wavespace);
224 
225  /// Set Modified Basis for the stability analysis
226  inline void SetModifiedBasis(const bool modbasis);
227 
228  /// Set the \a i th value of \a m_phys to value \a val
229  inline void SetPhys(int i, NekDouble val);
230 
231  /// This function returns the third direction expansion condition,
232  /// which can be in wave space (coefficient) or not
233  /// It is stored in the variable m_WaveSpace.
234  inline bool GetWaveSpace(void) const;
235 
236  /// Fills the array #m_phys
237  inline void SetPhys(const Array<OneD, const NekDouble> &inarray);
238 
239  /// Sets the array #m_phys
240  inline void SetPhysArray(Array<OneD, NekDouble> &inarray);
241 
242  /// This function manually sets whether the array of physical
243  /// values \f$\boldsymbol{u}_l\f$ (implemented as #m_phys) is
244  /// filled or not.
245  inline void SetPhysState(const bool physState);
246 
247  /// This function indicates whether the array of physical values
248  /// \f$\boldsymbol{u}_l\f$ (implemented as #m_phys) is filled or
249  /// not.
250  inline bool GetPhysState(void) const;
251 
252  /// multiply the metric jacobi and quadrature weights
254  const Array<OneD, const NekDouble> &inarray,
255  Array<OneD, NekDouble> &outarray);
256 
257  /// Divided by the metric jacobi and quadrature weights
259  const Array<OneD, const NekDouble> &inarray,
260  Array<OneD, NekDouble> &outarray);
261 
262  /// This function calculates the inner product of a function
263  /// \f$f(\boldsymbol{x})\f$ with respect to all \em local
264  /// expansion modes \f$\phi_n^e(\boldsymbol{x})\f$.
265  inline void IProductWRTBase(const Array<OneD, const NekDouble> &inarray,
266  Array<OneD, NekDouble> &outarray);
267 
268  /// This function calculates the inner product of a function
269  /// \f$f(\boldsymbol{x})\f$ with respect to the derivative (in
270  /// direction \param dir) of all \em local expansion modes
271  /// \f$\phi_n^e(\boldsymbol{x})\f$.
273  const int dir, const Array<OneD, const NekDouble> &inarray,
274  Array<OneD, NekDouble> &outarray);
275 
277  const Array<OneD, const NekDouble> &direction,
278  const Array<OneD, const NekDouble> &inarray,
279  Array<OneD, NekDouble> &outarray);
280 
281  /// This function calculates the inner product of a
282  /// function \f$f(\boldsymbol{x})\f$ with respect to the
283  /// derivative of all \em local expansion modes
284  /// \f$\phi_n^e(\boldsymbol{x})\f$.
286  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
287  Array<OneD, NekDouble> &outarray);
288 
289  /// This function elementally evaluates the forward transformation
290  /// of a function \f$u(\boldsymbol{x})\f$ onto the global
291  /// spectral/hp expansion.
292  inline void FwdTransLocalElmt(const Array<OneD, const NekDouble> &inarray,
293  Array<OneD, NekDouble> &outarray);
294 
295  ///
296  inline void FwdTrans(const Array<OneD, const NekDouble> &inarray,
297  Array<OneD, NekDouble> &outarray);
298 
300  const NekDouble alpha,
301  const NekDouble exponent,
302  const NekDouble cutoff);
303 
304  /// This function elementally mulplies the coefficient space of
305  /// Sin my the elemental inverse of the mass matrix.
307  const Array<OneD, const NekDouble> &inarray,
308  Array<OneD, NekDouble> &outarray);
309 
310  inline void MultiplyByInvMassMatrix(
311  const Array<OneD, const NekDouble> &inarray,
312  Array<OneD, NekDouble> &outarray);
313 
314  inline void MultiplyByMassMatrix(
315  const Array<OneD, const NekDouble> &inarray,
316  Array<OneD, NekDouble> &outarray)
317  {
319  BwdTrans(inarray, tmp);
320  IProductWRTBase(tmp, outarray);
321  }
322 
323  /// Smooth a field across elements
324  inline void SmoothField(Array<OneD, NekDouble> &field);
325 
326  /// Solve helmholtz problem
327  inline void HelmSolve(
328  const Array<OneD, const NekDouble> &inarray,
329  Array<OneD, NekDouble> &outarray,
330  const StdRegions::ConstFactorMap &factors,
332  const MultiRegions::VarFactorsMap &varfactors =
335  const bool PhysSpaceForcing = true);
336 
337  /// Solve Advection Diffusion Reaction
339  const Array<OneD, Array<OneD, NekDouble>> &velocity,
340  const Array<OneD, const NekDouble> &inarray,
341  Array<OneD, NekDouble> &outarray, const NekDouble lambda,
343 
344  /// Solve Advection Diffusion Reaction
345  inline void LinearAdvectionReactionSolve(
346  const Array<OneD, Array<OneD, NekDouble>> &velocity,
347  const Array<OneD, const NekDouble> &inarray,
348  Array<OneD, NekDouble> &outarray, const NekDouble lambda,
350 
351  ///
353  const Array<OneD, const NekDouble> &inarray,
354  Array<OneD, NekDouble> &outarray);
355 
356  /// This function elementally evaluates the backward transformation
357  /// of the global spectral/hp element expansion.
358  inline void BwdTrans(const Array<OneD, const NekDouble> &inarray,
359  Array<OneD, NekDouble> &outarray);
360 
361  /// This function calculates the coordinates of all the elemental
362  /// quadrature points \f$\boldsymbol{x}_i\f$.
363  inline void GetCoords(
364  Array<OneD, NekDouble> &coord_0,
367 
368  inline void GetCoords(
369  const int eid, Array<OneD, NekDouble> &coord_0,
372 
373  // Homogeneous transforms
374  inline void HomogeneousFwdTrans(const Array<OneD, const NekDouble> &inarray,
375  Array<OneD, NekDouble> &outarray,
376  bool Shuff = true, bool UnShuff = true);
377 
378  inline void HomogeneousBwdTrans(const Array<OneD, const NekDouble> &inarray,
379  Array<OneD, NekDouble> &outarray,
380  bool Shuff = true, bool UnShuff = true);
381 
382  inline void DealiasedProd(const Array<OneD, NekDouble> &inarray1,
383  const Array<OneD, NekDouble> &inarray2,
384  Array<OneD, NekDouble> &outarray);
385 
386  inline void DealiasedDotProd(
387  const Array<OneD, Array<OneD, NekDouble>> &inarray1,
388  const Array<OneD, Array<OneD, NekDouble>> &inarray2,
389  Array<OneD, Array<OneD, NekDouble>> &outarray);
390 
391  inline void GetBCValues(Array<OneD, NekDouble> &BndVals,
392  const Array<OneD, NekDouble> &TotField, int BndID);
393 
396  Array<OneD, NekDouble> &outarray,
397  int BndID);
398 
399  inline void NormVectorIProductWRTBase(
401  Array<OneD, NekDouble> &outarray);
402 
403  /// Apply geometry information to each expansion.
405 
406  /// Reset geometry information and reset matrices
408  {
409  v_Reset();
410  }
411 
412  void WriteTecplotHeader(std::ostream &outfile, std::string var = "")
413  {
414  v_WriteTecplotHeader(outfile, var);
415  }
416 
417  void WriteTecplotZone(std::ostream &outfile, int expansion = -1)
418  {
419  v_WriteTecplotZone(outfile, expansion);
420  }
421 
422  void WriteTecplotField(std::ostream &outfile, int expansion = -1)
423  {
424  v_WriteTecplotField(outfile, expansion);
425  }
426 
427  void WriteTecplotConnectivity(std::ostream &outfile, int expansion = -1)
428  {
429  v_WriteTecplotConnectivity(outfile, expansion);
430  }
431 
432  MULTI_REGIONS_EXPORT void WriteVtkHeader(std::ostream &outfile);
433  MULTI_REGIONS_EXPORT void WriteVtkFooter(std::ostream &outfile);
434 
435  void WriteVtkPieceHeader(std::ostream &outfile, int expansion,
436  int istrip = 0)
437  {
438  v_WriteVtkPieceHeader(outfile, expansion, istrip);
439  }
440 
441  MULTI_REGIONS_EXPORT void WriteVtkPieceFooter(std::ostream &outfile,
442  int expansion);
443 
444  void WriteVtkPieceData(std::ostream &outfile, int expansion,
445  std::string var = "v")
446  {
447  v_WriteVtkPieceData(outfile, expansion, var);
448  }
449 
450  /// This function returns the dimension of the coordinates of the
451  /// element \a eid.
452  // inline
453  MULTI_REGIONS_EXPORT int GetCoordim(int eid);
454 
455  /// Set the \a i th coefficiient in \a m_coeffs to value \a val
456  inline void SetCoeff(int i, NekDouble val);
457 
458  /// Set the \a i th coefficiient in #m_coeffs to value \a val
459  inline void SetCoeffs(int i, NekDouble val);
460 
461  /// Set the #m_coeffs array to inarray
462  inline void SetCoeffsArray(Array<OneD, NekDouble> &inarray);
463 
464  /// This function returns the dimension of the shape of the
465  /// element \a eid.
466  // inline
468 
469  /// This function returns (a reference to) the array
470  /// \f$\boldsymbol{\hat{u}}_l\f$ (implemented as #m_coeffs)
471  /// containing all local expansion coefficients.
472  inline const Array<OneD, const NekDouble> &GetCoeffs() const;
473 
474  /// Impose Dirichlet Boundary Conditions onto Array
475  inline void ImposeDirichletConditions(Array<OneD, NekDouble> &outarray);
476 
477  /// Fill Bnd Condition expansion from the values stored in expansion
478  inline void FillBndCondFromField(void);
479 
480  /// Fill Bnd Condition expansion in nreg from the values
481  /// stored in expansion
482  inline void FillBndCondFromField(const int nreg);
483 
484  /// Gathers the global coefficients \f$\boldsymbol{\hat{u}}_g\f$
485  /// from the local coefficients \f$\boldsymbol{\hat{u}}_l\f$.
486  // inline
487  MULTI_REGIONS_EXPORT inline void LocalToGlobal(bool useComm = true);
488 
490  const Array<OneD, const NekDouble> &inarray,
491  Array<OneD, NekDouble> &outarray, bool useComm = true);
492 
493  /// Scatters from the global coefficients
494  /// \f$\boldsymbol{\hat{u}}_g\f$ to the local coefficients
495  /// \f$\boldsymbol{\hat{u}}_l\f$.
496  // inline
497  MULTI_REGIONS_EXPORT inline void GlobalToLocal(void);
498 
499  /**
500  * This operation is evaluated as:
501  * \f{tabbing}
502  * \hspace{1cm} \= Do \= $e=$ $1, N_{\mathrm{el}}$ \\
503  * > > Do \= $i=$ $0,N_m^e-1$ \\
504  * > > > $\boldsymbol{\hat{u}}^{e}[i] = \mbox{sign}[e][i] \cdot
505  * \boldsymbol{\hat{u}}_g[\mbox{map}[e][i]]$ \\
506  * > > continue \\
507  * > continue
508  * \f}
509  * where \a map\f$[e][i]\f$ is the mapping array and \a
510  * sign\f$[e][i]\f$ is an array of similar dimensions ensuring the
511  * correct modal connectivity between the different elements (both
512  * these arrays are contained in the data member #m_locToGloMap). This
513  * operation is equivalent to the scatter operation
514  * \f$\boldsymbol{\hat{u}}_l=\mathcal{A}\boldsymbol{\hat{u}}_g\f$,
515  * where \f$\mathcal{A}\f$ is the
516  * \f$N_{\mathrm{eof}}\times N_{\mathrm{dof}}\f$ permutation matrix.
517  *
518  * @param inarray An array of size \f$N_\mathrm{dof}\f$
519  * containing the global degrees of freedom
520  * \f$\boldsymbol{x}_g\f$.
521  * @param outarray The resulting local degrees of freedom
522  * \f$\boldsymbol{x}_l\f$ will be stored in this
523  * array of size \f$N_\mathrm{eof}\f$.
524  */
526  const Array<OneD, const NekDouble> &inarray,
527  Array<OneD, NekDouble> &outarray);
528 
529  /// Get the \a i th value (coefficient) of #m_coeffs
530  inline NekDouble GetCoeff(int i);
531 
532  /// Get the \a i th value (coefficient) of #m_coeffs
533  inline NekDouble GetCoeffs(int i);
534 
535  /// This function returns (a reference to) the array
536  /// \f$\boldsymbol{u}_l\f$ (implemented as #m_phys) containing the
537  /// function \f$u^{\delta}(\boldsymbol{x})\f$ evaluated at the
538  /// quadrature points.
539  // inline
541 
542  /// This function calculates the \f$L_\infty\f$ error of the global
543  /// spectral/hp element approximation.
545  Linf(const Array<OneD, const NekDouble> &inarray,
547 
548  /// This function calculates the \f$L_2\f$ error with
549  /// respect to soln of the global
550  /// spectral/hp element approximation.
552  const Array<OneD, const NekDouble> &inarray,
554  {
555  return v_L2(inarray, soln);
556  }
557 
558  /// Calculates the \f$H^1\f$ error of the global spectral/hp
559  /// element approximation.
561  H1(const Array<OneD, const NekDouble> &inarray,
563 
564  /**
565  * The integration is evaluated locally, that is
566  * \f[\int
567  * f(\boldsymbol{x})d\boldsymbol{x}=\sum_{e=1}^{{N_{\mathrm{el}}}}
568  * \left\{\int_{\Omega_e}f(\boldsymbol{x})d\boldsymbol{x}\right\}, \f]
569  * where the integration over the separate elements is done by the
570  * function StdRegions#StdExpansion#Integral, which discretely
571  * evaluates the integral using Gaussian quadrature.
572  *
573  * Note that the array #m_phys should be filled with the values of the
574  * function \f$f(\boldsymbol{x})\f$ at the quadrature points
575  * \f$\boldsymbol{x}_i\f$.
576  *
577  * @return The value of the discretely evaluated integral
578  * \f$\int f(\boldsymbol{x})d\boldsymbol{x}\f$.
579  */
581  {
582  ASSERTL1(m_physState == true, "local physical space is not true ");
583 
584  return Integral(m_phys);
585  }
586 
587  /**
588  * The integration is evaluated locally, that is
589  * \f[\int
590  * f(\boldsymbol{x})d\boldsymbol{x}=\sum_{e=1}^{{N_{\mathrm{el}}}}
591  * \left\{\int_{\Omega_e}f(\boldsymbol{x})d\boldsymbol{x}\right\}, \f]
592  * where the integration over the separate elements is done by the
593  * function StdRegions#StdExpansion#Integral, which discretely
594  * evaluates the integral using Gaussian quadrature.
595  *
596  * @param inarray An array of size \f$Q_{\mathrm{tot}}\f$
597  * containing the values of the function
598  * \f$f(\boldsymbol{x})\f$ at the quadrature
599  * points \f$\boldsymbol{x}_i\f$.
600  * @return The value of the discretely evaluated integral
601  * \f$\int f(\boldsymbol{x})d\boldsymbol{x}\f$.
602  */
604  {
605  return v_Integral(inarray);
606  }
607 
609  {
610  return v_VectorFlux(inarray);
611  }
612 
613  /// This function calculates the energy associated with
614  /// each one of the modesof a 3D homogeneous nD expansion
616  {
617  return v_HomogeneousEnergy();
618  }
619 
620  /// This function sets the Spectral Vanishing Viscosity
621  /// in homogeneous1D expansion.
623  {
625  }
626 
627  /// This function returns a vector containing the wave
628  /// numbers in z-direction associated
629  /// with the 3D homogenous expansion. Required if a
630  /// parellelisation is applied in the Fourier direction
632  {
633  return v_GetZIDs();
634  }
635 
636  /// This function returns the transposition class
637  /// associated with the homogeneous expansion.
639  {
640  return v_GetTransposition();
641  }
642 
643  /// This function returns the Width of homogeneous direction
644  /// associated with the homogeneous expansion.
646  {
647  return v_GetHomoLen();
648  }
649 
650  /// This function sets the Width of homogeneous direction
651  /// associated with the homogeneous expansion.
652  void SetHomoLen(const NekDouble lhom)
653  {
654  return v_SetHomoLen(lhom);
655  }
656 
657  /// This function returns a vector containing the wave
658  /// numbers in y-direction associated
659  /// with the 3D homogenous expansion. Required if a
660  /// parellelisation is applied in the Fourier direction
662  {
663  return v_GetYIDs();
664  }
665 
666  /// This function interpolates the physical space points in
667  /// \a inarray to \a outarray using the same points defined in the
668  /// expansion but where the number of points are rescaled
669  /// by \a 1DScale
670  void PhysInterp1DScaled(const NekDouble scale,
671  const Array<OneD, NekDouble> &inarray,
672  Array<OneD, NekDouble> &outarray)
673  {
674  v_PhysInterp1DScaled(scale, inarray, outarray);
675  }
676 
677  /// This function Galerkin projects the physical space points in
678  /// \a inarray to \a outarray where inarray is assumed to
679  /// be defined in the expansion but where the number of
680  /// points are rescaled by \a 1DScale
682  const Array<OneD, NekDouble> &inarray,
683  Array<OneD, NekDouble> &outarray)
684  {
685  v_PhysGalerkinProjection1DScaled(scale, inarray, outarray);
686  }
687 
688  /// This function returns the number of elements in the expansion.
689  inline int GetExpSize(void);
690 
691  /// This function returns the number of elements in the
692  /// expansion which may be different for a homogeoenous extended
693  /// expansionp.
694  inline int GetNumElmts(void)
695  {
696  return v_GetNumElmts();
697  }
698 
699  /// This function returns the vector of elements in the expansion.
700  inline const std::shared_ptr<LocalRegions::ExpansionVector> GetExp() const;
701 
702  /// This function returns (a shared pointer to) the local elemental
703  /// expansion of the \f$n^{\mathrm{th}}\f$ element.
704  inline LocalRegions::ExpansionSharedPtr &GetExp(int n) const;
705 
706  /// This function returns (a shared pointer to) the local elemental
707  /// expansion containing the arbitrary point given by \a gloCoord.
709  const Array<OneD, const NekDouble> &gloCoord);
710 
711  /**
712  * @brief This function returns the index of the local elemental
713  * expansion containing the arbitrary point given by \a gloCoord,
714  * within a distance tolerance of tol.
715  *
716  * If returnNearestElmt is true and no element contains the point,
717  * this function returns the nearest element whose bounding box contains
718  * the point. The bounding box has a 10% margin in each direction.
719  *
720  * @param gloCoord (input) coordinate in physics space
721  * @param locCoords (output) local coordinate xi in the returned element
722  * @param tol distance tolerance to judge if a point lies in an
723  * element
724  * @param returnNearestElmt if true and no element contains this point, the
725  * nearest element whose bounding box contains this
726  * point is returned
727  * @param cachedId an initial guess of the most possible element index
728  * @param maxDistance if returnNearestElmt is set as true, the nearest
729  * element will be returned. But the distance of the
730  * nearest element and this point should be <=
731  * maxDistance.
732  *
733  * @return element index; if no element is found, -1 is returned.
734  */
736  const Array<OneD, const NekDouble> &gloCoord, NekDouble tol = 0.0,
737  bool returnNearestElmt = false, int cachedId = -1,
738  NekDouble maxDistance = 1e6);
739 
740  /** This function returns the index and the Local
741  * Cartesian Coordinates \a locCoords of the local
742  * elemental expansion containing the arbitrary point
743  * given by \a gloCoords.
744  **/
746  const Array<OneD, const NekDouble> &gloCoords,
747  Array<OneD, NekDouble> &locCoords, NekDouble tol = 0.0,
748  bool returnNearestElmt = false, int cachedId = -1,
749  NekDouble maxDistance = 1e6);
750 
751  /** This function return the expansion field value
752  * at the coordinates given as input.
753  **/
756  const Array<OneD, const NekDouble> &phys);
757 
758  /// Get the start offset position for a global list of #m_coeffs
759  /// correspoinding to element n.
760  inline int GetCoeff_Offset(int n) const;
761 
762  /// Get the start offset position for a global list of m_phys
763  /// correspoinding to element n.
764  inline int GetPhys_Offset(int n) const;
765 
766  /// This function returns (a reference to) the array
767  /// \f$\boldsymbol{\hat{u}}_l\f$ (implemented as #m_coeffs)
768  /// containing all local expansion coefficients.
770 
771  /// This function returns (a reference to) the array
772  /// \f$\boldsymbol{u}_l\f$ (implemented as #m_phys) containing the
773  /// function \f$u^{\delta}(\boldsymbol{x})\f$ evaluated at the
774  /// quadrature points.
776 
777  inline void PhysDeriv(Direction edir,
778  const Array<OneD, const NekDouble> &inarray,
779  Array<OneD, NekDouble> &out_d);
780 
781  /// This function discretely evaluates the derivative of a function
782  /// \f$f(\boldsymbol{x})\f$ on the domain consisting of all
783  /// elements of the expansion.
784  inline void PhysDeriv(
785  const Array<OneD, const NekDouble> &inarray,
786  Array<OneD, NekDouble> &out_d0,
789 
790  inline void PhysDeriv(const int dir,
791  const Array<OneD, const NekDouble> &inarray,
792  Array<OneD, NekDouble> &out_d);
793 
794  inline void CurlCurl(Array<OneD, Array<OneD, NekDouble>> &Vel,
796 
797  inline void PhysDirectionalDeriv(
798  const Array<OneD, const NekDouble> &direction,
799  const Array<OneD, const NekDouble> &inarray,
800  Array<OneD, NekDouble> &outarray);
801 
802  inline void GetMovingFrames(const SpatialDomains::GeomMMF MMFdir,
803  const Array<OneD, const NekDouble> &CircCentre,
804  Array<OneD, Array<OneD, NekDouble>> &outarray);
805 
806  // functions associated with DisContField
809 
810  /// Get the weight value for boundary conditions
812 
813  /// Set the weight value for boundary conditions
814  inline void SetBndCondBwdWeight(const int index, const NekDouble value);
815 
816  inline std::shared_ptr<ExpList> &UpdateBndCondExpansion(int i);
817 
818  inline void Upwind(const Array<OneD, const NekDouble> &Vn,
819  const Array<OneD, const NekDouble> &Fwd,
820  const Array<OneD, const NekDouble> &Bwd,
822 
823  inline void Upwind(const Array<OneD, const Array<OneD, NekDouble>> &Vec,
824  const Array<OneD, const NekDouble> &Fwd,
825  const Array<OneD, const NekDouble> &Bwd,
827 
828  /**
829  * Return a reference to the trace space associated with this
830  * expansion list.
831  */
832  inline std::shared_ptr<ExpList> &GetTrace();
833 
834  inline std::shared_ptr<AssemblyMapDG> &GetTraceMap(void);
835 
836  inline const Array<OneD, const int> &GetTraceBndMap(void);
837 
838  inline void GetNormals(Array<OneD, Array<OneD, NekDouble>> &normals);
839 
840  /// Get the length of elements in boundary normal direction
842  Array<OneD, NekDouble> &lengthsFwd, Array<OneD, NekDouble> &lengthsBwd);
843 
844  /// Get the weight value for boundary conditions
845  /// for boundary average and jump calculations
847  Array<OneD, NekDouble> &weightJump);
848 
849  inline void AddTraceIntegral(const Array<OneD, const NekDouble> &Fx,
851  Array<OneD, NekDouble> &outarray);
852 
853  inline void AddTraceIntegral(const Array<OneD, const NekDouble> &Fn,
854  Array<OneD, NekDouble> &outarray);
855 
857  const Array<OneD, const NekDouble> &Bwd,
858  Array<OneD, NekDouble> &outarray);
859 
862 
863  inline void GetFwdBwdTracePhys(const Array<OneD, const NekDouble> &field,
866  bool FillBnd = true,
867  bool PutFwdInBwdOnBCs = false,
868  bool DoExchange = true);
869 
870  inline void FillBwdWithBoundCond(const Array<OneD, NekDouble> &Fwd,
872  bool PutFwdInBwdOnBCs = false);
873 
874  /// Add Fwd and Bwd value to field,
875  /// a reverse procedure of GetFwdBwdTracePhys
877  const Array<OneD, const NekDouble> &Bwd,
878  Array<OneD, NekDouble> &field);
879 
881  const Array<OneD, const NekDouble> &Fwd,
883  {
884  v_AddTraceQuadPhysToOffDiag(Fwd, Bwd, field);
885  }
886 
888  const Array<OneD, const NekDouble> &Bwd,
889  Array<OneD, NekDouble> &locTraceFwd,
890  Array<OneD, NekDouble> &locTraceBwd);
891 
892  /// Fill Bwd with boundary conditions
893  inline void FillBwdWithBwdWeight(Array<OneD, NekDouble> &weightave,
894  Array<OneD, NekDouble> &weightjmp);
895 
896  /// Copy and fill the Periodic boundaries
897  inline void PeriodicBwdCopy(const Array<OneD, const NekDouble> &Fwd,
899 
900  inline const std::vector<bool> &GetLeftAdjacentFaces(void) const;
901 
902  inline void ExtractTracePhys(Array<OneD, NekDouble> &outarray);
903 
904  inline void ExtractTracePhys(const Array<OneD, const NekDouble> &inarray,
905  Array<OneD, NekDouble> &outarray);
906 
908  &GetBndConditions();
909 
912 
913  inline void EvaluateBoundaryConditions(
914  const NekDouble time = 0.0, const std::string varName = "",
917 
918  // Routines for continous matrix solution
919  /// This function calculates the result of the multiplication of a
920  /// matrix of type specified by \a mkey with a vector given by \a
921  /// inarray.
923  const GlobalMatrixKey &gkey,
924  const Array<OneD, const NekDouble> &inarray,
925  Array<OneD, NekDouble> &outarray);
926 
927  inline void SetUpPhysNormals();
928 
929  inline void GetBoundaryToElmtMap(Array<OneD, int> &ElmtID,
930  Array<OneD, int> &EdgeID);
931 
932  inline void GetBndElmtExpansion(int i, std::shared_ptr<ExpList> &result,
933  const bool DeclareCoeffPhysArrays = true);
934 
935  inline void ExtractElmtToBndPhys(int i, const Array<OneD, NekDouble> &elmt,
936  Array<OneD, NekDouble> &boundary);
937 
938  inline void ExtractPhysToBndElmt(int i,
939  const Array<OneD, const NekDouble> &phys,
940  Array<OneD, NekDouble> &bndElmt);
941 
942  inline void ExtractPhysToBnd(int i,
943  const Array<OneD, const NekDouble> &phys,
945 
946  inline void GetBoundaryNormals(
947  int i, Array<OneD, Array<OneD, NekDouble>> &normals);
948 
950  std::vector<LibUtilities::FieldDefinitionsSharedPtr> &fielddef,
951  int NumHomoDir = 0,
954  std::vector<NekDouble> &HomoLen = LibUtilities::NullNekDoubleVector,
955  bool homoStrips = false,
956  std::vector<unsigned int> &HomoSIDs =
958  std::vector<unsigned int> &HomoZIDs =
960  std::vector<unsigned int> &HomoYIDs =
962 
963  std::map<int, RobinBCInfoSharedPtr> GetRobinBCInfo()
964  {
965  return v_GetRobinBCInfo();
966  }
967 
968  void GetPeriodicEntities(PeriodicMap &periodicVerts,
969  PeriodicMap &periodicEdges,
970  PeriodicMap &periodicFaces = NullPeriodicMap)
971  {
972  v_GetPeriodicEntities(periodicVerts, periodicEdges, periodicFaces);
973  }
974 
975  std::vector<LibUtilities::FieldDefinitionsSharedPtr> GetFieldDefinitions()
976  {
977  return v_GetFieldDefinitions();
978  }
979 
981  std::vector<LibUtilities::FieldDefinitionsSharedPtr> &fielddef)
982  {
983  v_GetFieldDefinitions(fielddef);
984  }
985 
986  /// Append the element data listed in elements
987  /// fielddef->m_ElementIDs onto fielddata
989  std::vector<NekDouble> &fielddata)
990  {
991  v_AppendFieldData(fielddef, fielddata);
992  }
993 
994  /// Append the data in coeffs listed in elements
995  /// fielddef->m_ElementIDs onto fielddata
997  std::vector<NekDouble> &fielddata,
998  Array<OneD, NekDouble> &coeffs)
999  {
1000  v_AppendFieldData(fielddef, fielddata, coeffs);
1001  }
1002 
1003  /** \brief Extract the data in fielddata into the coeffs
1004  * using the basic ExpList Elemental expansions rather
1005  * than planes in homogeneous case
1006  */
1009  std::vector<NekDouble> &fielddata, std::string &field,
1010  Array<OneD, NekDouble> &coeffs);
1011 
1012  /** \brief Extract the data from fromField using
1013  * fromExpList the coeffs using the basic ExpList
1014  * Elemental expansions rather than planes in homogeneous
1015  * case
1016  */
1018  const std::shared_ptr<ExpList> &fromExpList,
1019  const Array<OneD, const NekDouble> &fromCoeffs,
1020  Array<OneD, NekDouble> &toCoeffs);
1021 
1022  // Extract data in fielddata into the m_coeffs_list for the 3D stability
1023  // analysis (base flow is 2D)
1026  std::vector<NekDouble> &fielddata, std::string &field,
1027  Array<OneD, NekDouble> &coeffs);
1028 
1030  const int ElementID, const NekDouble scalar1, const NekDouble scalar2,
1031  Array<OneD, NekDouble> &outarray);
1032 
1033  /// Returns a shared pointer to the current object.
1034  std::shared_ptr<ExpList> GetSharedThisPtr()
1035  {
1036  return shared_from_this();
1037  }
1038 
1039  /// Returns the session object
1040  std::shared_ptr<LibUtilities::SessionReader> GetSession() const
1041  {
1042  return m_session;
1043  }
1044 
1045  /// Returns the comm object
1046  std::shared_ptr<LibUtilities::Comm> GetComm() const
1047  {
1048  return m_comm;
1049  }
1050 
1052  {
1053  return m_graph;
1054  }
1055 
1056  // Wrapper functions for Homogeneous Expansions
1058  {
1059  return v_GetHomogeneousBasis();
1060  }
1061 
1062  std::shared_ptr<ExpList> &GetPlane(int n)
1063  {
1064  return v_GetPlane(n);
1065  }
1066 
1069 
1071 
1072  /// Get m_coeffs to elemental value map
1074  &GetCoeffsToElmt() const;
1075 
1079  Array<OneD, DNekMatSharedPtr> &fieldMat);
1080 
1082  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
1083  const int nDirctn, Array<OneD, DNekMatSharedPtr> &mtxPerVar);
1084 
1086  const TensorOfArray3D<NekDouble> &inarray,
1087  Array<OneD, DNekMatSharedPtr> &mtxPerVar);
1088 
1090  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
1091  Array<OneD, DNekMatSharedPtr> &mtxPerVar);
1092 
1094  const Array<OneD, const NekDouble> &FwdFlux,
1095  const Array<OneD, const NekDouble> &BwdFlux,
1096  Array<OneD, NekDouble> &outarray)
1097  {
1098  v_AddTraceIntegralToOffDiag(FwdFlux, BwdFlux, outarray);
1099  }
1100 
1102  const int dir, const Array<OneD, const DNekMatSharedPtr> ElmtJacQuad,
1103  Array<OneD, DNekMatSharedPtr> ElmtJacCoef);
1104 
1106  const Array<OneD, const DNekMatSharedPtr> ElmtJacQuad,
1107  Array<OneD, DNekMatSharedPtr> ElmtJacCoef);
1108 
1110  GetLocTraceToTraceMap() const;
1111 
1112  // Return the internal vector which identifieds if trace
1113  // is left adjacent definiing which trace the normal
1114  // points otwards from
1115  MULTI_REGIONS_EXPORT std::vector<bool> &GetLeftAdjacentTraces(void);
1116 
1117 protected:
1118  /// Exapnsion type
1120 
1121  std::shared_ptr<DNekMat> GenGlobalMatrixFull(
1122  const GlobalLinSysKey &mkey,
1123  const std::shared_ptr<AssemblyMapCG> &locToGloMap);
1124 
1125  /// Communicator
1127 
1128  /// Session
1130 
1131  /// Mesh associated with this expansion list.
1133 
1134  /// The total number of local degrees of freedom. #m_ncoeffs
1135  /// \f$=N_{\mathrm{eof}}=\sum_{e=1}^{{N_{\mathrm{el}}}}N^{e}_l\f$
1137 
1138  /** The total number of quadrature points. #m_npoints
1139  *\f$=Q_{\mathrm{tot}}=\sum_{e=1}^{{N_{\mathrm{el}}}}N^{e}_Q\f$
1140  **/
1142 
1143  /**
1144  * \brief Concatenation of all local expansion coefficients.
1145  *
1146  * The array of length #m_ncoeffs\f$=N_{\mathrm{eof}}\f$ which is
1147  * the concatenation of the local expansion coefficients
1148  * \f$\hat{u}_n^e\f$ over all \f$N_{\mathrm{el}}\f$ elements
1149  * \f[\mathrm{\texttt{m\_coeffs}}=\boldsymbol{\hat{u}}_{l} =
1150  * \underline{\boldsymbol{\hat{u}}}^e = \left [ \begin{array}{c}
1151  * \boldsymbol{\hat{u}}^{1} \ \
1152  * \boldsymbol{\hat{u}}^{2} \ \
1153  * \vdots \ \
1154  * \boldsymbol{\hat{u}}^{{{N_{\mathrm{el}}}}} \end{array} \right ],
1155  * \quad
1156  * \mathrm{where}\quad \boldsymbol{\hat{u}}^{e}[n]=\hat{u}_n^{e}\f]
1157  */
1159 
1160  /**
1161  * \brief The global expansion evaluated at the quadrature points
1162  *
1163  * The array of length #m_npoints\f$=Q_{\mathrm{tot}}\f$ containing
1164  * the evaluation of \f$u^{\delta}(\boldsymbol{x})\f$ at the
1165  * quadrature points \f$\boldsymbol{x}_i\f$.
1166  * \f[\mathrm{\texttt{m\_phys}}=\boldsymbol{u}_{l} =
1167  * \underline{\boldsymbol{u}}^e = \left [ \begin{array}{c}
1168  * \boldsymbol{u}^{1} \ \
1169  * \boldsymbol{u}^{2} \ \
1170  * \vdots \ \
1171  * \boldsymbol{u}^{{{N_{\mathrm{el}}}}} \end{array} \right ],\quad
1172  * \mathrm{where}\quad
1173  * \boldsymbol{u}^{e}[i]=u^{\delta}(\boldsymbol{x}_i)\f]
1174  */
1176 
1177  /**
1178  * \brief The state of the array #m_phys.
1179  *
1180  * Indicates whether the array #m_phys, created to contain the
1181  * evaluation of \f$u^{\delta}(\boldsymbol{x})\f$ at the quadrature
1182  * points \f$\boldsymbol{x}_i\f$, is filled with these values.
1183  */
1185 
1186  /**
1187  * \brief The list of local expansions.
1188  *
1189  * The (shared pointer to the) vector containing (shared pointers
1190  * to) all local expansions. The fact that the local expansions are
1191  * all stored as a (pointer to a) #StdExpansion, the abstract base
1192  * class for all local expansions, allows a general implementation
1193  * where most of the routines for the derived classes are defined
1194  * in the #ExpList base class.
1195  */
1196  std::shared_ptr<LocalRegions::ExpansionVector> m_exp;
1197 
1199 
1200  /// Vector of bools to act as an initialise on first call flag
1201  std::vector<bool> m_collectionsDoInit;
1202 
1203  /// Offset of elemental data into the array #m_coeffs
1204  std::vector<int> m_coll_coeff_offset;
1205 
1206  /// Offset of elemental data into the array #m_phys
1207  std::vector<int> m_coll_phys_offset;
1208 
1209  /// Offset of elemental data into the array #m_coeffs
1211 
1212  /// Offset of elemental data into the array #m_phys
1214 
1215  /// m_coeffs to elemental value map
1217 
1219 
1220  //@todo should this be in ExpList or
1221  // ExpListHomogeneous1D.cpp it's a bool which determine if
1222  // the expansion is in the wave space (coefficient space)
1223  // or not
1225 
1226  /// Mapping from geometry ID of element to index inside #m_exp
1227  std::unordered_map<int, int> m_elmtToExpId;
1228 
1229  /// This function assembles the block diagonal matrix of local
1230  /// matrices of the type \a mtype.
1232 
1234 
1235  void MultiplyByBlockMatrix(const GlobalMatrixKey &gkey,
1236  const Array<OneD, const NekDouble> &inarray,
1237  Array<OneD, NekDouble> &outarray);
1238 
1239  /// Generates a global matrix from the given key and map.
1240  std::shared_ptr<GlobalMatrix> GenGlobalMatrix(
1241  const GlobalMatrixKey &mkey,
1242  const std::shared_ptr<AssemblyMapCG> &locToGloMap);
1243 
1245  const std::shared_ptr<DNekMat> &Gmat,
1246  Array<OneD, NekDouble> &EigValsReal,
1247  Array<OneD, NekDouble> &EigValsImag,
1249 
1250  /// This operation constructs the global linear system of type \a
1251  /// mkey.
1252  std::shared_ptr<GlobalLinSys> GenGlobalLinSys(
1253  const GlobalLinSysKey &mkey,
1254  const std::shared_ptr<AssemblyMapCG> &locToGloMap);
1255 
1256  /// Generate a GlobalLinSys from information provided by the key
1257  /// "mkey" and the mapping provided in LocToGloBaseMap.
1258  std::shared_ptr<GlobalLinSys> GenGlobalBndLinSys(
1259  const GlobalLinSysKey &mkey, const AssemblyMapSharedPtr &locToGloMap);
1260 
1261  // Virtual prototypes
1262 
1263  virtual int v_GetNumElmts(void)
1264  {
1265  return (*m_exp).size();
1266  }
1267 
1269  &v_GetBndCondExpansions(void);
1270 
1272 
1273  virtual void v_SetBndCondBwdWeight(const int index, const NekDouble value);
1274 
1275  virtual std::shared_ptr<ExpList> &v_UpdateBndCondExpansion(int i);
1276 
1277  virtual void v_Upwind(const Array<OneD, const Array<OneD, NekDouble>> &Vec,
1278  const Array<OneD, const NekDouble> &Fwd,
1279  const Array<OneD, const NekDouble> &Bwd,
1281 
1282  virtual void v_Upwind(const Array<OneD, const NekDouble> &Vn,
1283  const Array<OneD, const NekDouble> &Fwd,
1284  const Array<OneD, const NekDouble> &Bwd,
1286 
1287  virtual std::shared_ptr<ExpList> &v_GetTrace();
1288 
1289  virtual std::shared_ptr<AssemblyMapDG> &v_GetTraceMap();
1290 
1291  virtual const Array<OneD, const int> &v_GetTraceBndMap();
1292 
1293  virtual const std::shared_ptr<LocTraceToTraceMap> &v_GetLocTraceToTraceMap(
1294  void) const;
1295 
1296  virtual std::vector<bool> &v_GetLeftAdjacentTraces(void);
1297 
1298  /// Populate \a normals with the normals of all expansions.
1299  virtual void v_GetNormals(Array<OneD, Array<OneD, NekDouble>> &normals);
1300 
1301  virtual void v_AddTraceIntegral(const Array<OneD, const NekDouble> &Fx,
1302  const Array<OneD, const NekDouble> &Fy,
1303  Array<OneD, NekDouble> &outarray);
1304 
1305  virtual void v_AddTraceIntegral(const Array<OneD, const NekDouble> &Fn,
1306  Array<OneD, NekDouble> &outarray);
1307 
1308  virtual void v_AddFwdBwdTraceIntegral(
1309  const Array<OneD, const NekDouble> &Fwd,
1310  const Array<OneD, const NekDouble> &Bwd,
1311  Array<OneD, NekDouble> &outarray);
1312 
1313  virtual void v_GetFwdBwdTracePhys(Array<OneD, NekDouble> &Fwd,
1314  Array<OneD, NekDouble> &Bwd);
1315 
1316  virtual void v_GetFwdBwdTracePhys(const Array<OneD, const NekDouble> &field,
1319  bool FillBnd = true,
1320  bool PutFwdInBwdOnBCs = false,
1321  bool DoExchange = true);
1322 
1323  virtual void v_FillBwdWithBoundCond(const Array<OneD, NekDouble> &Fwd,
1325  bool PutFwdInBwdOnBCs);
1326 
1327  virtual void v_AddTraceQuadPhysToField(
1328  const Array<OneD, const NekDouble> &Fwd,
1330 
1331  virtual void v_AddTraceQuadPhysToOffDiag(
1332  const Array<OneD, const NekDouble> &Fwd,
1334 
1335  virtual void v_GetLocTraceFromTracePts(
1336  const Array<OneD, const NekDouble> &Fwd,
1337  const Array<OneD, const NekDouble> &Bwd,
1338  Array<OneD, NekDouble> &locTraceFwd,
1339  Array<OneD, NekDouble> &locTraceBwd);
1340 
1341  virtual void v_FillBwdWithBwdWeight(Array<OneD, NekDouble> &weightave,
1342  Array<OneD, NekDouble> &weightjmp);
1343 
1344  virtual void v_PeriodicBwdCopy(const Array<OneD, const NekDouble> &Fwd,
1345  Array<OneD, NekDouble> &Bwd);
1346 
1347  virtual const std::vector<bool> &v_GetLeftAdjacentFaces(void) const;
1348 
1349  virtual void v_ExtractTracePhys(Array<OneD, NekDouble> &outarray);
1350 
1351  virtual void v_ExtractTracePhys(const Array<OneD, const NekDouble> &inarray,
1352  Array<OneD, NekDouble> &outarray);
1353 
1354  virtual void v_MultiplyByInvMassMatrix(
1355  const Array<OneD, const NekDouble> &inarray,
1356  Array<OneD, NekDouble> &outarray);
1357 
1358  virtual void v_HelmSolve(const Array<OneD, const NekDouble> &inarray,
1359  Array<OneD, NekDouble> &outarray,
1360  const StdRegions::ConstFactorMap &factors,
1361  const StdRegions::VarCoeffMap &varcoeff,
1362  const MultiRegions::VarFactorsMap &varfactors,
1363  const Array<OneD, const NekDouble> &dirForcing,
1364  const bool PhysSpaceForcing);
1365 
1367  const Array<OneD, Array<OneD, NekDouble>> &velocity,
1368  const Array<OneD, const NekDouble> &inarray,
1369  Array<OneD, NekDouble> &outarray, const NekDouble lambda,
1371 
1372  virtual void v_LinearAdvectionReactionSolve(
1373  const Array<OneD, Array<OneD, NekDouble>> &velocity,
1374  const Array<OneD, const NekDouble> &inarray,
1375  Array<OneD, NekDouble> &outarray, const NekDouble lambda,
1377 
1378  // wrapper functions about virtual functions
1379  virtual void v_ImposeDirichletConditions(Array<OneD, NekDouble> &outarray);
1380 
1381  virtual void v_FillBndCondFromField();
1382 
1383  virtual void v_FillBndCondFromField(const int nreg);
1384 
1385  virtual void v_Reset();
1386 
1387  virtual void v_LocalToGlobal(bool UseComm);
1388 
1389  virtual void v_LocalToGlobal(const Array<OneD, const NekDouble> &inarray,
1390  Array<OneD, NekDouble> &outarray,
1391  bool UseComm);
1392 
1393  virtual void v_GlobalToLocal(void);
1394 
1395  virtual void v_GlobalToLocal(const Array<OneD, const NekDouble> &inarray,
1396  Array<OneD, NekDouble> &outarray);
1397 
1398  virtual void v_BwdTrans(const Array<OneD, const NekDouble> &inarray,
1399  Array<OneD, NekDouble> &outarray);
1400 
1401  virtual void v_FwdTrans(const Array<OneD, const NekDouble> &inarray,
1402  Array<OneD, NekDouble> &outarray);
1403 
1404  virtual void v_FwdTransLocalElmt(
1405  const Array<OneD, const NekDouble> &inarray,
1406  Array<OneD, NekDouble> &outarray);
1407 
1408  virtual void v_FwdTransBndConstrained(
1409  const Array<OneD, const NekDouble> &inarray,
1410  Array<OneD, NekDouble> &outarray);
1411 
1412  virtual void v_SmoothField(Array<OneD, NekDouble> &field);
1413 
1414  virtual void v_IProductWRTBase(const Array<OneD, const NekDouble> &inarray,
1415  Array<OneD, NekDouble> &outarray);
1416 
1417  virtual void v_GetCoords(
1420 
1421  virtual void v_GetCoords(const int eid, Array<OneD, NekDouble> &xc0,
1423  Array<OneD, NekDouble> &xc2);
1424 
1425  virtual void v_PhysDeriv(const Array<OneD, const NekDouble> &inarray,
1426  Array<OneD, NekDouble> &out_d0,
1427  Array<OneD, NekDouble> &out_d1,
1428  Array<OneD, NekDouble> &out_d2);
1429 
1430  virtual void v_PhysDeriv(const int dir,
1431  const Array<OneD, const NekDouble> &inarray,
1432  Array<OneD, NekDouble> &out_d);
1433 
1434  virtual void v_PhysDeriv(Direction edir,
1435  const Array<OneD, const NekDouble> &inarray,
1436  Array<OneD, NekDouble> &out_d);
1437 
1438  virtual void v_CurlCurl(Array<OneD, Array<OneD, NekDouble>> &Vel,
1440 
1441  virtual void v_PhysDirectionalDeriv(
1442  const Array<OneD, const NekDouble> &direction,
1443  const Array<OneD, const NekDouble> &inarray,
1444  Array<OneD, NekDouble> &outarray);
1445 
1446  virtual void v_GetMovingFrames(
1447  const SpatialDomains::GeomMMF MMFdir,
1448  const Array<OneD, const NekDouble> &CircCentre,
1449  Array<OneD, Array<OneD, NekDouble>> &outarray);
1450 
1451  virtual void v_HomogeneousFwdTrans(
1452  const Array<OneD, const NekDouble> &inarray,
1453  Array<OneD, NekDouble> &outarray, bool Shuff = true,
1454  bool UnShuff = true);
1455 
1456  virtual void v_HomogeneousBwdTrans(
1457  const Array<OneD, const NekDouble> &inarray,
1458  Array<OneD, NekDouble> &outarray, bool Shuff = true,
1459  bool UnShuff = true);
1460 
1461  virtual void v_DealiasedProd(const Array<OneD, NekDouble> &inarray1,
1462  const Array<OneD, NekDouble> &inarray2,
1463  Array<OneD, NekDouble> &outarray);
1464 
1465  virtual void v_DealiasedDotProd(
1466  const Array<OneD, Array<OneD, NekDouble>> &inarray1,
1467  const Array<OneD, Array<OneD, NekDouble>> &inarray2,
1468  Array<OneD, Array<OneD, NekDouble>> &outarray);
1469 
1470  virtual void v_GetBCValues(Array<OneD, NekDouble> &BndVals,
1471  const Array<OneD, NekDouble> &TotField,
1472  int BndID);
1473 
1476  Array<OneD, NekDouble> &outarray,
1477  int BndID);
1478 
1479  virtual void v_NormVectorIProductWRTBase(
1481  Array<OneD, NekDouble> &outarray);
1482 
1483  virtual void v_SetUpPhysNormals();
1484 
1485  virtual void v_GetBoundaryToElmtMap(Array<OneD, int> &ElmtID,
1486  Array<OneD, int> &EdgeID);
1487 
1488  virtual void v_GetBndElmtExpansion(int i, std::shared_ptr<ExpList> &result,
1489  const bool DeclareCoeffPhysArrays);
1490 
1491  virtual void v_ExtractElmtToBndPhys(const int i,
1492  const Array<OneD, NekDouble> &elmt,
1493  Array<OneD, NekDouble> &boundary);
1494 
1495  virtual void v_ExtractPhysToBndElmt(
1496  const int i, const Array<OneD, const NekDouble> &phys,
1497  Array<OneD, NekDouble> &bndElmt);
1498 
1499  virtual void v_ExtractPhysToBnd(const int i,
1500  const Array<OneD, const NekDouble> &phys,
1501  Array<OneD, NekDouble> &bnd);
1502 
1503  virtual void v_GetBoundaryNormals(
1504  int i, Array<OneD, Array<OneD, NekDouble>> &normals);
1505 
1506  virtual std::vector<LibUtilities::FieldDefinitionsSharedPtr>
1507  v_GetFieldDefinitions(void);
1508 
1509  virtual void v_GetFieldDefinitions(
1510  std::vector<LibUtilities::FieldDefinitionsSharedPtr> &fielddef);
1511 
1512  virtual void v_AppendFieldData(
1514  std::vector<NekDouble> &fielddata);
1515 
1516  virtual void v_AppendFieldData(
1518  std::vector<NekDouble> &fielddata, Array<OneD, NekDouble> &coeffs);
1519 
1520  virtual void v_ExtractDataToCoeffs(
1522  std::vector<NekDouble> &fielddata, std::string &field,
1523  Array<OneD, NekDouble> &coeffs);
1524 
1525  virtual void v_ExtractCoeffsToCoeffs(
1526  const std::shared_ptr<ExpList> &fromExpList,
1527  const Array<OneD, const NekDouble> &fromCoeffs,
1528  Array<OneD, NekDouble> &toCoeffs);
1529 
1530  virtual void v_WriteTecplotHeader(std::ostream &outfile,
1531  std::string var = "");
1532  virtual void v_WriteTecplotZone(std::ostream &outfile, int expansion);
1533  virtual void v_WriteTecplotField(std::ostream &outfile, int expansion);
1534  virtual void v_WriteTecplotConnectivity(std::ostream &outfile,
1535  int expansion);
1536  virtual void v_WriteVtkPieceHeader(std::ostream &outfile, int expansion,
1537  int istrip);
1538 
1539  virtual void v_WriteVtkPieceData(std::ostream &outfile, int expansion,
1540  std::string var);
1541 
1542  virtual NekDouble v_L2(
1543  const Array<OneD, const NekDouble> &phys,
1545 
1546  virtual NekDouble v_Integral(const Array<OneD, const NekDouble> &inarray);
1547  virtual NekDouble v_VectorFlux(
1548  const Array<OneD, Array<OneD, NekDouble>> &inarray);
1549 
1552  virtual NekDouble v_GetHomoLen(void);
1553  virtual void v_SetHomoLen(const NekDouble lhom);
1556 
1557  // 1D Scaling and projection
1558  virtual void v_PhysInterp1DScaled(const NekDouble scale,
1559  const Array<OneD, NekDouble> &inarray,
1560  Array<OneD, NekDouble> &outarray);
1561 
1562  virtual void v_PhysGalerkinProjection1DScaled(
1563  const NekDouble scale, const Array<OneD, NekDouble> &inarray,
1564  Array<OneD, NekDouble> &outarray);
1565 
1566  virtual void v_ClearGlobalLinSysManager(void);
1567 
1568  void ExtractFileBCs(const std::string &fileName,
1570  const std::string &varName,
1571  const std::shared_ptr<ExpList> locExpList);
1572 
1573  // Utility function for a common case of retrieving a
1574  // BoundaryCondition from a boundary condition collection.
1578  unsigned int index, const std::string &variable);
1579 
1580 private:
1581  /// Definition of the total number of degrees of freedom and
1582  /// quadrature points and offsets to access data
1583  void SetupCoeffPhys(bool DeclareCoeffPhysArrays = true,
1584  bool SetupOffsets = true);
1585 
1586  /// Define a list of elements using the geometry and basis
1587  /// key information in expmap;
1589 
1591  &v_GetBndConditions();
1592 
1595 
1596  virtual void v_EvaluateBoundaryConditions(
1597  const NekDouble time = 0.0, const std::string varName = "",
1600 
1601  virtual std::map<int, RobinBCInfoSharedPtr> v_GetRobinBCInfo(void);
1602 
1603  virtual void v_GetPeriodicEntities(PeriodicMap &periodicVerts,
1604  PeriodicMap &periodicEdges,
1605  PeriodicMap &periodicFaces);
1606 
1607  // Homogeneous direction wrapper functions.
1609  {
1611  "This method is not defined or valid for this class type");
1613  }
1614 
1615  // wrapper function to set viscosity for Homo1D expansion
1617  {
1618  boost::ignore_unused(visc);
1620  "This method is not defined or valid for this class type");
1621  }
1622 
1623  virtual std::shared_ptr<ExpList> &v_GetPlane(int n);
1624 
1625  virtual void v_AddTraceIntegralToOffDiag(
1626  const Array<OneD, const NekDouble> &FwdFlux,
1627  const Array<OneD, const NekDouble> &BwdFlux,
1628  Array<OneD, NekDouble> &outarray);
1629 };
1630 
1631 /// An empty ExpList object.
1634 
1635 // Inline routines follow.
1636 
1637 /**
1638  * Returns the total number of local degrees of freedom
1639  * \f$N_{\mathrm{eof}}=\sum_{e=1}^{{N_{\mathrm{el}}}}N^{e}_m\f$.
1640  */
1641 inline int ExpList::GetNcoeffs() const
1642 {
1643  return m_ncoeffs;
1644 }
1645 
1646 inline int ExpList::GetNcoeffs(const int eid) const
1647 {
1648  return (*m_exp)[eid]->GetNcoeffs();
1649 }
1650 
1651 /**
1652  * Evaulates the maximum number of modes in the elemental basis
1653  * order over all elements
1654  */
1656 {
1657  unsigned int i;
1658  int returnval = 0;
1659 
1660  for (i = 0; i < (*m_exp).size(); ++i)
1661  {
1662  returnval = (std::max)(returnval, (*m_exp)[i]->EvalBasisNumModesMax());
1663  }
1664 
1665  return returnval;
1666 }
1667 
1668 /**
1669  *
1670  */
1672 {
1673  unsigned int i;
1674  Array<OneD, int> returnval((*m_exp).size(), 0);
1675 
1676  for (i = 0; i < (*m_exp).size(); ++i)
1677  {
1678  returnval[i] =
1679  (std::max)(returnval[i], (*m_exp)[i]->EvalBasisNumModesMax());
1680  }
1681 
1682  return returnval;
1683 }
1684 
1685 /**
1686  *
1687  */
1688 inline int ExpList::GetTotPoints() const
1689 {
1690  return m_npoints;
1691 }
1692 
1693 inline int ExpList::GetTotPoints(const int eid) const
1694 {
1695  return (*m_exp)[eid]->GetTotPoints();
1696 }
1697 
1698 inline int ExpList::Get1DScaledTotPoints(const NekDouble scale) const
1699 {
1700  int returnval = 0;
1701  int cnt;
1702  int nbase = (*m_exp)[0]->GetNumBases();
1703 
1704  for (int i = 0; i < (*m_exp).size(); ++i)
1705  {
1706  cnt = 1;
1707  for (int j = 0; j < nbase; ++j)
1708  {
1709  cnt *= (int)(scale * ((*m_exp)[i]->GetNumPoints(j)));
1710  }
1711  returnval += cnt;
1712  }
1713  return returnval;
1714 }
1715 
1716 /**
1717  *
1718  */
1719 inline int ExpList::GetNpoints() const
1720 {
1721  return m_npoints;
1722 }
1723 
1724 /**
1725  *
1726  */
1727 inline void ExpList::SetWaveSpace(const bool wavespace)
1728 {
1729  m_WaveSpace = wavespace;
1730 }
1731 
1732 /**
1733  *
1734  */
1735 inline bool ExpList::GetWaveSpace() const
1736 {
1737  return m_WaveSpace;
1738 }
1739 
1740 /// Set the \a i th value of\a m_phys to value \a val
1741 inline void ExpList::SetPhys(int i, NekDouble val)
1742 {
1743  m_phys[i] = val;
1744 }
1745 
1746 /**
1747  * This function fills the array \f$\boldsymbol{u}_l\f$, the evaluation
1748  * of the expansion at the quadrature points (implemented as #m_phys),
1749  * with the values of the array \a inarray.
1750  *
1751  * @param inarray The array containing the values where
1752  * #m_phys should be filled with.
1753  */
1755 {
1756  ASSERTL0(inarray.size() == m_npoints,
1757  "Input array does not have correct number of elements.");
1758 
1759  Vmath::Vcopy(m_npoints, &inarray[0], 1, &m_phys[0], 1);
1760  m_physState = true;
1761 }
1762 
1764 {
1765  m_phys = inarray;
1766 }
1767 
1768 /**
1769  * @param physState \a true (=filled) or \a false (=not filled).
1770  */
1771 inline void ExpList::SetPhysState(const bool physState)
1772 {
1773  m_physState = physState;
1774 }
1775 
1776 /**
1777  * @return physState \a true (=filled) or \a false (=not filled).
1778  */
1779 inline bool ExpList::GetPhysState() const
1780 {
1781  return m_physState;
1782 }
1783 
1784 /**
1785  *
1786  */
1788  const Array<OneD, const NekDouble> &inarray,
1789  Array<OneD, NekDouble> &outarray)
1790 {
1791  v_IProductWRTBase(inarray, outarray);
1792 }
1793 
1794 /**
1795  *
1796  */
1798  Array<OneD, NekDouble> &outarray)
1799 {
1800  v_FwdTrans(inarray, outarray);
1801 }
1802 
1803 /**
1804  *
1805  */
1807  const Array<OneD, const NekDouble> &inarray,
1808  Array<OneD, NekDouble> &outarray)
1809 {
1810  v_FwdTransLocalElmt(inarray, outarray);
1811 }
1812 
1813 /**
1814  *
1815  */
1817  const Array<OneD, const NekDouble> &inarray,
1818  Array<OneD, NekDouble> &outarray)
1819 {
1820  v_FwdTransBndConstrained(inarray, outarray);
1821 }
1822 
1823 /**
1824  *
1825  */
1827 {
1828  v_SmoothField(field);
1829 }
1830 
1831 /**
1832  *
1833  */
1834 
1835 /**
1836  * Given the coefficients of an expansion, this function evaluates the
1837  * spectral/hp expansion \f$u^{\delta}(\boldsymbol{x})\f$ at the
1838  * quadrature points \f$\boldsymbol{x}_i\f$. This operation is
1839  * evaluated locally by the function ExpList#BwdTrans.
1840  *
1841  * The coefficients of the expansion should be contained in the variable
1842  * #m_coeffs of the ExpList object \a In. The resulting physical values
1843  * at the quadrature points \f$u^{\delta}(\boldsymbol{x}_i)\f$ are
1844  * stored in the array #m_phys.
1845  *
1846  * @param In An ExpList, containing the local coefficients
1847  * \f$\hat{u}_n^e\f$ in its array #m_coeffs.
1848  */
1850  Array<OneD, NekDouble> &outarray)
1851 {
1852  v_BwdTrans(inarray, outarray);
1853 }
1854 
1855 /**
1856  *
1857  */
1859  const Array<OneD, const NekDouble> &inarray,
1860  Array<OneD, NekDouble> &outarray)
1861 {
1862  v_MultiplyByInvMassMatrix(inarray, outarray);
1863 }
1864 
1865 /**
1866  *
1867  */
1869  Array<OneD, NekDouble> &outarray,
1870  const StdRegions::ConstFactorMap &factors,
1871  const StdRegions::VarCoeffMap &varcoeff,
1872  const MultiRegions::VarFactorsMap &varfactors,
1873  const Array<OneD, const NekDouble> &dirForcing,
1874  const bool PhysSpaceForcing)
1875 
1876 {
1877  v_HelmSolve(inarray, outarray, factors, varcoeff, varfactors, dirForcing,
1878  PhysSpaceForcing);
1879 }
1880 
1881 /**
1882  *
1883  */
1885  const Array<OneD, Array<OneD, NekDouble>> &velocity,
1886  const Array<OneD, const NekDouble> &inarray,
1887  Array<OneD, NekDouble> &outarray, const NekDouble lambda,
1888  const Array<OneD, const NekDouble> &dirForcing)
1889 {
1890  v_LinearAdvectionDiffusionReactionSolve(velocity, inarray, outarray, lambda,
1891  dirForcing);
1892 }
1893 
1895  const Array<OneD, Array<OneD, NekDouble>> &velocity,
1896  const Array<OneD, const NekDouble> &inarray,
1897  Array<OneD, NekDouble> &outarray, const NekDouble lambda,
1898  const Array<OneD, const NekDouble> &dirForcing)
1899 {
1900  v_LinearAdvectionReactionSolve(velocity, inarray, outarray, lambda,
1901  dirForcing);
1902 }
1903 
1904 /**
1905  *
1906  */
1908  Array<OneD, NekDouble> &coord_1,
1909  Array<OneD, NekDouble> &coord_2)
1910 
1911 {
1912  v_GetCoords(coord_0, coord_1, coord_2);
1913 }
1914 
1915 inline void ExpList::GetCoords(const int eid, Array<OneD, NekDouble> &xc0,
1918 {
1919  v_GetCoords(eid, xc0, xc1, xc2);
1920 }
1921 
1922 /**
1923  *
1924  */
1926  const SpatialDomains::GeomMMF MMFdir,
1927  const Array<OneD, const NekDouble> &CircCentre,
1928  Array<OneD, Array<OneD, NekDouble>> &outarray)
1929 {
1930  v_GetMovingFrames(MMFdir, CircCentre, outarray);
1931 }
1932 
1933 /**
1934  *
1935  */
1937  Array<OneD, NekDouble> &out_d0,
1938  Array<OneD, NekDouble> &out_d1,
1939  Array<OneD, NekDouble> &out_d2)
1940 {
1941  v_PhysDeriv(inarray, out_d0, out_d1, out_d2);
1942 }
1943 
1944 /**
1945  *
1946  */
1947 inline void ExpList::PhysDeriv(const int dir,
1948  const Array<OneD, const NekDouble> &inarray,
1949  Array<OneD, NekDouble> &out_d)
1950 {
1951  v_PhysDeriv(dir, inarray, out_d);
1952 }
1953 
1955  const Array<OneD, const NekDouble> &inarray,
1956  Array<OneD, NekDouble> &out_d)
1957 {
1958  v_PhysDeriv(edir, inarray, out_d);
1959 }
1960 
1961 /**
1962  *
1963  */
1965  const Array<OneD, const NekDouble> &direction,
1966  const Array<OneD, const NekDouble> &inarray,
1967  Array<OneD, NekDouble> &outarray)
1968 {
1969  v_PhysDirectionalDeriv(direction, inarray, outarray);
1970 }
1971 
1972 /**
1973  *
1974  */
1975 
1978 {
1979  v_CurlCurl(Vel, Q);
1980 }
1981 
1982 /**
1983  *
1984  */
1986  const Array<OneD, const NekDouble> &inarray,
1987  Array<OneD, NekDouble> &outarray, bool Shuff, bool UnShuff)
1988 {
1989  v_HomogeneousFwdTrans(inarray, outarray, Shuff, UnShuff);
1990 }
1991 
1992 /**
1993  *
1994  */
1996  const Array<OneD, const NekDouble> &inarray,
1997  Array<OneD, NekDouble> &outarray, bool Shuff, bool UnShuff)
1998 {
1999  v_HomogeneousBwdTrans(inarray, outarray, Shuff, UnShuff);
2000 }
2001 
2002 /**
2003  *
2004  */
2006  const Array<OneD, NekDouble> &inarray2,
2007  Array<OneD, NekDouble> &outarray)
2008 {
2009  v_DealiasedProd(inarray1, inarray2, outarray);
2010 }
2011 
2012 /**
2013  *
2014  */
2016  const Array<OneD, Array<OneD, NekDouble>> &inarray1,
2017  const Array<OneD, Array<OneD, NekDouble>> &inarray2,
2018  Array<OneD, Array<OneD, NekDouble>> &outarray)
2019 {
2020  v_DealiasedDotProd(inarray1, inarray2, outarray);
2021 }
2022 
2023 /**
2024  *
2025  */
2027  const Array<OneD, NekDouble> &TotField,
2028  int BndID)
2029 {
2030  v_GetBCValues(BndVals, TotField, BndID);
2031 }
2032 
2033 /**
2034  *
2035  */
2038  Array<OneD, NekDouble> &outarray,
2039  int BndID)
2040 {
2041  v_NormVectorIProductWRTBase(V1, V2, outarray, BndID);
2042 }
2043 
2046 {
2047  v_NormVectorIProductWRTBase(V, outarray);
2048 }
2049 
2050 /**
2051  * @param eid The index of the element to be checked.
2052  * @return The dimension of the coordinates of the specific element.
2053  */
2054 inline int ExpList::GetCoordim(int eid)
2055 {
2056  ASSERTL2(eid <= (*m_exp).size(), "eid is larger than number of elements");
2057  return (*m_exp)[eid]->GetCoordim();
2058 }
2059 
2060 /**
2061  * @param eid The index of the element to be checked.
2062  * @return The dimension of the shape of the specific element.
2063  */
2065 {
2066  return (*m_exp)[0]->GetShapeDimension();
2067 }
2068 
2069 /**
2070  * @param i The index of m_coeffs to be set
2071  * @param val The value which m_coeffs[i] is to be set to.
2072  */
2073 inline void ExpList::SetCoeff(int i, NekDouble val)
2074 {
2075  m_coeffs[i] = val;
2076 }
2077 
2078 /**
2079  * @param i The index of #m_coeffs to be set.
2080  * @param val The value which #m_coeffs[i] is to be set to.
2081  */
2082 inline void ExpList::SetCoeffs(int i, NekDouble val)
2083 {
2084  m_coeffs[i] = val;
2085 }
2086 
2088 {
2089  m_coeffs = inarray;
2090 }
2091 
2092 /**
2093  * As the function returns a constant reference to a
2094  * <em>const Array</em>, it is not possible to modify the
2095  * underlying data of the array #m_coeffs. In order to do
2096  * so, use the function #UpdateCoeffs instead.
2097  *
2098  * @return (A constant reference to) the array #m_coeffs.
2099  */
2101 {
2102  return m_coeffs;
2103 }
2104 
2106 {
2107  v_ImposeDirichletConditions(outarray);
2108 }
2109 
2111 {
2113 }
2114 
2115 inline void ExpList::FillBndCondFromField(const int nreg)
2116 {
2117  v_FillBndCondFromField(nreg);
2118 }
2119 
2120 inline void ExpList::LocalToGlobal(bool useComm)
2121 {
2122  v_LocalToGlobal(useComm);
2123 }
2124 
2126  Array<OneD, NekDouble> &outarray,
2127  bool useComm)
2128 {
2129  v_LocalToGlobal(inarray, outarray, useComm);
2130 }
2131 
2132 inline void ExpList::GlobalToLocal(void)
2133 {
2134  v_GlobalToLocal();
2135 }
2136 
2137 /**
2138  * This operation is evaluated as:
2139  * \f{tabbing}
2140  * \hspace{1cm} \= Do \= $e=$ $1, N_{\mathrm{el}}$ \\
2141  * > > Do \= $i=$ $0,N_m^e-1$ \\
2142  * > > > $\boldsymbol{\hat{u}}^{e}[i] = \mbox{sign}[e][i] \cdot
2143  * \boldsymbol{\hat{u}}_g[\mbox{map}[e][i]]$ \\
2144  * > > continue \\
2145  * > continue
2146  * \f}
2147  * where \a map\f$[e][i]\f$ is the mapping array and \a
2148  * sign\f$[e][i]\f$ is an array of similar dimensions ensuring the
2149  * correct modal connectivity between the different elements (both
2150  * these arrays are contained in the data member #m_locToGloMap). This
2151  * operation is equivalent to the scatter operation
2152  * \f$\boldsymbol{\hat{u}}_l=\mathcal{A}\boldsymbol{\hat{u}}_g\f$, where
2153  * \f$\mathcal{A}\f$ is the
2154  * \f$N_{\mathrm{eof}}\times N_{\mathrm{dof}}\f$ permutation matrix.
2155  *
2156  * @param inarray An array of size \f$N_\mathrm{dof}\f$
2157  * containing the global degrees of freedom
2158  * \f$\boldsymbol{x}_g\f$.
2159  * @param outarray The resulting local degrees of freedom
2160  * \f$\boldsymbol{x}_l\f$ will be stored in this
2161  * array of size \f$N_\mathrm{eof}\f$.
2162  */
2164  Array<OneD, NekDouble> &outarray)
2165 {
2166  v_GlobalToLocal(inarray, outarray);
2167 }
2168 
2169 /**
2170  * @param i The index of #m_coeffs to be returned
2171  * @return The NekDouble held in #m_coeffs[i].
2172  */
2174 {
2175  return m_coeffs[i];
2176 }
2177 
2178 /**
2179  * @param i The index of #m_coeffs to be returned
2180  * @return The NekDouble held in #m_coeffs[i].
2181  */
2183 {
2184  return m_coeffs[i];
2185 }
2186 
2187 /**
2188  * As the function returns a constant reference to a
2189  * <em>const Array</em> it is not possible to modify the
2190  * underlying data of the array #m_phys. In order to do
2191  * so, use the function #UpdatePhys instead.
2192  *
2193  * @return (A constant reference to) the array #m_phys.
2194  */
2196 {
2197  return m_phys;
2198 }
2199 
2200 /**
2201  * @return \f$N_{\mathrm{el}}\f$, the number of elements in the
2202  * expansion.
2203  */
2204 inline int ExpList::GetExpSize(void)
2205 {
2206  return (*m_exp).size();
2207 }
2208 
2209 /**
2210  * @param n The index of the element concerned.
2211  *
2212  * @return (A shared pointer to) the local expansion of the
2213  * \f$n^{\mathrm{th}}\f$ element.
2214  */
2216 {
2217  return (*m_exp)[n];
2218 }
2219 
2220 /**
2221  * @return (A const shared pointer to) the local expansion vector #m_exp
2222  */
2223 inline const std::shared_ptr<LocalRegions::ExpansionVector> ExpList::GetExp(
2224  void) const
2225 {
2226  return m_exp;
2227 }
2228 
2229 /**
2230  *
2231  */
2232 inline int ExpList::GetCoeff_Offset(int n) const
2233 {
2234  return m_coeff_offset[n];
2235 }
2236 
2237 /**
2238  *
2239  */
2240 inline int ExpList::GetPhys_Offset(int n) const
2241 {
2242  return m_phys_offset[n];
2243 }
2244 
2245 /**
2246  * If one wants to get hold of the underlying data without modifying
2247  * them, rather use the function #GetCoeffs instead.
2248  *
2249  * @return (A reference to) the array #m_coeffs.
2250  */
2252 {
2253  return m_coeffs;
2254 }
2255 
2256 /**
2257  * If one wants to get hold of the underlying data without modifying
2258  * them, rather use the function #GetPhys instead.
2259  *
2260  * @return (A reference to) the array #m_phys.
2261  */
2263 {
2264  m_physState = true;
2265  return m_phys;
2266 }
2267 
2268 // functions associated with DisContField
2271 {
2272  return v_GetBndCondExpansions();
2273 }
2274 
2275 /// Get m_coeffs to elemental value map
2278 {
2279  return m_coeffsToElmt;
2280 }
2281 
2283  GetLocTraceToTraceMap() const
2284 {
2285  return v_GetLocTraceToTraceMap();
2286 }
2287 
2289 {
2290  return v_GetBndCondBwdWeight();
2291 }
2292 
2293 inline void ExpList::SetBndCondBwdWeight(const int index, const NekDouble value)
2294 {
2295  v_SetBndCondBwdWeight(index, value);
2296 }
2297 
2298 inline std::shared_ptr<ExpList> &ExpList::UpdateBndCondExpansion(int i)
2299 {
2300  return v_UpdateBndCondExpansion(i);
2301 }
2302 
2303 inline void ExpList::Upwind(
2304  const Array<OneD, const Array<OneD, NekDouble>> &Vec,
2305  const Array<OneD, const NekDouble> &Fwd,
2307 {
2308  v_Upwind(Vec, Fwd, Bwd, Upwind);
2309 }
2310 
2312  const Array<OneD, const NekDouble> &Fwd,
2313  const Array<OneD, const NekDouble> &Bwd,
2314  Array<OneD, NekDouble> &Upwind)
2315 {
2316  v_Upwind(Vn, Fwd, Bwd, Upwind);
2317 }
2318 
2319 inline std::shared_ptr<ExpList> &ExpList::GetTrace()
2320 {
2321  return v_GetTrace();
2322 }
2323 
2324 inline std::shared_ptr<AssemblyMapDG> &ExpList::GetTraceMap()
2325 {
2326  return v_GetTraceMap();
2327 }
2328 
2330 {
2331  return v_GetTraceBndMap();
2332 }
2333 
2335 {
2336  v_GetNormals(normals);
2337 }
2338 
2340  const Array<OneD, const NekDouble> &Fy,
2341  Array<OneD, NekDouble> &outarray)
2342 {
2343  v_AddTraceIntegral(Fx, Fy, outarray);
2344 }
2345 
2347  Array<OneD, NekDouble> &outarray)
2348 {
2349  v_AddTraceIntegral(Fn, outarray);
2350 }
2351 
2353  const Array<OneD, const NekDouble> &Fwd,
2355 {
2356  v_AddFwdBwdTraceIntegral(Fwd, Bwd, outarray);
2357 }
2358 
2361 {
2362  v_GetFwdBwdTracePhys(Fwd, Bwd);
2363 }
2364 
2367  Array<OneD, NekDouble> &Bwd, bool FillBnd, bool PutFwdInBwdOnBCs,
2368  bool DoExchange)
2369 {
2370  v_GetFwdBwdTracePhys(field, Fwd, Bwd, FillBnd, PutFwdInBwdOnBCs,
2371  DoExchange);
2372 }
2373 
2376  bool PutFwdInBwdOnBCs)
2377 {
2378  v_FillBwdWithBoundCond(Fwd, Bwd, PutFwdInBwdOnBCs);
2379 }
2380 
2382  const Array<OneD, const NekDouble> &Fwd,
2384 {
2385  v_AddTraceQuadPhysToField(Fwd, Bwd, field);
2386 }
2387 
2389  const Array<OneD, const NekDouble> &Fwd,
2390  const Array<OneD, const NekDouble> &Bwd,
2391  Array<OneD, NekDouble> &locTraceFwd, Array<OneD, NekDouble> &locTraceBwd)
2392 {
2393  v_GetLocTraceFromTracePts(Fwd, Bwd, locTraceFwd, locTraceBwd);
2394 }
2395 
2397  Array<OneD, NekDouble> &weightjmp)
2398 {
2399  v_FillBwdWithBwdWeight(weightave, weightjmp);
2400 }
2401 
2404 {
2405  v_PeriodicBwdCopy(Fwd, Bwd);
2406 }
2407 
2408 inline const std::vector<bool> &ExpList::GetLeftAdjacentFaces(void) const
2409 {
2410  return v_GetLeftAdjacentFaces();
2411 }
2412 
2414 {
2415  v_ExtractTracePhys(outarray);
2416 }
2417 
2419  const Array<OneD, const NekDouble> &inarray,
2420  Array<OneD, NekDouble> &outarray)
2421 {
2422  v_ExtractTracePhys(inarray, outarray);
2423 }
2424 
2427 {
2428  return v_GetBndConditions();
2429 }
2430 
2433 {
2434  return v_UpdateBndConditions();
2435 }
2436 
2438  const std::string varName,
2439  const NekDouble x2_in,
2440  const NekDouble x3_in)
2441 {
2442  v_EvaluateBoundaryConditions(time, varName, x2_in, x3_in);
2443 }
2444 
2446 {
2448 }
2449 
2451  Array<OneD, int> &EdgeID)
2452 {
2453  v_GetBoundaryToElmtMap(ElmtID, EdgeID);
2454 }
2455 
2457  std::shared_ptr<ExpList> &result,
2458  const bool DeclareCoeffPhysArrays)
2459 {
2460  v_GetBndElmtExpansion(i, result, DeclareCoeffPhysArrays);
2461 }
2462 
2464  const Array<OneD, NekDouble> &elmt,
2465  Array<OneD, NekDouble> &boundary)
2466 {
2467  v_ExtractElmtToBndPhys(i, elmt, boundary);
2468 }
2469 
2471  int i, const Array<OneD, const NekDouble> &phys,
2472  Array<OneD, NekDouble> &bndElmt)
2473 {
2474  v_ExtractPhysToBndElmt(i, phys, bndElmt);
2475 }
2476 
2477 inline void ExpList::ExtractPhysToBnd(int i,
2478  const Array<OneD, const NekDouble> &phys,
2480 {
2481  v_ExtractPhysToBnd(i, phys, bnd);
2482 }
2483 
2485  int i, Array<OneD, Array<OneD, NekDouble>> &normals)
2486 {
2487  v_GetBoundaryNormals(i, normals);
2488 }
2489 
2490 inline std::vector<bool> &ExpList::GetLeftAdjacentTraces(void)
2491 {
2492  return v_GetLeftAdjacentTraces();
2493 }
2494 
2496 } // namespace MultiRegions
2497 } // namespace Nektar
2498 
2499 #endif // EXPLIST_H
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:209
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
Definition: ErrorUtil.hpp:272
#define MULTI_REGIONS_EXPORT
Base class for all multi-elemental spectral/hp expansions.
Definition: ExpList.h:101
virtual void v_Reset()
Reset geometry information, metrics, matrix managers and geometry information.
Definition: ExpList.cpp:2670
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Definition: ExpList.h:1158
Array< OneD, NekDouble > & UpdateCoeffs()
This function returns (a reference to) the array (implemented as m_coeffs) containing all local expa...
Definition: ExpList.h:2251
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:5567
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:1787
const Array< OneD, const std::shared_ptr< ExpList > > & GetBndCondExpansions()
Definition: ExpList.h:2270
std::shared_ptr< ExpList > & GetTrace()
Definition: ExpList.h:2319
void WriteTecplotConnectivity(std::ostream &outfile, int expansion=-1)
Definition: ExpList.h:427
int GetNcoeffs(void) const
Returns the total number of local degrees of freedom .
Definition: ExpList.h:1641
virtual void v_ExtractCoeffsToCoeffs(const std::shared_ptr< ExpList > &fromExpList, const Array< OneD, const NekDouble > &fromCoeffs, Array< OneD, NekDouble > &toCoeffs)
Definition: ExpList.cpp:3649
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:3792
const DNekScalBlkMatSharedPtr & GetBlockMatrix(const GlobalMatrixKey &gkey)
Definition: ExpList.cpp:2010
virtual void v_AddFwdBwdTraceIntegral(const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:4371
ExpList(const ExpansionType Type=eNoType)
The default constructor using a type.
Definition: ExpList.cpp:104
void SetCoeffs(int i, NekDouble val)
Set the i th coefficiient in m_coeffs to value val.
Definition: ExpList.h:2082
void ExtractDataToCoeffs(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, std::string &field, Array< OneD, NekDouble > &coeffs)
Extract the data in fielddata into the coeffs.
Definition: ExpList.cpp:3538
void LinearAdvectionDiffusionReactionSolve(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:1884
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:5071
void CurlCurl(Array< OneD, Array< OneD, NekDouble >> &Vel, Array< OneD, Array< OneD, NekDouble >> &Q)
Definition: ExpList.h:1976
virtual void v_GetBoundaryToElmtMap(Array< OneD, int > &ElmtID, Array< OneD, int > &EdgeID)
Definition: ExpList.cpp:4961
void ExponentialFilter(Array< OneD, NekDouble > &array, const NekDouble alpha, const NekDouble exponent, const NekDouble cutoff)
Definition: ExpList.cpp:1743
virtual void v_HomogeneousBwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true)
Definition: ExpList.cpp:4527
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:2126
virtual void v_HomogeneousFwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true)
Definition: ExpList.cpp:4518
virtual void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:4686
std::shared_ptr< LibUtilities::Comm > GetComm() const
Returns the comm object.
Definition: ExpList.h:1046
void MultiplyByMassMatrix(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:314
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:2374
void AddRightIPTPhysDerivBase(const int dir, const Array< OneD, const DNekMatSharedPtr > ElmtJacQuad, Array< OneD, DNekMatSharedPtr > ElmtJacCoef)
Definition: ExpList.cpp:5803
virtual void v_ExtractTracePhys(Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:4459
void ExtractPhysToBndElmt(int i, const Array< OneD, const NekDouble > &phys, Array< OneD, NekDouble > &bndElmt)
Definition: ExpList.h:2470
void Reset()
Reset geometry information and reset matrices.
Definition: ExpList.h:407
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:2518
void GetFieldDefinitions(std::vector< LibUtilities::FieldDefinitionsSharedPtr > &fielddef)
Definition: ExpList.h:980
virtual void v_WriteTecplotConnectivity(std::ostream &outfile, int expansion)
Definition: ExpList.cpp:2828
BlockMatrixMapShPtr m_blockMat
Definition: ExpList.h:1218
virtual Array< OneD, const unsigned int > v_GetYIDs(void)
Definition: ExpList.cpp:3290
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
multiply the metric jacobi and quadrature weights
Definition: ExpList.cpp:1340
void SetPhys(int i, NekDouble val)
Set the i th value of m_phys to value val.
Definition: ExpList.h:1741
void ExtractTracePhys(Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:2413
void GetBoundaryToElmtMap(Array< OneD, int > &ElmtID, Array< OneD, int > &EdgeID)
Definition: ExpList.h:2450
void PeriodicBwdCopy(const Array< OneD, const NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
Copy and fill the Periodic boundaries.
Definition: ExpList.h:2402
void WriteVtkPieceFooter(std::ostream &outfile, int expansion)
Definition: ExpList.cpp:3089
virtual void v_ExtractPhysToBnd(const int i, const Array< OneD, const NekDouble > &phys, Array< OneD, NekDouble > &bnd)
Definition: ExpList.cpp:4880
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:3678
virtual void v_AddTraceIntegralToOffDiag(const Array< OneD, const NekDouble > &FwdFlux, const Array< OneD, const NekDouble > &BwdFlux, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:5344
virtual std::vector< bool > & v_GetLeftAdjacentTraces(void)
Definition: ExpList.cpp:3911
void InitialiseExpVector(const SpatialDomains::ExpansionInfoMap &expmap)
Define a list of elements using the geometry and basis key information in expmap;.
Definition: ExpList.cpp:1135
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:1092
std::shared_ptr< ExpList > & GetPlane(int n)
Definition: ExpList.h:1062
virtual LibUtilities::BasisSharedPtr v_GetHomogeneousBasis(void)
Definition: ExpList.h:1608
static SpatialDomains::BoundaryConditionShPtr GetBoundaryCondition(const SpatialDomains::BoundaryConditionCollection &collection, unsigned int index, const std::string &variable)
Definition: ExpList.cpp:5038
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:2100
void GetBndElmtExpansion(int i, std::shared_ptr< ExpList > &result, const bool DeclareCoeffPhysArrays=true)
Definition: ExpList.h:2456
virtual void v_DealiasedDotProd(const Array< OneD, Array< OneD, NekDouble >> &inarray1, const Array< OneD, Array< OneD, NekDouble >> &inarray2, Array< OneD, Array< OneD, NekDouble >> &outarray)
Definition: ExpList.cpp:4545
virtual void v_CurlCurl(Array< OneD, Array< OneD, NekDouble >> &Vel, Array< OneD, Array< OneD, NekDouble >> &Q)
Definition: ExpList.cpp:1645
void SetPhysArray(Array< OneD, NekDouble > &inarray)
Sets the array m_phys.
Definition: ExpList.h:1763
virtual void v_PhysDirectionalDeriv(const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:1713
virtual void v_GetBCValues(Array< OneD, NekDouble > &BndVals, const Array< OneD, NekDouble > &TotField, int BndID)
Definition: ExpList.cpp:4555
const Array< OneD, const std::pair< int, int > > & GetCoeffsToElmt() const
Get m_coeffs to elemental value map.
Definition: ExpList.h:2277
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:4426
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:2195
std::shared_ptr< ExpList > & UpdateBndCondExpansion(int i)
Definition: ExpList.h:2298
virtual void v_WriteTecplotHeader(std::ostream &outfile, std::string var="")
Definition: ExpList.cpp:2697
virtual std::shared_ptr< ExpList > & v_GetPlane(int n)
Definition: ExpList.cpp:5059
void NormVectorIProductWRTBase(Array< OneD, const NekDouble > &V1, Array< OneD, const NekDouble > &V2, Array< OneD, NekDouble > &outarray, int BndID)
Definition: ExpList.h:2036
int GetCoeff_Offset(int n) const
Get the start offset position for a global list of m_coeffs correspoinding to element n.
Definition: ExpList.h:2232
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:2204
virtual void v_AddTraceQuadPhysToOffDiag(const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &field)
Definition: ExpList.cpp:4417
virtual void v_ClearGlobalLinSysManager(void)
Definition: ExpList.cpp:3298
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:3379
virtual void v_GetPeriodicEntities(PeriodicMap &periodicVerts, PeriodicMap &periodicEdges, PeriodicMap &periodicFaces)
Definition: ExpList.cpp:5029
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.cpp:1384
void SetCoeffsArray(Array< OneD, NekDouble > &inarray)
Set the m_coeffs array to inarray.
Definition: ExpList.h:2087
LibUtilities::TranspositionSharedPtr GetTransposition(void)
This function returns the transposition class associated with the homogeneous expansion.
Definition: ExpList.h:638
virtual void v_LocalToGlobal(bool UseComm)
Definition: ExpList.cpp:4640
int GetNumElmts(void)
This function returns the number of elements in the expansion which may be different for a homogeoeno...
Definition: ExpList.h:694
void SetBndCondBwdWeight(const int index, const NekDouble value)
Set the weight value for boundary conditions.
Definition: ExpList.h:2293
virtual const Array< OneD, const SpatialDomains::BoundaryConditionShPtr > & v_GetBndConditions()
Definition: ExpList.cpp:4986
void SetHomoLen(const NekDouble lhom)
This function sets the Width of homogeneous direction associated with the homogeneous expansion.
Definition: ExpList.h:652
virtual int v_GetNumElmts(void)
Definition: ExpList.h:1263
int GetNpoints(void) const
Returns the total number of quadrature points m_npoints .
Definition: ExpList.h:1719
SpatialDomains::MeshGraphSharedPtr GetGraph()
Definition: ExpList.h:1051
virtual ~ExpList()
The default destructor.
Definition: ExpList.cpp:1308
void GetFwdBwdTracePhys(Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
Definition: ExpList.h:2359
const Array< OneD, const int > & GetTraceBndMap(void)
Definition: ExpList.h:2329
virtual void v_FillBndCondFromField()
Definition: ExpList.cpp:4625
virtual void v_FwdTransLocalElmt(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:1802
virtual void v_ImposeDirichletConditions(Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:4616
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:1806
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:5007
void FwdTransBndConstrained(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:1816
virtual const std::vector< bool > & v_GetLeftAdjacentFaces(void) const
Definition: ExpList.cpp:4451
virtual void v_GetNormals(Array< OneD, Array< OneD, NekDouble >> &normals)
Populate normals with the normals of all expansions.
Definition: ExpList.cpp:3989
void 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:1868
Array< OneD, int > m_coeff_offset
Offset of elemental data into the array m_coeffs.
Definition: ExpList.h:1210
int EvalBasisNumModesMax(void) const
Evaulates the maximum number of modes in the elemental basis order over all elements.
Definition: ExpList.h:1655
void AddFwdBwdTraceIntegral(const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:2352
virtual NekDouble v_VectorFlux(const Array< OneD, Array< OneD, NekDouble >> &inarray)
Definition: ExpList.cpp:3231
std::shared_ptr< AssemblyMapDG > & GetTraceMap(void)
Definition: ExpList.h:2324
void WriteTecplotField(std::ostream &outfile, int expansion=-1)
Definition: ExpList.h:422
virtual std::shared_ptr< AssemblyMapDG > & v_GetTraceMap()
Definition: ExpList.cpp:3898
Array< OneD, std::pair< int, int > > m_coeffsToElmt
m_coeffs to elemental value map
Definition: ExpList.h:1216
virtual void v_GetFwdBwdTracePhys(Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
Definition: ExpList.cpp:4380
NekDouble PhysEvaluate(const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &phys)
Definition: ExpList.cpp:2635
virtual NekDouble v_Integral(const Array< OneD, const NekDouble > &inarray)
Definition: ExpList.cpp:3217
void PhysDirectionalDeriv(const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:1964
virtual std::map< int, RobinBCInfoSharedPtr > v_GetRobinBCInfo(void)
Definition: ExpList.cpp:5019
virtual const Array< OneD, const int > & v_GetTraceBndMap()
Definition: ExpList.cpp:3906
virtual Array< OneD, const NekDouble > v_HomogeneousEnergy(void)
Definition: ExpList.cpp:3251
virtual NekDouble v_GetHomoLen(void)
Definition: ExpList.cpp:3267
void ExtractFileBCs(const std::string &fileName, LibUtilities::CommSharedPtr comm, const std::string &varName, const std::shared_ptr< ExpList > locExpList)
Definition: ExpList.cpp:3304
void GetNormals(Array< OneD, Array< OneD, NekDouble >> &normals)
Definition: ExpList.h:2334
const LocTraceToTraceMapSharedPtr & GetLocTraceToTraceMap() const
Definition: ExpList.h:2283
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:4507
virtual Array< OneD, const unsigned int > v_GetZIDs(void)
Definition: ExpList.cpp:3282
void SetHomo1DSpecVanVisc(Array< OneD, NekDouble > visc)
This function sets the Spectral Vanishing Viscosity in homogeneous1D expansion.
Definition: ExpList.h:622
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:3546
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:1907
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:1671
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:2431
void GetBoundaryNormals(int i, Array< OneD, Array< OneD, NekDouble >> &normals)
Definition: ExpList.h:2484
virtual void v_SetBndCondBwdWeight(const int index, const NekDouble value)
Definition: ExpList.cpp:4444
void GetMatIpwrtDeriveBase(const Array< OneD, const Array< OneD, NekDouble >> &inarray, const int nDirctn, Array< OneD, DNekMatSharedPtr > &mtxPerVar)
Definition: ExpList.cpp:5353
void ExtractPhysToBnd(int i, const Array< OneD, const NekDouble > &phys, Array< OneD, NekDouble > &bnd)
Definition: ExpList.h:2477
virtual void v_SmoothField(Array< OneD, NekDouble > &field)
Definition: ExpList.cpp:1834
virtual void v_AppendFieldData(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata)
Definition: ExpList.cpp:3506
std::vector< bool > m_collectionsDoInit
Vector of bools to act as an initialise on first call flag.
Definition: ExpList.h:1201
NekDouble GetCoeff(int i)
Get the i th value (coefficient) of m_coeffs.
Definition: ExpList.h:2173
bool m_physState
The state of the array m_phys.
Definition: ExpList.h:1184
void FillBndCondFromField(void)
Fill Bnd Condition expansion from the values stored in expansion.
Definition: ExpList.h:2110
void AddRightIPTBaseMatrix(const Array< OneD, const DNekMatSharedPtr > ElmtJacQuad, Array< OneD, DNekMatSharedPtr > ElmtJacCoef)
Definition: ExpList.cpp:5860
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:2388
virtual void v_WriteVtkPieceData(std::ostream &outfile, int expansion, std::string var)
Definition: ExpList.cpp:3096
virtual void v_WriteVtkPieceHeader(std::ostream &outfile, int expansion, int istrip)
Definition: ExpList.cpp:2955
int GetShapeDimension()
This function returns the dimension of the shape of the element eid.
Definition: ExpList.h:2064
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...
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:661
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:3361
std::vector< int > m_coll_phys_offset
Offset of elemental data into the array m_phys.
Definition: ExpList.h:1207
LibUtilities::CommSharedPtr m_comm
Communicator.
Definition: ExpList.h:1126
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:2413
void WriteVtkHeader(std::ostream &outfile)
Definition: ExpList.cpp:2941
bool GetPhysState(void) const
This function indicates whether the array of physical values (implemented as m_phys) is filled or no...
Definition: ExpList.h:1779
std::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:1196
virtual void v_GetCoords(Array< OneD, NekDouble > &coord_0, Array< OneD, NekDouble > &coord_1, Array< OneD, NekDouble > &coord_2=NullNekDouble1DArray)
Definition: ExpList.cpp:4720
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
Definition: ExpList.h:1136
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:3735
void AddTraceIntegralToOffDiag(const Array< OneD, const NekDouble > &FwdFlux, const Array< OneD, const NekDouble > &BwdFlux, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:1093
void HomogeneousFwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true)
Definition: ExpList.h:1985
std::shared_ptr< DNekMat > GenGlobalMatrixFull(const GlobalLinSysKey &mkey, const std::shared_ptr< AssemblyMapCG > &locToGloMap)
Definition: ExpList.cpp:2260
const Array< OneD, const SpatialDomains::BoundaryConditionShPtr > & GetBndConditions()
Definition: ExpList.h:2426
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:1771
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:681
virtual void v_NormVectorIProductWRTBase(Array< OneD, const NekDouble > &V1, Array< OneD, const NekDouble > &V2, Array< OneD, NekDouble > &outarray, int BndID)
Definition: ExpList.cpp:4563
virtual void v_ExtractElmtToBndPhys(const int i, const Array< OneD, NekDouble > &elmt, Array< OneD, NekDouble > &boundary)
Definition: ExpList.cpp:4803
virtual void v_ExtractPhysToBndElmt(const int i, const Array< OneD, const NekDouble > &phys, Array< OneD, NekDouble > &bndElmt)
Definition: ExpList.cpp:4841
std::vector< LibUtilities::FieldDefinitionsSharedPtr > GetFieldDefinitions()
Definition: ExpList.h:975
void EvaluateBoundaryConditions(const NekDouble time=0.0, const std::string varName="", const NekDouble=NekConstants::kNekUnsetDouble, const NekDouble=NekConstants::kNekUnsetDouble)
Definition: ExpList.h:2437
virtual void v_AddTraceIntegral(const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:4354
virtual void v_FillBwdWithBwdWeight(Array< OneD, NekDouble > &weightave, Array< OneD, NekDouble > &weightjmp)
Definition: ExpList.cpp:4969
Array< OneD, NekDouble > & UpdatePhys()
This function returns (a reference to) the array (implemented as m_phys) containing the function ev...
Definition: ExpList.h:2262
virtual void v_PhysGalerkinProjection1DScaled(const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:5917
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:2054
virtual void v_SetUpPhysNormals()
: Set up a normal along the trace elements between two elements at elemental level
Definition: ExpList.cpp:4780
virtual void v_GlobalToLocal(void)
Definition: ExpList.cpp:4655
virtual LibUtilities::TranspositionSharedPtr v_GetTransposition(void)
Definition: ExpList.cpp:3259
void ExtractElmtToBndPhys(int i, const Array< OneD, NekDouble > &elmt, Array< OneD, NekDouble > &boundary)
Definition: ExpList.h:2463
void SetWaveSpace(const bool wavespace)
Sets the wave space to the one of the possible configuration true or false.
Definition: ExpList.h:1727
const std::shared_ptr< LocalRegions::ExpansionVector > GetExp() const
This function returns the vector of elements in the expansion.
Definition: ExpList.h:2223
void DealiasedDotProd(const Array< OneD, Array< OneD, NekDouble >> &inarray1, const Array< OneD, Array< OneD, NekDouble >> &inarray2, Array< OneD, Array< OneD, NekDouble >> &outarray)
Definition: ExpList.h:2015
virtual void 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:4483
SpatialDomains::MeshGraphSharedPtr m_graph
Mesh associated with this expansion list.
Definition: ExpList.h:1132
virtual void v_MultiplyByInvMassMatrix(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:4474
virtual void v_AddTraceQuadPhysToField(const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &field)
Definition: ExpList.cpp:4408
virtual const Array< OneD, const NekDouble > & v_GetBndCondBwdWeight()
Definition: ExpList.cpp:4436
virtual void v_GetMovingFrames(const SpatialDomains::GeomMMF MMFdir, const Array< OneD, const NekDouble > &CircCentre, Array< OneD, Array< OneD, NekDouble >> &outarray)
Definition: ExpList.cpp:3691
std::shared_ptr< LibUtilities::SessionReader > GetSession() const
Returns the session object.
Definition: ExpList.h:1040
void GlobalToLocal(void)
Scatters from the global coefficients to the local coefficients .
Definition: ExpList.h:2132
virtual NekDouble v_L2(const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &soln=NullNekDouble1DArray)
Definition: ExpList.cpp:3172
const Array< OneD, const NekDouble > & GetBndCondBwdWeight()
Get the weight value for boundary conditions.
Definition: ExpList.h:2288
void WriteVtkFooter(std::ostream &outfile)
Definition: ExpList.cpp:2949
void MultiplyByBlockMatrix(const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:1320
virtual std::shared_ptr< ExpList > & v_GetTrace()
Definition: ExpList.cpp:3890
void GetElmtNormalLength(Array< OneD, NekDouble > &lengthsFwd, Array< OneD, NekDouble > &lengthsBwd)
Get the length of elements in boundary normal direction.
Definition: ExpList.cpp:4206
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:1894
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:1698
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:670
NekDouble L2(const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &soln=NullNekDouble1DArray)
This function calculates the error with respect to soln of the global spectral/hp element approximat...
Definition: ExpList.h:551
virtual std::vector< LibUtilities::FieldDefinitionsSharedPtr > v_GetFieldDefinitions(void)
Definition: ExpList.cpp:3491
std::vector< int > m_coll_coeff_offset
Offset of elemental data into the array m_coeffs.
Definition: ExpList.h:1204
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:1404
void WriteVtkPieceData(std::ostream &outfile, int expansion, std::string var="v")
Definition: ExpList.h:444
void FillBwdWithBwdWeight(Array< OneD, NekDouble > &weightave, Array< OneD, NekDouble > &weightjmp)
Fill Bwd with boundary conditions.
Definition: ExpList.h:2396
virtual void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:2467
virtual Array< OneD, SpatialDomains::BoundaryConditionShPtr > & v_UpdateBndConditions()
Definition: ExpList.cpp:4997
void SetCoeff(int i, NekDouble val)
Set the i th coefficiient in m_coeffs to value val.
Definition: ExpList.h:2073
void LocalToGlobal(bool useComm=true)
Gathers the global coefficients from the local coefficients .
Definition: ExpList.h:2120
void DealiasedProd(const Array< OneD, NekDouble > &inarray1, const Array< OneD, NekDouble > &inarray2, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:2005
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:615
void GetPeriodicEntities(PeriodicMap &periodicVerts, PeriodicMap &periodicEdges, PeriodicMap &periodicFaces=NullPeriodicMap)
Definition: ExpList.h:968
NekDouble Integral(const Array< OneD, const NekDouble > &inarray)
Definition: ExpList.h:603
virtual void v_LinearAdvectionDiffusionReactionSolve(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:4496
virtual void v_WriteTecplotField(std::ostream &outfile, int expansion)
Definition: ExpList.cpp:2908
std::shared_ptr< ExpList > GetSharedThisPtr()
Returns a shared pointer to the current object.
Definition: ExpList.h:1034
Collections::CollectionVector m_collections
Definition: ExpList.h:1198
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:2381
virtual void v_PhysInterp1DScaled(const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:5268
virtual void v_GetBoundaryNormals(int i, Array< OneD, Array< OneD, NekDouble >> &normals)
Definition: ExpList.cpp:4917
void ApplyGeomInfo()
Apply geometry information to each expansion.
Definition: ExpList.cpp:2658
void GetMovingFrames(const SpatialDomains::GeomMMF MMFdir, const Array< OneD, const NekDouble > &CircCentre, Array< OneD, Array< OneD, NekDouble >> &outarray)
Definition: ExpList.h:1925
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:631
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:1548
virtual void v_FwdTransBndConstrained(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:1811
std::unordered_map< int, int > m_elmtToExpId
Mapping from geometry ID of element to index inside m_exp.
Definition: ExpList.h:1227
void GetDiagMatIpwrtBase(const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, DNekMatSharedPtr > &mtxPerVar)
Definition: ExpList.cpp:5496
void FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:1797
void SmoothField(Array< OneD, NekDouble > &field)
Smooth a field across elements.
Definition: ExpList.h:1826
virtual void v_WriteTecplotZone(std::ostream &outfile, int expansion)
Definition: ExpList.cpp:2735
LibUtilities::SessionReaderSharedPtr m_session
Session.
Definition: ExpList.h:1129
void HomogeneousBwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true)
Definition: ExpList.h:1995
virtual void v_PeriodicBwdCopy(const Array< OneD, const NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
Definition: ExpList.cpp:4976
void MultiplyByInvMassMatrix(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:1858
std::map< int, RobinBCInfoSharedPtr > GetRobinBCInfo()
Definition: ExpList.h:963
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:996
void AddTraceIntegral(const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:2339
virtual void v_FillBwdWithBoundCond(const Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd, bool PutFwdInBwdOnBCs)
Definition: ExpList.cpp:4399
void DivideByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Divided by the metric jacobi and quadrature weights.
Definition: ExpList.cpp:1357
virtual const Array< OneD, const std::shared_ptr< ExpList > > & v_GetBndCondExpansions(void)
Definition: ExpList.cpp:3765
void ImposeDirichletConditions(Array< OneD, NekDouble > &outarray)
Impose Dirichlet Boundary Conditions onto Array.
Definition: ExpList.h:2105
Array< OneD, int > m_phys_offset
Offset of elemental data into the array m_phys.
Definition: ExpList.h:1213
virtual void v_SetHomoLen(const NekDouble lhom)
Definition: ExpList.cpp:3275
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:2311
Array< OneD, SpatialDomains::BoundaryConditionShPtr > & UpdateBndConditions()
Definition: ExpList.h:2432
const DNekScalBlkMatSharedPtr GenBlockMatrix(const GlobalMatrixKey &gkey)
This function assembles the block diagonal matrix of local matrices of the type mtype.
Definition: ExpList.cpp:1874
void PhysDeriv(Direction edir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d)
Definition: ExpList.h:1954
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:1764
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:1849
NekDouble GetHomoLen(void)
This function returns the Width of homogeneous direction associated with the homogeneous expansion.
Definition: ExpList.h:645
virtual void v_GetBndElmtExpansion(int i, std::shared_ptr< ExpList > &result, const bool DeclareCoeffPhysArrays)
Definition: ExpList.cpp:4793
virtual void v_FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:4669
bool GetWaveSpace(void) const
This function returns the third direction expansion condition, which can be in wave space (coefficien...
Definition: ExpList.h:1735
virtual std::shared_ptr< ExpList > & v_UpdateBndCondExpansion(int i)
Definition: ExpList.cpp:3773
int GetPhys_Offset(int n) const
Get the start offset position for a global list of m_phys correspoinding to element n.
Definition: ExpList.h:2240
virtual const std::shared_ptr< LocTraceToTraceMap > & v_GetLocTraceToTraceMap(void) const
Definition: ExpList.cpp:5992
void SetModifiedBasis(const bool modbasis)
Set Modified Basis for the stability analysis.
virtual void v_ExtractDataToCoeffs(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, std::string &field, Array< OneD, NekDouble > &coeffs)
Extract data from raw field data into expansion list.
Definition: ExpList.cpp:3562
ExpansionType m_expType
Exapnsion type.
Definition: ExpList.h:1119
void WriteTecplotHeader(std::ostream &outfile, std::string var="")
Definition: ExpList.h:412
NekDouble VectorFlux(const Array< OneD, Array< OneD, NekDouble >> &inarray)
Definition: ExpList.h:608
const std::vector< bool > & GetLeftAdjacentFaces(void) const
Definition: ExpList.h:2408
void AppendFieldData(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata)
Append the element data listed in elements fielddef->m_ElementIDs onto fielddata.
Definition: ExpList.h:988
Array< OneD, NekDouble > m_phys
The global expansion evaluated at the quadrature points.
Definition: ExpList.h:1175
std::vector< bool > & GetLeftAdjacentTraces(void)
Definition: ExpList.h:2490
ExpansionType GetExpType(void)
Returns the type of the expansion.
Definition: ExpList.cpp:1303
void AddTraceQuadPhysToOffDiag(const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &field)
Definition: ExpList.h:880
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:3134
virtual void v_DealiasedProd(const Array< OneD, NekDouble > &inarray1, const Array< OneD, NekDouble > &inarray2, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:4536
void WriteVtkPieceHeader(std::ostream &outfile, int expansion, int istrip=0)
Definition: ExpList.h:435
int GetTotPoints(void) const
Returns the total number of quadrature points m_npoints .
Definition: ExpList.h:1688
void SetExpType(ExpansionType Type)
Returns the type of the expansion.
LibUtilities::BasisSharedPtr GetHomogeneousBasis(void)
Definition: ExpList.h:1057
void WriteTecplotZone(std::ostream &outfile, int expansion=-1)
Definition: ExpList.h:417
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:1616
void GetBCValues(Array< OneD, NekDouble > &BndVals, const Array< OneD, NekDouble > &TotField, int BndID)
Definition: ExpList.h:2026
int GetCoordim(int eid)
This function returns the dimension of the coordinates of the element eid.
Definition: ExpList.h:2054
Describes a matrix with ordering defined by a local to global map.
std::vector< Collection > CollectionVector
Definition: Collection.h:110
static BasisSharedPtr NullBasisSharedPtr
Definition: Basis.h:369
static std::vector< unsigned int > NullUnsignedIntVector
std::shared_ptr< Basis > BasisSharedPtr
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< FieldDefinitions > FieldDefinitionsSharedPtr
Definition: FieldIO.h:177
static std::vector< NekDouble > NullNekDoubleVector
std::shared_ptr< Transposition > TranspositionSharedPtr
static Array< OneD, BasisSharedPtr > NullBasisSharedPtr1DArray
Definition: Basis.h:370
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:54
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:68
std::vector< ExpansionSharedPtr > ExpansionVector
Definition: Expansion.h:70
static ExpList NullExpList
An empty ExpList object.
Definition: ExpList.h:1632
static PeriodicMap NullPeriodicMap
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:89
std::shared_ptr< BlockMatrixMap > BlockMatrixMapShPtr
A shared pointer to a BlockMatrixMap.
Definition: ExpList.h:95
static VarFactorsMap NullVarFactorsMap
static ExpListSharedPtr NullExpListSharedPtr
Definition: ExpList.h:1633
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
std::shared_ptr< LocTraceToTraceMap > LocTraceToTraceMapSharedPtr
static const Array< OneD, ExpListSharedPtr > NullExpListSharedPtrArray
Definition: ExpList.h:2495
std::map< StdRegions::ConstFactorType, Array< OneD, NekDouble > > VarFactorsMap
std::shared_ptr< AssemblyMap > AssemblyMapSharedPtr
Definition: AssemblyMap.h:51
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:93
static const NekDouble kNekUnsetDouble
GeomMMF
Principle direction for MMF.
std::shared_ptr< BoundaryConditionBase > BoundaryConditionShPtr
Definition: Conditions.h:212
std::map< int, CompositeSharedPtr > CompositeMap
Definition: MeshGraph.h:138
std::map< int, BoundaryConditionMapShPtr > BoundaryConditionCollection
Definition: Conditions.h:219
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:172
std::shared_ptr< PointGeom > PointGeomSharedPtr
Definition: Geometry.h:59
std::map< int, ExpansionInfoShPtr > ExpansionInfoMap
Definition: MeshGraph.h:143
std::map< StdRegions::VarCoeffType, Array< OneD, NekDouble > > VarCoeffMap
Definition: StdRegions.hpp:240
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:282
static VarCoeffMap NullVarCoeffMap
Definition: StdRegions.hpp:241
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
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.cpp:1255