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 
44 #include <LocalRegions/Expansion.h>
45 #include <Collections/Collection.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 
65  class AssemblyMapCG;
66  class GlobalLinSysKey;
67  class GlobalMatrix;
68 
69  enum Direction
70  {
71  eX,
72  eY,
73  eZ,
74  eS,
75  eN
76  };
77 
79  {
80  e0D,
81  e1D,
82  e2D,
86  e3D,
87  eNoType
88  };
89 
91  {
92  eX,
93  eY,
94  eZ
95  };
96 
97  /// A map between global matrix keys and their associated block
98  /// matrices.
99  typedef std::map<GlobalMatrixKey,DNekScalBlkMatSharedPtr> BlockMatrixMap;
100  /// A shared pointer to a BlockMatrixMap.
101  typedef std::shared_ptr<BlockMatrixMap> BlockMatrixMapShPtr;
102  /// Shared pointer to an ExpList object.
103  typedef std::shared_ptr<ExpList> ExpListSharedPtr;
104 
105  /// Base class for all multi-elemental spectral/hp expansions.
106  class ExpList: public std::enable_shared_from_this<ExpList>
107  {
108  public:
109  /// The default constructor using a type
111  const ExpansionType Type = eNoType);
112 
113  /// The copy constructor.
115  const ExpList &in,
116  const bool DeclareCoeffPhysArrays = true);
117 
118  /// Constructor copying only elements defined in eIds.
120  const ExpList &in,
121  const std::vector<unsigned int> &eIDs,
122  const bool DeclareCoeffPhysArrays = true,
123  const Collections::ImplementationType ImpType
125 
126  /// Generate an ExpList from a meshgraph \a graph and session file
128  const LibUtilities::SessionReaderSharedPtr &pSession,
130  const bool DeclareCoeffPhysArrays = true,
131  const std::string &var = "DefaultVar",
132  const Collections::ImplementationType ImpType
134 
135  /// Sets up a list of local expansions based on an expansion Map
137  const LibUtilities::SessionReaderSharedPtr &pSession,
138  const SpatialDomains::ExpansionInfoMap &expansions,
139  const bool DeclareCoeffPhysArrays = true,
140  const Collections::ImplementationType ImpType
142 
143  //---------------------------------------------------------
144  // Specialised constructors in ExpListConstructor.cpp
145  //---------------------------------------------------------
146  /// Specialised constructors for 0D Expansions
147  /// Wrapper around LocalRegion::PointExp - used in PrePacing.cpp
150 
151  /// Generate expansions for the trace space expansions used in
152  /// DisContField.
154  const LibUtilities::SessionReaderSharedPtr &pSession,
155  const Array<OneD,const ExpListSharedPtr> &bndConstraint,
156  const Array<OneD, const SpatialDomains
157  ::BoundaryConditionShPtr> &bndCond,
158  const LocalRegions::ExpansionVector &locexp,
160  const bool DeclareCoeffPhysArrays = true,
161  const std::string variable = "DefaultVar",
162  const Collections::ImplementationType ImpType
164 
165  /// Constructor based on domain information only for 1D &
166  /// 2D boundary conditions
168  const LibUtilities::SessionReaderSharedPtr &pSession,
169  const SpatialDomains::CompositeMap &domain,
171  const bool DeclareCoeffPhysArrays = true,
172  const std::string variable = "DefaultVar",
173  bool SetToOneSpaceDimension = false,
174  const LibUtilities::CommSharedPtr comm
176  const Collections::ImplementationType ImpType
178 
179 
180  /// The default destructor.
181  MULTI_REGIONS_EXPORT virtual ~ExpList();
182 
183  /// Returns the total number of local degrees of freedom
184  /// \f$N_{\mathrm{eof}}=\sum_{e=1}^{{N_{\mathrm{el}}}}N^{e}_m\f$.
185  inline int GetNcoeffs(void) const;
186 
187  /// Returns the total number of local degrees of freedom
188  /// for element eid
189  MULTI_REGIONS_EXPORT int GetNcoeffs(const int eid) const;
190 
191  /// Returns the type of the expansion
193 
194  /// Returns the type of the expansion
196 
197  /// Evaulates the maximum number of modes in the elemental basis
198  /// order over all elements
199  inline int EvalBasisNumModesMax(void) const;
200 
201  /// Returns the vector of the number of modes in the elemental
202  /// basis order over all elements.
204  EvalBasisNumModesMaxPerExp(void) const;
205 
206  /// Returns the total number of quadrature points #m_npoints
207  /// \f$=Q_{\mathrm{tot}}\f$.
208  inline int GetTotPoints(void) const;
209 
210  /// Returns the total number of quadrature points for eid's element
211  /// \f$=Q_{\mathrm{tot}}\f$.
212  inline int GetTotPoints(const int eid) const;
213 
214  /// Returns the total number of quadrature points #m_npoints
215  /// \f$=Q_{\mathrm{tot}}\f$.
216  inline int GetNpoints(void) const;
217 
218 
219  /// Returns the total number of qudature points scaled by
220  /// the factor scale on each 1D direction
221  inline int Get1DScaledTotPoints(const NekDouble scale) const;
222 
223  /// Sets the wave space to the one of the possible configuration
224  /// true or false
225  inline void SetWaveSpace(const bool wavespace);
226 
227 
228  ///Set Modified Basis for the stability analysis
229  inline void SetModifiedBasis(const bool modbasis);
230 
231  /// Set the \a i th value of \a m_phys to value \a val
232  inline void SetPhys(int i, NekDouble val);
233 
234  /// This function returns the third direction expansion condition,
235  /// which can be in wave space (coefficient) or not
236  /// It is stored in the variable m_WaveSpace.
237  inline bool GetWaveSpace(void) const;
238 
239  /// Fills the array #m_phys
240  inline void SetPhys(const Array<OneD, const NekDouble> &inarray);
241 
242  /// Sets the array #m_phys
243  inline void SetPhysArray(Array<OneD, NekDouble> &inarray);
244 
245  /// This function manually sets whether the array of physical
246  /// values \f$\boldsymbol{u}_l\f$ (implemented as #m_phys) is
247  /// filled or not.
248  inline void SetPhysState(const bool physState);
249 
250  /// This function indicates whether the array of physical values
251  /// \f$\boldsymbol{u}_l\f$ (implemented as #m_phys) is filled or
252  /// not.
253  inline bool GetPhysState(void) const;
254 
255  /// multiply the metric jacobi and quadrature weights
257  const Array<OneD, const NekDouble> &inarray,
258  Array<OneD, NekDouble> &outarray);
259 
260  /// Divided by the metric jacobi and quadrature weights
262  const Array<OneD, const NekDouble> &inarray,
263  Array<OneD, NekDouble> &outarray);
264 
265  /// This function calculates the inner product of a function
266  /// \f$f(\boldsymbol{x})\f$ with respect to all \em local
267  /// expansion modes \f$\phi_n^e(\boldsymbol{x})\f$.
268  inline void IProductWRTBase_IterPerExp(
269  const Array<OneD, const NekDouble> &inarray,
270  Array<OneD, NekDouble> &outarray);
271 
272  ///
273  inline void IProductWRTBase(
274  const Array<OneD, const NekDouble> &inarray,
275  Array<OneD, NekDouble> &outarray);
276 
277  /// This function calculates the inner product of a function
278  /// \f$f(\boldsymbol{x})\f$ with respect to the derivative (in
279  /// direction \param dir) of all \em local expansion modes
280  /// \f$\phi_n^e(\boldsymbol{x})\f$.
282  const int dir,
283  const Array<OneD, const NekDouble> &inarray,
284  Array<OneD, NekDouble> &outarray);
285 
287  const Array<OneD, const NekDouble> &direction,
288  const Array<OneD, const NekDouble> &inarray,
289  Array<OneD, NekDouble> &outarray);
290 
291  /// This function calculates the inner product of a
292  /// function \f$f(\boldsymbol{x})\f$ with respect to the
293  /// derivative of all \em local expansion modes
294  /// \f$\phi_n^e(\boldsymbol{x})\f$.
296  (const Array<OneD, const Array<OneD, NekDouble> > &inarray,
297  Array<OneD, NekDouble> &outarray);
298 
299  /// This function elementally evaluates the forward transformation
300  /// of a function \f$u(\boldsymbol{x})\f$ onto the global
301  /// spectral/hp expansion.
302  inline void FwdTrans_IterPerExp (
303  const Array<OneD,
304  const NekDouble> &inarray,
305  Array<OneD,NekDouble> &outarray);
306 
307  ///
308  inline void FwdTrans(
309  const Array<OneD,
310  const NekDouble> &inarray,
311  Array<OneD, NekDouble> &outarray);
312 
314  Array<OneD, NekDouble> &array,
315  const NekDouble alpha,
316  const NekDouble exponent,
317  const NekDouble cutoff);
318 
319  /// This function elementally mulplies the coefficient space of
320  /// Sin my the elemental inverse of the mass matrix.
322  const Array<OneD,
323  const NekDouble> &inarray,
324  Array<OneD, NekDouble> &outarray);
325 
326  inline void MultiplyByInvMassMatrix(
327  const Array<OneD, const NekDouble> &inarray,
328  Array<OneD, NekDouble> &outarray);
329 
330  inline void MultiplyByMassMatrix(
331  const Array<OneD, const NekDouble> &inarray,
332  Array<OneD, NekDouble> &outarray)
333  {
335  BwdTrans(inarray,tmp);
336  IProductWRTBase(tmp,outarray);
337  }
338 
339  /// Smooth a field across elements
340  inline void SmoothField(Array<OneD,NekDouble> &field);
341 
342  /// Solve helmholtz problem
343  inline void HelmSolve(
344  const Array<OneD, const NekDouble> &inarray,
345  Array<OneD, NekDouble> &outarray,
346  const StdRegions::ConstFactorMap &factors,
347  const StdRegions::VarCoeffMap &varcoeff =
349  const MultiRegions::VarFactorsMap &varfactors =
351  const Array<OneD, const NekDouble> &dirForcing =
353  const bool PhysSpaceForcing = true);
354 
355  /// Solve Advection Diffusion Reaction
357  const Array<OneD, Array<OneD, NekDouble> > &velocity,
358  const Array<OneD, const NekDouble> &inarray,
359  Array<OneD, NekDouble> &outarray,
360  const NekDouble lambda,
362  dirForcing = NullNekDouble1DArray);
363 
364 
365  /// Solve Advection Diffusion Reaction
366  inline void LinearAdvectionReactionSolve(
367  const Array<OneD, Array<OneD, NekDouble> > &velocity,
368  const Array<OneD, const NekDouble> &inarray,
369  Array<OneD, NekDouble> &outarray,
370  const NekDouble lambda,
372  dirForcing = NullNekDouble1DArray);
373 
374  ///
376  const Array<OneD, const NekDouble> &inarray,
377  Array<OneD, NekDouble> &outarray);
378 
379 
380  /// This function elementally evaluates the backward transformation
381  /// of the global spectral/hp element expansion.
382  inline void BwdTrans_IterPerExp (
383  const Array<OneD, const NekDouble> &inarray,
384  Array<OneD,NekDouble> &outarray);
385 
386  ///
387  inline void BwdTrans (
388  const Array<OneD, const NekDouble> &inarray,
389  Array<OneD,NekDouble> &outarray);
390 
391  /// This function calculates the coordinates of all the elemental
392  /// quadrature points \f$\boldsymbol{x}_i\f$.
393  inline void GetCoords(
394  Array<OneD, NekDouble> &coord_0,
397 
398  // Homogeneous transforms
399  inline void HomogeneousFwdTrans(
400  const Array<OneD, const NekDouble> &inarray,
401  Array<OneD, NekDouble> &outarray,
402  bool Shuff = true,
403  bool UnShuff = true);
404 
405  inline void HomogeneousBwdTrans(
406  const Array<OneD, const NekDouble> &inarray,
407  Array<OneD, NekDouble> &outarray,
408  bool Shuff = true,
409  bool UnShuff = true);
410 
411  inline void DealiasedProd(
412  const Array<OneD, NekDouble> &inarray1,
413  const Array<OneD, NekDouble> &inarray2,
414  Array<OneD, NekDouble> &outarray);
415 
416  inline void DealiasedDotProd(
417  const Array<OneD, Array<OneD, NekDouble> > &inarray1,
418  const Array<OneD, Array<OneD, NekDouble> > &inarray2,
419  Array<OneD, Array<OneD, NekDouble> > &outarray);
420 
421  inline void GetBCValues(
422  Array<OneD, NekDouble> &BndVals,
423  const Array<OneD, NekDouble> &TotField,
424  int BndID);
425 
426  inline void NormVectorIProductWRTBase(
429  Array<OneD, NekDouble> &outarray,
430  int BndID);
431 
432  inline void NormVectorIProductWRTBase(
434  Array<OneD, NekDouble> &outarray);
435 
436  /// Apply geometry information to each expansion.
438 
439  /// Reset geometry information and reset matrices
441  {
442  v_Reset();
443  }
444 
445  void WriteTecplotHeader(std::ostream &outfile,
446  std::string var = "")
447  {
448  v_WriteTecplotHeader(outfile, var);
449  }
450 
452  std::ostream &outfile,
453  int expansion = -1)
454  {
455  v_WriteTecplotZone(outfile, expansion);
456  }
457 
458  void WriteTecplotField(std::ostream &outfile,
459  int expansion = -1)
460  {
461  v_WriteTecplotField(outfile, expansion);
462  }
463 
464  void WriteTecplotConnectivity(std::ostream &outfile,
465  int expansion = -1)
466  {
467  v_WriteTecplotConnectivity(outfile, expansion);
468  }
469 
470  MULTI_REGIONS_EXPORT void WriteVtkHeader(std::ostream &outfile);
471  MULTI_REGIONS_EXPORT void WriteVtkFooter(std::ostream &outfile);
472 
473  void WriteVtkPieceHeader(std::ostream &outfile, int expansion,
474  int istrip = 0)
475  {
476  v_WriteVtkPieceHeader(outfile, expansion, istrip);
477  }
478 
480  std::ostream &outfile,
481  int expansion);
482 
484  std::ostream &outfile,
485  int expansion,
486  std::string var = "v")
487  {
488  v_WriteVtkPieceData(outfile, expansion, var);
489  }
490 
491  /// This function returns the dimension of the coordinates of the
492  /// element \a eid.
493  // inline
494  MULTI_REGIONS_EXPORT int GetCoordim(int eid);
495 
496  /// Set the \a i th coefficiient in \a m_coeffs to value \a val
497  inline void SetCoeff(int i, NekDouble val);
498 
499  /// Set the \a i th coefficiient in #m_coeffs to value \a val
500  inline void SetCoeffs(int i, NekDouble val);
501 
502  /// Set the #m_coeffs array to inarray
503  inline void SetCoeffsArray(Array<OneD, NekDouble> &inarray);
504 
505  /// This function returns the dimension of the shape of the
506  /// element \a eid.
507  // inline
509 
510  /// This function returns (a reference to) the array
511  /// \f$\boldsymbol{\hat{u}}_l\f$ (implemented as #m_coeffs)
512  /// containing all local expansion coefficients.
513  inline const Array<OneD, const NekDouble> &GetCoeffs() const;
514 
515  /// Impose Dirichlet Boundary Conditions onto Array
516  inline void ImposeDirichletConditions(
517  Array<OneD,NekDouble>& outarray);
518 
519 
520  /// Fill Bnd Condition expansion from the values stored in expansion
521  inline void FillBndCondFromField(void);
522 
523  /// Fill Bnd Condition expansion in nreg from the values
524  /// stored in expansion
525  inline void FillBndCondFromField(const int nreg);
526 
527  /// Gathers the global coefficients \f$\boldsymbol{\hat{u}}_g\f$
528  /// from the local coefficients \f$\boldsymbol{\hat{u}}_l\f$.
529  // inline
530  MULTI_REGIONS_EXPORT inline void LocalToGlobal(bool useComm = true);
531 
533  const Array<OneD, const NekDouble> &inarray,
534  Array<OneD,NekDouble> &outarray,
535  bool useComm = true);
536 
537  /// Scatters from the global coefficients
538  /// \f$\boldsymbol{\hat{u}}_g\f$ to the local coefficients
539  /// \f$\boldsymbol{\hat{u}}_l\f$.
540  // inline
541  MULTI_REGIONS_EXPORT inline void GlobalToLocal(void);
542 
543  /**
544  * This operation is evaluated as:
545  * \f{tabbing}
546  * \hspace{1cm} \= Do \= $e=$ $1, N_{\mathrm{el}}$ \\
547  * > > Do \= $i=$ $0,N_m^e-1$ \\
548  * > > > $\boldsymbol{\hat{u}}^{e}[i] = \mbox{sign}[e][i] \cdot
549  * \boldsymbol{\hat{u}}_g[\mbox{map}[e][i]]$ \\
550  * > > continue \\
551  * > continue
552  * \f}
553  * where \a map\f$[e][i]\f$ is the mapping array and \a
554  * sign\f$[e][i]\f$ is an array of similar dimensions ensuring the
555  * correct modal connectivity between the different elements (both
556  * these arrays are contained in the data member #m_locToGloMap). This
557  * operation is equivalent to the scatter operation
558  * \f$\boldsymbol{\hat{u}}_l=\mathcal{A}\boldsymbol{\hat{u}}_g\f$,
559  * where \f$\mathcal{A}\f$ is the
560  * \f$N_{\mathrm{eof}}\times N_{\mathrm{dof}}\f$ permutation matrix.
561  *
562  * @param inarray An array of size \f$N_\mathrm{dof}\f$
563  * containing the global degrees of freedom
564  * \f$\boldsymbol{x}_g\f$.
565  * @param outarray The resulting local degrees of freedom
566  * \f$\boldsymbol{x}_l\f$ will be stored in this
567  * array of size \f$N_\mathrm{eof}\f$.
568  */
570  const Array<OneD, const NekDouble> &inarray,
571  Array<OneD,NekDouble> &outarray);
572 
573  /// Get the \a i th value (coefficient) of #m_coeffs
574  inline NekDouble GetCoeff(int i);
575 
576  /// Get the \a i th value (coefficient) of #m_coeffs
577  inline NekDouble GetCoeffs(int i);
578 
579  /// This function returns (a reference to) the array
580  /// \f$\boldsymbol{u}_l\f$ (implemented as #m_phys) containing the
581  /// function \f$u^{\delta}(\boldsymbol{x})\f$ evaluated at the
582  /// quadrature points.
583  // inline
585  &GetPhys() const;
586 
587  /// This function calculates the \f$L_\infty\f$ error of the global
588  /// spectral/hp element approximation.
590  const Array<OneD, const NekDouble> &inarray,
592 
593  /// This function calculates the \f$L_2\f$ error with
594  /// respect to soln of the global
595  /// spectral/hp element approximation.
597  const Array<OneD, const NekDouble> &inarray,
599  {
600  return v_L2(inarray, soln);
601  }
602 
603  /// Calculates the \f$H^1\f$ error of the global spectral/hp
604  /// element approximation.
606  const Array<OneD, const NekDouble> &inarray,
608 
609  /**
610  * The integration is evaluated locally, that is
611  * \f[\int
612  * f(\boldsymbol{x})d\boldsymbol{x}=\sum_{e=1}^{{N_{\mathrm{el}}}}
613  * \left\{\int_{\Omega_e}f(\boldsymbol{x})d\boldsymbol{x}\right\}, \f]
614  * where the integration over the separate elements is done by the
615  * function StdRegions#StdExpansion#Integral, which discretely
616  * evaluates the integral using Gaussian quadrature.
617  *
618  * Note that the array #m_phys should be filled with the values of the
619  * function \f$f(\boldsymbol{x})\f$ at the quadrature points
620  * \f$\boldsymbol{x}_i\f$.
621  *
622  * @return The value of the discretely evaluated integral
623  * \f$\int f(\boldsymbol{x})d\boldsymbol{x}\f$.
624  */
626  {
627  ASSERTL1(m_physState == true,
628  "local physical space is not true ");
629 
630  return Integral(m_phys);
631  }
632 
633  /**
634  * The integration is evaluated locally, that is
635  * \f[\int
636  * f(\boldsymbol{x})d\boldsymbol{x}=\sum_{e=1}^{{N_{\mathrm{el}}}}
637  * \left\{\int_{\Omega_e}f(\boldsymbol{x})d\boldsymbol{x}\right\}, \f]
638  * where the integration over the separate elements is done by the
639  * function StdRegions#StdExpansion#Integral, which discretely
640  * evaluates the integral using Gaussian quadrature.
641  *
642  * @param inarray An array of size \f$Q_{\mathrm{tot}}\f$
643  * containing the values of the function
644  * \f$f(\boldsymbol{x})\f$ at the quadrature
645  * points \f$\boldsymbol{x}_i\f$.
646  * @return The value of the discretely evaluated integral
647  * \f$\int f(\boldsymbol{x})d\boldsymbol{x}\f$.
648  */
650  {
651  return v_Integral(inarray);
652  }
653 
655  {
656  return v_VectorFlux(inarray);
657  }
658 
659  /// This function calculates the energy associated with
660  /// each one of the modesof a 3D homogeneous nD expansion
662  {
663  return v_HomogeneousEnergy();
664  }
665 
666  /// This function sets the Spectral Vanishing Viscosity
667  /// in homogeneous1D expansion.
669  {
671  }
672 
673  /// This function returns a vector containing the wave
674  /// numbers in z-direction associated
675  /// with the 3D homogenous expansion. Required if a
676  /// parellelisation is applied in the Fourier direction
678  {
679  return v_GetZIDs();
680  }
681 
682  /// This function returns the transposition class
683  /// associated with the homogeneous expansion.
685  {
686  return v_GetTransposition();
687  }
688 
689  /// This function returns the Width of homogeneous direction
690  /// associated with the homogeneous expansion.
692  {
693  return v_GetHomoLen();
694  }
695 
696  /// This function sets the Width of homogeneous direction
697  /// associated with the homogeneous expansion.
698  void SetHomoLen(const NekDouble lhom)
699  {
700  return v_SetHomoLen(lhom);
701  }
702 
703  /// This function returns a vector containing the wave
704  /// numbers in y-direction associated
705  /// with the 3D homogenous expansion. Required if a
706  /// parellelisation is applied in the Fourier direction
708  {
709  return v_GetYIDs();
710  }
711 
712  /// This function interpolates the physical space points in
713  /// \a inarray to \a outarray using the same points defined in the
714  /// expansion but where the number of points are rescaled
715  /// by \a 1DScale
717  const NekDouble scale,
718  const Array<OneD, NekDouble> &inarray,
719  Array<OneD, NekDouble> &outarray)
720  {
721  v_PhysInterp1DScaled(scale, inarray,outarray);
722  }
723 
724  /// This function Galerkin projects the physical space points in
725  /// \a inarray to \a outarray where inarray is assumed to
726  /// be defined in the expansion but where the number of
727  /// points are rescaled by \a 1DScale
729  const NekDouble scale,
730  const Array<OneD, NekDouble> &inarray,
731  Array<OneD, NekDouble> &outarray)
732  {
733  v_PhysGalerkinProjection1DScaled(scale, inarray, outarray);
734  }
735 
736  /// This function returns the number of elements in the expansion.
737  inline int GetExpSize(void);
738 
739 
740  /// This function returns the number of elements in the
741  /// expansion which may be different for a homogeoenous extended
742  /// expansionp.
743  inline int GetNumElmts(void)
744  {
745  return v_GetNumElmts();
746  }
747 
748  /// This function returns the vector of elements in the expansion.
749  inline const std::shared_ptr<LocalRegions::ExpansionVector>
750  GetExp() const;
751 
752  /// This function returns (a shared pointer to) the local elemental
753  /// expansion of the \f$n^{\mathrm{th}}\f$ element.
754  inline LocalRegions::ExpansionSharedPtr& GetExp(int n) const;
755 
756  /// This function returns (a shared pointer to) the local elemental
757  /// expansion containing the arbitrary point given by \a gloCoord.
759  const Array<OneD, const NekDouble> &gloCoord);
760 
761  /**
762  * @brief This function returns the index of the local elemental
763  * expansion containing the arbitrary point given by \a gloCoord,
764  * within a distance tolerance of tol.
765  *
766  * If returnNearestElmt is true and no element contains the point,
767  * this function returns the nearest element whose bounding box contains
768  * the point. The bounding box has a 10% margin in each direction.
769  *
770  * @param gloCoord (input) coordinate in physics space
771  * @param locCoords (output) local coordinate xi in the returned element
772  * @param tol distance tolerance to judge if a point lies in an element
773  * @param returnNearestElmt if true and no element contains this point, the
774  * nearest element whose bounding box contains this
775  * point is returned
776  * @param cachedId an initial guess of the most possible element index
777  * @param maxDistance if returnNearestElmt is set as true, the nearest
778  * element will be returned. But the distance of the
779  * nearest element and this point should be <= maxDistance.
780  *
781  * @return element index; if no element is found, -1 is returned.
782  */
784  const Array<OneD, const NekDouble> &gloCoord,
785  NekDouble tol = 0.0,
786  bool returnNearestElmt = false,
787  int cachedId = -1,
788  NekDouble maxDistance = 1e6);
789 
790  /** This function returns the index and the Local
791  * Cartesian Coordinates \a locCoords of the local
792  * elemental expansion containing the arbitrary point
793  * given by \a gloCoords.
794  **/
796  const Array<OneD, const NekDouble> &gloCoords,
797  Array<OneD, NekDouble> &locCoords,
798  NekDouble tol = 0.0,
799  bool returnNearestElmt = false,
800  int cachedId = -1,
801  NekDouble maxDistance = 1e6);
802 
803  /** This function return the expansion field value
804  * at the coordinates given as input.
805  **/
807  const Array<OneD, const NekDouble> &coords,
808  const Array<OneD, const NekDouble> &phys);
809 
810  /// Get the start offset position for a global list of #m_coeffs
811  /// correspoinding to element n.
812  inline int GetCoeff_Offset(int n) const;
813 
814  /// Get the start offset position for a global list of m_phys
815  /// correspoinding to element n.
816  inline int GetPhys_Offset(int n) const;
817 
818  /// This function returns (a reference to) the array
819  /// \f$\boldsymbol{\hat{u}}_l\f$ (implemented as #m_coeffs)
820  /// containing all local expansion coefficients.
822 
823  /// This function returns (a reference to) the array
824  /// \f$\boldsymbol{u}_l\f$ (implemented as #m_phys) containing the
825  /// function \f$u^{\delta}(\boldsymbol{x})\f$ evaluated at the
826  /// quadrature points.
828 
829  inline void PhysDeriv(
830  Direction edir,
831  const Array<OneD, const NekDouble> &inarray,
832  Array<OneD, NekDouble> &out_d);
833 
834  /// This function discretely evaluates the derivative of a function
835  /// \f$f(\boldsymbol{x})\f$ on the domain consisting of all
836  /// elements of the expansion.
837  inline void PhysDeriv(
838  const Array<OneD, const NekDouble> &inarray,
839  Array<OneD, NekDouble> &out_d0,
842 
843  inline void PhysDeriv(
844  const int dir,
845  const Array<OneD, const NekDouble> &inarray,
846  Array<OneD, NekDouble> &out_d);
847 
848  inline void CurlCurl(
851 
852  inline void PhysDirectionalDeriv(
853  const Array<OneD, const NekDouble> &direction,
854  const Array<OneD, const NekDouble> &inarray,
855  Array<OneD, NekDouble> &outarray);
856 
857  inline void GetMovingFrames(
858  const SpatialDomains::GeomMMF MMFdir,
859  const Array<OneD, const NekDouble> &CircCentre,
860  Array<OneD, Array<OneD, NekDouble> > &outarray);
861 
862  // functions associated with DisContField
865 
866  /// Get the weight value for boundary conditions
867  inline const Array<OneD, const NekDouble>
869 
870  /// Set the weight value for boundary conditions
871  inline void SetBndCondBwdWeight(
872  const int index,
873  const NekDouble value);
874 
875  inline std::shared_ptr<ExpList> &UpdateBndCondExpansion(int i);
876 
877  inline void Upwind(
879  const Array<OneD, const NekDouble> &Fwd,
880  const Array<OneD, const NekDouble> &Bwd,
882 
883  inline void Upwind(
884  const Array<OneD, const Array<OneD, NekDouble> > &Vec,
885  const Array<OneD, const NekDouble> &Fwd,
886  const Array<OneD, const NekDouble> &Bwd,
888 
889  /**
890  * Return a reference to the trace space associated with this
891  * expansion list.
892  */
893  inline std::shared_ptr<ExpList> &GetTrace();
894 
895  inline std::shared_ptr<AssemblyMapDG> &GetTraceMap(void);
896 
897  inline const Array<OneD, const int> &GetTraceBndMap(void);
898 
899  inline void GetNormals(Array<OneD, Array<OneD, NekDouble> >&normals);
900 
901  /// Get the length of elements in boundary normal direction
903  Array<OneD, NekDouble> &lengthsFwd,
904  Array<OneD, NekDouble> &lengthsBwd);
905 
906  /// Get the weight value for boundary conditions
907  /// for boundary average and jump calculations
909  Array<OneD, NekDouble> &weightAver,
910  Array<OneD, NekDouble> &weightJump);
911 
912  inline void AddTraceIntegral(
915  Array<OneD, NekDouble> &outarray);
916 
917  inline void AddTraceIntegral(
919  Array<OneD, NekDouble> &outarray);
920 
921  inline void AddFwdBwdTraceIntegral(
922  const Array<OneD, const NekDouble> &Fwd,
923  const Array<OneD, const NekDouble> &Bwd,
924  Array<OneD, NekDouble> &outarray);
925 
926  inline void GetFwdBwdTracePhys(
928  Array<OneD,NekDouble> &Bwd);
929 
930  inline void GetFwdBwdTracePhys(
931  const Array<OneD,const NekDouble> &field,
934  bool FillBnd = true,
935  bool PutFwdInBwdOnBCs = false,
936  bool DoExchange = true);
937 
938  inline void FillBwdWithBoundCond(
939  const Array<OneD, NekDouble> &Fwd,
941  bool PutFwdInBwdOnBCs = false);
942 
943  /// Add Fwd and Bwd value to field,
944  /// a reverse procedure of GetFwdBwdTracePhys
945  inline void AddTraceQuadPhysToField(
946  const Array<OneD, const NekDouble> &Fwd,
947  const Array<OneD, const NekDouble> &Bwd,
948  Array<OneD, NekDouble> &field);
949 
951  const Array<OneD, const NekDouble> &Fwd,
952  const Array<OneD, const NekDouble> &Bwd,
953  Array<OneD, NekDouble> &field)
954  {
955  v_AddTraceQuadPhysToOffDiag(Fwd,Bwd,field);
956  }
957 
958  inline void GetLocTraceFromTracePts(
959  const Array<OneD, const NekDouble> &Fwd,
960  const Array<OneD, const NekDouble> &Bwd,
961  Array<OneD, NekDouble> &locTraceFwd,
962  Array<OneD, NekDouble> &locTraceBwd);
963 
964  /// Fill Bwd with boundary conditions
965  inline void FillBwdWithBwdWeight(
966  Array<OneD, NekDouble> &weightave,
967  Array<OneD, NekDouble> &weightjmp);
968 
969  /// Copy and fill the Periodic boundaries
970  inline void PeriodicBwdCopy(
971  const Array<OneD, const NekDouble> &Fwd,
973 
974  inline const std::vector<bool> &GetLeftAdjacentFaces(void) const;
975 
976  inline void ExtractTracePhys(Array<OneD,NekDouble> &outarray);
977 
978  inline void ExtractTracePhys(
979  const Array<OneD, const NekDouble> &inarray,
980  Array<OneD,NekDouble> &outarray);
981 
982  inline const Array<OneD, const SpatialDomains::
983  BoundaryConditionShPtr>& GetBndConditions();
984 
985  inline Array<OneD, SpatialDomains::
986  BoundaryConditionShPtr>& UpdateBndConditions();
987 
988  inline void EvaluateBoundaryConditions(
989  const NekDouble time = 0.0,
990  const std::string varName = "",
993 
994  // Routines for continous matrix solution
995  /// This function calculates the result of the multiplication of a
996  /// matrix of type specified by \a mkey with a vector given by \a
997  /// inarray.
998  inline void GeneralMatrixOp(
999  const GlobalMatrixKey &gkey,
1000  const Array<OneD,const NekDouble> &inarray,
1001  Array<OneD, NekDouble> &outarray);
1002 
1004  const GlobalMatrixKey &gkey,
1005  const Array<OneD,const NekDouble> &inarray,
1006  Array<OneD, NekDouble> &outarray);
1007 
1008  inline void SetUpPhysNormals();
1009 
1010  inline void GetBoundaryToElmtMap(Array<OneD, int> &ElmtID,
1011  Array<OneD,int> &EdgeID);
1012 
1013  inline void GetBndElmtExpansion(int i,
1014  std::shared_ptr<ExpList> &result,
1015  const bool DeclareCoeffPhysArrays = true);
1016 
1017  inline void ExtractElmtToBndPhys(int i,
1018  const Array<OneD, NekDouble> &elmt,
1019  Array<OneD, NekDouble> &boundary);
1020 
1021  inline void ExtractPhysToBndElmt(int i,
1022  const Array<OneD, const NekDouble> &phys,
1023  Array<OneD, NekDouble> &bndElmt);
1024 
1025  inline void ExtractPhysToBnd(int i,
1026  const Array<OneD, const NekDouble> &phys,
1027  Array<OneD, NekDouble> &bnd);
1028 
1029  inline void GetBoundaryNormals(int i,
1030  Array<OneD, Array<OneD, NekDouble> > &normals);
1031 
1033  std::vector<LibUtilities::FieldDefinitionsSharedPtr> &fielddef,
1034  int NumHomoDir = 0,
1037  std::vector<NekDouble> &HomoLen =
1039  bool homoStrips = false,
1040  std::vector<unsigned int> &HomoSIDs =
1042  std::vector<unsigned int> &HomoZIDs =
1044  std::vector<unsigned int> &HomoYIDs =
1046 
1047 
1048  std::map<int, RobinBCInfoSharedPtr> GetRobinBCInfo()
1049  {
1050  return v_GetRobinBCInfo();
1051  }
1052 
1054  PeriodicMap &periodicVerts,
1055  PeriodicMap &periodicEdges,
1056  PeriodicMap &periodicFaces = NullPeriodicMap)
1057  {
1058  v_GetPeriodicEntities(periodicVerts, periodicEdges, periodicFaces);
1059  }
1060 
1061  std::vector<LibUtilities::FieldDefinitionsSharedPtr>
1063  {
1064  return v_GetFieldDefinitions();
1065  }
1066 
1067  void GetFieldDefinitions(std::vector<LibUtilities::FieldDefinitionsSharedPtr> &fielddef)
1068  {
1069  v_GetFieldDefinitions(fielddef);
1070  }
1071 
1072  /// Append the element data listed in elements
1073  /// fielddef->m_ElementIDs onto fielddata
1076  std::vector<NekDouble> &fielddata)
1077  {
1078  v_AppendFieldData(fielddef,fielddata);
1079  }
1080 
1081 
1082  /// Append the data in coeffs listed in elements
1083  /// fielddef->m_ElementIDs onto fielddata
1086  std::vector<NekDouble> &fielddata,
1087  Array<OneD, NekDouble> &coeffs)
1088  {
1089  v_AppendFieldData(fielddef,fielddata,coeffs);
1090  }
1091 
1092 
1093  /** \brief Extract the data in fielddata into the coeffs
1094  * using the basic ExpList Elemental expansions rather
1095  * than planes in homogeneous case
1096  */
1099  std::vector<NekDouble> &fielddata,
1100  std::string &field,
1101  Array<OneD, NekDouble> &coeffs);
1102 
1103 
1104  /** \brief Extract the data from fromField using
1105  * fromExpList the coeffs using the basic ExpList
1106  * Elemental expansions rather than planes in homogeneous
1107  * case
1108  */
1110  const std::shared_ptr<ExpList> &fromExpList,
1111  const Array<OneD, const NekDouble> &fromCoeffs,
1112  Array<OneD, NekDouble> &toCoeffs);
1113 
1114 
1115  //Extract data in fielddata into the m_coeffs_list for the 3D stability analysis (base flow is 2D)
1118  std::vector<NekDouble> &fielddata,
1119  std::string &field,
1120  Array<OneD, NekDouble> &coeffs);
1121 
1123  const int ElementID,
1124  const NekDouble scalar1,
1125  const NekDouble scalar2,
1126  Array<OneD, NekDouble> &outarray);
1127 
1128  /// Returns a shared pointer to the current object.
1129  std::shared_ptr<ExpList> GetSharedThisPtr()
1130  {
1131  return shared_from_this();
1132  }
1133 
1134  /// Returns the session object
1135  std::shared_ptr<LibUtilities::SessionReader> GetSession() const
1136  {
1137  return m_session;
1138  }
1139 
1140  /// Returns the comm object
1141  std::shared_ptr<LibUtilities::Comm> GetComm()
1142  {
1143  return m_comm;
1144  }
1145 
1147  {
1148  return m_graph;
1149  }
1150 
1151  // Wrapper functions for Homogeneous Expansions
1153  {
1154  return v_GetHomogeneousBasis();
1155  }
1156 
1157  std::shared_ptr<ExpList> &GetPlane(int n)
1158  {
1159  return v_GetPlane(n);
1160  }
1161 
1165 
1167 
1168  /// Get m_coeffs to elemental value map
1169  MULTI_REGIONS_EXPORT inline const
1171  &GetCoeffsToElmt() const;
1172 
1176  Array<OneD, DNekMatSharedPtr> &fieldMat);
1177 
1179  const Array<OneD, const Array<OneD, NekDouble> >&inarray,
1180  const int nDirctn, Array<OneD, DNekMatSharedPtr> &mtxPerVar);
1181 
1183  const TensorOfArray3D<NekDouble> &inarray,
1184  Array<OneD, DNekMatSharedPtr> &mtxPerVar);
1185 
1187  const Array<OneD, const Array<OneD, NekDouble> >&inarray,
1188  Array<OneD, DNekMatSharedPtr> &mtxPerVar);
1189 
1191  const Array<OneD, const NekDouble> &FwdFlux,
1192  const Array<OneD, const NekDouble> &BwdFlux,
1193  Array<OneD, NekDouble> &outarray)
1194  {
1195  v_AddTraceIntegralToOffDiag(FwdFlux,BwdFlux,outarray);
1196  }
1197 
1199  const int dir,
1200  const Array<OneD, const DNekMatSharedPtr> ElmtJacQuad,
1201  Array<OneD, DNekMatSharedPtr> ElmtJacCoef);
1202 
1204  const Array<OneD, const DNekMatSharedPtr> ElmtJacQuad,
1205  Array<OneD, DNekMatSharedPtr> ElmtJacCoef);
1206 
1208  &GetLocTraceToTraceMap() const;
1209 
1210  protected:
1211  /// Exapnsion type
1213 
1214 
1215  std::shared_ptr<DNekMat> GenGlobalMatrixFull(
1216  const GlobalLinSysKey &mkey,
1217  const std::shared_ptr<AssemblyMapCG> &locToGloMap);
1218 
1219  /// Communicator
1221 
1222  /// Session
1224 
1225  /// Mesh associated with this expansion list.
1227 
1228  /// The total number of local degrees of freedom. #m_ncoeffs
1229  /// \f$=N_{\mathrm{eof}}=\sum_{e=1}^{{N_{\mathrm{el}}}}N^{e}_l\f$
1231 
1232  /** The total number of quadrature points. #m_npoints
1233  *\f$=Q_{\mathrm{tot}}=\sum_{e=1}^{{N_{\mathrm{el}}}}N^{e}_Q\f$
1234  **/
1236 
1237  /**
1238  * \brief Concatenation of all local expansion coefficients.
1239  *
1240  * The array of length #m_ncoeffs\f$=N_{\mathrm{eof}}\f$ which is
1241  * the concatenation of the local expansion coefficients
1242  * \f$\hat{u}_n^e\f$ over all \f$N_{\mathrm{el}}\f$ elements
1243  * \f[\mathrm{\texttt{m\_coeffs}}=\boldsymbol{\hat{u}}_{l} =
1244  * \underline{\boldsymbol{\hat{u}}}^e = \left [ \begin{array}{c}
1245  * \boldsymbol{\hat{u}}^{1} \ \
1246  * \boldsymbol{\hat{u}}^{2} \ \
1247  * \vdots \ \
1248  * \boldsymbol{\hat{u}}^{{{N_{\mathrm{el}}}}} \end{array} \right ],
1249  * \quad
1250  * \mathrm{where}\quad \boldsymbol{\hat{u}}^{e}[n]=\hat{u}_n^{e}\f]
1251  */
1253 
1254  /**
1255  * \brief The global expansion evaluated at the quadrature points
1256  *
1257  * The array of length #m_npoints\f$=Q_{\mathrm{tot}}\f$ containing
1258  * the evaluation of \f$u^{\delta}(\boldsymbol{x})\f$ at the
1259  * quadrature points \f$\boldsymbol{x}_i\f$.
1260  * \f[\mathrm{\texttt{m\_phys}}=\boldsymbol{u}_{l} =
1261  * \underline{\boldsymbol{u}}^e = \left [ \begin{array}{c}
1262  * \boldsymbol{u}^{1} \ \
1263  * \boldsymbol{u}^{2} \ \
1264  * \vdots \ \
1265  * \boldsymbol{u}^{{{N_{\mathrm{el}}}}} \end{array} \right ],\quad
1266  * \mathrm{where}\quad
1267  * \boldsymbol{u}^{e}[i]=u^{\delta}(\boldsymbol{x}_i)\f]
1268  */
1270 
1271  /**
1272  * \brief The state of the array #m_phys.
1273  *
1274  * Indicates whether the array #m_phys, created to contain the
1275  * evaluation of \f$u^{\delta}(\boldsymbol{x})\f$ at the quadrature
1276  * points \f$\boldsymbol{x}_i\f$, is filled with these values.
1277  */
1279 
1280  /**
1281  * \brief The list of local expansions.
1282  *
1283  * The (shared pointer to the) vector containing (shared pointers
1284  * to) all local expansions. The fact that the local expansions are
1285  * all stored as a (pointer to a) #StdExpansion, the abstract base
1286  * class for all local expansions, allows a general implementation
1287  * where most of the routines for the derived classes are defined
1288  * in the #ExpList base class.
1289  */
1290  std::shared_ptr<LocalRegions::ExpansionVector> m_exp;
1291 
1293 
1294  /// Vector of bools to act as an initialise on first call flag
1295  std::vector<bool> m_collectionsDoInit;
1296 
1297  /// Offset of elemental data into the array #m_coeffs
1298  std::vector<int> m_coll_coeff_offset;
1299 
1300  /// Offset of elemental data into the array #m_phys
1301  std::vector<int> m_coll_phys_offset;
1302 
1303  /// Offset of elemental data into the array #m_coeffs
1305 
1306  /// Offset of elemental data into the array #m_phys
1308 
1309  /// m_coeffs to elemental value map
1311 
1313 
1314  //@todo should this be in ExpList or
1315  // ExpListHomogeneous1D.cpp it's a bool which determine if
1316  // the expansion is in the wave space (coefficient space)
1317  // or not
1319 
1320  /// Mapping from geometry ID of element to index inside #m_exp
1321  std::unordered_map<int, int> m_elmtToExpId;
1322 
1323  /// This function assembles the block diagonal matrix of local
1324  /// matrices of the type \a mtype.
1326  const GlobalMatrixKey &gkey);
1327 
1329  const GlobalMatrixKey &gkey);
1330 
1331  void MultiplyByBlockMatrix(
1332  const GlobalMatrixKey &gkey,
1333  const Array<OneD,const NekDouble> &inarray,
1334  Array<OneD, NekDouble> &outarray);
1335 
1336  /// Generates a global matrix from the given key and map.
1337  std::shared_ptr<GlobalMatrix> GenGlobalMatrix(
1338  const GlobalMatrixKey &mkey,
1339  const std::shared_ptr<AssemblyMapCG> &locToGloMap);
1340 
1341 
1343  const std::shared_ptr<DNekMat> &Gmat,
1344  Array<OneD, NekDouble> &EigValsReal,
1345  Array<OneD, NekDouble> &EigValsImag,
1346  Array<OneD, NekDouble> &EigVecs
1348 
1349 
1350  /// This operation constructs the global linear system of type \a
1351  /// mkey.
1352  std::shared_ptr<GlobalLinSys> GenGlobalLinSys(
1353  const GlobalLinSysKey &mkey,
1354  const std::shared_ptr<AssemblyMapCG> &locToGloMap);
1355 
1356  /// Generate a GlobalLinSys from information provided by the key
1357  /// "mkey" and the mapping provided in LocToGloBaseMap.
1358  std::shared_ptr<GlobalLinSys> GenGlobalBndLinSys(
1359  const GlobalLinSysKey &mkey,
1360  const AssemblyMapSharedPtr &locToGloMap);
1361 
1362  // Virtual prototypes
1363 
1364  virtual int v_GetNumElmts(void)
1365  {
1366  return (*m_exp).size();
1367  }
1368 
1370  &v_GetBndCondExpansions(void);
1371 
1372  virtual const Array<OneD, const NekDouble>
1374 
1375  virtual void v_SetBndCondBwdWeight(
1376  const int index,
1377  const NekDouble value);
1378 
1379  virtual std::shared_ptr<ExpList> &v_UpdateBndCondExpansion(int i);
1380 
1381  virtual void v_Upwind(
1382  const Array<OneD, const Array<OneD, NekDouble> > &Vec,
1383  const Array<OneD, const NekDouble> &Fwd,
1384  const Array<OneD, const NekDouble> &Bwd,
1386 
1387  virtual void v_Upwind(
1388  const Array<OneD, const NekDouble> &Vn,
1389  const Array<OneD, const NekDouble> &Fwd,
1390  const Array<OneD, const NekDouble> &Bwd,
1392 
1393  virtual std::shared_ptr<ExpList> &v_GetTrace();
1394 
1395  virtual std::shared_ptr<AssemblyMapDG> &v_GetTraceMap();
1396 
1397  virtual const Array<OneD, const int> &v_GetTraceBndMap();
1398 
1399  virtual const std::shared_ptr<LocTraceToTraceMap>
1400  &v_GetLocTraceToTraceMap(void) const;
1401 
1402  /// Populate \a normals with the normals of all expansions.
1403  virtual void v_GetNormals(
1404  Array<OneD, Array<OneD, NekDouble> > &normals);
1405 
1406  virtual void v_AddTraceIntegral(
1407  const Array<OneD, const NekDouble> &Fx,
1408  const Array<OneD, const NekDouble> &Fy,
1409  Array<OneD, NekDouble> &outarray);
1410 
1411  virtual void v_AddTraceIntegral(
1412  const Array<OneD, const NekDouble> &Fn,
1413  Array<OneD, NekDouble> &outarray);
1414 
1415  virtual void v_AddFwdBwdTraceIntegral(
1416  const Array<OneD, const NekDouble> &Fwd,
1417  const Array<OneD, const NekDouble> &Bwd,
1418  Array<OneD, NekDouble> &outarray);
1419 
1420  virtual void v_GetFwdBwdTracePhys(
1421  Array<OneD,NekDouble> &Fwd,
1422  Array<OneD,NekDouble> &Bwd);
1423 
1424  virtual void v_GetFwdBwdTracePhys(
1425  const Array<OneD,const NekDouble> &field,
1426  Array<OneD,NekDouble> &Fwd,
1427  Array<OneD,NekDouble> &Bwd,
1428  bool FillBnd = true,
1429  bool PutFwdInBwdOnBCs = false,
1430  bool DoExchange = true);
1431 
1432  virtual void v_FillBwdWithBoundCond(
1433  const Array<OneD, NekDouble> &Fwd,
1435  bool PutFwdInBwdOnBCs);
1436 
1437  virtual void v_AddTraceQuadPhysToField(
1438  const Array<OneD, const NekDouble> &Fwd,
1439  const Array<OneD, const NekDouble> &Bwd,
1440  Array<OneD, NekDouble> &field);
1441 
1442  virtual void v_AddTraceQuadPhysToOffDiag(
1443  const Array<OneD, const NekDouble> &Fwd,
1444  const Array<OneD, const NekDouble> &Bwd,
1445  Array<OneD, NekDouble> &field);
1446 
1447  virtual void v_GetLocTraceFromTracePts(
1448  const Array<OneD, const NekDouble> &Fwd,
1449  const Array<OneD, const NekDouble> &Bwd,
1450  Array<OneD, NekDouble> &locTraceFwd,
1451  Array<OneD, NekDouble> &locTraceBwd);
1452 
1453  virtual void v_FillBwdWithBwdWeight(
1454  Array<OneD, NekDouble> &weightave,
1455  Array<OneD, NekDouble> &weightjmp);
1456 
1457  virtual void v_PeriodicBwdCopy(
1458  const Array<OneD, const NekDouble> &Fwd,
1459  Array<OneD, NekDouble> &Bwd);
1460 
1461  virtual const std::vector<bool> &v_GetLeftAdjacentFaces(void) const;
1462 
1463  virtual void v_ExtractTracePhys(
1464  Array<OneD,NekDouble> &outarray);
1465 
1466  virtual void v_ExtractTracePhys(
1467  const Array<OneD, const NekDouble> &inarray,
1468  Array<OneD,NekDouble> &outarray);
1469 
1470  virtual void v_MultiplyByInvMassMatrix(
1471  const Array<OneD,const NekDouble> &inarray,
1472  Array<OneD, NekDouble> &outarray);
1473 
1474  virtual void v_HelmSolve(
1475  const Array<OneD, const NekDouble> &inarray,
1476  Array<OneD, NekDouble> &outarray,
1477  const StdRegions::ConstFactorMap &factors,
1478  const StdRegions::VarCoeffMap &varcoeff,
1479  const MultiRegions::VarFactorsMap &varfactors,
1480  const Array<OneD, const NekDouble> &dirForcing,
1481  const bool PhysSpaceForcing);
1482 
1484  const Array<OneD, Array<OneD, NekDouble> > &velocity,
1485  const Array<OneD, const NekDouble> &inarray,
1486  Array<OneD, NekDouble> &outarray,
1487  const NekDouble lambda,
1489  dirForcing = NullNekDouble1DArray);
1490 
1491  virtual void v_LinearAdvectionReactionSolve(
1492  const Array<OneD, Array<OneD, NekDouble> > &velocity,
1493  const Array<OneD, const NekDouble> &inarray,
1494  Array<OneD, NekDouble> &outarray,
1495  const NekDouble lambda,
1497  dirForcing = NullNekDouble1DArray);
1498 
1499  // wrapper functions about virtual functions
1500  virtual void v_ImposeDirichletConditions(Array<OneD,NekDouble>& outarray);
1501 
1502  virtual void v_FillBndCondFromField();
1503 
1504  virtual void v_FillBndCondFromField(const int nreg);
1505 
1506  virtual void v_Reset();
1507 
1508  virtual void v_LocalToGlobal(bool UseComm);
1509 
1510  virtual void v_LocalToGlobal(
1511  const Array<OneD, const NekDouble> &inarray,
1512  Array<OneD,NekDouble> &outarray,
1513  bool UseComm);
1514 
1515  virtual void v_GlobalToLocal(void);
1516 
1517  virtual void v_GlobalToLocal(
1518  const Array<OneD, const NekDouble> &inarray,
1519  Array<OneD,NekDouble> &outarray);
1520 
1521  virtual void v_BwdTrans(
1522  const Array<OneD,const NekDouble> &inarray,
1523  Array<OneD, NekDouble> &outarray);
1524 
1525  virtual void v_BwdTrans_IterPerExp(
1526  const Array<OneD,const NekDouble> &inarray,
1527  Array<OneD,NekDouble> &outarray);
1528 
1529  virtual void v_FwdTrans(
1530  const Array<OneD,const NekDouble> &inarray,
1531  Array<OneD, NekDouble> &outarray);
1532 
1533  virtual void v_FwdTrans_IterPerExp(
1534  const Array<OneD,const NekDouble> &inarray,
1535  Array<OneD,NekDouble> &outarray);
1536 
1537  virtual void v_FwdTrans_BndConstrained(
1538  const Array<OneD,const NekDouble> &inarray,
1539  Array<OneD,NekDouble> &outarray);
1540 
1541  virtual void v_SmoothField(Array<OneD,NekDouble> &field);
1542 
1543  virtual void v_IProductWRTBase(
1544  const Array<OneD, const NekDouble> &inarray,
1545  Array<OneD, NekDouble> &outarray);
1546 
1547  virtual void v_IProductWRTBase_IterPerExp(
1548  const Array<OneD,const NekDouble> &inarray,
1549  Array<OneD, NekDouble> &outarray);
1550 
1551  virtual void v_GeneralMatrixOp(
1552  const GlobalMatrixKey &gkey,
1553  const Array<OneD,const NekDouble> &inarray,
1554  Array<OneD, NekDouble> &outarray);
1555 
1556  virtual void v_GetCoords(
1557  Array<OneD, NekDouble> &coord_0,
1558  Array<OneD, NekDouble> &coord_1,
1560 
1561  virtual void v_PhysDeriv(
1562  const Array<OneD, const NekDouble> &inarray,
1563  Array<OneD, NekDouble> &out_d0,
1564  Array<OneD, NekDouble> &out_d1,
1565  Array<OneD, NekDouble> &out_d2);
1566 
1567  virtual void v_PhysDeriv(
1568  const int dir,
1569  const Array<OneD, const NekDouble> &inarray,
1570  Array<OneD, NekDouble> &out_d);
1571 
1572  virtual void v_PhysDeriv(
1573  Direction edir,
1574  const Array<OneD, const NekDouble> &inarray,
1575  Array<OneD, NekDouble> &out_d);
1576 
1577  virtual void v_CurlCurl(
1580 
1581  virtual void v_PhysDirectionalDeriv(
1582  const Array<OneD, const NekDouble> &direction,
1583  const Array<OneD, const NekDouble> &inarray,
1584  Array<OneD, NekDouble> &outarray);
1585 
1586  virtual void v_GetMovingFrames(
1587  const SpatialDomains::GeomMMF MMFdir,
1588  const Array<OneD, const NekDouble> &CircCentre,
1589  Array<OneD, Array<OneD, NekDouble> > &outarray);
1590 
1591  virtual void v_HomogeneousFwdTrans(
1592  const Array<OneD, const NekDouble> &inarray,
1593  Array<OneD, NekDouble> &outarray,
1594  bool Shuff = true,
1595  bool UnShuff = true);
1596 
1597  virtual void v_HomogeneousBwdTrans(
1598  const Array<OneD, const NekDouble> &inarray,
1599  Array<OneD, NekDouble> &outarray,
1600  bool Shuff = true,
1601  bool UnShuff = true);
1602 
1603  virtual void v_DealiasedProd(
1604  const Array<OneD, NekDouble> &inarray1,
1605  const Array<OneD, NekDouble> &inarray2,
1606  Array<OneD, NekDouble> &outarray);
1607 
1608  virtual void v_DealiasedDotProd(
1609  const Array<OneD, Array<OneD, NekDouble> > &inarray1,
1610  const Array<OneD, Array<OneD, NekDouble> > &inarray2,
1611  Array<OneD, Array<OneD, NekDouble> > &outarray);
1612 
1613  virtual void v_GetBCValues(
1614  Array<OneD, NekDouble> &BndVals,
1615  const Array<OneD, NekDouble> &TotField,
1616  int BndID);
1617 
1618  virtual void v_NormVectorIProductWRTBase(
1621  Array<OneD, NekDouble> &outarray,
1622  int BndID);
1623 
1624  virtual void v_NormVectorIProductWRTBase(
1626  Array<OneD, NekDouble> &outarray);
1627 
1628  virtual void v_SetUpPhysNormals();
1629 
1630  virtual void v_GetBoundaryToElmtMap(Array<OneD, int> &ElmtID,
1631  Array<OneD,int> &EdgeID);
1632 
1633  virtual void v_GetBndElmtExpansion(int i,
1634  std::shared_ptr<ExpList> &result,
1635  const bool DeclareCoeffPhysArrays);
1636 
1637  virtual void v_ExtractElmtToBndPhys( const int i,
1638  const Array<OneD, NekDouble> & elmt,
1639  Array<OneD, NekDouble> & boundary );
1640 
1641  virtual void v_ExtractPhysToBndElmt( const int i,
1642  const Array<OneD, const NekDouble> & phys,
1643  Array<OneD, NekDouble> & bndElmt );
1644 
1645  virtual void v_ExtractPhysToBnd( const int i,
1646  const Array<OneD, const NekDouble> & phys,
1647  Array<OneD, NekDouble> & bnd );
1648 
1649  virtual void v_GetBoundaryNormals(int i,
1650  Array<OneD, Array<OneD, NekDouble> > &normals);
1651 
1652  virtual std::vector<LibUtilities::FieldDefinitionsSharedPtr>
1653  v_GetFieldDefinitions(void);
1654 
1655  virtual void v_GetFieldDefinitions(
1656  std::vector<LibUtilities::FieldDefinitionsSharedPtr> &fielddef);
1657 
1658  virtual void v_AppendFieldData(
1660  std::vector<NekDouble> &fielddata);
1661 
1662  virtual void v_AppendFieldData(
1664  std::vector<NekDouble> &fielddata,
1665  Array<OneD, NekDouble> &coeffs);
1666 
1667  virtual void v_ExtractDataToCoeffs(
1669  std::vector<NekDouble> &fielddata, std::string &field,
1670  Array<OneD, NekDouble> &coeffs);
1671 
1672  virtual void v_ExtractCoeffsToCoeffs(const std::shared_ptr<ExpList> &fromExpList, const Array<OneD, const NekDouble> &fromCoeffs, Array<OneD, NekDouble> &toCoeffs);
1673 
1674  virtual void v_WriteTecplotHeader(std::ostream &outfile,
1675  std::string var = "");
1676  virtual void v_WriteTecplotZone(std::ostream &outfile,
1677  int expansion);
1678  virtual void v_WriteTecplotField(std::ostream &outfile,
1679  int expansion);
1680  virtual void v_WriteTecplotConnectivity(std::ostream &outfile,
1681  int expansion);
1682  virtual void v_WriteVtkPieceHeader(
1683  std::ostream &outfile,
1684  int expansion,
1685  int istrip);
1686 
1687  virtual void v_WriteVtkPieceData(
1688  std::ostream &outfile,
1689  int expansion,
1690  std::string var);
1691 
1692  virtual NekDouble v_L2(
1693  const Array<OneD, const NekDouble> &phys,
1695 
1696  virtual NekDouble v_Integral (
1697  const Array<OneD, const NekDouble> &inarray);
1698  virtual NekDouble v_VectorFlux (
1699  const Array<OneD, Array<OneD, NekDouble> > &inarray);
1700 
1703  virtual NekDouble v_GetHomoLen(void);
1704  virtual void v_SetHomoLen(const NekDouble lhom);
1707 
1708  // 1D Scaling and projection
1709  virtual void v_PhysInterp1DScaled(
1710  const NekDouble scale, const Array<OneD, NekDouble> &inarray,
1711  Array<OneD, NekDouble> &outarray);
1712 
1713  virtual void v_PhysGalerkinProjection1DScaled(
1714  const NekDouble scale,
1715  const Array<OneD, NekDouble> &inarray,
1716  Array<OneD, NekDouble> &outarray);
1717 
1718  virtual void v_ClearGlobalLinSysManager(void);
1719 
1720  void ExtractFileBCs(const std::string &fileName,
1722  const std::string &varName,
1723  const std::shared_ptr<ExpList> locExpList);
1724 
1725  // Utility function for a common case of retrieving a
1726  // BoundaryCondition from a boundary condition collection.
1729  GetBoundaryCondition(const SpatialDomains::
1730  BoundaryConditionCollection& collection,
1731  unsigned int index, const std::string& variable);
1732 
1733 
1734  private:
1735  /// Definition of the total number of degrees of freedom and
1736  /// quadrature points and offsets to access data
1737  void SetupCoeffPhys(bool DeclareCoeffPhysArrays = true,
1738  bool SetupOffsets = true);
1739 
1740  /// Define a list of elements using the geometry and basis
1741  /// key information in expmap;
1743 
1744  virtual const Array<OneD,
1746 
1749 
1750  virtual void v_EvaluateBoundaryConditions(
1751  const NekDouble time = 0.0,
1752  const std::string varName = "",
1755 
1756  virtual std::map<int, RobinBCInfoSharedPtr> v_GetRobinBCInfo(void);
1757 
1758 
1759  virtual void v_GetPeriodicEntities(
1760  PeriodicMap &periodicVerts,
1761  PeriodicMap &periodicEdges,
1762  PeriodicMap &periodicFaces);
1763 
1764  // Homogeneous direction wrapper functions.
1766  {
1768  "This method is not defined or valid for this class type");
1770  }
1771 
1772  // wrapper function to set viscosity for Homo1D expansion
1774  {
1775  boost::ignore_unused(visc);
1777  "This method is not defined or valid for this class type");
1778  }
1779 
1780 
1781  virtual std::shared_ptr<ExpList> &v_GetPlane(int n);
1782 
1783  virtual void v_AddTraceIntegralToOffDiag(
1784  const Array<OneD, const NekDouble> &FwdFlux,
1785  const Array<OneD, const NekDouble> &BwdFlux,
1786  Array<OneD, NekDouble> &outarray);
1787 
1788  };
1789 
1790  /// An empty ExpList object.
1793 
1794  // Inline routines follow.
1795 
1796  /**
1797  * Returns the total number of local degrees of freedom
1798  * \f$N_{\mathrm{eof}}=\sum_{e=1}^{{N_{\mathrm{el}}}}N^{e}_m\f$.
1799  */
1800  inline int ExpList::GetNcoeffs() const
1801  {
1802  return m_ncoeffs;
1803  }
1804 
1805  inline int ExpList::GetNcoeffs(const int eid) const
1806  {
1807  return (*m_exp)[eid]->GetNcoeffs();
1808  }
1809 
1810  /**
1811  * Evaulates the maximum number of modes in the elemental basis
1812  * order over all elements
1813  */
1815  {
1816  unsigned int i;
1817  int returnval = 0;
1818 
1819  for(i= 0; i < (*m_exp).size(); ++i)
1820  {
1821  returnval = (std::max)(returnval,
1822  (*m_exp)[i]->EvalBasisNumModesMax());
1823  }
1824 
1825  return returnval;
1826  }
1827 
1828  /**
1829  *
1830  */
1832  const
1833  {
1834  unsigned int i;
1835  Array<OneD,int> returnval((*m_exp).size(),0);
1836 
1837  for(i= 0; i < (*m_exp).size(); ++i)
1838  {
1839  returnval[i]
1840  = (std::max)(returnval[i],(*m_exp)[i]->EvalBasisNumModesMax());
1841  }
1842 
1843  return returnval;
1844  }
1845 
1846 
1847  /**
1848  *
1849  */
1850  inline int ExpList::GetTotPoints() const
1851  {
1852  return m_npoints;
1853  }
1854 
1855  inline int ExpList::GetTotPoints(const int eid) const
1856  {
1857  return (*m_exp)[eid]->GetTotPoints();
1858  }
1859 
1860 
1861  inline int ExpList::Get1DScaledTotPoints(const NekDouble scale) const
1862  {
1863  int returnval = 0;
1864  int cnt;
1865  int nbase = (*m_exp)[0]->GetNumBases();
1866 
1867  for(int i = 0; i < (*m_exp).size(); ++i)
1868  {
1869  cnt = 1;
1870  for(int j = 0; j < nbase; ++j)
1871  {
1872  cnt *= (int)(scale*((*m_exp)[i]->GetNumPoints(j)));
1873  }
1874  returnval += cnt;
1875  }
1876  return returnval;
1877  }
1878 
1879  /**
1880  *
1881  */
1882  inline int ExpList::GetNpoints() const
1883  {
1884  return m_npoints;
1885  }
1886 
1887 
1888  /**
1889  *
1890  */
1891  inline void ExpList::SetWaveSpace(const bool wavespace)
1892  {
1893  m_WaveSpace = wavespace;
1894  }
1895 
1896  /**
1897  *
1898  */
1899  inline bool ExpList::GetWaveSpace() const
1900  {
1901  return m_WaveSpace;
1902  }
1903 
1904  /// Set the \a i th value of\a m_phys to value \a val
1905  inline void ExpList::SetPhys(int i, NekDouble val)
1906  {
1907  m_phys[i] = val;
1908  }
1909 
1910 
1911  /**
1912  * This function fills the array \f$\boldsymbol{u}_l\f$, the evaluation
1913  * of the expansion at the quadrature points (implemented as #m_phys),
1914  * with the values of the array \a inarray.
1915  *
1916  * @param inarray The array containing the values where
1917  * #m_phys should be filled with.
1918  */
1919  inline void ExpList::SetPhys(
1920  const Array<OneD, const NekDouble> &inarray)
1921  {
1922  ASSERTL0(inarray.size() == m_npoints,
1923  "Input array does not have correct number of elements.");
1924 
1925  Vmath::Vcopy(m_npoints,&inarray[0],1,&m_phys[0],1);
1926  m_physState = true;
1927  }
1928 
1929 
1931  {
1932  m_phys = inarray;
1933  }
1934 
1935 
1936  /**
1937  * @param physState \a true (=filled) or \a false (=not filled).
1938  */
1939  inline void ExpList::SetPhysState(const bool physState)
1940  {
1941  m_physState = physState;
1942  }
1943 
1944 
1945  /**
1946  * @return physState \a true (=filled) or \a false (=not filled).
1947  */
1948  inline bool ExpList::GetPhysState() const
1949  {
1950  return m_physState;
1951  }
1952 
1953  /**
1954  *
1955  */
1957  const Array<OneD, const NekDouble> &inarray,
1958  Array<OneD, NekDouble> &outarray)
1959  {
1960  v_IProductWRTBase(inarray,outarray);
1961  }
1962 
1963  /**
1964  *
1965  */
1967  const Array<OneD, const NekDouble> &inarray,
1968  Array<OneD, NekDouble> &outarray)
1969  {
1970  v_IProductWRTBase_IterPerExp(inarray,outarray);
1971  }
1972 
1973  /**
1974  *
1975  */
1976  inline void ExpList::FwdTrans(
1977  const Array<OneD, const NekDouble> &inarray,
1978  Array<OneD, NekDouble> &outarray)
1979  {
1980  v_FwdTrans(inarray,outarray);
1981  }
1982 
1983  /**
1984  *
1985  */
1987  const Array<OneD, const NekDouble> &inarray,
1988  Array<OneD,NekDouble> &outarray)
1989  {
1990  v_FwdTrans_IterPerExp(inarray,outarray);
1991  }
1992 
1993  /**
1994  *
1995  */
1997  const Array<OneD, const NekDouble> &inarray,
1998  Array<OneD,NekDouble> &outarray)
1999  {
2000  v_FwdTrans_BndConstrained(inarray,outarray);
2001  }
2002 
2003 
2004  /**
2005  *
2006  */
2008  {
2009  v_SmoothField(field);
2010  }
2011 
2012  /**
2013  *
2014  */
2015  inline void ExpList::BwdTrans (
2016  const Array<OneD, const NekDouble> &inarray,
2017  Array<OneD, NekDouble> &outarray)
2018  {
2019  v_BwdTrans(inarray,outarray);
2020  }
2021 
2022  /**
2023  *
2024  */
2026  const Array<OneD, const NekDouble> &inarray,
2027  Array<OneD, NekDouble> &outarray)
2028  {
2029  v_BwdTrans_IterPerExp(inarray,outarray);
2030  }
2031 
2032 
2033  /**
2034  *
2035  */
2037  const Array<OneD,const NekDouble> &inarray,
2038  Array<OneD, NekDouble> &outarray)
2039  {
2040  v_MultiplyByInvMassMatrix(inarray,outarray);
2041  }
2042 
2043  /**
2044  *
2045  */
2046  inline void ExpList::HelmSolve(
2047  const Array<OneD, const NekDouble> &inarray,
2048  Array<OneD, NekDouble> &outarray,
2049  const StdRegions::ConstFactorMap &factors,
2050  const StdRegions::VarCoeffMap &varcoeff,
2051  const MultiRegions::VarFactorsMap &varfactors,
2052  const Array<OneD, const NekDouble> &dirForcing,
2053  const bool PhysSpaceForcing)
2054 
2055  {
2056  v_HelmSolve(inarray, outarray, factors, varcoeff,
2057  varfactors, dirForcing, PhysSpaceForcing);
2058  }
2059 
2060 
2061  /**
2062  *
2063  */
2065  const Array<OneD, Array<OneD, NekDouble> > &velocity,
2066  const Array<OneD, const NekDouble> &inarray,
2067  Array<OneD, NekDouble> &outarray,
2068  const NekDouble lambda,
2069  const Array<OneD, const NekDouble>& dirForcing)
2070  {
2071  v_LinearAdvectionDiffusionReactionSolve(velocity,inarray, outarray,
2072  lambda, dirForcing);
2073  }
2074 
2076  const Array<OneD, Array<OneD, NekDouble> > &velocity,
2077  const Array<OneD, const NekDouble> &inarray,
2078  Array<OneD, NekDouble> &outarray,
2079  const NekDouble lambda,
2080  const Array<OneD, const NekDouble>& dirForcing)
2081  {
2082  v_LinearAdvectionReactionSolve(velocity,inarray, outarray,
2083  lambda, dirForcing);
2084  }
2085 
2086  /**
2087  *
2088  */
2090  Array<OneD, NekDouble> &coord_1,
2091  Array<OneD, NekDouble> &coord_2)
2092 
2093  {
2094  v_GetCoords(coord_0,coord_1,coord_2);
2095  }
2096 
2097 
2098  /**
2099  *
2100  */
2102  const SpatialDomains::GeomMMF MMFdir,
2103  const Array<OneD, const NekDouble> &CircCentre,
2104  Array<OneD, Array<OneD, NekDouble> > &outarray)
2105  {
2106  v_GetMovingFrames(MMFdir,CircCentre,outarray);
2107  }
2108 
2109 
2110  /**
2111  *
2112  */
2114  Array<OneD, NekDouble> &out_d0,
2115  Array<OneD, NekDouble> &out_d1,
2116  Array<OneD, NekDouble> &out_d2)
2117  {
2118  v_PhysDeriv(inarray,out_d0,out_d1,out_d2);
2119  }
2120 
2121  /**
2122  *
2123  */
2124  inline void ExpList::PhysDeriv(
2125  const int dir,
2126  const Array<OneD, const NekDouble> &inarray,
2127  Array<OneD, NekDouble> &out_d)
2128  {
2129  v_PhysDeriv(dir,inarray,out_d);
2130  }
2131 
2132  inline void ExpList::PhysDeriv(
2133  Direction edir,
2134  const Array<OneD, const NekDouble> &inarray,
2135  Array<OneD, NekDouble> &out_d)
2136  {
2137  v_PhysDeriv(edir, inarray,out_d);
2138  }
2139 
2140 
2141  /**
2142  *
2143  */
2145  const Array<OneD, const NekDouble> &direction,
2146  const Array<OneD, const NekDouble> &inarray,
2147  Array<OneD, NekDouble> &outarray)
2148  {
2149  v_PhysDirectionalDeriv(direction, inarray, outarray);
2150  }
2151 
2152 
2153  /**
2154  *
2155  */
2156 
2157  inline void ExpList::CurlCurl(
2160  {
2161  v_CurlCurl(Vel, Q);
2162  }
2163 
2164  /**
2165  *
2166  */
2168  const Array<OneD, const NekDouble> &inarray,
2169  Array<OneD, NekDouble> &outarray,
2170  bool Shuff,
2171  bool UnShuff)
2172  {
2173  v_HomogeneousFwdTrans(inarray,outarray,Shuff,UnShuff);
2174  }
2175 
2176  /**
2177  *
2178  */
2180  const Array<OneD, const NekDouble> &inarray,
2181  Array<OneD, NekDouble> &outarray,
2182  bool Shuff,
2183  bool UnShuff)
2184  {
2185  v_HomogeneousBwdTrans(inarray,outarray,Shuff,UnShuff);
2186  }
2187 
2188  /**
2189  *
2190  */
2192  const Array<OneD, NekDouble> &inarray1,
2193  const Array<OneD, NekDouble> &inarray2,
2194  Array<OneD, NekDouble> &outarray)
2195  {
2196  v_DealiasedProd(inarray1,inarray2,outarray);
2197  }
2198 
2199  /**
2200  *
2201  */
2203  const Array<OneD, Array<OneD, NekDouble> > &inarray1,
2204  const Array<OneD, Array<OneD, NekDouble> > &inarray2,
2205  Array<OneD, Array<OneD, NekDouble> > &outarray)
2206  {
2207  v_DealiasedDotProd(inarray1,inarray2,outarray);
2208  }
2209 
2210  /**
2211  *
2212  */
2214  Array<OneD, NekDouble> &BndVals,
2215  const Array<OneD, NekDouble> &TotField,
2216  int BndID)
2217  {
2218  v_GetBCValues(BndVals,TotField,BndID);
2219  }
2220 
2221  /**
2222  *
2223  */
2227  Array<OneD, NekDouble> &outarray,
2228  int BndID)
2229  {
2230  v_NormVectorIProductWRTBase(V1,V2,outarray,BndID);
2231  }
2232 
2235  Array<OneD, NekDouble> &outarray)
2236  {
2237  v_NormVectorIProductWRTBase(V, outarray);
2238  }
2239 
2240  /**
2241  * @param eid The index of the element to be checked.
2242  * @return The dimension of the coordinates of the specific element.
2243  */
2244  inline int ExpList::GetCoordim(int eid)
2245  {
2246  ASSERTL2(eid <= (*m_exp).size(),
2247  "eid is larger than number of elements");
2248  return (*m_exp)[eid]->GetCoordim();
2249  }
2250 
2251  /**
2252  * @param eid The index of the element to be checked.
2253  * @return The dimension of the shape of the specific element.
2254  */
2256  {
2257  return (*m_exp)[0]->GetShapeDimension();
2258  }
2259 
2260  /**
2261  * @param i The index of m_coeffs to be set
2262  * @param val The value which m_coeffs[i] is to be set to.
2263  */
2264  inline void ExpList::SetCoeff(int i, NekDouble val)
2265  {
2266  m_coeffs[i] = val;
2267  }
2268 
2269 
2270  /**
2271  * @param i The index of #m_coeffs to be set.
2272  * @param val The value which #m_coeffs[i] is to be set to.
2273  */
2274  inline void ExpList::SetCoeffs(int i, NekDouble val)
2275  {
2276  m_coeffs[i] = val;
2277  }
2278 
2279 
2281  {
2282  m_coeffs = inarray;
2283  }
2284 
2285  /**
2286  * As the function returns a constant reference to a
2287  * <em>const Array</em>, it is not possible to modify the
2288  * underlying data of the array #m_coeffs. In order to do
2289  * so, use the function #UpdateCoeffs instead.
2290  *
2291  * @return (A constant reference to) the array #m_coeffs.
2292  */
2294  {
2295  return m_coeffs;
2296  }
2297 
2299  Array<OneD,NekDouble>& outarray)
2300  {
2301  v_ImposeDirichletConditions(outarray);
2302  }
2303 
2305  {
2307  }
2308 
2309  inline void ExpList::FillBndCondFromField(const int nreg)
2310  {
2311  v_FillBndCondFromField(nreg);
2312  }
2313 
2314  inline void ExpList::LocalToGlobal(bool useComm)
2315  {
2316  v_LocalToGlobal(useComm);
2317  }
2318 
2320  const Array<OneD, const NekDouble> &inarray,
2321  Array<OneD,NekDouble> &outarray,
2322  bool useComm)
2323  {
2324  v_LocalToGlobal(inarray, outarray,useComm);
2325  }
2326 
2327  inline void ExpList::GlobalToLocal(void)
2328  {
2329  v_GlobalToLocal();
2330  }
2331 
2332  /**
2333  * This operation is evaluated as:
2334  * \f{tabbing}
2335  * \hspace{1cm} \= Do \= $e=$ $1, N_{\mathrm{el}}$ \\
2336  * > > Do \= $i=$ $0,N_m^e-1$ \\
2337  * > > > $\boldsymbol{\hat{u}}^{e}[i] = \mbox{sign}[e][i] \cdot
2338  * \boldsymbol{\hat{u}}_g[\mbox{map}[e][i]]$ \\
2339  * > > continue \\
2340  * > continue
2341  * \f}
2342  * where \a map\f$[e][i]\f$ is the mapping array and \a
2343  * sign\f$[e][i]\f$ is an array of similar dimensions ensuring the
2344  * correct modal connectivity between the different elements (both
2345  * these arrays are contained in the data member #m_locToGloMap). This
2346  * operation is equivalent to the scatter operation
2347  * \f$\boldsymbol{\hat{u}}_l=\mathcal{A}\boldsymbol{\hat{u}}_g\f$, where
2348  * \f$\mathcal{A}\f$ is the
2349  * \f$N_{\mathrm{eof}}\times N_{\mathrm{dof}}\f$ permutation matrix.
2350  *
2351  * @param inarray An array of size \f$N_\mathrm{dof}\f$
2352  * containing the global degrees of freedom
2353  * \f$\boldsymbol{x}_g\f$.
2354  * @param outarray The resulting local degrees of freedom
2355  * \f$\boldsymbol{x}_l\f$ will be stored in this
2356  * array of size \f$N_\mathrm{eof}\f$.
2357  */
2359  const Array<OneD, const NekDouble> &inarray,
2360  Array<OneD,NekDouble> &outarray)
2361  {
2362  v_GlobalToLocal(inarray, outarray);
2363  }
2364 
2365 
2366  /**
2367  * @param i The index of #m_coeffs to be returned
2368  * @return The NekDouble held in #m_coeffs[i].
2369  */
2371  {
2372  return m_coeffs[i];
2373  }
2374 
2375  /**
2376  * @param i The index of #m_coeffs to be returned
2377  * @return The NekDouble held in #m_coeffs[i].
2378  */
2380  {
2381  return m_coeffs[i];
2382  }
2383 
2384  /**
2385  * As the function returns a constant reference to a
2386  * <em>const Array</em> it is not possible to modify the
2387  * underlying data of the array #m_phys. In order to do
2388  * so, use the function #UpdatePhys instead.
2389  *
2390  * @return (A constant reference to) the array #m_phys.
2391  */
2393  {
2394  return m_phys;
2395  }
2396 
2397  /**
2398  * @return \f$N_{\mathrm{el}}\f$, the number of elements in the
2399  * expansion.
2400  */
2401  inline int ExpList::GetExpSize(void)
2402  {
2403  return (*m_exp).size();
2404  }
2405 
2406 
2407  /**
2408  * @param n The index of the element concerned.
2409  *
2410  * @return (A shared pointer to) the local expansion of the
2411  * \f$n^{\mathrm{th}}\f$ element.
2412  */
2414  {
2415  return (*m_exp)[n];
2416  }
2417 
2418  /**
2419  * @return (A const shared pointer to) the local expansion vector #m_exp
2420  */
2421  inline const std::shared_ptr<LocalRegions::ExpansionVector>
2422  ExpList::GetExp(void) const
2423  {
2424  return m_exp;
2425  }
2426 
2427 
2428  /**
2429  *
2430  */
2431  inline int ExpList::GetCoeff_Offset(int n) const
2432  {
2433  return m_coeff_offset[n];
2434  }
2435 
2436  /**
2437  *
2438  */
2439  inline int ExpList::GetPhys_Offset(int n) const
2440  {
2441  return m_phys_offset[n];
2442  }
2443 
2444  /**
2445  * If one wants to get hold of the underlying data without modifying
2446  * them, rather use the function #GetCoeffs instead.
2447  *
2448  * @return (A reference to) the array #m_coeffs.
2449  */
2451  {
2452  return m_coeffs;
2453  }
2454 
2455  /**
2456  * If one wants to get hold of the underlying data without modifying
2457  * them, rather use the function #GetPhys instead.
2458  *
2459  * @return (A reference to) the array #m_phys.
2460  */
2462  {
2463  m_physState = true;
2464  return m_phys;
2465  }
2466 
2467 
2468  // functions associated with DisContField
2471  {
2472  return v_GetBndCondExpansions();
2473  }
2474 
2475  /// Get m_coeffs to elemental value map
2476  MULTI_REGIONS_EXPORT inline const
2479  {
2480  return m_coeffsToElmt;
2481  }
2482 
2485  {
2486  return v_GetLocTraceToTraceMap();
2487  }
2488 
2489  inline const Array<OneD, const NekDouble >
2491  {
2492  return v_GetBndCondBwdWeight();
2493  }
2494 
2496  const int index,
2497  const NekDouble value)
2498  {
2499  v_SetBndCondBwdWeight(index, value);
2500  }
2501 
2502  inline std::shared_ptr<ExpList> &ExpList::UpdateBndCondExpansion(int i)
2503  {
2504  return v_UpdateBndCondExpansion(i);
2505  }
2506 
2507  inline void ExpList::Upwind(
2508  const Array<OneD, const Array<OneD, NekDouble> > &Vec,
2509  const Array<OneD, const NekDouble> &Fwd,
2510  const Array<OneD, const NekDouble> &Bwd,
2511  Array<OneD, NekDouble> &Upwind)
2512  {
2513  v_Upwind(Vec, Fwd, Bwd, Upwind);
2514  }
2515 
2516  inline void ExpList::Upwind(
2517  const Array<OneD, const NekDouble> &Vn,
2518  const Array<OneD, const NekDouble> &Fwd,
2519  const Array<OneD, const NekDouble> &Bwd,
2520  Array<OneD, NekDouble> &Upwind)
2521  {
2522  v_Upwind(Vn, Fwd, Bwd, Upwind);
2523  }
2524 
2525  inline std::shared_ptr<ExpList> &ExpList::GetTrace()
2526  {
2527  return v_GetTrace();
2528  }
2529 
2530  inline std::shared_ptr<AssemblyMapDG> &ExpList::GetTraceMap()
2531  {
2532  return v_GetTraceMap();
2533  }
2534 
2536  {
2537  return v_GetTraceBndMap();
2538  }
2539 
2540  inline void ExpList::GetNormals(
2541  Array<OneD, Array<OneD, NekDouble> > &normals)
2542  {
2543  v_GetNormals(normals);
2544  }
2545 
2547  const Array<OneD, const NekDouble> &Fx,
2548  const Array<OneD, const NekDouble> &Fy,
2549  Array<OneD, NekDouble> &outarray)
2550  {
2551  v_AddTraceIntegral(Fx,Fy,outarray);
2552  }
2553 
2555  const Array<OneD, const NekDouble> &Fn,
2556  Array<OneD, NekDouble> &outarray)
2557  {
2558  v_AddTraceIntegral(Fn,outarray);
2559  }
2560 
2562  const Array<OneD, const NekDouble> &Fwd,
2563  const Array<OneD, const NekDouble> &Bwd,
2564  Array<OneD, NekDouble> &outarray)
2565  {
2566  v_AddFwdBwdTraceIntegral(Fwd,Bwd,outarray);
2567  }
2568 
2570  Array<OneD,NekDouble> &Fwd,
2571  Array<OneD,NekDouble> &Bwd)
2572  {
2573  v_GetFwdBwdTracePhys(Fwd,Bwd);
2574  }
2575 
2577  const Array<OneD,const NekDouble> &field,
2578  Array<OneD,NekDouble> &Fwd,
2579  Array<OneD,NekDouble> &Bwd,
2580  bool FillBnd,
2581  bool PutFwdInBwdOnBCs,
2582  bool DoExchange)
2583  {
2584  v_GetFwdBwdTracePhys(field,Fwd,Bwd,FillBnd,
2585  PutFwdInBwdOnBCs,DoExchange);
2586  }
2587 
2589  const Array<OneD,NekDouble> &Fwd,
2590  Array<OneD,NekDouble> &Bwd,
2591  bool PutFwdInBwdOnBCs)
2592  {
2593  v_FillBwdWithBoundCond(Fwd, Bwd, PutFwdInBwdOnBCs);
2594  }
2595 
2597  const Array<OneD, const NekDouble> &Fwd,
2598  const Array<OneD, const NekDouble> &Bwd,
2599  Array<OneD, NekDouble> &field)
2600  {
2601  v_AddTraceQuadPhysToField(Fwd, Bwd, field);
2602  }
2603 
2605  const Array<OneD, const NekDouble> &Fwd,
2606  const Array<OneD, const NekDouble> &Bwd,
2607  Array<OneD, NekDouble> &locTraceFwd,
2608  Array<OneD, NekDouble> &locTraceBwd)
2609  {
2610  v_GetLocTraceFromTracePts(Fwd,Bwd,locTraceFwd,locTraceBwd);
2611  }
2612 
2614  Array<OneD, NekDouble> &weightave,
2615  Array<OneD, NekDouble> &weightjmp)
2616  {
2617  v_FillBwdWithBwdWeight(weightave,weightjmp);
2618  }
2619 
2621  const Array<OneD, const NekDouble> &Fwd,
2623  {
2624  v_PeriodicBwdCopy(Fwd, Bwd);
2625  }
2626 
2627  inline const std::vector<bool> &ExpList::GetLeftAdjacentFaces(void) const
2628  {
2629  return v_GetLeftAdjacentFaces();
2630  }
2631 
2633  {
2634  v_ExtractTracePhys(outarray);
2635  }
2636 
2637 
2639  const Array<OneD, const NekDouble> &inarray,
2640  Array<OneD,NekDouble> &outarray)
2641  {
2642  v_ExtractTracePhys(inarray,outarray);
2643  }
2644 
2647  {
2648  return v_GetBndConditions();
2649  }
2650 
2653  {
2654  return v_UpdateBndConditions();
2655  }
2656 
2658  const NekDouble time,
2659  const std::string varName,
2660  const NekDouble x2_in,
2661  const NekDouble x3_in)
2662  {
2663  v_EvaluateBoundaryConditions(time, varName, x2_in, x3_in);
2664  }
2665 
2666  // Routines for continous matrix solution
2667  /**
2668  * This operation is equivalent to the evaluation of
2669  * \f$\underline{\boldsymbol{M}}^e\boldsymbol{\hat{u}}_l\f$, that is,
2670  * \f[ \left[
2671  * \begin{array}{cccc}
2672  * \boldsymbol{M}^1 & 0 & \hspace{3mm}0 \hspace{3mm}& 0 \\
2673  * 0 & \boldsymbol{M}^2 & 0 & 0 \\
2674  * 0 & 0 & \ddots & 0 \\
2675  * 0 & 0 & 0 & \boldsymbol{M}^{N_{\mathrm{el}}} \end{array} \right]
2676  *\left [ \begin{array}{c}
2677  * \boldsymbol{\hat{u}}^{1} \\
2678  * \boldsymbol{\hat{u}}^{2} \\
2679  * \vdots \\
2680  * \boldsymbol{\hat{u}}^{{{N_{\mathrm{el}}}}} \end{array} \right ]\f]
2681  * where \f$\boldsymbol{M}^e\f$ are the local matrices of type
2682  * specified by the key \a mkey. The decoupling of the local matrices
2683  * allows for a local evaluation of the operation. However, rather than
2684  * a local matrix-vector multiplication, the local operations are
2685  * evaluated as implemented in the function
2686  * StdRegions#StdExpansion#GeneralMatrixOp.
2687  *
2688  * @param mkey This key uniquely defines the type matrix
2689  * required for the operation.
2690  * @param inarray The vector \f$\boldsymbol{\hat{u}}_l\f$ of
2691  * size \f$N_{\mathrm{eof}}\f$.
2692  * @param outarray The resulting vector of size
2693  * \f$N_{\mathrm{eof}}\f$.
2694  */
2696  const GlobalMatrixKey &gkey,
2697  const Array<OneD,const NekDouble> &inarray,
2698  Array<OneD, NekDouble> &outarray)
2699  {
2700  v_GeneralMatrixOp(gkey,inarray,outarray);
2701  }
2702 
2703 
2705  {
2707  }
2708 
2710  Array<OneD,int> &EdgeID)
2711  {
2712  v_GetBoundaryToElmtMap(ElmtID,EdgeID);
2713  }
2714 
2715  inline void ExpList::GetBndElmtExpansion(int i,
2716  std::shared_ptr<ExpList> &result,
2717  const bool DeclareCoeffPhysArrays)
2718  {
2719  v_GetBndElmtExpansion(i, result, DeclareCoeffPhysArrays);
2720  }
2721 
2723  const Array<OneD, NekDouble> &elmt,
2724  Array<OneD, NekDouble> &boundary)
2725  {
2726  v_ExtractElmtToBndPhys(i, elmt, boundary);
2727  }
2728 
2730  const Array<OneD, const NekDouble> &phys,
2731  Array<OneD, NekDouble> &bndElmt)
2732  {
2733  v_ExtractPhysToBndElmt(i, phys, bndElmt);
2734  }
2735 
2736  inline void ExpList::ExtractPhysToBnd(int i,
2737  const Array<OneD, const NekDouble> &phys,
2739  {
2740  v_ExtractPhysToBnd(i, phys, bnd);
2741  }
2742 
2743  inline void ExpList::GetBoundaryNormals(int i,
2744  Array<OneD, Array<OneD, NekDouble> > &normals)
2745  {
2746  v_GetBoundaryNormals(i, normals);
2747  }
2748 
2750 
2751  } //end of namespace
2752 } //end of namespace
2753 
2754 #endif // EXPLIST_H
2755 
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
#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:250
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
Definition: ErrorUtil.hpp:274
#define MULTI_REGIONS_EXPORT
Base class for all multi-elemental spectral/hp expansions.
Definition: ExpList.h:107
virtual void v_Reset()
Reset geometry information, metrics, matrix managers and geometry information.
Definition: ExpList.cpp:2580
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Definition: ExpList.h:1252
Array< OneD, NekDouble > & UpdateCoeffs()
This function returns (a reference to) the array (implemented as m_coeffs) containing all local expa...
Definition: ExpList.h:2450
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:5477
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:2075
void IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:1956
const Array< OneD, const std::shared_ptr< ExpList > > & GetBndCondExpansions()
Definition: ExpList.h:2470
void IProductWRTBase_IterPerExp(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:1966
std::shared_ptr< ExpList > & GetTrace()
Definition: ExpList.h:2525
void WriteTecplotConnectivity(std::ostream &outfile, int expansion=-1)
Definition: ExpList.h:464
int GetNcoeffs(void) const
Returns the total number of local degrees of freedom .
Definition: ExpList.h:1800
virtual void v_ExtractCoeffsToCoeffs(const std::shared_ptr< ExpList > &fromExpList, const Array< OneD, const NekDouble > &fromCoeffs, Array< OneD, NekDouble > &toCoeffs)
Definition: ExpList.cpp:3546
const DNekScalBlkMatSharedPtr & GetBlockMatrix(const GlobalMatrixKey &gkey)
Definition: ExpList.cpp:1946
virtual void v_AddFwdBwdTraceIntegral(const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:4283
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:2274
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:3435
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:5017
virtual void v_FwdTrans_BndConstrained(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:1745
virtual void v_GetBoundaryToElmtMap(Array< OneD, int > &ElmtID, Array< OneD, int > &EdgeID)
Definition: ExpList.cpp:4902
void ExponentialFilter(Array< OneD, NekDouble > &array, const NekDouble alpha, const NekDouble exponent, const NekDouble cutoff)
Definition: ExpList.cpp:1670
virtual void v_HomogeneousBwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true)
Definition: ExpList.cpp:4456
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:2034
virtual void v_HomogeneousFwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true)
Definition: ExpList.cpp:4446
virtual void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:4620
void MultiplyByMassMatrix(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:330
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:2588
void AddRightIPTPhysDerivBase(const int dir, const Array< OneD, const DNekMatSharedPtr > ElmtJacQuad, Array< OneD, DNekMatSharedPtr > ElmtJacCoef)
Definition: ExpList.cpp:5718
virtual void v_ExtractTracePhys(Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:4383
void ExtractPhysToBndElmt(int i, const Array< OneD, const NekDouble > &phys, Array< OneD, NekDouble > &bndElmt)
Definition: ExpList.h:2729
void Reset()
Reset geometry information and reset matrices.
Definition: ExpList.h:440
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:2418
void GetFieldDefinitions(std::vector< LibUtilities::FieldDefinitionsSharedPtr > &fielddef)
Definition: ExpList.h:1067
virtual void v_WriteTecplotConnectivity(std::ostream &outfile, int expansion)
Definition: ExpList.cpp:2741
BlockMatrixMapShPtr m_blockMat
Definition: ExpList.h:1312
virtual Array< OneD, const unsigned int > v_GetYIDs(void)
Definition: ExpList.cpp:3192
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
multiply the metric jacobi and quadrature weights
Definition: ExpList.cpp:1220
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:4434
void SetPhys(int i, NekDouble val)
Set the i th value of m_phys to value val.
Definition: ExpList.h:1905
void ExtractTracePhys(Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:2632
void GetBoundaryToElmtMap(Array< OneD, int > &ElmtID, Array< OneD, int > &EdgeID)
Definition: ExpList.h:2709
void PeriodicBwdCopy(const Array< OneD, const NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
Copy and fill the Periodic boundaries.
Definition: ExpList.h:2620
void WriteVtkPieceFooter(std::ostream &outfile, int expansion)
Definition: ExpList.cpp:2991
virtual void v_ExtractPhysToBnd(const int i, const Array< OneD, const NekDouble > &phys, Array< OneD, NekDouble > &bnd)
Definition: ExpList.cpp:4819
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:3572
virtual void v_AddTraceIntegralToOffDiag(const Array< OneD, const NekDouble > &FwdFlux, const Array< OneD, const NekDouble > &BwdFlux, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:5244
void InitialiseExpVector(const SpatialDomains::ExpansionInfoMap &expmap)
Define a list of elements using the geometry and basis key information in expmap;.
Definition: ExpList.cpp:1022
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:978
std::shared_ptr< ExpList > & GetPlane(int n)
Definition: ExpList.h:1157
virtual LibUtilities::BasisSharedPtr v_GetHomogeneousBasis(void)
Definition: ExpList.h:1765
static SpatialDomains::BoundaryConditionShPtr GetBoundaryCondition(const SpatialDomains::BoundaryConditionCollection &collection, unsigned int index, const std::string &variable)
Definition: ExpList.cpp:4983
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:2293
void GetBndElmtExpansion(int i, std::shared_ptr< ExpList > &result, const bool DeclareCoeffPhysArrays=true)
Definition: ExpList.h:2715
void SetPhysArray(Array< OneD, NekDouble > &inarray)
Sets the array m_phys.
Definition: ExpList.h:1930
virtual void v_PhysDirectionalDeriv(const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:1638
virtual void v_GetBCValues(Array< OneD, NekDouble > &BndVals, const Array< OneD, NekDouble > &TotField, int BndID)
Definition: ExpList.cpp:4485
const Array< OneD, const std::pair< int, int > > & GetCoeffsToElmt() const
Get m_coeffs to elemental value map.
Definition: ExpList.h:2478
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:4345
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:2392
std::shared_ptr< ExpList > & UpdateBndCondExpansion(int i)
Definition: ExpList.h:2502
virtual void v_WriteTecplotHeader(std::ostream &outfile, std::string var="")
Definition: ExpList.cpp:2607
virtual std::shared_ptr< ExpList > & v_GetPlane(int n)
Definition: ExpList.cpp:5005
void NormVectorIProductWRTBase(Array< OneD, const NekDouble > &V1, Array< OneD, const NekDouble > &V2, Array< OneD, NekDouble > &outarray, int BndID)
Definition: ExpList.h:2224
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:2431
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:2401
virtual void v_AddTraceQuadPhysToOffDiag(const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &field)
Definition: ExpList.cpp:4335
virtual void v_ClearGlobalLinSysManager(void)
Definition: ExpList.cpp:3201
void FwdTrans_IterPerExp(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:1986
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:3286
virtual void v_GetPeriodicEntities(PeriodicMap &periodicVerts, PeriodicMap &periodicEdges, PeriodicMap &periodicFaces)
Definition: ExpList.cpp:4973
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:1299
void SetCoeffsArray(Array< OneD, NekDouble > &inarray)
Set the m_coeffs array to inarray.
Definition: ExpList.h:2280
virtual void v_CurlCurl(Array< OneD, Array< OneD, NekDouble > > &Vel, Array< OneD, Array< OneD, NekDouble > > &Q)
Definition: ExpList.cpp:1566
LibUtilities::TranspositionSharedPtr GetTransposition(void)
This function returns the transposition class associated with the homogeneous expansion.
Definition: ExpList.h:684
virtual void v_LocalToGlobal(bool UseComm)
Definition: ExpList.cpp:4574
int GetNumElmts(void)
This function returns the number of elements in the expansion which may be different for a homogeoeno...
Definition: ExpList.h:743
void SetBndCondBwdWeight(const int index, const NekDouble value)
Set the weight value for boundary conditions.
Definition: ExpList.h:2495
virtual const Array< OneD, const SpatialDomains::BoundaryConditionShPtr > & v_GetBndConditions()
Definition: ExpList.cpp:4929
void SetHomoLen(const NekDouble lhom)
This function sets the Width of homogeneous direction associated with the homogeneous expansion.
Definition: ExpList.h:698
virtual int v_GetNumElmts(void)
Definition: ExpList.h:1364
int GetNpoints(void) const
Returns the total number of quadrature points m_npoints .
Definition: ExpList.h:1882
SpatialDomains::MeshGraphSharedPtr GetGraph()
Definition: ExpList.h:1146
virtual ~ExpList()
The default destructor.
Definition: ExpList.cpp:1187
void GetFwdBwdTracePhys(Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
Definition: ExpList.h:2569
const Array< OneD, const int > & GetTraceBndMap(void)
Definition: ExpList.h:2535
virtual void v_FillBndCondFromField()
Definition: ExpList.cpp:4559
virtual void v_ImposeDirichletConditions(Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:4550
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:4950
virtual const std::vector< bool > & v_GetLeftAdjacentFaces(void) const
Definition: ExpList.cpp:4374
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:2046
Array< OneD, int > m_coeff_offset
Offset of elemental data into the array m_coeffs.
Definition: ExpList.h:1304
int EvalBasisNumModesMax(void) const
Evaulates the maximum number of modes in the elemental basis order over all elements.
Definition: ExpList.h:1814
void AddFwdBwdTraceIntegral(const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:2561
std::shared_ptr< AssemblyMapDG > & GetTraceMap(void)
Definition: ExpList.h:2530
void WriteTecplotField(std::ostream &outfile, int expansion=-1)
Definition: ExpList.h:458
virtual void v_GetMovingFrames(const SpatialDomains::GeomMMF MMFdir, const Array< OneD, const NekDouble > &CircCentre, Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: ExpList.cpp:3586
virtual std::shared_ptr< AssemblyMapDG > & v_GetTraceMap()
Definition: ExpList.cpp:3805
Array< OneD, std::pair< int, int > > m_coeffsToElmt
m_coeffs to elemental value map
Definition: ExpList.h:1310
virtual void v_GetFwdBwdTracePhys(Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
Definition: ExpList.cpp:4293
NekDouble PhysEvaluate(const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &phys)
Definition: ExpList.cpp:2542
virtual NekDouble v_Integral(const Array< OneD, const NekDouble > &inarray)
Definition: ExpList.cpp:3120
void PhysDirectionalDeriv(const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:2144
virtual std::map< int, RobinBCInfoSharedPtr > v_GetRobinBCInfo(void)
Definition: ExpList.cpp:4963
virtual const Array< OneD, const int > & v_GetTraceBndMap()
Definition: ExpList.cpp:3813
virtual void v_GetBoundaryNormals(int i, Array< OneD, Array< OneD, NekDouble > > &normals)
Definition: ExpList.cpp:4857
void GetBoundaryNormals(int i, Array< OneD, Array< OneD, NekDouble > > &normals)
Definition: ExpList.h:2743
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:4422
virtual Array< OneD, const NekDouble > v_HomogeneousEnergy(void)
Definition: ExpList.cpp:3153
virtual NekDouble v_GetHomoLen(void)
Definition: ExpList.cpp:3169
void ExtractFileBCs(const std::string &fileName, LibUtilities::CommSharedPtr comm, const std::string &varName, const std::shared_ptr< ExpList > locExpList)
Definition: ExpList.cpp:3207
void GetMovingFrames(const SpatialDomains::GeomMMF MMFdir, const Array< OneD, const NekDouble > &CircCentre, Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: ExpList.h:2101
const LocTraceToTraceMapSharedPtr & GetLocTraceToTraceMap() const
Definition: ExpList.h:2484
virtual Array< OneD, const unsigned int > v_GetZIDs(void)
Definition: ExpList.cpp:3184
void SetHomo1DSpecVanVisc(Array< OneD, NekDouble > visc)
This function sets the Spectral Vanishing Viscosity in homogeneous1D expansion.
Definition: ExpList.h:668
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:3444
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:2089
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:1831
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:2329
virtual void v_SetBndCondBwdWeight(const int index, const NekDouble value)
Definition: ExpList.cpp:4365
void GetNormals(Array< OneD, Array< OneD, NekDouble > > &normals)
Definition: ExpList.h:2540
void ExtractPhysToBnd(int i, const Array< OneD, const NekDouble > &phys, Array< OneD, NekDouble > &bnd)
Definition: ExpList.h:2736
virtual void v_SmoothField(Array< OneD, NekDouble > &field)
Definition: ExpList.cpp:1767
virtual void v_AppendFieldData(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata)
Definition: ExpList.cpp:3407
std::vector< bool > m_collectionsDoInit
Vector of bools to act as an initialise on first call flag.
Definition: ExpList.h:1295
virtual void v_FwdTrans_IterPerExp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:1735
NekDouble GetCoeff(int i)
Get the i th value (coefficient) of m_coeffs.
Definition: ExpList.h:2370
bool m_physState
The state of the array m_phys.
Definition: ExpList.h:1278
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:2202
void FillBndCondFromField(void)
Fill Bnd Condition expansion from the values stored in expansion.
Definition: ExpList.h:2304
void AddRightIPTBaseMatrix(const Array< OneD, const DNekMatSharedPtr > ElmtJacQuad, Array< OneD, DNekMatSharedPtr > ElmtJacCoef)
Definition: ExpList.cpp:5775
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:2604
virtual void v_GeneralMatrixOp(const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:4644
std::shared_ptr< LibUtilities::Comm > GetComm()
Returns the comm object.
Definition: ExpList.h:1141
virtual void v_WriteVtkPieceData(std::ostream &outfile, int expansion, std::string var)
Definition: ExpList.cpp:2998
virtual void v_WriteVtkPieceHeader(std::ostream &outfile, int expansion, int istrip)
Definition: ExpList.cpp:2866
virtual void v_BwdTrans_IterPerExp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:2366
int GetShapeDimension()
This function returns the dimension of the shape of the element eid.
Definition: ExpList.h:2255
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:707
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:3267
std::vector< int > m_coll_phys_offset
Offset of elemental data into the array m_phys.
Definition: ExpList.h:1301
LibUtilities::CommSharedPtr m_comm
Communicator.
Definition: ExpList.h:1220
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:2310
void WriteVtkHeader(std::ostream &outfile)
Definition: ExpList.cpp:2852
bool GetPhysState(void) const
This function indicates whether the array of physical values (implemented as m_phys) is filled or no...
Definition: ExpList.h:1948
std::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:1290
virtual void v_GetCoords(Array< OneD, NekDouble > &coord_0, Array< OneD, NekDouble > &coord_1, Array< OneD, NekDouble > &coord_2=NullNekDouble1DArray)
Definition: ExpList.cpp:4663
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
Definition: ExpList.h:1230
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:3638
NekDouble VectorFlux(const Array< OneD, Array< OneD, NekDouble > > &inarray)
Definition: ExpList.h:654
void AddTraceIntegralToOffDiag(const Array< OneD, const NekDouble > &FwdFlux, const Array< OneD, const NekDouble > &BwdFlux, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:1190
void HomogeneousFwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true)
Definition: ExpList.h:2167
std::shared_ptr< DNekMat > GenGlobalMatrixFull(const GlobalLinSysKey &mkey, const std::shared_ptr< AssemblyMapCG > &locToGloMap)
Definition: ExpList.cpp:2170
const Array< OneD, const SpatialDomains::BoundaryConditionShPtr > & GetBndConditions()
Definition: ExpList.h:2646
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:1939
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:728
virtual void v_NormVectorIProductWRTBase(Array< OneD, const NekDouble > &V1, Array< OneD, const NekDouble > &V2, Array< OneD, NekDouble > &outarray, int BndID)
Definition: ExpList.cpp:4494
virtual void v_ExtractElmtToBndPhys(const int i, const Array< OneD, NekDouble > &elmt, Array< OneD, NekDouble > &boundary)
Definition: ExpList.cpp:4740
virtual void v_ExtractPhysToBndElmt(const int i, const Array< OneD, const NekDouble > &phys, Array< OneD, NekDouble > &bndElmt)
Definition: ExpList.cpp:4779
void CurlCurl(Array< OneD, Array< OneD, NekDouble > > &Vel, Array< OneD, Array< OneD, NekDouble > > &Q)
Definition: ExpList.h:2157
std::vector< LibUtilities::FieldDefinitionsSharedPtr > GetFieldDefinitions()
Definition: ExpList.h:1062
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:3696
void EvaluateBoundaryConditions(const NekDouble time=0.0, const std::string varName="", const NekDouble=NekConstants::kNekUnsetDouble, const NekDouble=NekConstants::kNekUnsetDouble)
Definition: ExpList.h:2657
virtual void v_AddTraceIntegral(const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:4264
virtual void v_FillBwdWithBwdWeight(Array< OneD, NekDouble > &weightave, Array< OneD, NekDouble > &weightjmp)
Definition: ExpList.cpp:4910
Array< OneD, NekDouble > & UpdatePhys()
This function returns (a reference to) the array (implemented as m_phys) containing the function ev...
Definition: ExpList.h:2461
virtual void v_PhysGalerkinProjection1DScaled(const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:5833
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.h:2695
virtual void v_SetUpPhysNormals()
Definition: ExpList.cpp:4716
virtual void v_GlobalToLocal(void)
Definition: ExpList.cpp:4592
virtual LibUtilities::TranspositionSharedPtr v_GetTransposition(void)
Definition: ExpList.cpp:3161
void ExtractElmtToBndPhys(int i, const Array< OneD, NekDouble > &elmt, Array< OneD, NekDouble > &boundary)
Definition: ExpList.h:2722
void SetWaveSpace(const bool wavespace)
Sets the wave space to the one of the possible configuration true or false.
Definition: ExpList.h:1891
void GeneralMatrixOp_IterPerExp(const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:1961
const std::shared_ptr< LocalRegions::ExpansionVector > GetExp() const
This function returns the vector of elements in the expansion.
Definition: ExpList.h:2422
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:4408
SpatialDomains::MeshGraphSharedPtr m_graph
Mesh associated with this expansion list.
Definition: ExpList.h:1226
virtual void v_MultiplyByInvMassMatrix(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:4399
virtual void v_AddTraceQuadPhysToField(const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &field)
Definition: ExpList.cpp:4325
virtual const Array< OneD, const NekDouble > & v_GetBndCondBwdWeight()
Definition: ExpList.cpp:4357
std::shared_ptr< LibUtilities::SessionReader > GetSession() const
Returns the session object.
Definition: ExpList.h:1135
void GlobalToLocal(void)
Scatters from the global coefficients to the local coefficients .
Definition: ExpList.h:2327
virtual NekDouble v_L2(const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &soln=NullNekDouble1DArray)
Definition: ExpList.cpp:3074
const Array< OneD, const NekDouble > & GetBndCondBwdWeight()
Get the weight value for boundary conditions.
Definition: ExpList.h:2490
void WriteVtkFooter(std::ostream &outfile)
Definition: ExpList.cpp:2860
void MultiplyByBlockMatrix(const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:1199
virtual std::shared_ptr< ExpList > & v_GetTrace()
Definition: ExpList.cpp:3797
void GetElmtNormalLength(Array< OneD, NekDouble > &lengthsFwd, Array< OneD, NekDouble > &lengthsBwd)
Get the length of elements in boundary normal direction.
Definition: ExpList.cpp:4116
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:1861
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:716
virtual NekDouble v_VectorFlux(const Array< OneD, Array< OneD, NekDouble > > &inarray)
Definition: ExpList.cpp:3134
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:596
virtual std::vector< LibUtilities::FieldDefinitionsSharedPtr > v_GetFieldDefinitions(void)
Definition: ExpList.cpp:3393
std::vector< int > m_coll_coeff_offset
Offset of elemental data into the array m_coeffs.
Definition: ExpList.h:1298
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:1319
void WriteVtkPieceData(std::ostream &outfile, int expansion, std::string var="v")
Definition: ExpList.h:483
void FillBwdWithBwdWeight(Array< OneD, NekDouble > &weightave, Array< OneD, NekDouble > &weightjmp)
Fill Bwd with boundary conditions.
Definition: ExpList.h:2613
virtual void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:4608
virtual Array< OneD, SpatialDomains::BoundaryConditionShPtr > & v_UpdateBndConditions()
Definition: ExpList.cpp:4940
void SetCoeff(int i, NekDouble val)
Set the i th coefficiient in m_coeffs to value val.
Definition: ExpList.h:2264
void GetMatIpwrtDeriveBase(const Array< OneD, const Array< OneD, NekDouble > > &inarray, const int nDirctn, Array< OneD, DNekMatSharedPtr > &mtxPerVar)
Definition: ExpList.cpp:5253
void LocalToGlobal(bool useComm=true)
Gathers the global coefficients from the local coefficients .
Definition: ExpList.h:2314
void DealiasedProd(const Array< OneD, NekDouble > &inarray1, const Array< OneD, NekDouble > &inarray2, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:2191
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:661
void GetPeriodicEntities(PeriodicMap &periodicVerts, PeriodicMap &periodicEdges, PeriodicMap &periodicFaces=NullPeriodicMap)
Definition: ExpList.h:1053
NekDouble Integral(const Array< OneD, const NekDouble > &inarray)
Definition: ExpList.h:649
virtual void v_WriteTecplotField(std::ostream &outfile, int expansion)
Definition: ExpList.cpp:2818
std::shared_ptr< ExpList > GetSharedThisPtr()
Returns a shared pointer to the current object.
Definition: ExpList.h:1129
Collections::CollectionVector m_collections
Definition: ExpList.h:1292
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:2596
virtual void v_PhysInterp1DScaled(const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:5170
void ApplyGeomInfo()
Apply geometry information to each expansion.
Definition: ExpList.cpp:2567
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:677
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:1465
std::unordered_map< int, int > m_elmtToExpId
Mapping from geometry ID of element to index inside m_exp.
Definition: ExpList.h:1321
void FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:1976
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:2064
void SmoothField(Array< OneD, NekDouble > &field)
Smooth a field across elements.
Definition: ExpList.h:2007
virtual void v_WriteTecplotZone(std::ostream &outfile, int expansion)
Definition: ExpList.cpp:2646
LibUtilities::SessionReaderSharedPtr m_session
Session.
Definition: ExpList.h:1223
void HomogeneousBwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true)
Definition: ExpList.h:2179
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:4475
virtual void v_PeriodicBwdCopy(const Array< OneD, const NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
Definition: ExpList.cpp:4918
void MultiplyByInvMassMatrix(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:2036
std::map< int, RobinBCInfoSharedPtr > GetRobinBCInfo()
Definition: ExpList.h:1048
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:1084
void AddTraceIntegral(const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:2546
virtual void v_FillBwdWithBoundCond(const Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd, bool PutFwdInBwdOnBCs)
Definition: ExpList.cpp:4315
void DivideByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Divided by the metric jacobi and quadrature weights.
Definition: ExpList.cpp:1237
virtual const Array< OneD, const std::shared_ptr< ExpList > > & v_GetBndCondExpansions(void)
Definition: ExpList.cpp:3669
void ImposeDirichletConditions(Array< OneD, NekDouble > &outarray)
Impose Dirichlet Boundary Conditions onto Array.
Definition: ExpList.h:2298
Array< OneD, int > m_phys_offset
Offset of elemental data into the array m_phys.
Definition: ExpList.h:1307
virtual void v_SetHomoLen(const NekDouble lhom)
Definition: ExpList.cpp:3177
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:2516
virtual void v_GetNormals(Array< OneD, Array< OneD, NekDouble > > &normals)
Populate normals with the normals of all expansions.
Definition: ExpList.cpp:3893
Array< OneD, SpatialDomains::BoundaryConditionShPtr > & UpdateBndConditions()
Definition: ExpList.h:2652
const DNekScalBlkMatSharedPtr GenBlockMatrix(const GlobalMatrixKey &gkey)
This function assembles the block diagonal matrix of local matrices of the type mtype.
Definition: ExpList.cpp:1810
virtual void v_IProductWRTBase_IterPerExp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:1262
void PhysDeriv(Direction edir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d)
Definition: ExpList.h:2132
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:1696
void BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:2015
NekDouble GetHomoLen(void)
This function returns the Width of homogeneous direction associated with the homogeneous expansion.
Definition: ExpList.h:691
virtual void v_GetBndElmtExpansion(int i, std::shared_ptr< ExpList > &result, const bool DeclareCoeffPhysArrays)
Definition: ExpList.cpp:4729
virtual void v_FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:4614
bool GetWaveSpace(void) const
This function returns the third direction expansion condition, which can be in wave space (coefficien...
Definition: ExpList.h:1899
virtual std::shared_ptr< ExpList > & v_UpdateBndCondExpansion(int i)
Definition: ExpList.cpp:3677
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:2439
virtual const std::shared_ptr< LocTraceToTraceMap > & v_GetLocTraceToTraceMap(void) const
Definition: ExpList.cpp:5915
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:3457
ExpansionType m_expType
Exapnsion type.
Definition: ExpList.h:1212
void WriteTecplotHeader(std::ostream &outfile, std::string var="")
Definition: ExpList.h:445
const std::vector< bool > & GetLeftAdjacentFaces(void) const
Definition: ExpList.h:2627
void AppendFieldData(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata)
Append the element data listed in elements fielddef->m_ElementIDs onto fielddata.
Definition: ExpList.h:1074
Array< OneD, NekDouble > m_phys
The global expansion evaluated at the quadrature points.
Definition: ExpList.h:1269
ExpansionType GetExpType(void)
Returns the type of the expansion.
Definition: ExpList.cpp:1182
void AddTraceQuadPhysToOffDiag(const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &field)
Definition: ExpList.h:950
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:3035
virtual void v_DealiasedProd(const Array< OneD, NekDouble > &inarray1, const Array< OneD, NekDouble > &inarray2, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:4466
void WriteVtkPieceHeader(std::ostream &outfile, int expansion, int istrip=0)
Definition: ExpList.h:473
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 bool DeclareCoeffPhysArrays=true, const std::string variable="DefaultVar", const Collections::ImplementationType ImpType=Collections::eNoImpType)
Generate expansions for the trace space expansions used in DisContField.
void BwdTrans_IterPerExp(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:2025
void FwdTrans_BndConstrained(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:1996
int GetTotPoints(void) const
Returns the total number of quadrature points m_npoints .
Definition: ExpList.h:1850
void SetExpType(ExpansionType Type)
Returns the type of the expansion.
LibUtilities::BasisSharedPtr GetHomogeneousBasis(void)
Definition: ExpList.h:1152
void WriteTecplotZone(std::ostream &outfile, int expansion=-1)
Definition: ExpList.h:451
virtual void v_SetHomo1DSpecVanVisc(Array< OneD, NekDouble > visc)
Definition: ExpList.h:1773
void GetBCValues(Array< OneD, NekDouble > &BndVals, const Array< OneD, NekDouble > &TotField, int BndID)
Definition: ExpList.h:2213
int GetCoordim(int eid)
This function returns the dimension of the coordinates of the element eid.
Definition: ExpList.h:2244
void GetDiagMatIpwrtBase(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, DNekMatSharedPtr > &mtxPerVar)
Definition: ExpList.cpp:5405
Describes a matrix with ordering defined by a local to global map.
std::vector< Collection > CollectionVector
Definition: Collection.h:114
static BasisSharedPtr NullBasisSharedPtr
Definition: Basis.h:370
static std::vector< unsigned int > NullUnsignedIntVector
std::shared_ptr< Basis > BasisSharedPtr
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< FieldDefinitions > FieldDefinitionsSharedPtr
Definition: FieldIO.h:179
static std::vector< NekDouble > NullNekDoubleVector
std::shared_ptr< Transposition > TranspositionSharedPtr
static Array< OneD, BasisSharedPtr > NullBasisSharedPtr1DArray
Definition: Basis.h:371
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:1791
static PeriodicMap NullPeriodicMap
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:90
std::shared_ptr< BlockMatrixMap > BlockMatrixMapShPtr
A shared pointer to a BlockMatrixMap.
Definition: ExpList.h:101
static VarFactorsMap NullVarFactorsMap
static ExpListSharedPtr NullExpListSharedPtr
Definition: ExpList.h:1792
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:2749
std::map< StdRegions::ConstFactorType, Array< OneD, NekDouble > > VarFactorsMap
std::shared_ptr< AssemblyMap > AssemblyMapSharedPtr
Definition: AssemblyMap.h:52
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:99
static const NekDouble kNekUnsetDouble
GeomMMF
Principle direction for MMF.
std::shared_ptr< BoundaryConditionBase > BoundaryConditionShPtr
Definition: Conditions.h:219
std::map< int, CompositeSharedPtr > CompositeMap
Definition: MeshGraph.h:138
std::map< int, BoundaryConditionMapShPtr > BoundaryConditionCollection
Definition: Conditions.h:226
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174
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:272
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:314
static VarCoeffMap NullVarCoeffMap
Definition: StdRegions.hpp:273
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
std::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
Definition: NekTypeDefs.hpp:73
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:1199