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>
48 #include <boost/enable_shared_from_this.hpp>
50 
52 
53 #include <tinyxml.h>
54 
55 namespace Nektar
56 {
57  namespace MultiRegions
58  {
59  // Forward declarations
60  class GlobalLinSys;
61  class AssemblyMapDG;
62 
63  class AssemblyMapCG;
64  class GlobalLinSysKey;
65  class GlobalMatrix;
66 
67  enum Direction
68  {
69  eX,
70  eY,
71  eZ,
72  eS,
74  };
75 
77  {
78  e0D,
79  e1D,
80  e2D,
83  e3D,
85  };
86 
88  {
89  eX,
90  eY,
91  eZ
92  };
93 
94  /// A map between global matrix keys and their associated block
95  /// matrices.
96  typedef map<GlobalMatrixKey,DNekScalBlkMatSharedPtr> BlockMatrixMap;
97  /// A shared pointer to a BlockMatrixMap.
98  typedef boost::shared_ptr<BlockMatrixMap> BlockMatrixMapShPtr;
99 
100 
101  /// Base class for all multi-elemental spectral/hp expansions.
102  class ExpList: public boost::enable_shared_from_this<ExpList>
103  {
104  public:
105  /// The default constructor.
107 
108  /// The default constructor.
110  const LibUtilities::SessionReaderSharedPtr &pSession);
111 
112  /// The default constructor.
114  const LibUtilities::SessionReaderSharedPtr &pSession,
115  const SpatialDomains::MeshGraphSharedPtr &pGraph);
116 
117  /// The copy constructor.
119  const ExpList &in,
120  const bool DeclareCoeffPhysArrays = true);
121 
122  /// The default destructor.
123  MULTI_REGIONS_EXPORT virtual ~ExpList();
124 
125  /// Returns the total number of local degrees of freedom
126  /// \f$N_{\mathrm{eof}}=\sum_{e=1}^{{N_{\mathrm{el}}}}N^{e}_m\f$.
127  inline int GetNcoeffs(void) const;
128 
129  /// Returns the total number of local degrees of freedom
130  /// for element eid
131  MULTI_REGIONS_EXPORT int GetNcoeffs(const int eid) const;
132 
133  /// Returns the type of the expansion
135 
136  /// Returns the type of the expansion
138 
139  /// Evaulates the maximum number of modes in the elemental basis
140  /// order over all elements
141  inline int EvalBasisNumModesMax(void) const;
142 
143  /// Returns the vector of the number of modes in the elemental
144  /// basis order over all elements.
145  MULTI_REGIONS_EXPORT const Array<OneD,int>
146  EvalBasisNumModesMaxPerExp(void) const;
147 
148  /// Returns the total number of quadrature points #m_npoints
149  /// \f$=Q_{\mathrm{tot}}\f$.
150  inline int GetTotPoints(void) const;
151 
152  /// Returns the total number of quadrature points for eid's element
153  /// \f$=Q_{\mathrm{tot}}\f$.
154  inline int GetTotPoints(const int eid) const;
155 
156  /// Returns the total number of quadrature points #m_npoints
157  /// \f$=Q_{\mathrm{tot}}\f$.
158  inline int GetNpoints(void) const;
159 
160 
161  /// Returns the total number of qudature points scaled by
162  /// the factor scale on each 1D direction
163  inline int Get1DScaledTotPoints(const NekDouble scale) const;
164 
165  /// Sets the wave space to the one of the possible configuration
166  /// true or false
167  inline void SetWaveSpace(const bool wavespace);
168 
169 
170  ///Set Modified Basis for the stability analysis
171  inline void SetModifiedBasis(const bool modbasis);
172 
173  /// Set the \a i th value of \a m_phys to value \a val
174  inline void SetPhys(int i, NekDouble val);
175 
176  /// This function returns the third direction expansion condition,
177  /// which can be in wave space (coefficient) or not
178  /// It is stored in the variable m_WaveSpace.
179  inline bool GetWaveSpace(void) const;
180 
181  /// Fills the array #m_phys
182  inline void SetPhys(const Array<OneD, const NekDouble> &inarray);
183 
184  /// Sets the array #m_phys
185  inline void SetPhysArray(Array<OneD, NekDouble> &inarray);
186 
187  /// This function manually sets whether the array of physical
188  /// values \f$\boldsymbol{u}_l\f$ (implemented as #m_phys) is
189  /// filled or not.
190  inline void SetPhysState(const bool physState);
191 
192  /// This function indicates whether the array of physical values
193  /// \f$\boldsymbol{u}_l\f$ (implemented as #m_phys) is filled or
194  /// not.
195  inline bool GetPhysState(void) const;
196 
197  /// This function integrates a function \f$f(\boldsymbol{x})\f$
198  /// over the domain consisting of all the elements of the expansion.
200 
201  /// This function integrates a function \f$f(\boldsymbol{x})\f$
202  /// over the domain consisting of all the elements of the expansion.
204  const Array<OneD,
205  const NekDouble> &inarray);
206 
207  /// This function calculates the inner product of a function
208  /// \f$f(\boldsymbol{x})\f$ with respect to all \emph{local}
209  /// expansion modes \f$\phi_n^e(\boldsymbol{x})\f$.
210  inline void IProductWRTBase_IterPerExp(
211  const Array<OneD, const NekDouble> &inarray,
212  Array<OneD, NekDouble> &outarray);
213 
214  ///
215  inline void IProductWRTBase(
216  const Array<OneD, const NekDouble> &inarray,
217  Array<OneD, NekDouble> &outarray,
218  CoeffState coeffstate = eLocal);
219 
220  /// This function calculates the inner product of a function
221  /// \f$f(\boldsymbol{x})\f$ with respect to the derivative (in
222  /// direction \param dir) of all \emph{local} expansion modes
223  /// \f$\phi_n^e(\boldsymbol{x})\f$.
225  (const int dir,
226  const Array<OneD, const NekDouble> &inarray,
227  Array<OneD, NekDouble> &outarray);
228 
229  /// This function elementally evaluates the forward transformation
230  /// of a function \f$u(\boldsymbol{x})\f$ onto the global
231  /// spectral/hp expansion.
232  inline void FwdTrans_IterPerExp (
233  const Array<OneD,
234  const NekDouble> &inarray,
235  Array<OneD,NekDouble> &outarray);
236 
237  ///
238  inline void FwdTrans(
239  const Array<OneD,
240  const NekDouble> &inarray,
241  Array<OneD, NekDouble> &outarray,
242  CoeffState coeffstate = eLocal);
243 
244  /// This function elementally mulplies the coefficient space of
245  /// Sin my the elemental inverse of the mass matrix.
247  const Array<OneD,
248  const NekDouble> &inarray,
249  Array<OneD, NekDouble> &outarray);
250 
251  ///
252  inline void MultiplyByInvMassMatrix(
253  const Array<OneD,const NekDouble> &inarray,
254  Array<OneD, NekDouble> &outarray,
255  CoeffState coeffstate = eLocal);
256 
257  /// Smooth a field across elements
258  inline void SmoothField(Array<OneD,NekDouble> &field);
259 
260  /// Solve helmholtz problem
261  inline void HelmSolve(
262  const Array<OneD, const NekDouble> &inarray,
263  Array<OneD, NekDouble> &outarray,
264  const FlagList &flags,
265  const StdRegions::ConstFactorMap &factors,
266  const StdRegions::VarCoeffMap &varcoeff =
268  const Array<OneD, const NekDouble> &dirForcing =
270 
271  /// Solve Advection Diffusion Reaction
273  const Array<OneD, Array<OneD, NekDouble> > &velocity,
274  const Array<OneD, const NekDouble> &inarray,
275  Array<OneD, NekDouble> &outarray,
276  const NekDouble lambda,
277  CoeffState coeffstate = eLocal,
278  const Array<OneD, const NekDouble>&
279  dirForcing = NullNekDouble1DArray);
280 
281 
282  /// Solve Advection Diffusion Reaction
283  inline void LinearAdvectionReactionSolve(
284  const Array<OneD, Array<OneD, NekDouble> > &velocity,
285  const Array<OneD, const NekDouble> &inarray,
286  Array<OneD, NekDouble> &outarray,
287  const NekDouble lambda,
288  CoeffState coeffstate = eLocal,
289  const Array<OneD, const NekDouble>&
290  dirForcing = NullNekDouble1DArray);
291 
292  ///
294  const Array<OneD, const NekDouble> &inarray,
295  Array<OneD, NekDouble> &outarray);
296 
297 
298  /// This function elementally evaluates the backward transformation
299  /// of the global spectral/hp element expansion.
300  inline void BwdTrans_IterPerExp (
301  const Array<OneD, const NekDouble> &inarray,
302  Array<OneD,NekDouble> &outarray);
303 
304  ///
305  inline void BwdTrans (
306  const Array<OneD,
307  const NekDouble> &inarray,
308  Array<OneD,NekDouble> &outarray,
309  CoeffState coeffstate = eLocal);
310 
311  /// This function calculates the coordinates of all the elemental
312  /// quadrature points \f$\boldsymbol{x}_i\f$.
313  inline void GetCoords(
314  Array<OneD, NekDouble> &coord_0,
315  Array<OneD, NekDouble> &coord_1 = NullNekDouble1DArray,
316  Array<OneD, NekDouble> &coord_2 = NullNekDouble1DArray);
317 
318  // Homogeneous transforms
319  inline void HomogeneousFwdTrans(
320  const Array<OneD, const NekDouble> &inarray,
321  Array<OneD, NekDouble> &outarray,
322  CoeffState coeffstate = eLocal,
323  bool Shuff = true,
324  bool UnShuff = true);
325 
326  inline void HomogeneousBwdTrans(
327  const Array<OneD, const NekDouble> &inarray,
328  Array<OneD, NekDouble> &outarray,
329  CoeffState coeffstate = eLocal,
330  bool Shuff = true,
331  bool UnShuff = true);
332 
333  inline void DealiasedProd(
334  const Array<OneD, NekDouble> &inarray1,
335  const Array<OneD, NekDouble> &inarray2,
336  Array<OneD, NekDouble> &outarray,
337  CoeffState coeffstate = eLocal);
338 
339  inline void GetBCValues(
340  Array<OneD, NekDouble> &BndVals,
341  const Array<OneD, NekDouble> &TotField,
342  int BndID);
343 
344  inline void NormVectorIProductWRTBase(
345  Array<OneD, const NekDouble> &V1,
346  Array<OneD, const NekDouble> &V2,
347  Array<OneD, NekDouble> &outarray,
348  int BndID);
349 
350  /// Apply geometry information to each expansion.
352 
353  void WriteTecplotHeader(std::ofstream &outfile,
354  std::string var = "")
355  {
356  v_WriteTecplotHeader(outfile, var);
357  }
358 
360  std::ofstream &outfile,
361  int expansion = -1)
362  {
363  v_WriteTecplotZone(outfile, expansion);
364  }
365 
366  void WriteTecplotField(std::ofstream &outfile,
367  int expansion = -1)
368  {
369  v_WriteTecplotField(outfile, expansion);
370  }
371 
372  void WriteTecplotConnectivity(std::ofstream &outfile,
373  int expansion = -1)
374  {
375  v_WriteTecplotConnectivity(outfile, expansion);
376  }
377 
378  MULTI_REGIONS_EXPORT void WriteVtkHeader(std::ofstream &outfile);
379  MULTI_REGIONS_EXPORT void WriteVtkFooter(std::ofstream &outfile);
380 
381  void WriteVtkPieceHeader(std::ofstream &outfile, int expansion)
382  {
383  v_WriteVtkPieceHeader(outfile, expansion);
384  }
385 
387  std::ofstream &outfile,
388  int expansion);
389 
391  std::ofstream &outfile,
392  int expansion,
393  std::string var = "v")
394  {
395  v_WriteVtkPieceData(outfile, expansion, var);
396  }
397 
398  /// This function returns the dimension of the coordinates of the
399  /// element \a eid.
400  // inline
401  MULTI_REGIONS_EXPORT int GetCoordim(int eid);
402 
403  /// Set the \a i th coefficiient in \a m_coeffs to value \a val
404  inline void SetCoeff(int i, NekDouble val);
405 
406  /// Set the \a i th coefficiient in #m_coeffs to value \a val
407  inline void SetCoeffs(int i, NekDouble val);
408 
409  /// Set the #m_coeffs array to inarray
410  inline void SetCoeffsArray(Array<OneD, NekDouble> &inarray);
411 
412  /// This function returns (a reference to) the array
413  /// \f$\boldsymbol{\hat{u}}_l\f$ (implemented as #m_coeffs)
414  /// containing all local expansion coefficients.
415  inline const Array<OneD, const NekDouble> &GetCoeffs() const;
416 
417  /// Impose Dirichlet Boundary Conditions onto Array
418  inline void ImposeDirichletConditions(
419  Array<OneD,NekDouble>& outarray);
420 
421  /// Fill Bnd Condition expansion from the values stored in expansion
422  inline void FillBndCondFromField(void);
423 
424  /// Put the coefficients into global ordering using m_coeffs
425  inline void LocalToGlobal(void);
426 
427  /// Put the coefficients into local ordering and place in m_coeffs
428  inline void GlobalToLocal(void);
429 
430  /// Get the \a i th value (coefficient) of #m_coeffs
431  inline NekDouble GetCoeff(int i);
432 
433  /// Get the \a i th value (coefficient) of #m_coeffs
434  inline NekDouble GetCoeffs(int i);
435 
436  /// This function returns (a reference to) the array
437  /// \f$\boldsymbol{u}_l\f$ (implemented as #m_phys) containing the
438  /// function \f$u^{\delta}(\boldsymbol{x})\f$ evaluated at the
439  /// quadrature points.
440  // inline
441  MULTI_REGIONS_EXPORT const Array<OneD, const NekDouble>
442  &GetPhys() const;
443 
444  /// This function calculates the \f$L_\infty\f$ error of the global
445  /// spectral/hp element approximation.
447  const Array<OneD, const NekDouble> &inarray,
448  const Array<OneD, const NekDouble> &soln = NullNekDouble1DArray);
449 
450  /// This function calculates the \f$L_2\f$ error with
451  /// respect to soln of the global
452  /// spectral/hp element approximation.
454  const Array<OneD, const NekDouble> &inarray,
455  const Array<OneD, const NekDouble> &soln = NullNekDouble1DArray)
456  {
457  return v_L2(inarray, soln);
458  }
459 
460  /// Calculates the \f$H^1\f$ error of the global spectral/hp
461  /// element approximation.
463  const Array<OneD, const NekDouble> &inarray,
464  const Array<OneD, const NekDouble> &soln = NullNekDouble1DArray);
465 
466  NekDouble Integral (const Array<OneD, const NekDouble> &inarray)
467  {
468  return v_Integral(inarray);
469  }
470 
471  /// This function calculates the energy associated with
472  /// each one of the modesof a 3D homogeneous nD expansion
473  Array<OneD, const NekDouble> HomogeneousEnergy (void)
474  {
475  return v_HomogeneousEnergy();
476  }
477 
478  /// This function sets the Spectral Vanishing Viscosity
479  /// in homogeneous1D expansion.
480  void SetHomo1DSpecVanVisc(Array<OneD, NekDouble> visc)
481  {
483  }
484 
485  /// This function returns a vector containing the wave
486  /// numbers in z-direction associated
487  /// with the 3D homogenous expansion. Required if a
488  /// parellelisation is applied in the Fourier direction
489  Array<OneD, const unsigned int> GetZIDs(void)
490  {
491  return v_GetZIDs();
492  }
493 
494  /// This function returns the transposition class
495  /// associaed with the homogeneous expansion.
497  {
498  return v_GetTransposition();
499  }
500 
501  /// This function returns the Width of homogeneous direction
502  /// associaed with the homogeneous expansion.
504  {
505  return v_GetHomoLen();
506  }
507 
508  /// This function returns a vector containing the wave
509  /// numbers in y-direction associated
510  /// with the 3D homogenous expansion. Required if a
511  /// parellelisation is applied in the Fourier direction
512  Array<OneD, const unsigned int> GetYIDs(void)
513  {
514  return v_GetYIDs();
515  }
516 
517  /// This function interpolates the physical space points in
518  /// \a inarray to \a outarray using the same points defined in the
519  /// expansion but where the number of points are rescaled
520  /// by \a 1DScale
522  const NekDouble scale,
523  const Array<OneD, NekDouble> &inarray,
524  Array<OneD, NekDouble> &outarray)
525  {
526  v_PhysInterp1DScaled(scale, inarray,outarray);
527  }
528 
529  /// This function Galerkin projects the physical space points in
530  /// \a inarray to \a outarray where inarray is assumed to
531  /// be defined in the expansion but where the number of
532  /// points are rescaled by \a 1DScale
534  const NekDouble scale,
535  const Array<OneD, NekDouble> &inarray,
536  Array<OneD, NekDouble> &outarray)
537  {
538  v_PhysGalerkinProjection1DScaled(scale, inarray, outarray);
539  }
540 
541  /// This function returns the number of elements in the expansion.
542  inline int GetExpSize(void);
543 
544 
545  /// This function returns the number of elements in the
546  /// expansion which may be different for a homogeoenous extended
547  /// expansionp.
548  inline int GetNumElmts(void)
549  {
550  return v_GetNumElmts();
551  }
552 
553  /// This function returns the vector of elements in the expansion.
554  inline const boost::shared_ptr<LocalRegions::ExpansionVector>
555  GetExp() const;
556 
557  /// This function returns (a shared pointer to) the local elemental
558  /// expansion of the \f$n^{\mathrm{th}}\f$ element.
559  inline LocalRegions::ExpansionSharedPtr& GetExp(int n) const;
560 
561  /// This function returns (a shared pointer to) the local elemental
562  /// expansion containing the arbitrary point given by \a gloCoord.
564  const Array<OneD, const NekDouble> &gloCoord);
565 
566  /** This function returns the index of the local elemental
567  * expansion containing the arbitrary point given by \a gloCoord.
568  **/
570  const Array<OneD, const NekDouble> &gloCoord,
571  NekDouble tol = 0.0,
572  bool returnNearestElmt = false);
573 
574  /** This function returns the index and the Local
575  * Cartesian Coordinates \a locCoords of the local
576  * elemental expansion containing the arbitrary point
577  * given by \a gloCoords.
578  **/
580  const Array<OneD, const NekDouble> &gloCoords,
581  Array<OneD, NekDouble> &locCoords,
582  NekDouble tol = 0.0,
583  bool returnNearestElmt = false);
584 
585  /// Get the start offset position for a global list of #m_coeffs
586  /// correspoinding to element n.
587  inline int GetCoeff_Offset(int n) const;
588 
589  /// Get the start offset position for a global list of m_phys
590  /// correspoinding to element n.
591  inline int GetPhys_Offset(int n) const;
592 
593  /// Get the element id associated with the n th
594  /// consecutive block of data in #m_phys and #m_coeffs
595  inline int GetOffset_Elmt_Id(int n) const;
596 
597  /// This function returns (a reference to) the array
598  /// \f$\boldsymbol{\hat{u}}_l\f$ (implemented as #m_coeffs)
599  /// containing all local expansion coefficients.
600  inline Array<OneD, NekDouble> &UpdateCoeffs();
601 
602  /// This function returns (a reference to) the array
603  /// \f$\boldsymbol{u}_l\f$ (implemented as #m_phys) containing the
604  /// function \f$u^{\delta}(\boldsymbol{x})\f$ evaluated at the
605  /// quadrature points.
606  inline Array<OneD, NekDouble> &UpdatePhys();
607 
608  inline void PhysDeriv(
609  Direction edir,
610  const Array<OneD, const NekDouble> &inarray,
611  Array<OneD, NekDouble> &out_d);
612 
613  /// This function discretely evaluates the derivative of a function
614  /// \f$f(\boldsymbol{x})\f$ on the domain consisting of all
615  /// elements of the expansion.
616  inline void PhysDeriv(
617  const Array<OneD, const NekDouble> &inarray,
618  Array<OneD, NekDouble> &out_d0,
619  Array<OneD, NekDouble> &out_d1 = NullNekDouble1DArray,
620  Array<OneD, NekDouble> &out_d2 = NullNekDouble1DArray);
621 
622  inline void PhysDeriv(
623  const int dir,
624  const Array<OneD, const NekDouble> &inarray,
625  Array<OneD, NekDouble> &out_d);
626 
627 
628  // functions associated with DisContField
629  inline const Array<OneD, const boost::shared_ptr<ExpList> >
631 
632  inline boost::shared_ptr<ExpList> &UpdateBndCondExpansion(int i);
633 
634  inline void Upwind(
635  const Array<OneD, const Array<OneD, NekDouble> > &Vec,
636  const Array<OneD, const NekDouble> &Fwd,
637  const Array<OneD, const NekDouble> &Bwd,
638  Array<OneD, NekDouble> &Upwind);
639 
640  inline void Upwind(
641  const Array<OneD, const NekDouble> &Vn,
642  const Array<OneD, const NekDouble> &Fwd,
643  const Array<OneD, const NekDouble> &Bwd,
644  Array<OneD, NekDouble> &Upwind);
645 
646  /**
647  * Return a reference to the trace space associated with this
648  * expansion list.
649  */
650  inline boost::shared_ptr<ExpList> &GetTrace();
651 
652  inline boost::shared_ptr<AssemblyMapDG> &GetTraceMap(void);
653 
654  inline const Array<OneD, const int> &GetTraceBndMap(void);
655 
656  inline void GetNormals(Array<OneD, Array<OneD, NekDouble> > &normals);
657 
658  inline void AddTraceIntegral(
659  const Array<OneD, const NekDouble> &Fx,
660  const Array<OneD, const NekDouble> &Fy,
661  Array<OneD, NekDouble> &outarray);
662 
663  inline void AddTraceIntegral(
664  const Array<OneD, const NekDouble> &Fn,
665  Array<OneD, NekDouble> &outarray);
666 
667  inline void AddFwdBwdTraceIntegral(
668  const Array<OneD, const NekDouble> &Fwd,
669  const Array<OneD, const NekDouble> &Bwd,
670  Array<OneD, NekDouble> &outarray);
671 
672  inline void GetFwdBwdTracePhys(
673  Array<OneD,NekDouble> &Fwd,
674  Array<OneD,NekDouble> &Bwd);
675 
676  inline void GetFwdBwdTracePhys(
677  const Array<OneD,const NekDouble> &field,
678  Array<OneD,NekDouble> &Fwd,
679  Array<OneD,NekDouble> &Bwd);
680 
681  inline void ExtractTracePhys(Array<OneD,NekDouble> &outarray);
682 
683  inline void ExtractTracePhys(
684  const Array<OneD, const NekDouble> &inarray,
685  Array<OneD,NekDouble> &outarray);
686 
687  inline const Array<OneD, const SpatialDomains::
689 
690  inline Array<OneD, SpatialDomains::
692 
693  inline void EvaluateBoundaryConditions(
694  const NekDouble time = 0.0,
695  const std::string varName = "",
698 
699  // Routines for continous matrix solution
700  /// This function calculates the result of the multiplication of a
701  /// matrix of type specified by \a mkey with a vector given by \a
702  /// inarray.
703  inline void GeneralMatrixOp(
704  const GlobalMatrixKey &gkey,
705  const Array<OneD,const NekDouble> &inarray,
706  Array<OneD, NekDouble> &outarray,
707  CoeffState coeffstate = eLocal);
708 
710  const GlobalMatrixKey &gkey,
711  const Array<OneD,const NekDouble> &inarray,
712  Array<OneD, NekDouble> &outarray);
713 
714  inline void SetUpPhysNormals();
715 
716  inline void GetBoundaryToElmtMap(Array<OneD, int> &ElmtID,
717  Array<OneD,int> &EdgeID);
718 
720  std::vector<LibUtilities::FieldDefinitionsSharedPtr> &fielddef,
721  int NumHomoDir = 0,
722  Array<OneD, LibUtilities::BasisSharedPtr> &HomoBasis =
724  std::vector<NekDouble> &HomoLen =
726  std::vector<unsigned int> &HomoZIDs =
728  std::vector<unsigned int> &HomoYIDs =
730 
731 
733  {
734  return m_globalOptParam;
735  }
736 
737  map<int, RobinBCInfoSharedPtr> GetRobinBCInfo()
738  {
739  return v_GetRobinBCInfo();
740  }
741 
743  PeriodicMap &periodicVerts,
744  PeriodicMap &periodicEdges,
745  PeriodicMap &periodicFaces = NullPeriodicMap)
746  {
747  v_GetPeriodicEntities(periodicVerts, periodicEdges, periodicFaces);
748  }
749 
750  std::vector<LibUtilities::FieldDefinitionsSharedPtr>
752  {
753  return v_GetFieldDefinitions();
754  }
755 
756 
757  void GetFieldDefinitions(std::vector<LibUtilities::FieldDefinitionsSharedPtr> &fielddef)
758  {
759  v_GetFieldDefinitions(fielddef);
760  }
761 
762 
763 
764  /// Append the element data listed in elements
765  /// fielddef->m_ElementIDs onto fielddata
768  std::vector<NekDouble> &fielddata)
769  {
770  v_AppendFieldData(fielddef,fielddata);
771  }
772 
773 
774  /// Append the data in coeffs listed in elements
775  /// fielddef->m_ElementIDs onto fielddata
778  std::vector<NekDouble> &fielddata,
779  Array<OneD, NekDouble> &coeffs)
780  {
781  v_AppendFieldData(fielddef,fielddata,coeffs);
782  }
783 
784 
785  /** \brief Extract the data in fielddata into the coeffs
786  * using the basic ExpList Elemental expansions rather
787  * than planes in homogeneous case
788  */
791  std::vector<NekDouble> &fielddata,
792  std::string &field,
793  Array<OneD, NekDouble> &coeffs);
794 
795 
796  /** \brief Extract the data from fromField using
797  * fromExpList the coeffs using the basic ExpList
798  * Elemental expansions rather than planes in homogeneous
799  * case
800  */
802  const boost::shared_ptr<ExpList> &fromExpList,
803  const Array<OneD, const NekDouble> &fromCoeffs,
804  Array<OneD, NekDouble> &toCoeffs);
805 
806 
807  //Extract data in fielddata into the m_coeffs_list for the 3D stability analysis (base flow is 2D)
810  std::vector<NekDouble> &fielddata,
811  std::string &field,
812  Array<OneD, NekDouble> &coeffs);
813 
814 
815  /// Returns a shared pointer to the current object.
816  boost::shared_ptr<ExpList> GetSharedThisPtr()
817  {
818  return shared_from_this();
819  }
820 
821  /// Returns the session object
822  boost::shared_ptr<LibUtilities::SessionReader> GetSession()
823  {
824  return m_session;
825  }
826 
827  /// Returns the comm object
828  boost::shared_ptr<LibUtilities::Comm> GetComm()
829  {
830  return m_comm;
831  }
832 
834  {
835  return m_graph;
836  }
837 
838  // Wrapper functions for Homogeneous Expansions
840  {
841  return v_GetHomogeneousBasis();
842  }
843 
844  boost::shared_ptr<ExpList> &GetPlane(int n)
845  {
846  return v_GetPlane(n);
847  }
848 
849  //expansion type
851 
852  protected:
853  boost::shared_ptr<DNekMat> GenGlobalMatrixFull(
854  const GlobalLinSysKey &mkey,
855  const boost::shared_ptr<AssemblyMapCG> &locToGloMap);
856 
857  /// Communicator
859 
860  /// Session
862 
863  /// Mesh associated with this expansion list.
865 
866  /// The total number of local degrees of freedom. #m_ncoeffs
867  /// \f$=N_{\mathrm{eof}}=\sum_{e=1}^{{N_{\mathrm{el}}}}N^{e}_l\f$
869 
870  /** The total number of quadrature points. #m_npoints
871  *\f$=Q_{\mathrm{tot}}=\sum_{e=1}^{{N_{\mathrm{el}}}}N^{e}_Q\f$
872  **/
874 
875  /**
876  * \brief Concatenation of all local expansion coefficients.
877  *
878  * The array of length #m_ncoeffs\f$=N_{\mathrm{eof}}\f$ which is
879  * the concatenation of the local expansion coefficients
880  * \f$\hat{u}_n^e\f$ over all \f$N_{\mathrm{el}}\f$ elements
881  * \f[\mathrm{\texttt{m\_coeffs}}=\boldsymbol{\hat{u}}_{l} =
882  * \underline{\boldsymbol{\hat{u}}}^e = \left [ \begin{array}{c}
883  * \boldsymbol{\hat{u}}^{1} \ \
884  * \boldsymbol{\hat{u}}^{2} \ \
885  * \vdots \ \
886  * \boldsymbol{\hat{u}}^{{{N_{\mathrm{el}}}}} \end{array} \right ],
887  * \quad
888  * \mathrm{where}\quad \boldsymbol{\hat{u}}^{e}[n]=\hat{u}_n^{e}\f]
889  */
890  Array<OneD, NekDouble> m_coeffs;
891 
892  /**
893  * \brief The global expansion evaluated at the quadrature points
894  *
895  * The array of length #m_npoints\f$=Q_{\mathrm{tot}}\f$ containing
896  * the evaluation of \f$u^{\delta}(\boldsymbol{x})\f$ at the
897  * quadrature points \f$\boldsymbol{x}_i\f$.
898  * \f[\mathrm{\texttt{m\_phys}}=\boldsymbol{u}_{l} =
899  * \underline{\boldsymbol{u}}^e = \left [ \begin{array}{c}
900  * \boldsymbol{u}^{1} \ \
901  * \boldsymbol{u}^{2} \ \
902  * \vdots \ \
903  * \boldsymbol{u}^{{{N_{\mathrm{el}}}}} \end{array} \right ],\quad
904  * \mathrm{where}\quad
905  * \boldsymbol{u}^{e}[i]=u^{\delta}(\boldsymbol{x}_i)\f]
906  */
907  Array<OneD, NekDouble> m_phys;
908 
909  /**
910  * \brief The state of the array #m_phys.
911  *
912  * Indicates whether the array #m_phys, created to contain the
913  * evaluation of \f$u^{\delta}(\boldsymbol{x})\f$ at the quadrature
914  * points \f$\boldsymbol{x}_i\f$, is filled with these values.
915  */
917 
918  /**
919  * \brief The list of local expansions.
920  *
921  * The (shared pointer to the) vector containing (shared pointers
922  * to) all local expansions. The fact that the local expansions are
923  * all stored as a (pointer to a) #StdExpansion, the abstract base
924  * class for all local expansions, allows a general implementation
925  * where most of the routines for the derived classes are defined
926  * in the #ExpList base class.
927  */
928  boost::shared_ptr<LocalRegions::ExpansionVector> m_exp;
929 
930  /// Offset of elemental data into the array #m_coeffs
931  Array<OneD, int> m_coeff_offset;
932 
933  /// Offset of elemental data into the array #m_phys
934  Array<OneD, int> m_phys_offset;
935 
936  /// Array containing the element id #m_offset_elmt_id[n]
937  /// that the n^th consecutive block of data in #m_coeffs
938  /// and #m_phys is associated, i.e. for an array of
939  /// constant expansion size and single shape elements
940  /// m_phys[n*m_npoints] is the data related to
941  /// m_exp[m_offset_elmt_id[n]];
942  Array<OneD, int> m_offset_elmt_id;
943 
945 
947 
948  //@todo should this be in ExpList or ExpListHomogeneous1D.cpp
949  // it's a bool which determine if the expansion is in the wave space (coefficient space)
950  // or not
952 
953  /// This function assembles the block diagonal matrix of local
954  /// matrices of the type \a mtype.
956  const GlobalMatrixKey &gkey);
957 
959  const GlobalMatrixKey &gkey);
960 
962  const GlobalMatrixKey &gkey,
963  const Array<OneD,const NekDouble> &inarray,
964  Array<OneD, NekDouble> &outarray);
965 
966  /// Generates a global matrix from the given key and map.
967  boost::shared_ptr<GlobalMatrix> GenGlobalMatrix(
968  const GlobalMatrixKey &mkey,
969  const boost::shared_ptr<AssemblyMapCG> &locToGloMap);
970 
971 
972  void GlobalEigenSystem(
973  const boost::shared_ptr<DNekMat> &Gmat,
974  Array<OneD, NekDouble> &EigValsReal,
975  Array<OneD, NekDouble> &EigValsImag,
976  Array<OneD, NekDouble> &EigVecs
978 
979 
980  /// This operation constructs the global linear system of type \a
981  /// mkey.
982  boost::shared_ptr<GlobalLinSys> GenGlobalLinSys(
983  const GlobalLinSysKey &mkey,
984  const boost::shared_ptr<AssemblyMapCG> &locToGloMap);
985 
986  /// Generate a GlobalLinSys from information provided by the key
987  /// "mkey" and the mapping provided in LocToGloBaseMap.
988  boost::shared_ptr<GlobalLinSys> GenGlobalBndLinSys(
989  const GlobalLinSysKey &mkey,
990  const AssemblyMapSharedPtr &locToGloMap);
991 
993  {
995  }
996 
997  // Virtual prototypes
998 
999  virtual int v_GetNumElmts(void)
1000  {
1001  return (*m_exp).size();
1002  }
1003 
1004  virtual const Array<OneD,const boost::shared_ptr<ExpList> >
1005  &v_GetBndCondExpansions(void);
1006 
1007  virtual boost::shared_ptr<ExpList> &v_UpdateBndCondExpansion(int i);
1008 
1009  virtual void v_Upwind(
1010  const Array<OneD, const Array<OneD, NekDouble> > &Vec,
1011  const Array<OneD, const NekDouble> &Fwd,
1012  const Array<OneD, const NekDouble> &Bwd,
1013  Array<OneD, NekDouble> &Upwind);
1014 
1015  virtual void v_Upwind(
1016  const Array<OneD, const NekDouble> &Vn,
1017  const Array<OneD, const NekDouble> &Fwd,
1018  const Array<OneD, const NekDouble> &Bwd,
1019  Array<OneD, NekDouble> &Upwind);
1020 
1021  virtual boost::shared_ptr<ExpList> &v_GetTrace();
1022 
1023  virtual boost::shared_ptr<AssemblyMapDG> &v_GetTraceMap();
1024 
1025  virtual const Array<OneD, const int> &v_GetTraceBndMap();
1026 
1027  virtual void v_GetNormals(
1028  Array<OneD, Array<OneD, NekDouble> > &normals);
1029 
1030  virtual void v_AddTraceIntegral(
1031  const Array<OneD, const NekDouble> &Fx,
1032  const Array<OneD, const NekDouble> &Fy,
1033  Array<OneD, NekDouble> &outarray);
1034 
1035  virtual void v_AddTraceIntegral(
1036  const Array<OneD, const NekDouble> &Fn,
1037  Array<OneD, NekDouble> &outarray);
1038 
1039  virtual void v_AddFwdBwdTraceIntegral(
1040  const Array<OneD, const NekDouble> &Fwd,
1041  const Array<OneD, const NekDouble> &Bwd,
1042  Array<OneD, NekDouble> &outarray);
1043 
1044  virtual void v_GetFwdBwdTracePhys(
1045  Array<OneD,NekDouble> &Fwd,
1046  Array<OneD,NekDouble> &Bwd);
1047 
1048  virtual void v_GetFwdBwdTracePhys(
1049  const Array<OneD,const NekDouble> &field,
1050  Array<OneD,NekDouble> &Fwd,
1051  Array<OneD,NekDouble> &Bwd);
1052 
1053  virtual void v_ExtractTracePhys(
1054  Array<OneD,NekDouble> &outarray);
1055 
1056  virtual void v_ExtractTracePhys(
1057  const Array<OneD, const NekDouble> &inarray,
1058  Array<OneD,NekDouble> &outarray);
1059 
1060  virtual void v_MultiplyByInvMassMatrix(
1061  const Array<OneD,const NekDouble> &inarray,
1062  Array<OneD, NekDouble> &outarray,
1063  CoeffState coeffstate);
1064 
1065  virtual void v_HelmSolve(
1066  const Array<OneD, const NekDouble> &inarray,
1067  Array<OneD, NekDouble> &outarray,
1068  const FlagList &flags,
1069  const StdRegions::ConstFactorMap &factors,
1070  const StdRegions::VarCoeffMap &varcoeff,
1071  const Array<OneD, const NekDouble> &dirForcing);
1072 
1074  const Array<OneD, Array<OneD, NekDouble> > &velocity,
1075  const Array<OneD, const NekDouble> &inarray,
1076  Array<OneD, NekDouble> &outarray,
1077  const NekDouble lambda,
1078  CoeffState coeffstate = eLocal,
1079  const Array<OneD, const NekDouble>&
1080  dirForcing = NullNekDouble1DArray);
1081 
1082  virtual void v_LinearAdvectionReactionSolve(
1083  const Array<OneD, Array<OneD, NekDouble> > &velocity,
1084  const Array<OneD, const NekDouble> &inarray,
1085  Array<OneD, NekDouble> &outarray,
1086  const NekDouble lambda,
1087  CoeffState coeffstate = eLocal,
1088  const Array<OneD, const NekDouble>&
1089  dirForcing = NullNekDouble1DArray);
1090 
1091  // wrapper functions about virtual functions
1092  virtual void v_ImposeDirichletConditions(Array<OneD,NekDouble>& outarray);
1093 
1094  virtual void v_FillBndCondFromField();
1095 
1096 
1097  virtual void v_LocalToGlobal(void);
1098 
1099  virtual void v_GlobalToLocal(void);
1100 
1101  virtual void v_BwdTrans(
1102  const Array<OneD,const NekDouble> &inarray,
1103  Array<OneD, NekDouble> &outarray,
1104  CoeffState coeffstate);
1105 
1106  virtual void v_BwdTrans_IterPerExp(
1107  const Array<OneD,const NekDouble> &inarray,
1108  Array<OneD,NekDouble> &outarray);
1109 
1110  virtual void v_FwdTrans(
1111  const Array<OneD,const NekDouble> &inarray,
1112  Array<OneD, NekDouble> &outarray,
1113  CoeffState coeffstate);
1114 
1115  virtual void v_FwdTrans_IterPerExp(
1116  const Array<OneD,const NekDouble> &inarray,
1117  Array<OneD,NekDouble> &outarray);
1118 
1119  virtual void v_SmoothField(Array<OneD,NekDouble> &field);
1120 
1121  virtual void v_IProductWRTBase(
1122  const Array<OneD, const NekDouble> &inarray,
1123  Array<OneD, NekDouble> &outarray,
1124  CoeffState coeffstate);
1125 
1126  virtual void v_IProductWRTBase_IterPerExp(
1127  const Array<OneD,const NekDouble> &inarray,
1128  Array<OneD, NekDouble> &outarray);
1129 
1130  virtual void v_GeneralMatrixOp(
1131  const GlobalMatrixKey &gkey,
1132  const Array<OneD,const NekDouble> &inarray,
1133  Array<OneD, NekDouble> &outarray,
1134  CoeffState coeffstate);
1135 
1136  virtual void v_GetCoords(
1137  Array<OneD, NekDouble> &coord_0,
1138  Array<OneD, NekDouble> &coord_1,
1139  Array<OneD, NekDouble> &coord_2 = NullNekDouble1DArray);
1140 
1141  virtual void v_PhysDeriv(
1142  const Array<OneD, const NekDouble> &inarray,
1143  Array<OneD, NekDouble> &out_d0,
1144  Array<OneD, NekDouble> &out_d1,
1145  Array<OneD, NekDouble> &out_d2);
1146 
1147  virtual void v_PhysDeriv(
1148  const int dir,
1149  const Array<OneD, const NekDouble> &inarray,
1150  Array<OneD, NekDouble> &out_d);
1151 
1152  virtual void v_PhysDeriv(
1153  Direction edir,
1154  const Array<OneD, const NekDouble> &inarray,
1155  Array<OneD, NekDouble> &out_d);
1156 
1157  virtual void v_HomogeneousFwdTrans(
1158  const Array<OneD, const NekDouble> &inarray,
1159  Array<OneD, NekDouble> &outarray,
1160  CoeffState coeffstate = eLocal,
1161  bool Shuff = true,
1162  bool UnShuff = true);
1163 
1164  virtual void v_HomogeneousBwdTrans(
1165  const Array<OneD, const NekDouble> &inarray,
1166  Array<OneD, NekDouble> &outarray,
1167  CoeffState coeffstate = eLocal,
1168  bool Shuff = true,
1169  bool UnShuff = true);
1170 
1171  virtual void v_DealiasedProd(
1172  const Array<OneD, NekDouble> &inarray1,
1173  const Array<OneD, NekDouble> &inarray2,
1174  Array<OneD, NekDouble> &outarray,
1175  CoeffState coeffstate = eLocal);
1176 
1177  virtual void v_GetBCValues(
1178  Array<OneD, NekDouble> &BndVals,
1179  const Array<OneD, NekDouble> &TotField,
1180  int BndID);
1181 
1182  virtual void v_NormVectorIProductWRTBase(
1183  Array<OneD, const NekDouble> &V1,
1184  Array<OneD, const NekDouble> &V2,
1185  Array<OneD, NekDouble> &outarray,
1186  int BndID);
1187 
1188  virtual void v_SetUpPhysNormals();
1189 
1190  virtual void v_GetBoundaryToElmtMap(Array<OneD, int> &ElmtID,
1191  Array<OneD,int> &EdgeID);
1192 
1193  virtual void v_ReadGlobalOptimizationParameters();
1194 
1195  virtual std::vector<LibUtilities::FieldDefinitionsSharedPtr>
1196  v_GetFieldDefinitions(void);
1197 
1198  virtual void v_GetFieldDefinitions(
1199  std::vector<LibUtilities::FieldDefinitionsSharedPtr> &fielddef);
1200 
1201  virtual void v_AppendFieldData(
1203  std::vector<NekDouble> &fielddata);
1204 
1205  virtual void v_AppendFieldData(
1207  std::vector<NekDouble> &fielddata,
1208  Array<OneD, NekDouble> &coeffs);
1209 
1210  virtual void v_ExtractDataToCoeffs(
1212  std::vector<NekDouble> &fielddata, std::string &field,
1213  Array<OneD, NekDouble> &coeffs);
1214 
1215  virtual void v_ExtractCoeffsToCoeffs(const boost::shared_ptr<ExpList> &fromExpList, const Array<OneD, const NekDouble> &fromCoeffs, Array<OneD, NekDouble> &toCoeffs);
1216 
1217  virtual void v_WriteTecplotHeader(std::ofstream &outfile,
1218  std::string var = "");
1219  virtual void v_WriteTecplotZone(std::ofstream &outfile,
1220  int expansion);
1221  virtual void v_WriteTecplotField(std::ofstream &outfile,
1222  int expansion);
1223  virtual void v_WriteTecplotConnectivity(std::ofstream &outfile,
1224  int expansion);
1225  virtual void v_WriteVtkPieceHeader(std::ofstream &outfile, int expansion);
1226  virtual void v_WriteVtkPieceData(std::ofstream &outfile, int expansion,
1227  std::string var);
1228 
1229  virtual NekDouble v_L2(
1230  const Array<OneD, const NekDouble> &phys,
1231  const Array<OneD, const NekDouble> &soln = NullNekDouble1DArray);
1232 
1233  virtual NekDouble v_Integral (
1234  const Array<OneD, const NekDouble> &inarray);
1235 
1236  virtual Array<OneD, const NekDouble> v_HomogeneousEnergy(void);
1238  virtual NekDouble v_GetHomoLen(void);
1239  virtual Array<OneD, const unsigned int> v_GetZIDs(void);
1240  virtual Array<OneD, const unsigned int> v_GetYIDs(void);
1241 
1242  // 1D Scaling and projection
1243  virtual void v_PhysInterp1DScaled(
1244  const NekDouble scale, const Array<OneD, NekDouble> &inarray,
1245  Array<OneD, NekDouble> &outarray);
1246 
1247  virtual void v_PhysGalerkinProjection1DScaled(
1248  const NekDouble scale,
1249  const Array<OneD, NekDouble> &inarray,
1250  Array<OneD, NekDouble> &outarray);
1251 
1252  void ExtractFileBCs(const std::string &fileName,
1253  const std::string &varName,
1254  const boost::shared_ptr<ExpList> locExpList);
1255 
1256  // Utility function for a common case of retrieving a
1257  // BoundaryCondition from a boundary condition collection.
1260  GetBoundaryCondition(const SpatialDomains::
1261  BoundaryConditionCollection& collection,
1262  unsigned int index, const std::string& variable);
1263 
1264  private:
1265  virtual const Array<OneD,const SpatialDomains::BoundaryConditionShPtr> &v_GetBndConditions();
1266 
1267  virtual Array<OneD, SpatialDomains::BoundaryConditionShPtr>
1269 
1270  virtual void v_EvaluateBoundaryConditions(
1271  const NekDouble time = 0.0,
1272  const std::string varName = "",
1275 
1276  virtual map<int, RobinBCInfoSharedPtr> v_GetRobinBCInfo(void);
1277 
1278 
1279  virtual void v_GetPeriodicEntities(
1280  PeriodicMap &periodicVerts,
1281  PeriodicMap &periodicEdges,
1282  PeriodicMap &periodicFaces);
1283 
1284  // Homogeneous direction wrapper functions.
1286  {
1287  ASSERTL0(false,
1288  "This method is not defined or valid for this class type");
1290  }
1291 
1292  // wrapper function to set viscosity for Homo1D expansion
1293  virtual void v_SetHomo1DSpecVanVisc(Array<OneD, NekDouble> visc)
1294  {
1295  ASSERTL0(false,
1296  "This method is not defined or valid for this class type");
1297  }
1298 
1299 
1300  virtual boost::shared_ptr<ExpList> &v_GetPlane(int n);
1301  };
1302 
1303 
1304  /// Shared pointer to an ExpList object.
1305  typedef boost::shared_ptr<ExpList> ExpListSharedPtr;
1306  /// An empty ExpList object.
1309 
1310  // Inline routines follow.
1311 
1312  /**
1313  * Returns the total number of local degrees of freedom
1314  * \f$N_{\mathrm{eof}}=\sum_{e=1}^{{N_{\mathrm{el}}}}N^{e}_m\f$.
1315  */
1316  inline int ExpList::GetNcoeffs() const
1317  {
1318  return m_ncoeffs;
1319  }
1320 
1321  inline int ExpList::GetNcoeffs(const int eid) const
1322  {
1323  return (*m_exp)[eid]->GetNcoeffs();
1324  }
1325 
1326  /**
1327  * Evaulates the maximum number of modes in the elemental basis
1328  * order over all elements
1329  */
1331  {
1332  unsigned int i;
1333  int returnval = 0;
1334 
1335  for(i= 0; i < (*m_exp).size(); ++i)
1336  {
1337  returnval = max(returnval,
1338  (*m_exp)[i]->EvalBasisNumModesMax());
1339  }
1340 
1341  return returnval;
1342  }
1343 
1344  /**
1345  *
1346  */
1347  inline const Array<OneD,int> ExpList::EvalBasisNumModesMaxPerExp(void)
1348  const
1349  {
1350  unsigned int i;
1351  Array<OneD,int> returnval((*m_exp).size(),0);
1352 
1353  for(i= 0; i < (*m_exp).size(); ++i)
1354  {
1355  returnval[i]
1356  = max(returnval[i],(*m_exp)[i]->EvalBasisNumModesMax());
1357  }
1358 
1359  return returnval;
1360  }
1361 
1362 
1363  /**
1364  *
1365  */
1366  inline int ExpList::GetTotPoints() const
1367  {
1368  return m_npoints;
1369  }
1370 
1371  inline int ExpList::GetTotPoints(const int eid) const
1372  {
1373  return (*m_exp)[eid]->GetTotPoints();
1374  }
1375 
1376 
1377  inline int ExpList::Get1DScaledTotPoints(const NekDouble scale) const
1378  {
1379  int returnval = 0;
1380  int cnt;
1381  int nbase = (*m_exp)[0]->GetNumBases();
1382 
1383  for(int i = 0; i < (*m_exp).size(); ++i)
1384  {
1385  cnt = 1;
1386  for(int j = 0; j < nbase; ++j)
1387  {
1388  cnt *= (int)(scale*((*m_exp)[i]->GetNumPoints(j)));
1389  }
1390  returnval += cnt;
1391  }
1392  return returnval;
1393  }
1394 
1395  /**
1396  *
1397  */
1398  inline int ExpList::GetNpoints() const
1399  {
1400  return m_npoints;
1401  }
1402 
1403 
1404  /**
1405  *
1406  */
1407  inline void ExpList::SetWaveSpace(const bool wavespace)
1408  {
1409  m_WaveSpace = wavespace;
1410  }
1411 
1412  /**
1413  *
1414  */
1415  inline bool ExpList::GetWaveSpace() const
1416  {
1417  return m_WaveSpace;
1418  }
1419 
1420  /// Set the \a i th value of\a m_phys to value \a val
1421  inline void ExpList::SetPhys(int i, NekDouble val)
1422  {
1423  m_phys[i] = val;
1424  }
1425 
1426 
1427  /**
1428  * This function fills the array \f$\boldsymbol{u}_l\f$, the evaluation
1429  * of the expansion at the quadrature points (implemented as #m_phys),
1430  * with the values of the array \a inarray.
1431  *
1432  * @param inarray The array containing the values where
1433  * #m_phys should be filled with.
1434  */
1435  inline void ExpList::SetPhys(
1436  const Array<OneD, const NekDouble> &inarray)
1437  {
1438  ASSERTL0(inarray.num_elements() == m_npoints,
1439  "Input array does not have correct number of elements.");
1440 
1441  Vmath::Vcopy(m_npoints,&inarray[0],1,&m_phys[0],1);
1442  m_physState = true;
1443  }
1444 
1445 
1446  inline void ExpList::SetPhysArray(Array<OneD, NekDouble> &inarray)
1447  {
1448  m_phys = inarray;
1449  }
1450 
1451 
1452  /**
1453  * @param physState \a true (=filled) or \a false (=not filled).
1454  */
1455  inline void ExpList::SetPhysState(const bool physState)
1456  {
1457  m_physState = physState;
1458  }
1459 
1460 
1461  /**
1462  * @return physState \a true (=filled) or \a false (=not filled).
1463  */
1464  inline bool ExpList::GetPhysState() const
1465  {
1466  return m_physState;
1467  }
1468 
1469  /**
1470  *
1471  */
1473  const Array<OneD, const NekDouble> &inarray,
1474  Array<OneD, NekDouble> &outarray,
1475  CoeffState coeffstate)
1476  {
1477  v_IProductWRTBase(inarray,outarray, coeffstate);
1478  }
1479 
1480  /**
1481  *
1482  */
1484  const Array<OneD, const NekDouble> &inarray,
1485  Array<OneD, NekDouble> &outarray)
1486  {
1487  v_IProductWRTBase_IterPerExp(inarray,outarray);
1488  }
1489 
1490  /**
1491  *
1492  */
1493  inline void ExpList::FwdTrans(
1494  const Array<OneD, const NekDouble> &inarray,
1495  Array<OneD, NekDouble> &outarray,
1496  CoeffState coeffstate)
1497  {
1498  v_FwdTrans(inarray,outarray,coeffstate);
1499  }
1500 
1501  /**
1502  *
1503  */
1505  const Array<OneD, const NekDouble> &inarray,
1506  Array<OneD,NekDouble> &outarray)
1507  {
1508  v_FwdTrans_IterPerExp(inarray,outarray);
1509  }
1510 
1511  /**
1512  *
1513  */
1514  inline void ExpList::SmoothField(Array<OneD,NekDouble> &field)
1515  {
1516  v_SmoothField(field);
1517  }
1518 
1519  /**
1520  *
1521  */
1522  inline void ExpList::BwdTrans (
1523  const Array<OneD, const NekDouble> &inarray,
1524  Array<OneD, NekDouble> &outarray,
1525  CoeffState coeffstate)
1526  {
1527  v_BwdTrans(inarray,outarray,coeffstate);
1528  }
1529 
1530  /**
1531  *
1532  */
1534  const Array<OneD, const NekDouble> &inarray,
1535  Array<OneD, NekDouble> &outarray)
1536  {
1537  v_BwdTrans_IterPerExp(inarray,outarray);
1538  }
1539 
1540 
1541  /**
1542  *
1543  */
1545  const Array<OneD,const NekDouble> &inarray,
1546  Array<OneD, NekDouble> &outarray,
1547  CoeffState coeffstate)
1548  {
1549  v_MultiplyByInvMassMatrix(inarray,outarray,coeffstate);
1550  }
1551 
1552  /**
1553  *
1554  */
1555  inline void ExpList::HelmSolve(
1556  const Array<OneD, const NekDouble> &inarray,
1557  Array<OneD, NekDouble> &outarray,
1558  const FlagList &flags,
1559  const StdRegions::ConstFactorMap &factors,
1560  const StdRegions::VarCoeffMap &varcoeff,
1561  const Array<OneD, const NekDouble> &dirForcing)
1562  {
1563  v_HelmSolve(inarray, outarray, flags, factors, varcoeff, dirForcing);
1564  }
1565 
1566 
1567  /**
1568  *
1569  */
1571  const Array<OneD, Array<OneD, NekDouble> > &velocity,
1572  const Array<OneD, const NekDouble> &inarray,
1573  Array<OneD, NekDouble> &outarray,
1574  const NekDouble lambda,
1575  CoeffState coeffstate,
1576  const Array<OneD, const NekDouble>& dirForcing)
1577  {
1578  v_LinearAdvectionDiffusionReactionSolve(velocity,inarray, outarray,
1579  lambda, coeffstate,dirForcing);
1580  }
1581 
1583  const Array<OneD, Array<OneD, NekDouble> > &velocity,
1584  const Array<OneD, const NekDouble> &inarray,
1585  Array<OneD, NekDouble> &outarray,
1586  const NekDouble lambda,
1587  CoeffState coeffstate,
1588  const Array<OneD, const NekDouble>& dirForcing)
1589  {
1590  v_LinearAdvectionReactionSolve(velocity,inarray, outarray,
1591  lambda, coeffstate,dirForcing);
1592  }
1593 
1594  /**
1595  *
1596  */
1597  inline void ExpList::GetCoords(Array<OneD, NekDouble> &coord_0,
1598  Array<OneD, NekDouble> &coord_1,
1599  Array<OneD, NekDouble> &coord_2)
1600 
1601  {
1602  v_GetCoords(coord_0,coord_1,coord_2);
1603  }
1604 
1605  /**
1606  *
1607  */
1608  inline void ExpList::PhysDeriv(const Array<OneD, const NekDouble> &inarray,
1609  Array<OneD, NekDouble> &out_d0,
1610  Array<OneD, NekDouble> &out_d1,
1611  Array<OneD, NekDouble> &out_d2)
1612  {
1613  v_PhysDeriv(inarray,out_d0,out_d1,out_d2);
1614  }
1615 
1616  /**
1617  *
1618  */
1619  inline void ExpList::PhysDeriv(
1620  const int dir,
1621  const Array<OneD, const NekDouble> &inarray,
1622  Array<OneD, NekDouble> &out_d)
1623  {
1624  v_PhysDeriv(dir,inarray,out_d);
1625  }
1626 
1627  inline void ExpList::PhysDeriv(
1628  Direction edir,
1629  const Array<OneD, const NekDouble> &inarray,
1630  Array<OneD, NekDouble> &out_d)
1631  {
1632  v_PhysDeriv(edir, inarray,out_d);
1633  }
1634 
1635  /**
1636  *
1637  */
1639  const Array<OneD, const NekDouble> &inarray,
1640  Array<OneD, NekDouble> &outarray,
1641  CoeffState coeffstate,
1642  bool Shuff,
1643  bool UnShuff)
1644  {
1645  v_HomogeneousFwdTrans(inarray,outarray,coeffstate,Shuff,UnShuff);
1646  }
1647 
1648  /**
1649  *
1650  */
1652  const Array<OneD, const NekDouble> &inarray,
1653  Array<OneD, NekDouble> &outarray,
1654  CoeffState coeffstate,
1655  bool Shuff,
1656  bool UnShuff)
1657  {
1658  v_HomogeneousBwdTrans(inarray,outarray,coeffstate,Shuff,UnShuff);
1659  }
1660 
1661  /**
1662  *
1663  */
1665  const Array<OneD, NekDouble> &inarray1,
1666  const Array<OneD, NekDouble> &inarray2,
1667  Array<OneD, NekDouble> &outarray,
1668  CoeffState coeffstate)
1669  {
1670  v_DealiasedProd(inarray1,inarray2,outarray,coeffstate);
1671  }
1672 
1673  /**
1674  *
1675  */
1677  Array<OneD, NekDouble> &BndVals,
1678  const Array<OneD, NekDouble> &TotField,
1679  int BndID)
1680  {
1681  v_GetBCValues(BndVals,TotField,BndID);
1682  }
1683 
1684  /**
1685  *
1686  */
1688  Array<OneD, const NekDouble> &V1,
1689  Array<OneD, const NekDouble> &V2,
1690  Array<OneD, NekDouble> &outarray,
1691  int BndID)
1692  {
1693  v_NormVectorIProductWRTBase(V1,V2,outarray,BndID);
1694  }
1695 
1696  /**
1697  * @param eid The index of the element to be checked.
1698  * @return The dimension of the coordinates of the specific element.
1699  */
1700  inline int ExpList::GetCoordim(int eid)
1701  {
1702  ASSERTL2(eid <= (*m_exp).size(),
1703  "eid is larger than number of elements");
1704  return (*m_exp)[eid]->GetCoordim();
1705  }
1706 
1707  /**
1708  * @param i The index of m_coeffs to be set
1709  * @param val The value which m_coeffs[i] is to be set to.
1710  */
1711  inline void ExpList::SetCoeff(int i, NekDouble val)
1712  {
1713  m_coeffs[i] = val;
1714  }
1715 
1716 
1717  /**
1718  * @param i The index of #m_coeffs to be set.
1719  * @param val The value which #m_coeffs[i] is to be set to.
1720  */
1721  inline void ExpList::SetCoeffs(int i, NekDouble val)
1722  {
1723  m_coeffs[i] = val;
1724  }
1725 
1726 
1727  inline void ExpList::SetCoeffsArray(Array<OneD, NekDouble> &inarray)
1728  {
1729  m_coeffs = inarray;
1730  }
1731 
1732  /**
1733  * As the function returns a constant reference to a
1734  * <em>const Array</em>, it is not possible to modify the
1735  * underlying data of the array #m_coeffs. In order to do
1736  * so, use the function #UpdateCoeffs instead.
1737  *
1738  * @return (A constant reference to) the array #m_coeffs.
1739  */
1740  inline const Array<OneD, const NekDouble> &ExpList::GetCoeffs() const
1741  {
1742  return m_coeffs;
1743  }
1744 
1746  Array<OneD,NekDouble>& outarray)
1747  {
1748  v_ImposeDirichletConditions(outarray);
1749  }
1750 
1752  {
1754  }
1755 
1756  inline void ExpList::LocalToGlobal(void)
1757  {
1758  v_LocalToGlobal();
1759  }
1760 
1761  inline void ExpList::GlobalToLocal(void)
1762  {
1763  v_GlobalToLocal();
1764  }
1765 
1766 
1767  /**
1768  * @param i The index of #m_coeffs to be returned
1769  * @return The NekDouble held in #m_coeffs[i].
1770  */
1772  {
1773  return m_coeffs[i];
1774  }
1775 
1776  /**
1777  * @param i The index of #m_coeffs to be returned
1778  * @return The NekDouble held in #m_coeffs[i].
1779  */
1781  {
1782  return m_coeffs[i];
1783  }
1784 
1785  /**
1786  * As the function returns a constant reference to a
1787  * <em>const Array</em> it is not possible to modify the
1788  * underlying data of the array #m_phys. In order to do
1789  * so, use the function #UpdatePhys instead.
1790  *
1791  * @return (A constant reference to) the array #m_phys.
1792  */
1793  inline const Array<OneD, const NekDouble> &ExpList::GetPhys() const
1794  {
1795  return m_phys;
1796  }
1797 
1798  /**
1799  * @return \f$N_{\mathrm{el}}\f$, the number of elements in the
1800  * expansion.
1801  */
1802  inline int ExpList::GetExpSize(void)
1803  {
1804  return (*m_exp).size();
1805  }
1806 
1807 
1808  /**
1809  * @param n The index of the element concerned.
1810  *
1811  * @return (A shared pointer to) the local expansion of the
1812  * \f$n^{\mathrm{th}}\f$ element.
1813  */
1815  {
1816  return (*m_exp)[n];
1817  }
1818 
1819  /**
1820  * @return (A const shared pointer to) the local expansion vector #m_exp
1821  */
1822  inline const boost::shared_ptr<LocalRegions::ExpansionVector>
1823  ExpList::GetExp(void) const
1824  {
1825  return m_exp;
1826  }
1827 
1828 
1829  /**
1830  *
1831  */
1832  inline int ExpList::GetCoeff_Offset(int n) const
1833  {
1834  return m_coeff_offset[n];
1835  }
1836 
1837  /**
1838  *
1839  */
1840  inline int ExpList::GetPhys_Offset(int n) const
1841  {
1842  return m_phys_offset[n];
1843  }
1844 
1845  /**
1846  *
1847  */
1848  inline int ExpList::GetOffset_Elmt_Id(int n) const
1849  {
1850  return m_offset_elmt_id[n];
1851  }
1852 
1853  /**
1854  * If one wants to get hold of the underlying data without modifying
1855  * them, rather use the function #GetCoeffs instead.
1856  *
1857  * @return (A reference to) the array #m_coeffs.
1858  */
1859  inline Array<OneD, NekDouble> &ExpList::UpdateCoeffs()
1860  {
1861  return m_coeffs;
1862  }
1863 
1864  /**
1865  * If one wants to get hold of the underlying data without modifying
1866  * them, rather use the function #GetPhys instead.
1867  *
1868  * @return (A reference to) the array #m_phys.
1869  */
1870  inline Array<OneD, NekDouble> &ExpList::UpdatePhys()
1871  {
1872  m_physState = true;
1873  return m_phys;
1874  }
1875 
1876 
1877  // functions associated with DisContField
1878  inline const Array<OneD, const boost::shared_ptr<ExpList> >
1880  {
1881  return v_GetBndCondExpansions();
1882  }
1883 
1884  inline boost::shared_ptr<ExpList> &ExpList::UpdateBndCondExpansion(int i)
1885  {
1886  return v_UpdateBndCondExpansion(i);
1887  }
1888 
1889  inline void ExpList::Upwind(
1890  const Array<OneD, const Array<OneD, NekDouble> > &Vec,
1891  const Array<OneD, const NekDouble> &Fwd,
1892  const Array<OneD, const NekDouble> &Bwd,
1893  Array<OneD, NekDouble> &Upwind)
1894  {
1895  v_Upwind(Vec, Fwd, Bwd, Upwind);
1896  }
1897 
1898  inline void ExpList::Upwind(
1899  const Array<OneD, const NekDouble> &Vn,
1900  const Array<OneD, const NekDouble> &Fwd,
1901  const Array<OneD, const NekDouble> &Bwd,
1902  Array<OneD, NekDouble> &Upwind)
1903  {
1904  v_Upwind(Vn, Fwd, Bwd, Upwind);
1905  }
1906 
1907  inline boost::shared_ptr<ExpList> &ExpList::GetTrace()
1908  {
1909  return v_GetTrace();
1910  }
1911 
1912  inline boost::shared_ptr<AssemblyMapDG> &ExpList::GetTraceMap()
1913  {
1914  return v_GetTraceMap();
1915  }
1916 
1917  inline const Array<OneD, const int> &ExpList::GetTraceBndMap()
1918  {
1919  return v_GetTraceBndMap();
1920  }
1921 
1922  inline void ExpList::GetNormals(
1923  Array<OneD, Array<OneD, NekDouble> > &normals)
1924  {
1925  v_GetNormals(normals);
1926  }
1927 
1929  const Array<OneD, const NekDouble> &Fx,
1930  const Array<OneD, const NekDouble> &Fy,
1931  Array<OneD, NekDouble> &outarray)
1932  {
1933  v_AddTraceIntegral(Fx,Fy,outarray);
1934  }
1935 
1937  const Array<OneD, const NekDouble> &Fn,
1938  Array<OneD, NekDouble> &outarray)
1939  {
1940  v_AddTraceIntegral(Fn,outarray);
1941  }
1942 
1944  const Array<OneD, const NekDouble> &Fwd,
1945  const Array<OneD, const NekDouble> &Bwd,
1946  Array<OneD, NekDouble> &outarray)
1947  {
1948  v_AddFwdBwdTraceIntegral(Fwd,Bwd,outarray);
1949  }
1950 
1952  Array<OneD,NekDouble> &Fwd,
1953  Array<OneD,NekDouble> &Bwd)
1954  {
1955  v_GetFwdBwdTracePhys(Fwd,Bwd);
1956  }
1957 
1959  const Array<OneD,const NekDouble> &field,
1960  Array<OneD,NekDouble> &Fwd,
1961  Array<OneD,NekDouble> &Bwd)
1962  {
1963  v_GetFwdBwdTracePhys(field,Fwd,Bwd);
1964  }
1965 
1966  inline void ExpList::ExtractTracePhys(Array<OneD,NekDouble> &outarray)
1967  {
1968  v_ExtractTracePhys(outarray);
1969  }
1970 
1971 
1973  const Array<OneD, const NekDouble> &inarray,
1974  Array<OneD,NekDouble> &outarray)
1975  {
1976  v_ExtractTracePhys(inarray,outarray);
1977  }
1978 
1979  inline const Array<OneD,const SpatialDomains::BoundaryConditionShPtr>
1981  {
1982  return v_GetBndConditions();
1983  }
1984 
1985 
1986  inline Array<OneD, SpatialDomains::BoundaryConditionShPtr>
1988  {
1989  return v_UpdateBndConditions();
1990  }
1991 
1993  const NekDouble time,
1994  const std::string varName,
1995  const NekDouble x2_in,
1996  const NekDouble x3_in)
1997  {
1998  v_EvaluateBoundaryConditions(time, varName, x2_in, x3_in);
1999  }
2000 
2001  // Routines for continous matrix solution
2002  /**
2003  * This operation is equivalent to the evaluation of
2004  * \f$\underline{\boldsymbol{M}}^e\boldsymbol{\hat{u}}_l\f$, that is,
2005  * \f[ \left[
2006  * \begin{array}{cccc}
2007  * \boldsymbol{M}^1 & 0 & \hspace{3mm}0 \hspace{3mm}& 0 \\
2008  * 0 & \boldsymbol{M}^2 & 0 & 0 \\
2009  * 0 & 0 & \ddots & 0 \\
2010  * 0 & 0 & 0 & \boldsymbol{M}^{N_{\mathrm{el}}} \end{array} \right]
2011  *\left [ \begin{array}{c}
2012  * \boldsymbol{\hat{u}}^{1} \\
2013  * \boldsymbol{\hat{u}}^{2} \\
2014  * \vdots \\
2015  * \boldsymbol{\hat{u}}^{{{N_{\mathrm{el}}}}} \end{array} \right ]\f]
2016  * where \f$\boldsymbol{M}^e\f$ are the local matrices of type
2017  * specified by the key \a mkey. The decoupling of the local matrices
2018  * allows for a local evaluation of the operation. However, rather than
2019  * a local matrix-vector multiplication, the local operations are
2020  * evaluated as implemented in the function
2021  * StdRegions#StdExpansion#GeneralMatrixOp.
2022  *
2023  * @param mkey This key uniquely defines the type matrix
2024  * required for the operation.
2025  * @param inarray The vector \f$\boldsymbol{\hat{u}}_l\f$ of
2026  * size \f$N_{\mathrm{eof}}\f$.
2027  * @param outarray The resulting vector of size
2028  * \f$N_{\mathrm{eof}}\f$.
2029  */
2031  const GlobalMatrixKey &gkey,
2032  const Array<OneD,const NekDouble> &inarray,
2033  Array<OneD, NekDouble> &outarray,
2034  CoeffState coeffstate)
2035  {
2036  v_GeneralMatrixOp(gkey,inarray,outarray,coeffstate);
2037  }
2038 
2039 
2041  {
2043  }
2044 
2045  inline void ExpList::GetBoundaryToElmtMap( Array<OneD, int> &ElmtID,
2046  Array<OneD,int> &EdgeID)
2047  {
2048  v_GetBoundaryToElmtMap(ElmtID,EdgeID);
2049  }
2050 
2051  const static Array<OneD, ExpListSharedPtr> NullExpListSharedPtrArray;
2052 
2053  } //end of namespace
2054 } //end of namespace
2055 
2056 #endif // EXPLIST_H
2057