Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 // License for the specific language governing rights and limitations under
14 // Permission is hereby granted, free of charge, to any person obtaining a
15 // copy of this software and associated documentation files (the "Software"),
16 // to deal in the Software without restriction, including without limitation
17 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
18 // and/or sell copies of the Software, and to permit persons to whom the
19 // Software is furnished to do so, subject to the following conditions:
20 //
21 // The above copyright notice and this permission notice shall be included
22 // in all copies or substantial portions of the Software.
23 //
24 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30 // DEALINGS IN THE SOFTWARE.
31 //
32 // Description: Expansion list top class definition
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
36 #ifndef NEKTAR_LIBS_MULTIREGIONS_EXPLIST_H
37 #define NEKTAR_LIBS_MULTIREGIONS_EXPLIST_H
38 
43 #include <LocalRegions/Expansion.h>
44 #include <Collections/Collection.h>
51 #include <boost/enable_shared_from_this.hpp>
52 #include <tinyxml.h>
53 
54 namespace Nektar
55 {
56  namespace MultiRegions
57  {
58  // Forward declarations
59  class GlobalLinSys;
60  class AssemblyMapDG;
61 
62  class AssemblyMapCG;
63  class GlobalLinSysKey;
64  class GlobalMatrix;
65 
66  enum Direction
67  {
68  eX,
69  eY,
70  eZ,
71  eS,
73  };
74 
76  {
77  e0D,
78  e1D,
79  e2D,
82  e3D,
84  };
85 
87  {
88  eX,
89  eY,
90  eZ
91  };
92 
93  /// A map between global matrix keys and their associated block
94  /// matrices.
95  typedef map<GlobalMatrixKey,DNekScalBlkMatSharedPtr> BlockMatrixMap;
96  /// A shared pointer to a BlockMatrixMap.
97  typedef boost::shared_ptr<BlockMatrixMap> BlockMatrixMapShPtr;
98 
99 
100  /// Base class for all multi-elemental spectral/hp expansions.
101  class ExpList: public boost::enable_shared_from_this<ExpList>
102  {
103  public:
104  /// The default constructor.
106 
107  /// The default constructor.
109  const LibUtilities::SessionReaderSharedPtr &pSession);
110 
111  /// The default constructor.
113  const LibUtilities::SessionReaderSharedPtr &pSession,
114  const SpatialDomains::MeshGraphSharedPtr &pGraph);
115 
116  /// Constructor copying only elements defined in eIds.
118  const ExpList &in,
119  const std::vector<unsigned int> &eIDs,
120  const bool DeclareCoeffPhysArrays = true);
121 
122  /// The copy constructor.
124  const ExpList &in,
125  const bool DeclareCoeffPhysArrays = true);
126 
127  /// The default destructor.
128  MULTI_REGIONS_EXPORT virtual ~ExpList();
129 
130  /// Returns the total number of local degrees of freedom
131  /// \f$N_{\mathrm{eof}}=\sum_{e=1}^{{N_{\mathrm{el}}}}N^{e}_m\f$.
132  inline int GetNcoeffs(void) const;
133 
134  /// Returns the total number of local degrees of freedom
135  /// for element eid
136  MULTI_REGIONS_EXPORT int GetNcoeffs(const int eid) const;
137 
138  /// Returns the type of the expansion
140 
141  /// Returns the type of the expansion
143 
144  /// Evaulates the maximum number of modes in the elemental basis
145  /// order over all elements
146  inline int EvalBasisNumModesMax(void) const;
147 
148  /// Returns the vector of the number of modes in the elemental
149  /// basis order over all elements.
151  EvalBasisNumModesMaxPerExp(void) const;
152 
153  /// Returns the total number of quadrature points #m_npoints
154  /// \f$=Q_{\mathrm{tot}}\f$.
155  inline int GetTotPoints(void) const;
156 
157  /// Returns the total number of quadrature points for eid's element
158  /// \f$=Q_{\mathrm{tot}}\f$.
159  inline int GetTotPoints(const int eid) const;
160 
161  /// Returns the total number of quadrature points #m_npoints
162  /// \f$=Q_{\mathrm{tot}}\f$.
163  inline int GetNpoints(void) const;
164 
165 
166  /// Returns the total number of qudature points scaled by
167  /// the factor scale on each 1D direction
168  inline int Get1DScaledTotPoints(const NekDouble scale) const;
169 
170  /// Sets the wave space to the one of the possible configuration
171  /// true or false
172  inline void SetWaveSpace(const bool wavespace);
173 
174 
175  ///Set Modified Basis for the stability analysis
176  inline void SetModifiedBasis(const bool modbasis);
177 
178  /// Set the \a i th value of \a m_phys to value \a val
179  inline void SetPhys(int i, NekDouble val);
180 
181  /// This function returns the third direction expansion condition,
182  /// which can be in wave space (coefficient) or not
183  /// It is stored in the variable m_WaveSpace.
184  inline bool GetWaveSpace(void) const;
185 
186  /// Fills the array #m_phys
187  inline void SetPhys(const Array<OneD, const NekDouble> &inarray);
188 
189  /// Sets the array #m_phys
190  inline void SetPhysArray(Array<OneD, NekDouble> &inarray);
191 
192  /// This function manually sets whether the array of physical
193  /// values \f$\boldsymbol{u}_l\f$ (implemented as #m_phys) is
194  /// filled or not.
195  inline void SetPhysState(const bool physState);
196 
197  /// This function indicates whether the array of physical values
198  /// \f$\boldsymbol{u}_l\f$ (implemented as #m_phys) is filled or
199  /// not.
200  inline bool GetPhysState(void) const;
201 
202  /// This function integrates a function \f$f(\boldsymbol{x})\f$
203  /// over the domain consisting of all the elements of the expansion.
205 
206  /// This function integrates a function \f$f(\boldsymbol{x})\f$
207  /// over the domain consisting of all the elements of the expansion.
209  const Array<OneD,
210  const NekDouble> &inarray);
211 
212  /// This function calculates the inner product of a function
213  /// \f$f(\boldsymbol{x})\f$ with respect to all \emph{local}
214  /// expansion modes \f$\phi_n^e(\boldsymbol{x})\f$.
215  inline void IProductWRTBase_IterPerExp(
216  const Array<OneD, const NekDouble> &inarray,
217  Array<OneD, NekDouble> &outarray);
218 
219  ///
220  inline void IProductWRTBase(
221  const Array<OneD, const NekDouble> &inarray,
222  Array<OneD, NekDouble> &outarray,
223  CoeffState coeffstate = eLocal);
224 
225  /// This function calculates the inner product of a function
226  /// \f$f(\boldsymbol{x})\f$ with respect to the derivative (in
227  /// direction \param dir) of all \emph{local} expansion modes
228  /// \f$\phi_n^e(\boldsymbol{x})\f$.
230  (const int dir,
231  const Array<OneD, const NekDouble> &inarray,
232  Array<OneD, NekDouble> &outarray);
233 
234  /// This function calculates the inner product of a function
235  /// \f$f(\boldsymbol{x})\f$ with respect to the derivative (in
236  /// direction \param dir) of all \emph{local} expansion modes
237  /// \f$\phi_n^e(\boldsymbol{x})\f$.
239  (const Array<OneD, const Array<OneD, NekDouble> > &inarray,
240  Array<OneD, NekDouble> &outarray);
241 
242  /// This function elementally evaluates the forward transformation
243  /// of a function \f$u(\boldsymbol{x})\f$ onto the global
244  /// spectral/hp expansion.
245  inline void FwdTrans_IterPerExp (
246  const Array<OneD,
247  const NekDouble> &inarray,
248  Array<OneD,NekDouble> &outarray);
249 
250  ///
251  inline void FwdTrans(
252  const Array<OneD,
253  const NekDouble> &inarray,
254  Array<OneD, NekDouble> &outarray,
255  CoeffState coeffstate = eLocal);
256 
257  /// This function elementally mulplies the coefficient space of
258  /// Sin my the elemental inverse of the mass matrix.
260  const Array<OneD,
261  const NekDouble> &inarray,
262  Array<OneD, NekDouble> &outarray);
263 
264  ///
265  inline void MultiplyByInvMassMatrix(
266  const Array<OneD,const NekDouble> &inarray,
267  Array<OneD, NekDouble> &outarray,
268  CoeffState coeffstate = eLocal);
269 
270  /// Smooth a field across elements
271  inline void SmoothField(Array<OneD,NekDouble> &field);
272 
273  /// Solve helmholtz problem
274  inline void HelmSolve(
275  const Array<OneD, const NekDouble> &inarray,
276  Array<OneD, NekDouble> &outarray,
277  const FlagList &flags,
278  const StdRegions::ConstFactorMap &factors,
279  const StdRegions::VarCoeffMap &varcoeff =
281  const Array<OneD, const NekDouble> &dirForcing =
283 
284  /// Solve Advection Diffusion Reaction
286  const Array<OneD, Array<OneD, NekDouble> > &velocity,
287  const Array<OneD, const NekDouble> &inarray,
288  Array<OneD, NekDouble> &outarray,
289  const NekDouble lambda,
290  CoeffState coeffstate = eLocal,
292  dirForcing = NullNekDouble1DArray);
293 
294 
295  /// Solve Advection Diffusion Reaction
296  inline void LinearAdvectionReactionSolve(
297  const Array<OneD, Array<OneD, NekDouble> > &velocity,
298  const Array<OneD, const NekDouble> &inarray,
299  Array<OneD, NekDouble> &outarray,
300  const NekDouble lambda,
301  CoeffState coeffstate = eLocal,
303  dirForcing = NullNekDouble1DArray);
304 
305  ///
307  const Array<OneD, const NekDouble> &inarray,
308  Array<OneD, NekDouble> &outarray);
309 
310 
311  /// This function elementally evaluates the backward transformation
312  /// of the global spectral/hp element expansion.
313  inline void BwdTrans_IterPerExp (
314  const Array<OneD, const NekDouble> &inarray,
315  Array<OneD,NekDouble> &outarray);
316 
317  ///
318  inline void BwdTrans (
319  const Array<OneD,
320  const NekDouble> &inarray,
321  Array<OneD,NekDouble> &outarray,
322  CoeffState coeffstate = eLocal);
323 
324  /// This function calculates the coordinates of all the elemental
325  /// quadrature points \f$\boldsymbol{x}_i\f$.
326  inline void GetCoords(
327  Array<OneD, NekDouble> &coord_0,
330 
331  // Homogeneous transforms
332  inline void HomogeneousFwdTrans(
333  const Array<OneD, const NekDouble> &inarray,
334  Array<OneD, NekDouble> &outarray,
335  CoeffState coeffstate = eLocal,
336  bool Shuff = true,
337  bool UnShuff = true);
338 
339  inline void HomogeneousBwdTrans(
340  const Array<OneD, const NekDouble> &inarray,
341  Array<OneD, NekDouble> &outarray,
342  CoeffState coeffstate = eLocal,
343  bool Shuff = true,
344  bool UnShuff = true);
345 
346  inline void DealiasedProd(
347  const Array<OneD, NekDouble> &inarray1,
348  const Array<OneD, NekDouble> &inarray2,
349  Array<OneD, NekDouble> &outarray,
350  CoeffState coeffstate = eLocal);
351 
352  inline void GetBCValues(
353  Array<OneD, NekDouble> &BndVals,
354  const Array<OneD, NekDouble> &TotField,
355  int BndID);
356 
357  inline void NormVectorIProductWRTBase(
360  Array<OneD, NekDouble> &outarray,
361  int BndID);
362 
363  inline void NormVectorIProductWRTBase(
365  Array<OneD, NekDouble> &outarray);
366 
367  /// Apply geometry information to each expansion.
369 
370  /// Reset geometry information and reset matrices
372  {
373  v_Reset();
374  }
375 
376  void WriteTecplotHeader(std::ostream &outfile,
377  std::string var = "")
378  {
379  v_WriteTecplotHeader(outfile, var);
380  }
381 
383  std::ostream &outfile,
384  int expansion = -1)
385  {
386  v_WriteTecplotZone(outfile, expansion);
387  }
388 
389  void WriteTecplotField(std::ostream &outfile,
390  int expansion = -1)
391  {
392  v_WriteTecplotField(outfile, expansion);
393  }
394 
395  void WriteTecplotConnectivity(std::ostream &outfile,
396  int expansion = -1)
397  {
398  v_WriteTecplotConnectivity(outfile, expansion);
399  }
400 
401  MULTI_REGIONS_EXPORT void WriteVtkHeader(std::ostream &outfile);
402  MULTI_REGIONS_EXPORT void WriteVtkFooter(std::ostream &outfile);
403 
404  void WriteVtkPieceHeader(std::ostream &outfile, int expansion,
405  int istrip = 0)
406  {
407  v_WriteVtkPieceHeader(outfile, expansion, istrip);
408  }
409 
411  std::ostream &outfile,
412  int expansion);
413 
415  std::ostream &outfile,
416  int expansion,
417  std::string var = "v")
418  {
419  v_WriteVtkPieceData(outfile, expansion, var);
420  }
421 
422  /// This function returns the dimension of the coordinates of the
423  /// element \a eid.
424  // inline
425  MULTI_REGIONS_EXPORT int GetCoordim(int eid);
426 
427  /// Set the \a i th coefficiient in \a m_coeffs to value \a val
428  inline void SetCoeff(int i, NekDouble val);
429 
430  /// Set the \a i th coefficiient in #m_coeffs to value \a val
431  inline void SetCoeffs(int i, NekDouble val);
432 
433  /// Set the #m_coeffs array to inarray
434  inline void SetCoeffsArray(Array<OneD, NekDouble> &inarray);
435 
436  /// This function returns (a reference to) the array
437  /// \f$\boldsymbol{\hat{u}}_l\f$ (implemented as #m_coeffs)
438  /// containing all local expansion coefficients.
439  inline const Array<OneD, const NekDouble> &GetCoeffs() const;
440 
441  /// Impose Dirichlet Boundary Conditions onto Array
442  inline void ImposeDirichletConditions(
443  Array<OneD,NekDouble>& outarray);
444 
445  /// Fill Bnd Condition expansion from the values stored in expansion
446  inline void FillBndCondFromField(void);
447 
448  /// Put the coefficients into global ordering using m_coeffs
449  inline void LocalToGlobal(void);
450 
451  /// Put the coefficients into local ordering and place in m_coeffs
452  inline void GlobalToLocal(void);
453 
454  /// Get the \a i th value (coefficient) of #m_coeffs
455  inline NekDouble GetCoeff(int i);
456 
457  /// Get the \a i th value (coefficient) of #m_coeffs
458  inline NekDouble GetCoeffs(int i);
459 
460  /// This function returns (a reference to) the array
461  /// \f$\boldsymbol{u}_l\f$ (implemented as #m_phys) containing the
462  /// function \f$u^{\delta}(\boldsymbol{x})\f$ evaluated at the
463  /// quadrature points.
464  // inline
466  &GetPhys() const;
467 
468  /// This function calculates the \f$L_\infty\f$ error of the global
469  /// spectral/hp element approximation.
471  const Array<OneD, const NekDouble> &inarray,
473 
474  /// This function calculates the \f$L_2\f$ error with
475  /// respect to soln of the global
476  /// spectral/hp element approximation.
478  const Array<OneD, const NekDouble> &inarray,
480  {
481  return v_L2(inarray, soln);
482  }
483 
484  /// Calculates the \f$H^1\f$ error of the global spectral/hp
485  /// element approximation.
487  const Array<OneD, const NekDouble> &inarray,
489 
491  {
492  return v_Integral(inarray);
493  }
494 
495  /// This function calculates the energy associated with
496  /// each one of the modesof a 3D homogeneous nD expansion
498  {
499  return v_HomogeneousEnergy();
500  }
501 
502  /// This function sets the Spectral Vanishing Viscosity
503  /// in homogeneous1D expansion.
505  {
507  }
508 
509  /// This function returns a vector containing the wave
510  /// numbers in z-direction associated
511  /// with the 3D homogenous expansion. Required if a
512  /// parellelisation is applied in the Fourier direction
514  {
515  return v_GetZIDs();
516  }
517 
518  /// This function returns the transposition class
519  /// associaed with the homogeneous expansion.
521  {
522  return v_GetTransposition();
523  }
524 
525  /// This function returns the Width of homogeneous direction
526  /// associaed with the homogeneous expansion.
528  {
529  return v_GetHomoLen();
530  }
531 
532  /// This function returns a vector containing the wave
533  /// numbers in y-direction associated
534  /// with the 3D homogenous expansion. Required if a
535  /// parellelisation is applied in the Fourier direction
537  {
538  return v_GetYIDs();
539  }
540 
541  /// This function interpolates the physical space points in
542  /// \a inarray to \a outarray using the same points defined in the
543  /// expansion but where the number of points are rescaled
544  /// by \a 1DScale
546  const NekDouble scale,
547  const Array<OneD, NekDouble> &inarray,
548  Array<OneD, NekDouble> &outarray)
549  {
550  v_PhysInterp1DScaled(scale, inarray,outarray);
551  }
552 
553  /// This function Galerkin projects the physical space points in
554  /// \a inarray to \a outarray where inarray is assumed to
555  /// be defined in the expansion but where the number of
556  /// points are rescaled by \a 1DScale
558  const NekDouble scale,
559  const Array<OneD, NekDouble> &inarray,
560  Array<OneD, NekDouble> &outarray)
561  {
562  v_PhysGalerkinProjection1DScaled(scale, inarray, outarray);
563  }
564 
565  /// This function returns the number of elements in the expansion.
566  inline int GetExpSize(void);
567 
568 
569  /// This function returns the number of elements in the
570  /// expansion which may be different for a homogeoenous extended
571  /// expansionp.
572  inline int GetNumElmts(void)
573  {
574  return v_GetNumElmts();
575  }
576 
577  /// This function returns the vector of elements in the expansion.
578  inline const boost::shared_ptr<LocalRegions::ExpansionVector>
579  GetExp() const;
580 
581  /// This function returns (a shared pointer to) the local elemental
582  /// expansion of the \f$n^{\mathrm{th}}\f$ element.
583  inline LocalRegions::ExpansionSharedPtr& GetExp(int n) const;
584 
585  /// This function returns (a shared pointer to) the local elemental
586  /// expansion containing the arbitrary point given by \a gloCoord.
588  const Array<OneD, const NekDouble> &gloCoord);
589 
590  /** This function returns the index of the local elemental
591  * expansion containing the arbitrary point given by \a gloCoord.
592  **/
594  const Array<OneD, const NekDouble> &gloCoord,
595  NekDouble tol = 0.0,
596  bool returnNearestElmt = false);
597 
598  /** This function returns the index and the Local
599  * Cartesian Coordinates \a locCoords of the local
600  * elemental expansion containing the arbitrary point
601  * given by \a gloCoords.
602  **/
604  const Array<OneD, const NekDouble> &gloCoords,
605  Array<OneD, NekDouble> &locCoords,
606  NekDouble tol = 0.0,
607  bool returnNearestElmt = false);
608 
609  /// Get the start offset position for a global list of #m_coeffs
610  /// correspoinding to element n.
611  inline int GetCoeff_Offset(int n) const;
612 
613  /// Get the start offset position for a global list of m_phys
614  /// correspoinding to element n.
615  inline int GetPhys_Offset(int n) const;
616 
617  /// Get the element id associated with the n th
618  /// consecutive block of data in #m_phys and #m_coeffs
619  inline int GetOffset_Elmt_Id(int n) const;
620 
621  /// This function returns (a reference to) the array
622  /// \f$\boldsymbol{\hat{u}}_l\f$ (implemented as #m_coeffs)
623  /// containing all local expansion coefficients.
625 
626  /// This function returns (a reference to) the array
627  /// \f$\boldsymbol{u}_l\f$ (implemented as #m_phys) containing the
628  /// function \f$u^{\delta}(\boldsymbol{x})\f$ evaluated at the
629  /// quadrature points.
631 
632  inline void PhysDeriv(
633  Direction edir,
634  const Array<OneD, const NekDouble> &inarray,
635  Array<OneD, NekDouble> &out_d);
636 
637  /// This function discretely evaluates the derivative of a function
638  /// \f$f(\boldsymbol{x})\f$ on the domain consisting of all
639  /// elements of the expansion.
640  inline void PhysDeriv(
641  const Array<OneD, const NekDouble> &inarray,
642  Array<OneD, NekDouble> &out_d0,
645 
646  inline void PhysDeriv(
647  const int dir,
648  const Array<OneD, const NekDouble> &inarray,
649  Array<OneD, NekDouble> &out_d);
650 
651 
652  // functions associated with DisContField
655 
656  inline boost::shared_ptr<ExpList> &UpdateBndCondExpansion(int i);
657 
658  inline void Upwind(
659  const Array<OneD, const Array<OneD, NekDouble> > &Vec,
660  const Array<OneD, const NekDouble> &Fwd,
661  const Array<OneD, const NekDouble> &Bwd,
663 
664  inline void Upwind(
666  const Array<OneD, const NekDouble> &Fwd,
667  const Array<OneD, const NekDouble> &Bwd,
669 
670  /**
671  * Return a reference to the trace space associated with this
672  * expansion list.
673  */
674  inline boost::shared_ptr<ExpList> &GetTrace();
675 
676  inline boost::shared_ptr<AssemblyMapDG> &GetTraceMap(void);
677 
678  inline const Array<OneD, const int> &GetTraceBndMap(void);
679 
680  inline void GetNormals(Array<OneD, Array<OneD, NekDouble> > &normals);
681 
682  inline void AddTraceIntegral(
685  Array<OneD, NekDouble> &outarray);
686 
687  inline void AddTraceIntegral(
689  Array<OneD, NekDouble> &outarray);
690 
691  inline void AddFwdBwdTraceIntegral(
692  const Array<OneD, const NekDouble> &Fwd,
693  const Array<OneD, const NekDouble> &Bwd,
694  Array<OneD, NekDouble> &outarray);
695 
696  inline void GetFwdBwdTracePhys(
698  Array<OneD,NekDouble> &Bwd);
699 
700  inline void GetFwdBwdTracePhys(
701  const Array<OneD,const NekDouble> &field,
703  Array<OneD,NekDouble> &Bwd);
704 
705  inline const vector<bool> &GetLeftAdjacentFaces(void) const;
706 
707  inline void ExtractTracePhys(Array<OneD,NekDouble> &outarray);
708 
709  inline void ExtractTracePhys(
710  const Array<OneD, const NekDouble> &inarray,
711  Array<OneD,NekDouble> &outarray);
712 
713  inline const Array<OneD, const SpatialDomains::
715 
716  inline Array<OneD, SpatialDomains::
718 
719  inline void EvaluateBoundaryConditions(
720  const NekDouble time = 0.0,
721  const std::string varName = "",
724 
725  // Routines for continous matrix solution
726  /// This function calculates the result of the multiplication of a
727  /// matrix of type specified by \a mkey with a vector given by \a
728  /// inarray.
729  inline void GeneralMatrixOp(
730  const GlobalMatrixKey &gkey,
731  const Array<OneD,const NekDouble> &inarray,
732  Array<OneD, NekDouble> &outarray,
733  CoeffState coeffstate = eLocal);
734 
736  const GlobalMatrixKey &gkey,
737  const Array<OneD,const NekDouble> &inarray,
738  Array<OneD, NekDouble> &outarray);
739 
740  inline void SetUpPhysNormals();
741 
742  inline void GetBoundaryToElmtMap(Array<OneD, int> &ElmtID,
743  Array<OneD,int> &EdgeID);
744 
745  inline void GetBndElmtExpansion(int i,
746  boost::shared_ptr<ExpList> &result);
747 
748  inline void ExtractElmtToBndPhys(int i,
750  Array<OneD, NekDouble> &boundary);
751 
752  inline void ExtractPhysToBndElmt(int i,
753  const Array<OneD, const NekDouble> &phys,
754  Array<OneD, NekDouble> &bndElmt);
755 
756  inline void GetBoundaryNormals(int i,
757  Array<OneD, Array<OneD, NekDouble> > &normals);
758 
760  std::vector<LibUtilities::FieldDefinitionsSharedPtr> &fielddef,
761  int NumHomoDir = 0,
764  std::vector<NekDouble> &HomoLen =
766  bool homoStrips = false,
767  std::vector<unsigned int> &HomoSIDs =
769  std::vector<unsigned int> &HomoZIDs =
771  std::vector<unsigned int> &HomoYIDs =
773 
774 
776  {
777  return m_globalOptParam;
778  }
779 
780  map<int, RobinBCInfoSharedPtr> GetRobinBCInfo()
781  {
782  return v_GetRobinBCInfo();
783  }
784 
786  PeriodicMap &periodicVerts,
787  PeriodicMap &periodicEdges,
788  PeriodicMap &periodicFaces = NullPeriodicMap)
789  {
790  v_GetPeriodicEntities(periodicVerts, periodicEdges, periodicFaces);
791  }
792 
793  std::vector<LibUtilities::FieldDefinitionsSharedPtr>
795  {
796  return v_GetFieldDefinitions();
797  }
798 
799 
800  void GetFieldDefinitions(std::vector<LibUtilities::FieldDefinitionsSharedPtr> &fielddef)
801  {
802  v_GetFieldDefinitions(fielddef);
803  }
804 
805 
806 
807  /// Append the element data listed in elements
808  /// fielddef->m_ElementIDs onto fielddata
811  std::vector<NekDouble> &fielddata)
812  {
813  v_AppendFieldData(fielddef,fielddata);
814  }
815 
816 
817  /// Append the data in coeffs listed in elements
818  /// fielddef->m_ElementIDs onto fielddata
821  std::vector<NekDouble> &fielddata,
822  Array<OneD, NekDouble> &coeffs)
823  {
824  v_AppendFieldData(fielddef,fielddata,coeffs);
825  }
826 
827 
828  /** \brief Extract the data in fielddata into the coeffs
829  * using the basic ExpList Elemental expansions rather
830  * than planes in homogeneous case
831  */
834  std::vector<NekDouble> &fielddata,
835  std::string &field,
836  Array<OneD, NekDouble> &coeffs);
837 
838 
839  /** \brief Extract the data from fromField using
840  * fromExpList the coeffs using the basic ExpList
841  * Elemental expansions rather than planes in homogeneous
842  * case
843  */
845  const boost::shared_ptr<ExpList> &fromExpList,
846  const Array<OneD, const NekDouble> &fromCoeffs,
847  Array<OneD, NekDouble> &toCoeffs);
848 
849 
850  //Extract data in fielddata into the m_coeffs_list for the 3D stability analysis (base flow is 2D)
853  std::vector<NekDouble> &fielddata,
854  std::string &field,
855  Array<OneD, NekDouble> &coeffs);
856 
857 
858  /// Returns a shared pointer to the current object.
859  boost::shared_ptr<ExpList> GetSharedThisPtr()
860  {
861  return shared_from_this();
862  }
863 
864  /// Returns the session object
865  boost::shared_ptr<LibUtilities::SessionReader> GetSession()
866  {
867  return m_session;
868  }
869 
870  /// Returns the comm object
871  boost::shared_ptr<LibUtilities::Comm> GetComm()
872  {
873  return m_comm;
874  }
875 
877  {
878  return m_graph;
879  }
880 
881  // Wrapper functions for Homogeneous Expansions
883  {
884  return v_GetHomogeneousBasis();
885  }
886 
887  boost::shared_ptr<ExpList> &GetPlane(int n)
888  {
889  return v_GetPlane(n);
890  }
891 
892  //expansion type
894 
896  Collections::ImplementationType ImpType
897  = Collections::eNoImpType);
898 
900 
901  protected:
902  boost::shared_ptr<DNekMat> GenGlobalMatrixFull(
903  const GlobalLinSysKey &mkey,
904  const boost::shared_ptr<AssemblyMapCG> &locToGloMap);
905 
906  /// Communicator
908 
909  /// Session
911 
912  /// Mesh associated with this expansion list.
914 
915  /// The total number of local degrees of freedom. #m_ncoeffs
916  /// \f$=N_{\mathrm{eof}}=\sum_{e=1}^{{N_{\mathrm{el}}}}N^{e}_l\f$
918 
919  /** The total number of quadrature points. #m_npoints
920  *\f$=Q_{\mathrm{tot}}=\sum_{e=1}^{{N_{\mathrm{el}}}}N^{e}_Q\f$
921  **/
923 
924  /**
925  * \brief Concatenation of all local expansion coefficients.
926  *
927  * The array of length #m_ncoeffs\f$=N_{\mathrm{eof}}\f$ which is
928  * the concatenation of the local expansion coefficients
929  * \f$\hat{u}_n^e\f$ over all \f$N_{\mathrm{el}}\f$ elements
930  * \f[\mathrm{\texttt{m\_coeffs}}=\boldsymbol{\hat{u}}_{l} =
931  * \underline{\boldsymbol{\hat{u}}}^e = \left [ \begin{array}{c}
932  * \boldsymbol{\hat{u}}^{1} \ \
933  * \boldsymbol{\hat{u}}^{2} \ \
934  * \vdots \ \
935  * \boldsymbol{\hat{u}}^{{{N_{\mathrm{el}}}}} \end{array} \right ],
936  * \quad
937  * \mathrm{where}\quad \boldsymbol{\hat{u}}^{e}[n]=\hat{u}_n^{e}\f]
938  */
940 
941  /**
942  * \brief The global expansion evaluated at the quadrature points
943  *
944  * The array of length #m_npoints\f$=Q_{\mathrm{tot}}\f$ containing
945  * the evaluation of \f$u^{\delta}(\boldsymbol{x})\f$ at the
946  * quadrature points \f$\boldsymbol{x}_i\f$.
947  * \f[\mathrm{\texttt{m\_phys}}=\boldsymbol{u}_{l} =
948  * \underline{\boldsymbol{u}}^e = \left [ \begin{array}{c}
949  * \boldsymbol{u}^{1} \ \
950  * \boldsymbol{u}^{2} \ \
951  * \vdots \ \
952  * \boldsymbol{u}^{{{N_{\mathrm{el}}}}} \end{array} \right ],\quad
953  * \mathrm{where}\quad
954  * \boldsymbol{u}^{e}[i]=u^{\delta}(\boldsymbol{x}_i)\f]
955  */
957 
958  /**
959  * \brief The state of the array #m_phys.
960  *
961  * Indicates whether the array #m_phys, created to contain the
962  * evaluation of \f$u^{\delta}(\boldsymbol{x})\f$ at the quadrature
963  * points \f$\boldsymbol{x}_i\f$, is filled with these values.
964  */
966 
967  /**
968  * \brief The list of local expansions.
969  *
970  * The (shared pointer to the) vector containing (shared pointers
971  * to) all local expansions. The fact that the local expansions are
972  * all stored as a (pointer to a) #StdExpansion, the abstract base
973  * class for all local expansions, allows a general implementation
974  * where most of the routines for the derived classes are defined
975  * in the #ExpList base class.
976  */
977  boost::shared_ptr<LocalRegions::ExpansionVector> m_exp;
978 
979  Collections::CollectionVector m_collections;
980 
981  /// Offset of elemental data into the array #m_coeffs
982  std::vector<int> m_coll_coeff_offset;
983 
984  /// Offset of elemental data into the array #m_phys
985  std::vector<int> m_coll_phys_offset;
986 
987  /// Offset of elemental data into the array #m_coeffs
989 
990  /// Offset of elemental data into the array #m_phys
992 
993  /// Array containing the element id #m_offset_elmt_id[n]
994  /// that the n^th consecutive block of data in #m_coeffs
995  /// and #m_phys is associated, i.e. for an array of
996  /// constant expansion size and single shape elements
997  /// m_phys[n*m_npoints] is the data related to
998  /// m_exp[m_offset_elmt_id[n]];
1000 
1002 
1003  BlockMatrixMapShPtr m_blockMat;
1004 
1005  //@todo should this be in ExpList or ExpListHomogeneous1D.cpp
1006  // it's a bool which determine if the expansion is in the wave space (coefficient space)
1007  // or not
1009 
1010  /// This function assembles the block diagonal matrix of local
1011  /// matrices of the type \a mtype.
1013  const GlobalMatrixKey &gkey);
1014 
1016  const GlobalMatrixKey &gkey);
1017 
1018  void MultiplyByBlockMatrix(
1019  const GlobalMatrixKey &gkey,
1020  const Array<OneD,const NekDouble> &inarray,
1021  Array<OneD, NekDouble> &outarray);
1022 
1023  /// Generates a global matrix from the given key and map.
1024  boost::shared_ptr<GlobalMatrix> GenGlobalMatrix(
1025  const GlobalMatrixKey &mkey,
1026  const boost::shared_ptr<AssemblyMapCG> &locToGloMap);
1027 
1028 
1029  void GlobalEigenSystem(
1030  const boost::shared_ptr<DNekMat> &Gmat,
1031  Array<OneD, NekDouble> &EigValsReal,
1032  Array<OneD, NekDouble> &EigValsImag,
1033  Array<OneD, NekDouble> &EigVecs
1035 
1036 
1037  /// This operation constructs the global linear system of type \a
1038  /// mkey.
1039  boost::shared_ptr<GlobalLinSys> GenGlobalLinSys(
1040  const GlobalLinSysKey &mkey,
1041  const boost::shared_ptr<AssemblyMapCG> &locToGloMap);
1042 
1043  /// Generate a GlobalLinSys from information provided by the key
1044  /// "mkey" and the mapping provided in LocToGloBaseMap.
1045  boost::shared_ptr<GlobalLinSys> GenGlobalBndLinSys(
1046  const GlobalLinSysKey &mkey,
1047  const AssemblyMapSharedPtr &locToGloMap);
1048 
1050  {
1052  }
1053 
1054  // Virtual prototypes
1055 
1056  virtual int v_GetNumElmts(void)
1057  {
1058  return (*m_exp).size();
1059  }
1060 
1062  &v_GetBndCondExpansions(void);
1063 
1064  virtual boost::shared_ptr<ExpList> &v_UpdateBndCondExpansion(int i);
1065 
1066  virtual void v_Upwind(
1067  const Array<OneD, const Array<OneD, NekDouble> > &Vec,
1068  const Array<OneD, const NekDouble> &Fwd,
1069  const Array<OneD, const NekDouble> &Bwd,
1071 
1072  virtual void v_Upwind(
1073  const Array<OneD, const NekDouble> &Vn,
1074  const Array<OneD, const NekDouble> &Fwd,
1075  const Array<OneD, const NekDouble> &Bwd,
1077 
1078  virtual boost::shared_ptr<ExpList> &v_GetTrace();
1079 
1080  virtual boost::shared_ptr<AssemblyMapDG> &v_GetTraceMap();
1081 
1082  virtual const Array<OneD, const int> &v_GetTraceBndMap();
1083 
1084  virtual void v_GetNormals(
1085  Array<OneD, Array<OneD, NekDouble> > &normals);
1086 
1087  virtual void v_AddTraceIntegral(
1088  const Array<OneD, const NekDouble> &Fx,
1089  const Array<OneD, const NekDouble> &Fy,
1090  Array<OneD, NekDouble> &outarray);
1091 
1092  virtual void v_AddTraceIntegral(
1093  const Array<OneD, const NekDouble> &Fn,
1094  Array<OneD, NekDouble> &outarray);
1095 
1096  virtual void v_AddFwdBwdTraceIntegral(
1097  const Array<OneD, const NekDouble> &Fwd,
1098  const Array<OneD, const NekDouble> &Bwd,
1099  Array<OneD, NekDouble> &outarray);
1100 
1101  virtual void v_GetFwdBwdTracePhys(
1102  Array<OneD,NekDouble> &Fwd,
1103  Array<OneD,NekDouble> &Bwd);
1104 
1105  virtual void v_GetFwdBwdTracePhys(
1106  const Array<OneD,const NekDouble> &field,
1107  Array<OneD,NekDouble> &Fwd,
1108  Array<OneD,NekDouble> &Bwd);
1109 
1110  virtual const vector<bool> &v_GetLeftAdjacentFaces(void) const;
1111 
1112  virtual void v_ExtractTracePhys(
1113  Array<OneD,NekDouble> &outarray);
1114 
1115  virtual void v_ExtractTracePhys(
1116  const Array<OneD, const NekDouble> &inarray,
1117  Array<OneD,NekDouble> &outarray);
1118 
1119  virtual void v_MultiplyByInvMassMatrix(
1120  const Array<OneD,const NekDouble> &inarray,
1121  Array<OneD, NekDouble> &outarray,
1122  CoeffState coeffstate);
1123 
1124  virtual void v_HelmSolve(
1125  const Array<OneD, const NekDouble> &inarray,
1126  Array<OneD, NekDouble> &outarray,
1127  const FlagList &flags,
1128  const StdRegions::ConstFactorMap &factors,
1129  const StdRegions::VarCoeffMap &varcoeff,
1130  const Array<OneD, const NekDouble> &dirForcing);
1131 
1133  const Array<OneD, Array<OneD, NekDouble> > &velocity,
1134  const Array<OneD, const NekDouble> &inarray,
1135  Array<OneD, NekDouble> &outarray,
1136  const NekDouble lambda,
1137  CoeffState coeffstate = eLocal,
1139  dirForcing = NullNekDouble1DArray);
1140 
1141  virtual void v_LinearAdvectionReactionSolve(
1142  const Array<OneD, Array<OneD, NekDouble> > &velocity,
1143  const Array<OneD, const NekDouble> &inarray,
1144  Array<OneD, NekDouble> &outarray,
1145  const NekDouble lambda,
1146  CoeffState coeffstate = eLocal,
1148  dirForcing = NullNekDouble1DArray);
1149 
1150  // wrapper functions about virtual functions
1151  virtual void v_ImposeDirichletConditions(Array<OneD,NekDouble>& outarray);
1152 
1153  virtual void v_FillBndCondFromField();
1154 
1155  virtual void v_Reset();
1156 
1157  virtual void v_LocalToGlobal(void);
1158 
1159  virtual void v_GlobalToLocal(void);
1160 
1161  virtual void v_BwdTrans(
1162  const Array<OneD,const NekDouble> &inarray,
1163  Array<OneD, NekDouble> &outarray,
1164  CoeffState coeffstate);
1165 
1166  virtual void v_BwdTrans_IterPerExp(
1167  const Array<OneD,const NekDouble> &inarray,
1168  Array<OneD,NekDouble> &outarray);
1169 
1170  virtual void v_FwdTrans(
1171  const Array<OneD,const NekDouble> &inarray,
1172  Array<OneD, NekDouble> &outarray,
1173  CoeffState coeffstate);
1174 
1175  virtual void v_FwdTrans_IterPerExp(
1176  const Array<OneD,const NekDouble> &inarray,
1177  Array<OneD,NekDouble> &outarray);
1178 
1179  virtual void v_SmoothField(Array<OneD,NekDouble> &field);
1180 
1181  virtual void v_IProductWRTBase(
1182  const Array<OneD, const NekDouble> &inarray,
1183  Array<OneD, NekDouble> &outarray,
1184  CoeffState coeffstate);
1185 
1186  virtual void v_IProductWRTBase_IterPerExp(
1187  const Array<OneD,const NekDouble> &inarray,
1188  Array<OneD, NekDouble> &outarray);
1189 
1190  virtual void v_GeneralMatrixOp(
1191  const GlobalMatrixKey &gkey,
1192  const Array<OneD,const NekDouble> &inarray,
1193  Array<OneD, NekDouble> &outarray,
1194  CoeffState coeffstate);
1195 
1196  virtual void v_GetCoords(
1197  Array<OneD, NekDouble> &coord_0,
1198  Array<OneD, NekDouble> &coord_1,
1200 
1201  virtual void v_PhysDeriv(
1202  const Array<OneD, const NekDouble> &inarray,
1203  Array<OneD, NekDouble> &out_d0,
1204  Array<OneD, NekDouble> &out_d1,
1205  Array<OneD, NekDouble> &out_d2);
1206 
1207  virtual void v_PhysDeriv(
1208  const int dir,
1209  const Array<OneD, const NekDouble> &inarray,
1210  Array<OneD, NekDouble> &out_d);
1211 
1212  virtual void v_PhysDeriv(
1213  Direction edir,
1214  const Array<OneD, const NekDouble> &inarray,
1215  Array<OneD, NekDouble> &out_d);
1216 
1217  virtual void v_HomogeneousFwdTrans(
1218  const Array<OneD, const NekDouble> &inarray,
1219  Array<OneD, NekDouble> &outarray,
1220  CoeffState coeffstate = eLocal,
1221  bool Shuff = true,
1222  bool UnShuff = true);
1223 
1224  virtual void v_HomogeneousBwdTrans(
1225  const Array<OneD, const NekDouble> &inarray,
1226  Array<OneD, NekDouble> &outarray,
1227  CoeffState coeffstate = eLocal,
1228  bool Shuff = true,
1229  bool UnShuff = true);
1230 
1231  virtual void v_DealiasedProd(
1232  const Array<OneD, NekDouble> &inarray1,
1233  const Array<OneD, NekDouble> &inarray2,
1234  Array<OneD, NekDouble> &outarray,
1235  CoeffState coeffstate = eLocal);
1236 
1237  virtual void v_GetBCValues(
1238  Array<OneD, NekDouble> &BndVals,
1239  const Array<OneD, NekDouble> &TotField,
1240  int BndID);
1241 
1242  virtual void v_NormVectorIProductWRTBase(
1245  Array<OneD, NekDouble> &outarray,
1246  int BndID);
1247 
1248  virtual void v_NormVectorIProductWRTBase(
1249  Array<OneD, Array<OneD, NekDouble> > &V,
1250  Array<OneD, NekDouble> &outarray);
1251 
1252  virtual void v_SetUpPhysNormals();
1253 
1254  virtual void v_GetBoundaryToElmtMap(Array<OneD, int> &ElmtID,
1255  Array<OneD,int> &EdgeID);
1256 
1257  virtual void v_GetBndElmtExpansion(int i,
1258  boost::shared_ptr<ExpList> &result);
1259 
1260  virtual void v_ExtractElmtToBndPhys(int i,
1261  Array<OneD, NekDouble> &elmt,
1262  Array<OneD, NekDouble> &boundary);
1263 
1264  virtual void v_ExtractPhysToBndElmt(int i,
1265  const Array<OneD, const NekDouble> &phys,
1266  Array<OneD, NekDouble> &bndElmt);
1267 
1268  virtual void v_GetBoundaryNormals(int i,
1269  Array<OneD, Array<OneD, NekDouble> > &normals);
1270 
1271  virtual void v_ReadGlobalOptimizationParameters();
1272 
1273  virtual std::vector<LibUtilities::FieldDefinitionsSharedPtr>
1274  v_GetFieldDefinitions(void);
1275 
1276  virtual void v_GetFieldDefinitions(
1277  std::vector<LibUtilities::FieldDefinitionsSharedPtr> &fielddef);
1278 
1279  virtual void v_AppendFieldData(
1281  std::vector<NekDouble> &fielddata);
1282 
1283  virtual void v_AppendFieldData(
1285  std::vector<NekDouble> &fielddata,
1286  Array<OneD, NekDouble> &coeffs);
1287 
1288  virtual void v_ExtractDataToCoeffs(
1290  std::vector<NekDouble> &fielddata, std::string &field,
1291  Array<OneD, NekDouble> &coeffs);
1292 
1293  virtual void v_ExtractCoeffsToCoeffs(const boost::shared_ptr<ExpList> &fromExpList, const Array<OneD, const NekDouble> &fromCoeffs, Array<OneD, NekDouble> &toCoeffs);
1294 
1295  virtual void v_WriteTecplotHeader(std::ostream &outfile,
1296  std::string var = "");
1297  virtual void v_WriteTecplotZone(std::ostream &outfile,
1298  int expansion);
1299  virtual void v_WriteTecplotField(std::ostream &outfile,
1300  int expansion);
1301  virtual void v_WriteTecplotConnectivity(std::ostream &outfile,
1302  int expansion);
1303  virtual void v_WriteVtkPieceHeader(
1304  std::ostream &outfile,
1305  int expansion,
1306  int istrip);
1307 
1308  virtual void v_WriteVtkPieceData(
1309  std::ostream &outfile,
1310  int expansion,
1311  std::string var);
1312 
1313  virtual NekDouble v_L2(
1314  const Array<OneD, const NekDouble> &phys,
1316 
1317  virtual NekDouble v_Integral (
1318  const Array<OneD, const NekDouble> &inarray);
1319 
1322  virtual NekDouble v_GetHomoLen(void);
1325 
1326  // 1D Scaling and projection
1327  virtual void v_PhysInterp1DScaled(
1328  const NekDouble scale, const Array<OneD, NekDouble> &inarray,
1329  Array<OneD, NekDouble> &outarray);
1330 
1331  virtual void v_PhysGalerkinProjection1DScaled(
1332  const NekDouble scale,
1333  const Array<OneD, NekDouble> &inarray,
1334  Array<OneD, NekDouble> &outarray);
1335 
1336  virtual void v_ClearGlobalLinSysManager(void);
1337 
1338  void ExtractFileBCs(const std::string &fileName,
1339  const std::string &varName,
1340  const boost::shared_ptr<ExpList> locExpList);
1341 
1342  // Utility function for a common case of retrieving a
1343  // BoundaryCondition from a boundary condition collection.
1346  GetBoundaryCondition(const SpatialDomains::
1347  BoundaryConditionCollection& collection,
1348  unsigned int index, const std::string& variable);
1349 
1350 
1351  private:
1353 
1356 
1357  virtual void v_EvaluateBoundaryConditions(
1358  const NekDouble time = 0.0,
1359  const std::string varName = "",
1362 
1363  virtual map<int, RobinBCInfoSharedPtr> v_GetRobinBCInfo(void);
1364 
1365 
1366  virtual void v_GetPeriodicEntities(
1367  PeriodicMap &periodicVerts,
1368  PeriodicMap &periodicEdges,
1369  PeriodicMap &periodicFaces);
1370 
1371  // Homogeneous direction wrapper functions.
1373  {
1374  ASSERTL0(false,
1375  "This method is not defined or valid for this class type");
1377  }
1378 
1379  // wrapper function to set viscosity for Homo1D expansion
1381  {
1382  ASSERTL0(false,
1383  "This method is not defined or valid for this class type");
1384  }
1385 
1386 
1387  virtual boost::shared_ptr<ExpList> &v_GetPlane(int n);
1388  };
1389 
1390 
1391  /// Shared pointer to an ExpList object.
1392  typedef boost::shared_ptr<ExpList> ExpListSharedPtr;
1393  /// An empty ExpList object.
1395  static ExpListSharedPtr NullExpListSharedPtr;
1396 
1397  // Inline routines follow.
1398 
1399  /**
1400  * Returns the total number of local degrees of freedom
1401  * \f$N_{\mathrm{eof}}=\sum_{e=1}^{{N_{\mathrm{el}}}}N^{e}_m\f$.
1402  */
1403  inline int ExpList::GetNcoeffs() const
1404  {
1405  return m_ncoeffs;
1406  }
1407 
1408  inline int ExpList::GetNcoeffs(const int eid) const
1409  {
1410  return (*m_exp)[eid]->GetNcoeffs();
1411  }
1412 
1413  /**
1414  * Evaulates the maximum number of modes in the elemental basis
1415  * order over all elements
1416  */
1418  {
1419  unsigned int i;
1420  int returnval = 0;
1421 
1422  for(i= 0; i < (*m_exp).size(); ++i)
1423  {
1424  returnval = max(returnval,
1425  (*m_exp)[i]->EvalBasisNumModesMax());
1426  }
1427 
1428  return returnval;
1429  }
1430 
1431  /**
1432  *
1433  */
1435  const
1436  {
1437  unsigned int i;
1438  Array<OneD,int> returnval((*m_exp).size(),0);
1439 
1440  for(i= 0; i < (*m_exp).size(); ++i)
1441  {
1442  returnval[i]
1443  = max(returnval[i],(*m_exp)[i]->EvalBasisNumModesMax());
1444  }
1445 
1446  return returnval;
1447  }
1448 
1449 
1450  /**
1451  *
1452  */
1453  inline int ExpList::GetTotPoints() const
1454  {
1455  return m_npoints;
1456  }
1457 
1458  inline int ExpList::GetTotPoints(const int eid) const
1459  {
1460  return (*m_exp)[eid]->GetTotPoints();
1461  }
1462 
1463 
1464  inline int ExpList::Get1DScaledTotPoints(const NekDouble scale) const
1465  {
1466  int returnval = 0;
1467  int cnt;
1468  int nbase = (*m_exp)[0]->GetNumBases();
1469 
1470  for(int i = 0; i < (*m_exp).size(); ++i)
1471  {
1472  cnt = 1;
1473  for(int j = 0; j < nbase; ++j)
1474  {
1475  cnt *= (int)(scale*((*m_exp)[i]->GetNumPoints(j)));
1476  }
1477  returnval += cnt;
1478  }
1479  return returnval;
1480  }
1481 
1482  /**
1483  *
1484  */
1485  inline int ExpList::GetNpoints() const
1486  {
1487  return m_npoints;
1488  }
1489 
1490 
1491  /**
1492  *
1493  */
1494  inline void ExpList::SetWaveSpace(const bool wavespace)
1495  {
1496  m_WaveSpace = wavespace;
1497  }
1498 
1499  /**
1500  *
1501  */
1502  inline bool ExpList::GetWaveSpace() const
1503  {
1504  return m_WaveSpace;
1505  }
1506 
1507  /// Set the \a i th value of\a m_phys to value \a val
1508  inline void ExpList::SetPhys(int i, NekDouble val)
1509  {
1510  m_phys[i] = val;
1511  }
1512 
1513 
1514  /**
1515  * This function fills the array \f$\boldsymbol{u}_l\f$, the evaluation
1516  * of the expansion at the quadrature points (implemented as #m_phys),
1517  * with the values of the array \a inarray.
1518  *
1519  * @param inarray The array containing the values where
1520  * #m_phys should be filled with.
1521  */
1522  inline void ExpList::SetPhys(
1523  const Array<OneD, const NekDouble> &inarray)
1524  {
1525  ASSERTL0(inarray.num_elements() == m_npoints,
1526  "Input array does not have correct number of elements.");
1527 
1528  Vmath::Vcopy(m_npoints,&inarray[0],1,&m_phys[0],1);
1529  m_physState = true;
1530  }
1531 
1532 
1534  {
1535  m_phys = inarray;
1536  }
1537 
1538 
1539  /**
1540  * @param physState \a true (=filled) or \a false (=not filled).
1541  */
1542  inline void ExpList::SetPhysState(const bool physState)
1543  {
1544  m_physState = physState;
1545  }
1546 
1547 
1548  /**
1549  * @return physState \a true (=filled) or \a false (=not filled).
1550  */
1551  inline bool ExpList::GetPhysState() const
1552  {
1553  return m_physState;
1554  }
1555 
1556  /**
1557  *
1558  */
1560  const Array<OneD, const NekDouble> &inarray,
1561  Array<OneD, NekDouble> &outarray,
1562  CoeffState coeffstate)
1563  {
1564  v_IProductWRTBase(inarray,outarray, coeffstate);
1565  }
1566 
1567  /**
1568  *
1569  */
1571  const Array<OneD, const NekDouble> &inarray,
1572  Array<OneD, NekDouble> &outarray)
1573  {
1574  v_IProductWRTBase_IterPerExp(inarray,outarray);
1575  }
1576 
1577  /**
1578  *
1579  */
1580  inline void ExpList::FwdTrans(
1581  const Array<OneD, const NekDouble> &inarray,
1582  Array<OneD, NekDouble> &outarray,
1583  CoeffState coeffstate)
1584  {
1585  v_FwdTrans(inarray,outarray,coeffstate);
1586  }
1587 
1588  /**
1589  *
1590  */
1592  const Array<OneD, const NekDouble> &inarray,
1593  Array<OneD,NekDouble> &outarray)
1594  {
1595  v_FwdTrans_IterPerExp(inarray,outarray);
1596  }
1597 
1598  /**
1599  *
1600  */
1602  {
1603  v_SmoothField(field);
1604  }
1605 
1606  /**
1607  *
1608  */
1609  inline void ExpList::BwdTrans (
1610  const Array<OneD, const NekDouble> &inarray,
1611  Array<OneD, NekDouble> &outarray,
1612  CoeffState coeffstate)
1613  {
1614  v_BwdTrans(inarray,outarray,coeffstate);
1615  }
1616 
1617  /**
1618  *
1619  */
1621  const Array<OneD, const NekDouble> &inarray,
1622  Array<OneD, NekDouble> &outarray)
1623  {
1624  v_BwdTrans_IterPerExp(inarray,outarray);
1625  }
1626 
1627 
1628  /**
1629  *
1630  */
1632  const Array<OneD,const NekDouble> &inarray,
1633  Array<OneD, NekDouble> &outarray,
1634  CoeffState coeffstate)
1635  {
1636  v_MultiplyByInvMassMatrix(inarray,outarray,coeffstate);
1637  }
1638 
1639  /**
1640  *
1641  */
1642  inline void ExpList::HelmSolve(
1643  const Array<OneD, const NekDouble> &inarray,
1644  Array<OneD, NekDouble> &outarray,
1645  const FlagList &flags,
1646  const StdRegions::ConstFactorMap &factors,
1647  const StdRegions::VarCoeffMap &varcoeff,
1648  const Array<OneD, const NekDouble> &dirForcing)
1649  {
1650  v_HelmSolve(inarray, outarray, flags, factors, varcoeff, dirForcing);
1651  }
1652 
1653 
1654  /**
1655  *
1656  */
1658  const Array<OneD, Array<OneD, NekDouble> > &velocity,
1659  const Array<OneD, const NekDouble> &inarray,
1660  Array<OneD, NekDouble> &outarray,
1661  const NekDouble lambda,
1662  CoeffState coeffstate,
1663  const Array<OneD, const NekDouble>& dirForcing)
1664  {
1665  v_LinearAdvectionDiffusionReactionSolve(velocity,inarray, outarray,
1666  lambda, coeffstate,dirForcing);
1667  }
1668 
1670  const Array<OneD, Array<OneD, NekDouble> > &velocity,
1671  const Array<OneD, const NekDouble> &inarray,
1672  Array<OneD, NekDouble> &outarray,
1673  const NekDouble lambda,
1674  CoeffState coeffstate,
1675  const Array<OneD, const NekDouble>& dirForcing)
1676  {
1677  v_LinearAdvectionReactionSolve(velocity,inarray, outarray,
1678  lambda, coeffstate,dirForcing);
1679  }
1680 
1681  /**
1682  *
1683  */
1685  Array<OneD, NekDouble> &coord_1,
1686  Array<OneD, NekDouble> &coord_2)
1687 
1688  {
1689  v_GetCoords(coord_0,coord_1,coord_2);
1690  }
1691 
1692  /**
1693  *
1694  */
1696  Array<OneD, NekDouble> &out_d0,
1697  Array<OneD, NekDouble> &out_d1,
1698  Array<OneD, NekDouble> &out_d2)
1699  {
1700  v_PhysDeriv(inarray,out_d0,out_d1,out_d2);
1701  }
1702 
1703  /**
1704  *
1705  */
1706  inline void ExpList::PhysDeriv(
1707  const int dir,
1708  const Array<OneD, const NekDouble> &inarray,
1709  Array<OneD, NekDouble> &out_d)
1710  {
1711  v_PhysDeriv(dir,inarray,out_d);
1712  }
1713 
1714  inline void ExpList::PhysDeriv(
1715  Direction edir,
1716  const Array<OneD, const NekDouble> &inarray,
1717  Array<OneD, NekDouble> &out_d)
1718  {
1719  v_PhysDeriv(edir, inarray,out_d);
1720  }
1721 
1722  /**
1723  *
1724  */
1726  const Array<OneD, const NekDouble> &inarray,
1727  Array<OneD, NekDouble> &outarray,
1728  CoeffState coeffstate,
1729  bool Shuff,
1730  bool UnShuff)
1731  {
1732  v_HomogeneousFwdTrans(inarray,outarray,coeffstate,Shuff,UnShuff);
1733  }
1734 
1735  /**
1736  *
1737  */
1739  const Array<OneD, const NekDouble> &inarray,
1740  Array<OneD, NekDouble> &outarray,
1741  CoeffState coeffstate,
1742  bool Shuff,
1743  bool UnShuff)
1744  {
1745  v_HomogeneousBwdTrans(inarray,outarray,coeffstate,Shuff,UnShuff);
1746  }
1747 
1748  /**
1749  *
1750  */
1752  const Array<OneD, NekDouble> &inarray1,
1753  const Array<OneD, NekDouble> &inarray2,
1754  Array<OneD, NekDouble> &outarray,
1755  CoeffState coeffstate)
1756  {
1757  v_DealiasedProd(inarray1,inarray2,outarray,coeffstate);
1758  }
1759 
1760  /**
1761  *
1762  */
1764  Array<OneD, NekDouble> &BndVals,
1765  const Array<OneD, NekDouble> &TotField,
1766  int BndID)
1767  {
1768  v_GetBCValues(BndVals,TotField,BndID);
1769  }
1770 
1771  /**
1772  *
1773  */
1777  Array<OneD, NekDouble> &outarray,
1778  int BndID)
1779  {
1780  v_NormVectorIProductWRTBase(V1,V2,outarray,BndID);
1781  }
1782 
1785  Array<OneD, NekDouble> &outarray)
1786  {
1787  v_NormVectorIProductWRTBase(V, outarray);
1788  }
1789 
1790  /**
1791  * @param eid The index of the element to be checked.
1792  * @return The dimension of the coordinates of the specific element.
1793  */
1794  inline int ExpList::GetCoordim(int eid)
1795  {
1796  ASSERTL2(eid <= (*m_exp).size(),
1797  "eid is larger than number of elements");
1798  return (*m_exp)[eid]->GetCoordim();
1799  }
1800 
1801  /**
1802  * @param i The index of m_coeffs to be set
1803  * @param val The value which m_coeffs[i] is to be set to.
1804  */
1805  inline void ExpList::SetCoeff(int i, NekDouble val)
1806  {
1807  m_coeffs[i] = val;
1808  }
1809 
1810 
1811  /**
1812  * @param i The index of #m_coeffs to be set.
1813  * @param val The value which #m_coeffs[i] is to be set to.
1814  */
1815  inline void ExpList::SetCoeffs(int i, NekDouble val)
1816  {
1817  m_coeffs[i] = val;
1818  }
1819 
1820 
1822  {
1823  m_coeffs = inarray;
1824  }
1825 
1826  /**
1827  * As the function returns a constant reference to a
1828  * <em>const Array</em>, it is not possible to modify the
1829  * underlying data of the array #m_coeffs. In order to do
1830  * so, use the function #UpdateCoeffs instead.
1831  *
1832  * @return (A constant reference to) the array #m_coeffs.
1833  */
1835  {
1836  return m_coeffs;
1837  }
1838 
1840  Array<OneD,NekDouble>& outarray)
1841  {
1842  v_ImposeDirichletConditions(outarray);
1843  }
1844 
1846  {
1848  }
1849 
1850  inline void ExpList::LocalToGlobal(void)
1851  {
1852  v_LocalToGlobal();
1853  }
1854 
1855  inline void ExpList::GlobalToLocal(void)
1856  {
1857  v_GlobalToLocal();
1858  }
1859 
1860 
1861  /**
1862  * @param i The index of #m_coeffs to be returned
1863  * @return The NekDouble held in #m_coeffs[i].
1864  */
1866  {
1867  return m_coeffs[i];
1868  }
1869 
1870  /**
1871  * @param i The index of #m_coeffs to be returned
1872  * @return The NekDouble held in #m_coeffs[i].
1873  */
1875  {
1876  return m_coeffs[i];
1877  }
1878 
1879  /**
1880  * As the function returns a constant reference to a
1881  * <em>const Array</em> it is not possible to modify the
1882  * underlying data of the array #m_phys. In order to do
1883  * so, use the function #UpdatePhys instead.
1884  *
1885  * @return (A constant reference to) the array #m_phys.
1886  */
1888  {
1889  return m_phys;
1890  }
1891 
1892  /**
1893  * @return \f$N_{\mathrm{el}}\f$, the number of elements in the
1894  * expansion.
1895  */
1896  inline int ExpList::GetExpSize(void)
1897  {
1898  return (*m_exp).size();
1899  }
1900 
1901 
1902  /**
1903  * @param n The index of the element concerned.
1904  *
1905  * @return (A shared pointer to) the local expansion of the
1906  * \f$n^{\mathrm{th}}\f$ element.
1907  */
1909  {
1910  return (*m_exp)[n];
1911  }
1912 
1913  /**
1914  * @return (A const shared pointer to) the local expansion vector #m_exp
1915  */
1916  inline const boost::shared_ptr<LocalRegions::ExpansionVector>
1917  ExpList::GetExp(void) const
1918  {
1919  return m_exp;
1920  }
1921 
1922 
1923  /**
1924  *
1925  */
1926  inline int ExpList::GetCoeff_Offset(int n) const
1927  {
1928  return m_coeff_offset[n];
1929  }
1930 
1931  /**
1932  *
1933  */
1934  inline int ExpList::GetPhys_Offset(int n) const
1935  {
1936  return m_phys_offset[n];
1937  }
1938 
1939  /**
1940  *
1941  */
1942  inline int ExpList::GetOffset_Elmt_Id(int n) const
1943  {
1944  return m_offset_elmt_id[n];
1945  }
1946 
1947  /**
1948  * If one wants to get hold of the underlying data without modifying
1949  * them, rather use the function #GetCoeffs instead.
1950  *
1951  * @return (A reference to) the array #m_coeffs.
1952  */
1954  {
1955  return m_coeffs;
1956  }
1957 
1958  /**
1959  * If one wants to get hold of the underlying data without modifying
1960  * them, rather use the function #GetPhys instead.
1961  *
1962  * @return (A reference to) the array #m_phys.
1963  */
1965  {
1966  m_physState = true;
1967  return m_phys;
1968  }
1969 
1970 
1971  // functions associated with DisContField
1974  {
1975  return v_GetBndCondExpansions();
1976  }
1977 
1978  inline boost::shared_ptr<ExpList> &ExpList::UpdateBndCondExpansion(int i)
1979  {
1980  return v_UpdateBndCondExpansion(i);
1981  }
1982 
1983  inline void ExpList::Upwind(
1984  const Array<OneD, const Array<OneD, NekDouble> > &Vec,
1985  const Array<OneD, const NekDouble> &Fwd,
1986  const Array<OneD, const NekDouble> &Bwd,
1987  Array<OneD, NekDouble> &Upwind)
1988  {
1989  v_Upwind(Vec, Fwd, Bwd, Upwind);
1990  }
1991 
1992  inline void ExpList::Upwind(
1993  const Array<OneD, const NekDouble> &Vn,
1994  const Array<OneD, const NekDouble> &Fwd,
1995  const Array<OneD, const NekDouble> &Bwd,
1996  Array<OneD, NekDouble> &Upwind)
1997  {
1998  v_Upwind(Vn, Fwd, Bwd, Upwind);
1999  }
2000 
2001  inline boost::shared_ptr<ExpList> &ExpList::GetTrace()
2002  {
2003  return v_GetTrace();
2004  }
2005 
2006  inline boost::shared_ptr<AssemblyMapDG> &ExpList::GetTraceMap()
2007  {
2008  return v_GetTraceMap();
2009  }
2010 
2012  {
2013  return v_GetTraceBndMap();
2014  }
2015 
2016  inline void ExpList::GetNormals(
2017  Array<OneD, Array<OneD, NekDouble> > &normals)
2018  {
2019  v_GetNormals(normals);
2020  }
2021 
2023  const Array<OneD, const NekDouble> &Fx,
2024  const Array<OneD, const NekDouble> &Fy,
2025  Array<OneD, NekDouble> &outarray)
2026  {
2027  v_AddTraceIntegral(Fx,Fy,outarray);
2028  }
2029 
2031  const Array<OneD, const NekDouble> &Fn,
2032  Array<OneD, NekDouble> &outarray)
2033  {
2034  v_AddTraceIntegral(Fn,outarray);
2035  }
2036 
2038  const Array<OneD, const NekDouble> &Fwd,
2039  const Array<OneD, const NekDouble> &Bwd,
2040  Array<OneD, NekDouble> &outarray)
2041  {
2042  v_AddFwdBwdTraceIntegral(Fwd,Bwd,outarray);
2043  }
2044 
2046  Array<OneD,NekDouble> &Fwd,
2047  Array<OneD,NekDouble> &Bwd)
2048  {
2049  v_GetFwdBwdTracePhys(Fwd,Bwd);
2050  }
2051 
2053  const Array<OneD,const NekDouble> &field,
2054  Array<OneD,NekDouble> &Fwd,
2055  Array<OneD,NekDouble> &Bwd)
2056  {
2057  v_GetFwdBwdTracePhys(field,Fwd,Bwd);
2058  }
2059 
2060  inline const vector<bool> &ExpList::GetLeftAdjacentFaces(void) const
2061  {
2062  return v_GetLeftAdjacentFaces();
2063  }
2064 
2066  {
2067  v_ExtractTracePhys(outarray);
2068  }
2069 
2070 
2072  const Array<OneD, const NekDouble> &inarray,
2073  Array<OneD,NekDouble> &outarray)
2074  {
2075  v_ExtractTracePhys(inarray,outarray);
2076  }
2077 
2080  {
2081  return v_GetBndConditions();
2082  }
2083 
2084 
2087  {
2088  return v_UpdateBndConditions();
2089  }
2090 
2092  const NekDouble time,
2093  const std::string varName,
2094  const NekDouble x2_in,
2095  const NekDouble x3_in)
2096  {
2097  v_EvaluateBoundaryConditions(time, varName, x2_in, x3_in);
2098  }
2099 
2100  // Routines for continous matrix solution
2101  /**
2102  * This operation is equivalent to the evaluation of
2103  * \f$\underline{\boldsymbol{M}}^e\boldsymbol{\hat{u}}_l\f$, that is,
2104  * \f[ \left[
2105  * \begin{array}{cccc}
2106  * \boldsymbol{M}^1 & 0 & \hspace{3mm}0 \hspace{3mm}& 0 \\
2107  * 0 & \boldsymbol{M}^2 & 0 & 0 \\
2108  * 0 & 0 & \ddots & 0 \\
2109  * 0 & 0 & 0 & \boldsymbol{M}^{N_{\mathrm{el}}} \end{array} \right]
2110  *\left [ \begin{array}{c}
2111  * \boldsymbol{\hat{u}}^{1} \\
2112  * \boldsymbol{\hat{u}}^{2} \\
2113  * \vdots \\
2114  * \boldsymbol{\hat{u}}^{{{N_{\mathrm{el}}}}} \end{array} \right ]\f]
2115  * where \f$\boldsymbol{M}^e\f$ are the local matrices of type
2116  * specified by the key \a mkey. The decoupling of the local matrices
2117  * allows for a local evaluation of the operation. However, rather than
2118  * a local matrix-vector multiplication, the local operations are
2119  * evaluated as implemented in the function
2120  * StdRegions#StdExpansion#GeneralMatrixOp.
2121  *
2122  * @param mkey This key uniquely defines the type matrix
2123  * required for the operation.
2124  * @param inarray The vector \f$\boldsymbol{\hat{u}}_l\f$ of
2125  * size \f$N_{\mathrm{eof}}\f$.
2126  * @param outarray The resulting vector of size
2127  * \f$N_{\mathrm{eof}}\f$.
2128  */
2130  const GlobalMatrixKey &gkey,
2131  const Array<OneD,const NekDouble> &inarray,
2132  Array<OneD, NekDouble> &outarray,
2133  CoeffState coeffstate)
2134  {
2135  v_GeneralMatrixOp(gkey,inarray,outarray,coeffstate);
2136  }
2137 
2138 
2140  {
2142  }
2143 
2145  Array<OneD,int> &EdgeID)
2146  {
2147  v_GetBoundaryToElmtMap(ElmtID,EdgeID);
2148  }
2149 
2150  inline void ExpList::GetBndElmtExpansion(int i,
2151  boost::shared_ptr<ExpList> &result)
2152  {
2153  v_GetBndElmtExpansion(i, result);
2154  }
2155 
2157  Array<OneD, NekDouble> &elmt,
2158  Array<OneD, NekDouble> &boundary)
2159  {
2160  v_ExtractElmtToBndPhys(i, elmt, boundary);
2161  }
2162 
2164  const Array<OneD, const NekDouble> &phys,
2165  Array<OneD, NekDouble> &bndElmt)
2166  {
2167  v_ExtractPhysToBndElmt(i, phys, bndElmt);
2168  }
2169 
2170  inline void ExpList::GetBoundaryNormals(int i,
2171  Array<OneD, Array<OneD, NekDouble> > &normals)
2172  {
2173  v_GetBoundaryNormals(i, normals);
2174  }
2175 
2177 
2178  } //end of namespace
2179 } //end of namespace
2180 
2181 #endif // EXPLIST_H
2182 
virtual NekDouble v_L2(const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &soln=NullNekDouble1DArray)
Definition: ExpList.cpp:1807
virtual void v_GlobalToLocal(void)
Definition: ExpList.cpp:2537
virtual void v_AddFwdBwdTraceIntegral(const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:2336
virtual void v_SmoothField(Array< OneD, NekDouble > &field)
Definition: ExpList.cpp:634
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:488
void AddTraceIntegral(const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:2022
virtual void v_GeneralMatrixOp(const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
Definition: ExpList.cpp:2573
const DNekScalBlkMatSharedPtr & GetBlockMatrix(const GlobalMatrixKey &gkey)
Definition: ExpList.cpp:812
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:2156
virtual void v_GetBndElmtExpansion(int i, boost::shared_ptr< ExpList > &result)
Definition: ExpList.cpp:2654
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:1834
void SmoothField(Array< OneD, NekDouble > &field)
Smooth a field across elements.
Definition: ExpList.h:1601
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
int GetExpIndex(const Array< OneD, const NekDouble > &gloCoord, NekDouble tol=0.0, bool returnNearestElmt=false)
Definition: ExpList.cpp:1262
virtual void v_Reset()
Reset geometry information, metrics, matrix managers and geometry information.
Definition: ExpList.cpp:1435
boost::shared_ptr< GlobalMatrix > GenGlobalMatrix(const GlobalMatrixKey &mkey, const boost::shared_ptr< AssemblyMapCG > &locToGloMap)
Generates a global matrix from the given key and map.
Definition: ExpList.cpp:895
virtual void v_SetHomo1DSpecVanVisc(Array< OneD, NekDouble > visc)
Definition: ExpList.h:1380
virtual void v_GetBoundaryNormals(int i, Array< OneD, Array< OneD, NekDouble > > &normals)
Definition: ExpList.cpp:2743
void GeneralMatrixOp_IterPerExp(const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:827
virtual const Array< OneD, const boost::shared_ptr< ExpList > > & v_GetBndCondExpansions(void)
Definition: ExpList.cpp:2255
boost::shared_ptr< Transposition > TranspositionSharedPtr
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:1926
ExpList()
The default constructor.
Definition: ExpList.cpp:93
Local coefficients.
static ExpListSharedPtr NullExpListSharedPtr
Definition: ExpList.h:1395
static Array< OneD, NekDouble > NullNekDouble1DArray
void EvaluateBoundaryConditions(const NekDouble time=0.0, const std::string varName="", const NekDouble=NekConstants::kNekUnsetDouble, const NekDouble=NekConstants::kNekUnsetDouble)
Definition: ExpList.h:2091
std::map< int, vector< PeriodicEntity > > PeriodicMap
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:1464
NekOptimize::GlobalOptParamSharedPtr m_globalOptParam
Definition: ExpList.h:1001
boost::shared_ptr< AssemblyMap > AssemblyMapSharedPtr
Definition: AssemblyMap.h:53
boost::shared_ptr< ExpList > & GetPlane(int n)
Definition: ExpList.h:887
boost::shared_ptr< BlockMatrixMap > BlockMatrixMapShPtr
A shared pointer to a BlockMatrixMap.
Definition: ExpList.h:97
LibUtilities::TranspositionSharedPtr GetTransposition(void)
This function returns the transposition class associaed with the homogeneous expansion.
Definition: ExpList.h:520
void MultiplyByInvMassMatrix(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
Definition: ExpList.h:1631
virtual boost::shared_ptr< AssemblyMapDG > & v_GetTraceMap()
Definition: ExpList.cpp:2299
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:1934
Array< OneD, NekDouble > & UpdateCoeffs()
This function returns (a reference to) the array (implemented as m_coeffs) containing all local expa...
Definition: ExpList.h:1953
virtual void v_FillBndCondFromField()
Definition: ExpList.cpp:2525
void ExtractFileBCs(const std::string &fileName, const std::string &varName, const boost::shared_ptr< ExpList > locExpList)
Definition: ExpList.cpp:1909
boost::shared_ptr< DNekMat > GenGlobalMatrixFull(const GlobalLinSysKey &mkey, const boost::shared_ptr< AssemblyMapCG > &locToGloMap)
Definition: ExpList.cpp:1032
NekDouble PhysIntegral(void)
This function integrates a function over the domain consisting of all the elements of the expansion...
Definition: ExpList.cpp:276
virtual void v_WriteTecplotZone(std::ostream &outfile, int expansion)
Definition: ExpList.cpp:1501
void GetBndElmtExpansion(int i, boost::shared_ptr< ExpList > &result)
Definition: ExpList.h:2150
ExpansionType GetExpType(void)
Returns the type of the expansion.
Definition: ExpList.cpp:243
const Array< OneD, const int > & GetTraceBndMap(void)
Definition: ExpList.h:2011
boost::shared_ptr< FieldDefinitions > FieldDefinitionsSharedPtr
Definition: FieldIO.h:118
virtual void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
Definition: ExpList.cpp:2558
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:1434
const boost::shared_ptr< LocalRegions::ExpansionVector > GetExp() const
This function returns the vector of elements in the expansion.
Definition: ExpList.h:1917
std::vector< LibUtilities::FieldDefinitionsSharedPtr > GetFieldDefinitions()
Definition: ExpList.h:794
virtual void v_BwdTrans_IterPerExp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:1229
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:1984
static BasisSharedPtr NullBasisSharedPtr
Definition: Basis.h:358
#define MULTI_REGIONS_EXPORT
boost::shared_ptr< ExpList > & UpdateBndCondExpansion(int i)
Definition: ExpList.h:1978
void BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
Definition: ExpList.h:1609
virtual void v_GetNormals(Array< OneD, Array< OneD, NekDouble > > &normals)
Definition: ExpList.cpp:2312
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:513
map< int, RobinBCInfoSharedPtr > GetRobinBCInfo()
Definition: ExpList.h:780
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:251
boost::shared_ptr< LibUtilities::SessionReader > GetSession()
Returns the session object.
Definition: ExpList.h:865
void ImposeDirichletConditions(Array< OneD, NekDouble > &outarray)
Impose Dirichlet Boundary Conditions onto Array.
Definition: ExpList.h:1839
void GetBoundaryToElmtMap(Array< OneD, int > &ElmtID, Array< OneD, int > &EdgeID)
Definition: ExpList.h:2144
BlockMatrixMapShPtr m_blockMat
Definition: ExpList.h:1003
Array< OneD, NekDouble > m_phys
The global expansion evaluated at the quadrature points.
Definition: ExpList.h:956
void LinearAdvectionDiffusionReactionSolve(const Array< OneD, Array< OneD, NekDouble > > &velocity, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const NekDouble lambda, CoeffState coeffstate=eLocal, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
Solve Advection Diffusion Reaction.
Definition: ExpList.h:1657
virtual const vector< bool > & v_GetLeftAdjacentFaces(void) const
Definition: ExpList.cpp:2361
void LinearAdvectionReactionSolve(const Array< OneD, Array< OneD, NekDouble > > &velocity, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const NekDouble lambda, CoeffState coeffstate=eLocal, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
Solve Advection Diffusion Reaction.
Definition: ExpList.h:1669
void ApplyGeomInfo()
Apply geometry information to each expansion.
Definition: ExpList.cpp:1422
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Definition: ExpList.h:939
virtual void v_WriteVtkPieceHeader(std::ostream &outfile, int expansion, int istrip)
Definition: ExpList.cpp:1720
NekDouble Integral(const Array< OneD, const NekDouble > &inarray)
Definition: ExpList.h:490
virtual void v_NormVectorIProductWRTBase(Array< OneD, const NekDouble > &V1, Array< OneD, const NekDouble > &V2, Array< OneD, NekDouble > &outarray, int BndID)
Definition: ExpList.cpp:2462
virtual LibUtilities::TranspositionSharedPtr v_GetTransposition(void)
Definition: ExpList.cpp:1859
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:1896
void PhysDeriv(Direction edir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d)
Definition: ExpList.h:1714
boost::shared_ptr< SessionReader > SessionReaderSharedPtr
Definition: MeshPartition.h:51
virtual boost::shared_ptr< ExpList > & v_UpdateBndCondExpansion(int i)
Definition: ExpList.cpp:2263
boost::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:1192
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...
virtual void v_ExtractElmtToBndPhys(int i, Array< OneD, NekDouble > &elmt, Array< OneD, NekDouble > &boundary)
Definition: ExpList.cpp:2663
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:557
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:563
virtual void v_WriteTecplotHeader(std::ostream &outfile, std::string var="")
Definition: ExpList.cpp:1462
virtual NekDouble v_GetHomoLen(void)
Definition: ExpList.cpp:1867
void WriteVtkFooter(std::ostream &outfile)
Definition: ExpList.cpp:1714
boost::shared_ptr< GlobalOptParam > GlobalOptParamSharedPtr
Pointer to a GlobalOptParam object.
void SetCoeffs(int i, NekDouble val)
Set the i th coefficiient in m_coeffs to value val.
Definition: ExpList.h:1815
virtual void v_GetCoords(Array< OneD, NekDouble > &coord_0, Array< OneD, NekDouble > &coord_1, Array< OneD, NekDouble > &coord_2=NullNekDouble1DArray)
Definition: ExpList.cpp:2593
const Array< OneD, const boost::shared_ptr< ExpList > > & GetBndCondExpansions()
Definition: ExpList.h:1973
void WriteTecplotField(std::ostream &outfile, int expansion=-1)
Definition: ExpList.h:389
map< GlobalMatrixKey, DNekScalBlkMatSharedPtr > BlockMatrixMap
A map between global matrix keys and their associated block matrices.
Definition: ExpList.h:95
boost::shared_ptr< ExpList > GetSharedThisPtr()
Returns a shared pointer to the current object.
Definition: ExpList.h:859
void 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.h:1983
virtual Array< OneD, const NekDouble > v_HomogeneousEnergy(void)
Definition: ExpList.cpp:1851
virtual const Array< OneD, const int > & v_GetTraceBndMap()
Definition: ExpList.cpp:2307
virtual void v_WriteVtkPieceData(std::ostream &outfile, int expansion, std::string var)
Definition: ExpList.cpp:1731
static std::vector< NekDouble > NullNekDoubleVector
Definition: FieldIO.h:49
void ExtractPhysToBndElmt(int i, const Array< OneD, const NekDouble > &phys, Array< OneD, NekDouble > &bndElmt)
Definition: ExpList.h:2163
virtual void v_ExtractPhysToBndElmt(int i, const Array< OneD, const NekDouble > &phys, Array< OneD, NekDouble > &bndElmt)
Definition: ExpList.cpp:2702
void GeneralMatrixOp(const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
This function calculates the result of the multiplication of a matrix of type specified by mkey with ...
Definition: ExpList.h:2129
void NormVectorIProductWRTBase(Array< OneD, const NekDouble > &V1, Array< OneD, const NekDouble > &V2, Array< OneD, NekDouble > &outarray, int BndID)
Definition: ExpList.h:1774
virtual void v_GetBCValues(Array< OneD, NekDouble > &BndVals, const Array< OneD, NekDouble > &TotField, int BndID)
Definition: ExpList.cpp:2454
void GetNormals(Array< OneD, Array< OneD, NekDouble > > &normals)
Definition: ExpList.h:2016
void Reset()
Reset geometry information and reset matrices.
Definition: ExpList.h:371
boost::shared_ptr< LibUtilities::Comm > GetComm()
Returns the comm object.
Definition: ExpList.h:871
Array< OneD, int > m_coeff_offset
Offset of elemental data into the array m_coeffs.
Definition: ExpList.h:988
LibUtilities::BasisSharedPtr GetHomogeneousBasis(void)
Definition: ExpList.h:882
const NekOptimize::GlobalOptParamSharedPtr & GetGlobalOptParam(void)
Definition: ExpList.h:775
int GetNumElmts(void)
This function returns the number of elements in the expansion which may be different for a homogeoeno...
Definition: ExpList.h:572
Base class for all multi-elemental spectral/hp expansions.
Definition: ExpList.h:101
boost::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:53
boost::shared_ptr< ExpList > & GetTrace()
Definition: ExpList.h:2001
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:1684
virtual void v_ImposeDirichletConditions(Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:2517
virtual LibUtilities::BasisSharedPtr v_GetHomogeneousBasis(void)
Definition: ExpList.h:1372
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:2271
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:1768
boost::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:977
int GetTotPoints(void) const
Returns the total number of quadrature points m_npoints .
Definition: ExpList.h:1453
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:1620
virtual void v_FwdTrans_IterPerExp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:602
virtual ~ExpList()
The default destructor.
Definition: ExpList.cpp:256
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:381
std::map< StdRegions::VarCoeffType, Array< OneD, NekDouble > > VarCoeffMap
Definition: StdRegions.hpp:226
boost::shared_ptr< GlobalLinSys > GenGlobalLinSys(const GlobalLinSysKey &mkey, const boost::shared_ptr< AssemblyMapCG > &locToGloMap)
This operation constructs the global linear system of type mkey.
Definition: ExpList.cpp:1173
void WriteVtkPieceData(std::ostream &outfile, int expansion, std::string var="v")
Definition: ExpList.h:414
void AddFwdBwdTraceIntegral(const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:2037
Array< OneD, int > m_phys_offset
Offset of elemental data into the array m_phys.
Definition: ExpList.h:991
void FillBndCondFromField(void)
Fill Bnd Condition expansion from the values stored in expansion.
Definition: ExpList.h:1845
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:1570
NekDouble GetHomoLen(void)
This function returns the Width of homogeneous direction associaed with the homogeneous expansion...
Definition: ExpList.h:527
virtual NekDouble v_Integral(const Array< OneD, const NekDouble > &inarray)
Definition: ExpList.cpp:1837
virtual std::vector< LibUtilities::FieldDefinitionsSharedPtr > v_GetFieldDefinitions(void)
Definition: ExpList.cpp:2092
bool m_physState
The state of the array m_phys.
Definition: ExpList.h:965
virtual void v_WriteTecplotConnectivity(std::ostream &outfile, int expansion)
Definition: ExpList.cpp:1595
static PeriodicMap NullPeriodicMap
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
Definition: ExpList.h:917
void WriteTecplotConnectivity(std::ostream &outfile, int expansion=-1)
Definition: ExpList.h:395
Array< OneD, int > m_offset_elmt_id
Array containing the element id m_offset_elmt_id[n] that the n^th consecutive block of data in m_coef...
Definition: ExpList.h:999
boost::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
Definition: NekTypeDefs.hpp:74
void ExtractTracePhys(Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:2065
boost::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
static std::vector< unsigned int > NullUnsignedIntVector
Definition: FieldIO.h:51
virtual void v_HomogeneousFwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
Definition: ExpList.cpp:2428
virtual void v_LocalToGlobal(void)
Definition: ExpList.cpp:2531
void SetWaveSpace(const bool wavespace)
Sets the wave space to the one of the possible configuration true or false.
Definition: ExpList.h:1494
void SetHomo1DSpecVanVisc(Array< OneD, NekDouble > visc)
This function sets the Spectral Vanishing Viscosity in homogeneous1D expansion.
Definition: ExpList.h:504
const vector< bool > & GetLeftAdjacentFaces(void) const
Definition: ExpList.h:2060
SpatialDomains::MeshGraphSharedPtr m_graph
Mesh associated with this expansion list.
Definition: ExpList.h:913
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:2827
LibUtilities::SessionReaderSharedPtr m_session
Session.
Definition: ExpList.h:910
SpatialDomains::MeshGraphSharedPtr GetGraph()
Definition: ExpList.h:876
double NekDouble
void WriteTecplotHeader(std::ostream &outfile, std::string var="")
Definition: ExpList.h:376
const DNekScalBlkMatSharedPtr GenBlockMatrix(const GlobalMatrixKey &gkey)
This function assembles the block diagonal matrix of local matrices of the type mtype.
Definition: ExpList.cpp:676
boost::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:68
Defines a list of flags.
std::map< int, BoundaryConditionMapShPtr > BoundaryConditionCollection
Definition: Conditions.h:215
void HelmSolve(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const FlagList &flags, const StdRegions::ConstFactorMap &factors, const StdRegions::VarCoeffMap &varcoeff=StdRegions::NullVarCoeffMap, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
Solve helmholtz problem.
Definition: ExpList.h:1642
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:1591
void GetPeriodicEntities(PeriodicMap &periodicVerts, PeriodicMap &periodicEdges, PeriodicMap &periodicFaces=NullPeriodicMap)
Definition: ExpList.h:785
virtual Array< OneD, SpatialDomains::BoundaryConditionShPtr > & v_UpdateBndConditions()
Definition: ExpList.cpp:2817
virtual void v_ClearGlobalLinSysManager(void)
Definition: ExpList.cpp:1903
static const NekDouble kNekUnsetDouble
Describe a linear system.
void FwdTrans_BndConstrained(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:612
virtual boost::shared_ptr< ExpList > & v_GetTrace()
Definition: ExpList.cpp:2291
void HomogeneousFwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
Definition: ExpList.h:1725
virtual void v_DealiasedProd(const Array< OneD, NekDouble > &inarray1, const Array< OneD, NekDouble > &inarray2, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
Definition: ExpList.cpp:2448
virtual void v_ExtractCoeffsToCoeffs(const boost::shared_ptr< ExpList > &fromExpList, const Array< OneD, const NekDouble > &fromCoeffs, Array< OneD, NekDouble > &toCoeffs)
Definition: ExpList.cpp:2232
virtual boost::shared_ptr< ExpList > & v_GetPlane(int n)
Definition: ExpList.cpp:2872
Describes a matrix with ordering defined by a local to global map.
static ExpList NullExpList
An empty ExpList object.
Definition: ExpList.h:1394
void LocalToGlobal(void)
Put the coefficients into global ordering using m_coeffs.
Definition: ExpList.h:1850
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:1887
int GetNpoints(void) const
Returns the total number of quadrature points m_npoints .
Definition: ExpList.h:1485
Array< OneD, SpatialDomains::BoundaryConditionShPtr > & UpdateBndConditions()
Definition: ExpList.h:2086
virtual void v_GetPeriodicEntities(PeriodicMap &periodicVerts, PeriodicMap &periodicEdges, PeriodicMap &periodicFaces)
Definition: ExpList.cpp:2849
virtual void v_PhysInterp1DScaled(const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:1892
void GlobalEigenSystem(const boost::shared_ptr< DNekMat > &Gmat, Array< OneD, NekDouble > &EigValsReal, Array< OneD, NekDouble > &EigValsImag, Array< OneD, NekDouble > &EigVecs=NullNekDouble1DArray)
Array< OneD, NekDouble > & UpdatePhys()
This function returns (a reference to) the array (implemented as m_phys) containing the function ev...
Definition: ExpList.h:1964
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:86
LibUtilities::CommSharedPtr m_comm
Communicator.
Definition: ExpList.h:907
virtual void v_AddTraceIntegral(const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:2319
void GetFieldDefinitions(std::vector< LibUtilities::FieldDefinitionsSharedPtr > &fielddef)
Definition: ExpList.h:800
void ExtractCoeffsToCoeffs(const boost::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:2143
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:497
virtual void v_LinearAdvectionReactionSolve(const Array< OneD, Array< OneD, NekDouble > > &velocity, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const NekDouble lambda, CoeffState coeffstate=eLocal, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
Definition: ExpList.cpp:2416
virtual Array< OneD, const unsigned int > v_GetZIDs(void)
Definition: ExpList.cpp:1875
void ExtractElmtToBndPhys(int i, Array< OneD, NekDouble > &elmt, Array< OneD, NekDouble > &boundary)
Definition: ExpList.h:2156
void DealiasedProd(const Array< OneD, NekDouble > &inarray1, const Array< OneD, NekDouble > &inarray2, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
Definition: ExpList.h:1751
bool GetPhysState(void) const
This function indicates whether the array of physical values (implemented as m_phys) is filled or no...
Definition: ExpList.h:1551
static const Array< OneD, ExpListSharedPtr > NullExpListSharedPtrArray
Definition: ExpList.h:2176
void FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
Definition: ExpList.h:1580
virtual map< int, RobinBCInfoSharedPtr > v_GetRobinBCInfo(void)
Definition: ExpList.cpp:2839
static Array< OneD, BasisSharedPtr > NullBasisSharedPtr1DArray
Definition: Basis.h:359
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:213
NekDouble GetCoeff(int i)
Get the i th value (coefficient) of m_coeffs.
Definition: ExpList.h:1865
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:1965
void SetCoeff(int i, NekDouble val)
Set the i th coefficiient in m_coeffs to value val.
Definition: ExpList.h:1805
int GetOffset_Elmt_Id(int n) const
Get the element id associated with the n th consecutive block of data in m_phys and m_coeffs...
Definition: ExpList.h:1942
virtual void v_IProductWRTBase_IterPerExp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:354
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:536
void IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
Definition: ExpList.h:1559
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:1542
void WriteVtkPieceFooter(std::ostream &outfile, int expansion)
Definition: ExpList.cpp:1725
void AppendFieldData(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata)
Append the element data listed in elements fielddef->m_ElementIDs onto fielddata. ...
Definition: ExpList.h:809
int GetNcoeffs(void) const
Returns the total number of local degrees of freedom .
Definition: ExpList.h:1403
virtual void v_LinearAdvectionDiffusionReactionSolve(const Array< OneD, Array< OneD, NekDouble > > &velocity, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const NekDouble lambda, CoeffState coeffstate=eLocal, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
Definition: ExpList.cpp:2404
virtual void v_HelmSolve(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const FlagList &flags, const StdRegions::ConstFactorMap &factors, const StdRegions::VarCoeffMap &varcoeff, const Array< OneD, const NekDouble > &dirForcing)
Definition: ExpList.cpp:2393
std::vector< int > m_coll_phys_offset
Offset of elemental data into the array m_phys.
Definition: ExpList.h:985
virtual void v_AppendFieldData(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata)
Definition: ExpList.cpp:2106
static SpatialDomains::BoundaryConditionShPtr GetBoundaryCondition(const SpatialDomains::BoundaryConditionCollection &collection, unsigned int index, const std::string &variable)
Definition: ExpList.cpp:2858
boost::shared_ptr< AssemblyMapDG > & GetTraceMap(void)
Definition: ExpList.h:2006
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:2950
const Array< OneD, const SpatialDomains::BoundaryConditionShPtr > & GetBndConditions()
Definition: ExpList.h:2079
void WriteVtkPieceHeader(std::ostream &outfile, int expansion, int istrip=0)
Definition: ExpList.h:404
void MultiplyByBlockMatrix(const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:324
void SetCoeffsArray(Array< OneD, NekDouble > &inarray)
Set the m_coeffs array to inarray.
Definition: ExpList.h:1821
virtual void v_MultiplyByInvMassMatrix(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
Definition: ExpList.cpp:2384
int EvalBasisNumModesMax(void) const
Evaulates the maximum number of modes in the elemental basis order over all elements.
Definition: ExpList.h:1417
void GetBoundaryNormals(int i, Array< OneD, Array< OneD, NekDouble > > &normals)
Definition: ExpList.h:2170
bool GetWaveSpace(void) const
This function returns the third direction expansion condition, which can be in wave space (coefficien...
Definition: ExpList.h:1502
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:545
virtual void v_GetBoundaryToElmtMap(Array< OneD, int > &ElmtID, Array< OneD, int > &EdgeID)
Definition: ExpList.cpp:2788
virtual Array< OneD, const unsigned int > v_GetYIDs(void)
Definition: ExpList.cpp:1883
virtual void v_SetUpPhysNormals()
Definition: ExpList.cpp:2646
void GetBCValues(Array< OneD, NekDouble > &BndVals, const Array< OneD, NekDouble > &TotField, int BndID)
Definition: ExpList.h:1763
int GetCoordim(int eid)
This function returns the dimension of the coordinates of the element eid.
Definition: ExpList.h:1794
boost::shared_ptr< BoundaryConditionBase > BoundaryConditionShPtr
Definition: Conditions.h:208
boost::shared_ptr< Basis > BasisSharedPtr
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:819
virtual const Array< OneD, const SpatialDomains::BoundaryConditionShPtr > & v_GetBndConditions()
Definition: ExpList.cpp:2806
void SetExpType(ExpansionType Type)
Returns the type of the expansion.
Definition: ExpList.cpp:251
void SetPhysArray(Array< OneD, NekDouble > &inarray)
Sets the array m_phys.
Definition: ExpList.h:1533
void WriteTecplotZone(std::ostream &outfile, int expansion=-1)
Definition: ExpList.h:382
std::vector< int > m_coll_coeff_offset
Offset of elemental data into the array m_coeffs.
Definition: ExpList.h:982
void WriteVtkHeader(std::ostream &outfile)
Definition: ExpList.cpp:1706
virtual void v_ExtractTracePhys(Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:2370
virtual void v_PhysGalerkinProjection1DScaled(const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:1898
boost::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:442
void GlobalToLocal(void)
Put the coefficients into local ordering and place in m_coeffs.
Definition: ExpList.h:1855
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:477
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
virtual void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
Definition: ExpList.cpp:2544
void GetFwdBwdTracePhys(Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
Definition: ExpList.h:2045
virtual void v_GetFwdBwdTracePhys(Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
Definition: ExpList.cpp:2345
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:2134
void HomogeneousBwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
Definition: ExpList.h:1738
void SetPhys(int i, NekDouble val)
Set the i th value of m_phys to value val.
Definition: ExpList.h:1508
virtual void v_WriteTecplotField(std::ostream &outfile, int expansion)
Definition: ExpList.cpp:1672
virtual int v_GetNumElmts(void)
Definition: ExpList.h:1056
Collections::CollectionVector m_collections
Definition: ExpList.h:979
virtual void v_ReadGlobalOptimizationParameters()
Definition: ExpList.cpp:2797
static VarCoeffMap NullVarCoeffMap
Definition: StdRegions.hpp:227
virtual void v_HomogeneousBwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
Definition: ExpList.cpp:2438
virtual void v_FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
Definition: ExpList.cpp:2551