Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties 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 std::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  const bool PhysSpaceForcing = true);
284 
285  /// Solve Advection Diffusion Reaction
287  const Array<OneD, Array<OneD, NekDouble> > &velocity,
288  const Array<OneD, const NekDouble> &inarray,
289  Array<OneD, NekDouble> &outarray,
290  const NekDouble lambda,
291  CoeffState coeffstate = eLocal,
293  dirForcing = NullNekDouble1DArray);
294 
295 
296  /// Solve Advection Diffusion Reaction
297  inline void LinearAdvectionReactionSolve(
298  const Array<OneD, Array<OneD, NekDouble> > &velocity,
299  const Array<OneD, const NekDouble> &inarray,
300  Array<OneD, NekDouble> &outarray,
301  const NekDouble lambda,
302  CoeffState coeffstate = eLocal,
304  dirForcing = NullNekDouble1DArray);
305 
306  ///
308  const Array<OneD, const NekDouble> &inarray,
309  Array<OneD, NekDouble> &outarray);
310 
311 
312  /// This function elementally evaluates the backward transformation
313  /// of the global spectral/hp element expansion.
314  inline void BwdTrans_IterPerExp (
315  const Array<OneD, const NekDouble> &inarray,
316  Array<OneD,NekDouble> &outarray);
317 
318  ///
319  inline void BwdTrans (
320  const Array<OneD,
321  const NekDouble> &inarray,
322  Array<OneD,NekDouble> &outarray,
323  CoeffState coeffstate = eLocal);
324 
325  /// This function calculates the coordinates of all the elemental
326  /// quadrature points \f$\boldsymbol{x}_i\f$.
327  inline void GetCoords(
328  Array<OneD, NekDouble> &coord_0,
331 
332  // Homogeneous transforms
333  inline void HomogeneousFwdTrans(
334  const Array<OneD, const NekDouble> &inarray,
335  Array<OneD, NekDouble> &outarray,
336  CoeffState coeffstate = eLocal,
337  bool Shuff = true,
338  bool UnShuff = true);
339 
340  inline void HomogeneousBwdTrans(
341  const Array<OneD, const NekDouble> &inarray,
342  Array<OneD, NekDouble> &outarray,
343  CoeffState coeffstate = eLocal,
344  bool Shuff = true,
345  bool UnShuff = true);
346 
347  inline void DealiasedProd(
348  const Array<OneD, NekDouble> &inarray1,
349  const Array<OneD, NekDouble> &inarray2,
350  Array<OneD, NekDouble> &outarray,
351  CoeffState coeffstate = eLocal);
352 
353  inline void DealiasedDotProd(
354  const Array<OneD, Array<OneD, NekDouble> > &inarray1,
355  const Array<OneD, Array<OneD, NekDouble> > &inarray2,
356  Array<OneD, Array<OneD, NekDouble> > &outarray,
357  CoeffState coeffstate = eLocal);
358 
359  inline void GetBCValues(
360  Array<OneD, NekDouble> &BndVals,
361  const Array<OneD, NekDouble> &TotField,
362  int BndID);
363 
364  inline void NormVectorIProductWRTBase(
367  Array<OneD, NekDouble> &outarray,
368  int BndID);
369 
370  inline void NormVectorIProductWRTBase(
372  Array<OneD, NekDouble> &outarray);
373 
374  /// Apply geometry information to each expansion.
376 
377  /// Reset geometry information and reset matrices
379  {
380  v_Reset();
381  }
382 
383  void WriteTecplotHeader(std::ostream &outfile,
384  std::string var = "")
385  {
386  v_WriteTecplotHeader(outfile, var);
387  }
388 
390  std::ostream &outfile,
391  int expansion = -1)
392  {
393  v_WriteTecplotZone(outfile, expansion);
394  }
395 
396  void WriteTecplotField(std::ostream &outfile,
397  int expansion = -1)
398  {
399  v_WriteTecplotField(outfile, expansion);
400  }
401 
402  void WriteTecplotConnectivity(std::ostream &outfile,
403  int expansion = -1)
404  {
405  v_WriteTecplotConnectivity(outfile, expansion);
406  }
407 
408  MULTI_REGIONS_EXPORT void WriteVtkHeader(std::ostream &outfile);
409  MULTI_REGIONS_EXPORT void WriteVtkFooter(std::ostream &outfile);
410 
411  void WriteVtkPieceHeader(std::ostream &outfile, int expansion,
412  int istrip = 0)
413  {
414  v_WriteVtkPieceHeader(outfile, expansion, istrip);
415  }
416 
418  std::ostream &outfile,
419  int expansion);
420 
422  std::ostream &outfile,
423  int expansion,
424  std::string var = "v")
425  {
426  v_WriteVtkPieceData(outfile, expansion, var);
427  }
428 
429  /// This function returns the dimension of the coordinates of the
430  /// element \a eid.
431  // inline
432  MULTI_REGIONS_EXPORT int GetCoordim(int eid);
433 
434  /// Set the \a i th coefficiient in \a m_coeffs to value \a val
435  inline void SetCoeff(int i, NekDouble val);
436 
437  /// Set the \a i th coefficiient in #m_coeffs to value \a val
438  inline void SetCoeffs(int i, NekDouble val);
439 
440  /// Set the #m_coeffs array to inarray
441  inline void SetCoeffsArray(Array<OneD, NekDouble> &inarray);
442 
443  /// This function returns (a reference to) the array
444  /// \f$\boldsymbol{\hat{u}}_l\f$ (implemented as #m_coeffs)
445  /// containing all local expansion coefficients.
446  inline const Array<OneD, const NekDouble> &GetCoeffs() const;
447 
448  /// Impose Dirichlet Boundary Conditions onto Array
449  inline void ImposeDirichletConditions(
450  Array<OneD,NekDouble>& outarray);
451 
452 
453  /// Fill Bnd Condition expansion from the values stored in expansion
454  inline void FillBndCondFromField(void);
455 
456  /// Fill Bnd Condition expansion in nreg from the values stored in expansion
457  inline void FillBndCondFromField(const int nreg);
458 
459  /// Gathers the global coefficients \f$\boldsymbol{\hat{u}}_g\f$
460  /// from the local coefficients \f$\boldsymbol{\hat{u}}_l\f$.
461  // inline
462  MULTI_REGIONS_EXPORT inline void LocalToGlobal(bool useComm = true);
463 
465  const Array<OneD, const NekDouble> &inarray,
466  Array<OneD,NekDouble> &outarray,
467  bool useComm = true);
468 
469  /// Scatters from the global coefficients
470  /// \f$\boldsymbol{\hat{u}}_g\f$ to the local coefficients
471  /// \f$\boldsymbol{\hat{u}}_l\f$.
472  // inline
473  MULTI_REGIONS_EXPORT inline void GlobalToLocal(void);
474 
475  /**
476  * This operation is evaluated as:
477  * \f{tabbing}
478  * \hspace{1cm} \= Do \= $e=$ $1, N_{\mathrm{el}}$ \\
479  * > > Do \= $i=$ $0,N_m^e-1$ \\
480  * > > > $\boldsymbol{\hat{u}}^{e}[i] = \mbox{sign}[e][i] \cdot
481  * \boldsymbol{\hat{u}}_g[\mbox{map}[e][i]]$ \\
482  * > > continue \\
483  * > continue
484  * \f}
485  * where \a map\f$[e][i]\f$ is the mapping array and \a
486  * sign\f$[e][i]\f$ is an array of similar dimensions ensuring the
487  * correct modal connectivity between the different elements (both
488  * these arrays are contained in the data member #m_locToGloMap). This
489  * operation is equivalent to the scatter operation
490  * \f$\boldsymbol{\hat{u}}_l=\mathcal{A}\boldsymbol{\hat{u}}_g\f$,
491  * where \f$\mathcal{A}\f$ is the
492  * \f$N_{\mathrm{eof}}\times N_{\mathrm{dof}}\f$ permutation matrix.
493  *
494  * @param inarray An array of size \f$N_\mathrm{dof}\f$
495  * containing the global degrees of freedom
496  * \f$\boldsymbol{x}_g\f$.
497  * @param outarray The resulting local degrees of freedom
498  * \f$\boldsymbol{x}_l\f$ will be stored in this
499  * array of size \f$N_\mathrm{eof}\f$.
500  */
502  const Array<OneD, const NekDouble> &inarray,
503  Array<OneD,NekDouble> &outarray);
504 
505  /// Get the \a i th value (coefficient) of #m_coeffs
506  inline NekDouble GetCoeff(int i);
507 
508  /// Get the \a i th value (coefficient) of #m_coeffs
509  inline NekDouble GetCoeffs(int i);
510 
511  /// This function returns (a reference to) the array
512  /// \f$\boldsymbol{u}_l\f$ (implemented as #m_phys) containing the
513  /// function \f$u^{\delta}(\boldsymbol{x})\f$ evaluated at the
514  /// quadrature points.
515  // inline
517  &GetPhys() const;
518 
519  /// This function calculates the \f$L_\infty\f$ error of the global
520  /// spectral/hp element approximation.
522  const Array<OneD, const NekDouble> &inarray,
524 
525  /// This function calculates the \f$L_2\f$ error with
526  /// respect to soln of the global
527  /// spectral/hp element approximation.
529  const Array<OneD, const NekDouble> &inarray,
531  {
532  return v_L2(inarray, soln);
533  }
534 
535  /// Calculates the \f$H^1\f$ error of the global spectral/hp
536  /// element approximation.
538  const Array<OneD, const NekDouble> &inarray,
540 
542  {
543  return v_Integral(inarray);
544  }
545 
546  /// This function calculates the energy associated with
547  /// each one of the modesof a 3D homogeneous nD expansion
549  {
550  return v_HomogeneousEnergy();
551  }
552 
553  /// This function sets the Spectral Vanishing Viscosity
554  /// in homogeneous1D expansion.
556  {
558  }
559 
560  /// This function returns a vector containing the wave
561  /// numbers in z-direction associated
562  /// with the 3D homogenous expansion. Required if a
563  /// parellelisation is applied in the Fourier direction
565  {
566  return v_GetZIDs();
567  }
568 
569  /// This function returns the transposition class
570  /// associaed with the homogeneous expansion.
572  {
573  return v_GetTransposition();
574  }
575 
576  /// This function returns the Width of homogeneous direction
577  /// associaed with the homogeneous expansion.
579  {
580  return v_GetHomoLen();
581  }
582 
583  /// This function returns a vector containing the wave
584  /// numbers in y-direction associated
585  /// with the 3D homogenous expansion. Required if a
586  /// parellelisation is applied in the Fourier direction
588  {
589  return v_GetYIDs();
590  }
591 
592  /// This function interpolates the physical space points in
593  /// \a inarray to \a outarray using the same points defined in the
594  /// expansion but where the number of points are rescaled
595  /// by \a 1DScale
597  const NekDouble scale,
598  const Array<OneD, NekDouble> &inarray,
599  Array<OneD, NekDouble> &outarray)
600  {
601  v_PhysInterp1DScaled(scale, inarray,outarray);
602  }
603 
604  /// This function Galerkin projects the physical space points in
605  /// \a inarray to \a outarray where inarray is assumed to
606  /// be defined in the expansion but where the number of
607  /// points are rescaled by \a 1DScale
609  const NekDouble scale,
610  const Array<OneD, NekDouble> &inarray,
611  Array<OneD, NekDouble> &outarray)
612  {
613  v_PhysGalerkinProjection1DScaled(scale, inarray, outarray);
614  }
615 
616  /// This function returns the number of elements in the expansion.
617  inline int GetExpSize(void);
618 
619 
620  /// This function returns the number of elements in the
621  /// expansion which may be different for a homogeoenous extended
622  /// expansionp.
623  inline int GetNumElmts(void)
624  {
625  return v_GetNumElmts();
626  }
627 
628  /// This function returns the vector of elements in the expansion.
629  inline const boost::shared_ptr<LocalRegions::ExpansionVector>
630  GetExp() const;
631 
632  /// This function returns (a shared pointer to) the local elemental
633  /// expansion of the \f$n^{\mathrm{th}}\f$ element.
634  inline LocalRegions::ExpansionSharedPtr& GetExp(int n) const;
635 
636  /// This function returns (a shared pointer to) the local elemental
637  /// expansion containing the arbitrary point given by \a gloCoord.
639  const Array<OneD, const NekDouble> &gloCoord);
640 
641  /** This function returns the index of the local elemental
642  * expansion containing the arbitrary point given by \a gloCoord.
643  **/
645  const Array<OneD, const NekDouble> &gloCoord,
646  NekDouble tol = 0.0,
647  bool returnNearestElmt = false);
648 
649  /** This function returns the index and the Local
650  * Cartesian Coordinates \a locCoords of the local
651  * elemental expansion containing the arbitrary point
652  * given by \a gloCoords.
653  **/
655  const Array<OneD, const NekDouble> &gloCoords,
656  Array<OneD, NekDouble> &locCoords,
657  NekDouble tol = 0.0,
658  bool returnNearestElmt = false);
659 
660  /// Get the start offset position for a global list of #m_coeffs
661  /// correspoinding to element n.
662  inline int GetCoeff_Offset(int n) const;
663 
664  /// Get the start offset position for a global list of m_phys
665  /// correspoinding to element n.
666  inline int GetPhys_Offset(int n) const;
667 
668  /// Get the element id associated with the n th
669  /// consecutive block of data in #m_phys and #m_coeffs
670  inline int GetOffset_Elmt_Id(int n) const;
671 
672  /// This function returns (a reference to) the array
673  /// \f$\boldsymbol{\hat{u}}_l\f$ (implemented as #m_coeffs)
674  /// containing all local expansion coefficients.
676 
677  /// This function returns (a reference to) the array
678  /// \f$\boldsymbol{u}_l\f$ (implemented as #m_phys) containing the
679  /// function \f$u^{\delta}(\boldsymbol{x})\f$ evaluated at the
680  /// quadrature points.
682 
683  inline void PhysDeriv(
684  Direction edir,
685  const Array<OneD, const NekDouble> &inarray,
686  Array<OneD, NekDouble> &out_d);
687 
688  /// This function discretely evaluates the derivative of a function
689  /// \f$f(\boldsymbol{x})\f$ on the domain consisting of all
690  /// elements of the expansion.
691  inline void PhysDeriv(
692  const Array<OneD, const NekDouble> &inarray,
693  Array<OneD, NekDouble> &out_d0,
696 
697  inline void PhysDeriv(
698  const int dir,
699  const Array<OneD, const NekDouble> &inarray,
700  Array<OneD, NekDouble> &out_d);
701 
702  inline void CurlCurl(
705 
706  // functions associated with DisContField
709 
710  inline boost::shared_ptr<ExpList> &UpdateBndCondExpansion(int i);
711 
712  inline void Upwind(
713  const Array<OneD, const Array<OneD, NekDouble> > &Vec,
714  const Array<OneD, const NekDouble> &Fwd,
715  const Array<OneD, const NekDouble> &Bwd,
717 
718  inline void Upwind(
720  const Array<OneD, const NekDouble> &Fwd,
721  const Array<OneD, const NekDouble> &Bwd,
723 
724  /**
725  * Return a reference to the trace space associated with this
726  * expansion list.
727  */
728  inline boost::shared_ptr<ExpList> &GetTrace();
729 
730  inline boost::shared_ptr<AssemblyMapDG> &GetTraceMap(void);
731 
732  inline const Array<OneD, const int> &GetTraceBndMap(void);
733 
734  inline void GetNormals(Array<OneD, Array<OneD, NekDouble> > &normals);
735 
736  inline void AddTraceIntegral(
739  Array<OneD, NekDouble> &outarray);
740 
741  inline void AddTraceIntegral(
743  Array<OneD, NekDouble> &outarray);
744 
745  inline void AddFwdBwdTraceIntegral(
746  const Array<OneD, const NekDouble> &Fwd,
747  const Array<OneD, const NekDouble> &Bwd,
748  Array<OneD, NekDouble> &outarray);
749 
750  inline void GetFwdBwdTracePhys(
752  Array<OneD,NekDouble> &Bwd);
753 
754  inline void GetFwdBwdTracePhys(
755  const Array<OneD,const NekDouble> &field,
757  Array<OneD,NekDouble> &Bwd);
758 
759  inline const std::vector<bool> &GetLeftAdjacentFaces(void) const;
760 
761  inline void ExtractTracePhys(Array<OneD,NekDouble> &outarray);
762 
763  inline void ExtractTracePhys(
764  const Array<OneD, const NekDouble> &inarray,
765  Array<OneD,NekDouble> &outarray);
766 
767  inline const Array<OneD, const SpatialDomains::
769 
770  inline Array<OneD, SpatialDomains::
772 
773  inline void EvaluateBoundaryConditions(
774  const NekDouble time = 0.0,
775  const std::string varName = "",
778 
779  // Routines for continous matrix solution
780  /// This function calculates the result of the multiplication of a
781  /// matrix of type specified by \a mkey with a vector given by \a
782  /// inarray.
783  inline void GeneralMatrixOp(
784  const GlobalMatrixKey &gkey,
785  const Array<OneD,const NekDouble> &inarray,
786  Array<OneD, NekDouble> &outarray,
787  CoeffState coeffstate = eLocal);
788 
790  const GlobalMatrixKey &gkey,
791  const Array<OneD,const NekDouble> &inarray,
792  Array<OneD, NekDouble> &outarray);
793 
794  inline void SetUpPhysNormals();
795 
796  inline void GetBoundaryToElmtMap(Array<OneD, int> &ElmtID,
797  Array<OneD,int> &EdgeID);
798 
799  inline void GetBndElmtExpansion(int i,
800  boost::shared_ptr<ExpList> &result,
801  const bool DeclareCoeffPhysArrays = true);
802 
803  inline void ExtractElmtToBndPhys(int i,
805  Array<OneD, NekDouble> &boundary);
806 
807  inline void ExtractPhysToBndElmt(int i,
808  const Array<OneD, const NekDouble> &phys,
809  Array<OneD, NekDouble> &bndElmt);
810 
811  inline void ExtractPhysToBnd(int i,
812  const Array<OneD, const NekDouble> &phys,
814 
815  inline void GetBoundaryNormals(int i,
816  Array<OneD, Array<OneD, NekDouble> > &normals);
817 
819  std::vector<LibUtilities::FieldDefinitionsSharedPtr> &fielddef,
820  int NumHomoDir = 0,
823  std::vector<NekDouble> &HomoLen =
825  bool homoStrips = false,
826  std::vector<unsigned int> &HomoSIDs =
828  std::vector<unsigned int> &HomoZIDs =
830  std::vector<unsigned int> &HomoYIDs =
832 
833 
835  {
836  return m_globalOptParam;
837  }
838 
839  std::map<int, RobinBCInfoSharedPtr> GetRobinBCInfo()
840  {
841  return v_GetRobinBCInfo();
842  }
843 
845  PeriodicMap &periodicVerts,
846  PeriodicMap &periodicEdges,
847  PeriodicMap &periodicFaces = NullPeriodicMap)
848  {
849  v_GetPeriodicEntities(periodicVerts, periodicEdges, periodicFaces);
850  }
851 
852  std::vector<LibUtilities::FieldDefinitionsSharedPtr>
854  {
855  return v_GetFieldDefinitions();
856  }
857 
858 
859  void GetFieldDefinitions(std::vector<LibUtilities::FieldDefinitionsSharedPtr> &fielddef)
860  {
861  v_GetFieldDefinitions(fielddef);
862  }
863 
864 
865 
866  /// Append the element data listed in elements
867  /// fielddef->m_ElementIDs onto fielddata
870  std::vector<NekDouble> &fielddata)
871  {
872  v_AppendFieldData(fielddef,fielddata);
873  }
874 
875 
876  /// Append the data in coeffs listed in elements
877  /// fielddef->m_ElementIDs onto fielddata
880  std::vector<NekDouble> &fielddata,
881  Array<OneD, NekDouble> &coeffs)
882  {
883  v_AppendFieldData(fielddef,fielddata,coeffs);
884  }
885 
886 
887  /** \brief Extract the data in fielddata into the coeffs
888  * using the basic ExpList Elemental expansions rather
889  * than planes in homogeneous case
890  */
893  std::vector<NekDouble> &fielddata,
894  std::string &field,
895  Array<OneD, NekDouble> &coeffs);
896 
897 
898  /** \brief Extract the data from fromField using
899  * fromExpList the coeffs using the basic ExpList
900  * Elemental expansions rather than planes in homogeneous
901  * case
902  */
904  const boost::shared_ptr<ExpList> &fromExpList,
905  const Array<OneD, const NekDouble> &fromCoeffs,
906  Array<OneD, NekDouble> &toCoeffs);
907 
908 
909  //Extract data in fielddata into the m_coeffs_list for the 3D stability analysis (base flow is 2D)
912  std::vector<NekDouble> &fielddata,
913  std::string &field,
914  Array<OneD, NekDouble> &coeffs);
915 
916 
917  /// Returns a shared pointer to the current object.
918  boost::shared_ptr<ExpList> GetSharedThisPtr()
919  {
920  return shared_from_this();
921  }
922 
923  /// Returns the session object
924  boost::shared_ptr<LibUtilities::SessionReader> GetSession() const
925  {
926  return m_session;
927  }
928 
929  /// Returns the comm object
930  boost::shared_ptr<LibUtilities::Comm> GetComm()
931  {
932  return m_comm;
933  }
934 
936  {
937  return m_graph;
938  }
939 
940  // Wrapper functions for Homogeneous Expansions
942  {
943  return v_GetHomogeneousBasis();
944  }
945 
946  boost::shared_ptr<ExpList> &GetPlane(int n)
947  {
948  return v_GetPlane(n);
949  }
950 
951  //expansion type
953 
955  Collections::ImplementationType ImpType
956  = Collections::eNoImpType);
957 
959 
960  protected:
961  boost::shared_ptr<DNekMat> GenGlobalMatrixFull(
962  const GlobalLinSysKey &mkey,
963  const boost::shared_ptr<AssemblyMapCG> &locToGloMap);
964 
965  /// Communicator
967 
968  /// Session
970 
971  /// Mesh associated with this expansion list.
973 
974  /// The total number of local degrees of freedom. #m_ncoeffs
975  /// \f$=N_{\mathrm{eof}}=\sum_{e=1}^{{N_{\mathrm{el}}}}N^{e}_l\f$
977 
978  /** The total number of quadrature points. #m_npoints
979  *\f$=Q_{\mathrm{tot}}=\sum_{e=1}^{{N_{\mathrm{el}}}}N^{e}_Q\f$
980  **/
982 
983  /**
984  * \brief Concatenation of all local expansion coefficients.
985  *
986  * The array of length #m_ncoeffs\f$=N_{\mathrm{eof}}\f$ which is
987  * the concatenation of the local expansion coefficients
988  * \f$\hat{u}_n^e\f$ over all \f$N_{\mathrm{el}}\f$ elements
989  * \f[\mathrm{\texttt{m\_coeffs}}=\boldsymbol{\hat{u}}_{l} =
990  * \underline{\boldsymbol{\hat{u}}}^e = \left [ \begin{array}{c}
991  * \boldsymbol{\hat{u}}^{1} \ \
992  * \boldsymbol{\hat{u}}^{2} \ \
993  * \vdots \ \
994  * \boldsymbol{\hat{u}}^{{{N_{\mathrm{el}}}}} \end{array} \right ],
995  * \quad
996  * \mathrm{where}\quad \boldsymbol{\hat{u}}^{e}[n]=\hat{u}_n^{e}\f]
997  */
999 
1000  /**
1001  * \brief The global expansion evaluated at the quadrature points
1002  *
1003  * The array of length #m_npoints\f$=Q_{\mathrm{tot}}\f$ containing
1004  * the evaluation of \f$u^{\delta}(\boldsymbol{x})\f$ at the
1005  * quadrature points \f$\boldsymbol{x}_i\f$.
1006  * \f[\mathrm{\texttt{m\_phys}}=\boldsymbol{u}_{l} =
1007  * \underline{\boldsymbol{u}}^e = \left [ \begin{array}{c}
1008  * \boldsymbol{u}^{1} \ \
1009  * \boldsymbol{u}^{2} \ \
1010  * \vdots \ \
1011  * \boldsymbol{u}^{{{N_{\mathrm{el}}}}} \end{array} \right ],\quad
1012  * \mathrm{where}\quad
1013  * \boldsymbol{u}^{e}[i]=u^{\delta}(\boldsymbol{x}_i)\f]
1014  */
1016 
1017  /**
1018  * \brief The state of the array #m_phys.
1019  *
1020  * Indicates whether the array #m_phys, created to contain the
1021  * evaluation of \f$u^{\delta}(\boldsymbol{x})\f$ at the quadrature
1022  * points \f$\boldsymbol{x}_i\f$, is filled with these values.
1023  */
1025 
1026  /**
1027  * \brief The list of local expansions.
1028  *
1029  * The (shared pointer to the) vector containing (shared pointers
1030  * to) all local expansions. The fact that the local expansions are
1031  * all stored as a (pointer to a) #StdExpansion, the abstract base
1032  * class for all local expansions, allows a general implementation
1033  * where most of the routines for the derived classes are defined
1034  * in the #ExpList base class.
1035  */
1036  boost::shared_ptr<LocalRegions::ExpansionVector> m_exp;
1037 
1038  Collections::CollectionVector m_collections;
1039 
1040  /// Offset of elemental data into the array #m_coeffs
1041  std::vector<int> m_coll_coeff_offset;
1042 
1043  /// Offset of elemental data into the array #m_phys
1044  std::vector<int> m_coll_phys_offset;
1045 
1046  /// Offset of elemental data into the array #m_coeffs
1048 
1049  /// Offset of elemental data into the array #m_phys
1051 
1052  /// Array containing the element id #m_offset_elmt_id[n]
1053  /// that the n^th consecutive block of data in #m_coeffs
1054  /// and #m_phys is associated, i.e. for an array of
1055  /// constant expansion size and single shape elements
1056  /// m_phys[n*m_npoints] is the data related to
1057  /// m_exp[m_offset_elmt_id[n]];
1059 
1061 
1062  BlockMatrixMapShPtr m_blockMat;
1063 
1064  //@todo should this be in ExpList or ExpListHomogeneous1D.cpp
1065  // it's a bool which determine if the expansion is in the wave space (coefficient space)
1066  // or not
1068 
1069  /// Mapping from geometry ID of element to index inside #m_exp
1070  boost::unordered_map<int, int> m_elmtToExpId;
1071 
1072  /// This function assembles the block diagonal matrix of local
1073  /// matrices of the type \a mtype.
1075  const GlobalMatrixKey &gkey);
1076 
1078  const GlobalMatrixKey &gkey);
1079 
1080  void MultiplyByBlockMatrix(
1081  const GlobalMatrixKey &gkey,
1082  const Array<OneD,const NekDouble> &inarray,
1083  Array<OneD, NekDouble> &outarray);
1084 
1085  /// Generates a global matrix from the given key and map.
1086  boost::shared_ptr<GlobalMatrix> GenGlobalMatrix(
1087  const GlobalMatrixKey &mkey,
1088  const boost::shared_ptr<AssemblyMapCG> &locToGloMap);
1089 
1090 
1091  void GlobalEigenSystem(
1092  const boost::shared_ptr<DNekMat> &Gmat,
1093  Array<OneD, NekDouble> &EigValsReal,
1094  Array<OneD, NekDouble> &EigValsImag,
1095  Array<OneD, NekDouble> &EigVecs
1097 
1098 
1099  /// This operation constructs the global linear system of type \a
1100  /// mkey.
1101  boost::shared_ptr<GlobalLinSys> GenGlobalLinSys(
1102  const GlobalLinSysKey &mkey,
1103  const boost::shared_ptr<AssemblyMapCG> &locToGloMap);
1104 
1105  /// Generate a GlobalLinSys from information provided by the key
1106  /// "mkey" and the mapping provided in LocToGloBaseMap.
1107  boost::shared_ptr<GlobalLinSys> GenGlobalBndLinSys(
1108  const GlobalLinSysKey &mkey,
1109  const AssemblyMapSharedPtr &locToGloMap);
1110 
1112  {
1114  }
1115 
1116  // Virtual prototypes
1117 
1118  virtual int v_GetNumElmts(void)
1119  {
1120  return (*m_exp).size();
1121  }
1122 
1124  &v_GetBndCondExpansions(void);
1125 
1126  virtual boost::shared_ptr<ExpList> &v_UpdateBndCondExpansion(int i);
1127 
1128  virtual void v_Upwind(
1129  const Array<OneD, const Array<OneD, NekDouble> > &Vec,
1130  const Array<OneD, const NekDouble> &Fwd,
1131  const Array<OneD, const NekDouble> &Bwd,
1133 
1134  virtual void v_Upwind(
1135  const Array<OneD, const NekDouble> &Vn,
1136  const Array<OneD, const NekDouble> &Fwd,
1137  const Array<OneD, const NekDouble> &Bwd,
1139 
1140  virtual boost::shared_ptr<ExpList> &v_GetTrace();
1141 
1142  virtual boost::shared_ptr<AssemblyMapDG> &v_GetTraceMap();
1143 
1144  virtual const Array<OneD, const int> &v_GetTraceBndMap();
1145 
1146  virtual void v_GetNormals(
1147  Array<OneD, Array<OneD, NekDouble> > &normals);
1148 
1149  virtual void v_AddTraceIntegral(
1150  const Array<OneD, const NekDouble> &Fx,
1151  const Array<OneD, const NekDouble> &Fy,
1152  Array<OneD, NekDouble> &outarray);
1153 
1154  virtual void v_AddTraceIntegral(
1155  const Array<OneD, const NekDouble> &Fn,
1156  Array<OneD, NekDouble> &outarray);
1157 
1158  virtual void v_AddFwdBwdTraceIntegral(
1159  const Array<OneD, const NekDouble> &Fwd,
1160  const Array<OneD, const NekDouble> &Bwd,
1161  Array<OneD, NekDouble> &outarray);
1162 
1163  virtual void v_GetFwdBwdTracePhys(
1164  Array<OneD,NekDouble> &Fwd,
1165  Array<OneD,NekDouble> &Bwd);
1166 
1167  virtual void v_GetFwdBwdTracePhys(
1168  const Array<OneD,const NekDouble> &field,
1169  Array<OneD,NekDouble> &Fwd,
1170  Array<OneD,NekDouble> &Bwd);
1171 
1172  virtual const std::vector<bool> &v_GetLeftAdjacentFaces(void) const;
1173 
1174  virtual void v_ExtractTracePhys(
1175  Array<OneD,NekDouble> &outarray);
1176 
1177  virtual void v_ExtractTracePhys(
1178  const Array<OneD, const NekDouble> &inarray,
1179  Array<OneD,NekDouble> &outarray);
1180 
1181  virtual void v_MultiplyByInvMassMatrix(
1182  const Array<OneD,const NekDouble> &inarray,
1183  Array<OneD, NekDouble> &outarray,
1184  CoeffState coeffstate);
1185 
1186  virtual void v_HelmSolve(
1187  const Array<OneD, const NekDouble> &inarray,
1188  Array<OneD, NekDouble> &outarray,
1189  const FlagList &flags,
1190  const StdRegions::ConstFactorMap &factors,
1191  const StdRegions::VarCoeffMap &varcoeff,
1192  const Array<OneD, const NekDouble> &dirForcing,
1193  const bool PhysSpaceForcing);
1194 
1196  const Array<OneD, Array<OneD, NekDouble> > &velocity,
1197  const Array<OneD, const NekDouble> &inarray,
1198  Array<OneD, NekDouble> &outarray,
1199  const NekDouble lambda,
1200  CoeffState coeffstate = eLocal,
1202  dirForcing = NullNekDouble1DArray);
1203 
1204  virtual void v_LinearAdvectionReactionSolve(
1205  const Array<OneD, Array<OneD, NekDouble> > &velocity,
1206  const Array<OneD, const NekDouble> &inarray,
1207  Array<OneD, NekDouble> &outarray,
1208  const NekDouble lambda,
1209  CoeffState coeffstate = eLocal,
1211  dirForcing = NullNekDouble1DArray);
1212 
1213  // wrapper functions about virtual functions
1214  virtual void v_ImposeDirichletConditions(Array<OneD,NekDouble>& outarray);
1215 
1216  virtual void v_FillBndCondFromField();
1217 
1218  virtual void v_FillBndCondFromField(const int nreg);
1219 
1220  virtual void v_Reset();
1221 
1222  virtual void v_LocalToGlobal(bool UseComm);
1223 
1224  virtual void v_LocalToGlobal(
1225  const Array<OneD, const NekDouble> &inarray,
1226  Array<OneD,NekDouble> &outarray,
1227  bool UseComm);
1228 
1229  virtual void v_GlobalToLocal(void);
1230 
1231  virtual void v_GlobalToLocal(
1232  const Array<OneD, const NekDouble> &inarray,
1233  Array<OneD,NekDouble> &outarray);
1234 
1235  virtual void v_BwdTrans(
1236  const Array<OneD,const NekDouble> &inarray,
1237  Array<OneD, NekDouble> &outarray,
1238  CoeffState coeffstate);
1239 
1240  virtual void v_BwdTrans_IterPerExp(
1241  const Array<OneD,const NekDouble> &inarray,
1242  Array<OneD,NekDouble> &outarray);
1243 
1244  virtual void v_FwdTrans(
1245  const Array<OneD,const NekDouble> &inarray,
1246  Array<OneD, NekDouble> &outarray,
1247  CoeffState coeffstate);
1248 
1249  virtual void v_FwdTrans_IterPerExp(
1250  const Array<OneD,const NekDouble> &inarray,
1251  Array<OneD,NekDouble> &outarray);
1252 
1253  virtual void v_SmoothField(Array<OneD,NekDouble> &field);
1254 
1255  virtual void v_IProductWRTBase(
1256  const Array<OneD, const NekDouble> &inarray,
1257  Array<OneD, NekDouble> &outarray,
1258  CoeffState coeffstate);
1259 
1260  virtual void v_IProductWRTBase_IterPerExp(
1261  const Array<OneD,const NekDouble> &inarray,
1262  Array<OneD, NekDouble> &outarray);
1263 
1264  virtual void v_GeneralMatrixOp(
1265  const GlobalMatrixKey &gkey,
1266  const Array<OneD,const NekDouble> &inarray,
1267  Array<OneD, NekDouble> &outarray,
1268  CoeffState coeffstate);
1269 
1270  virtual void v_GetCoords(
1271  Array<OneD, NekDouble> &coord_0,
1272  Array<OneD, NekDouble> &coord_1,
1274 
1275  virtual void v_PhysDeriv(
1276  const Array<OneD, const NekDouble> &inarray,
1277  Array<OneD, NekDouble> &out_d0,
1278  Array<OneD, NekDouble> &out_d1,
1279  Array<OneD, NekDouble> &out_d2);
1280 
1281  virtual void v_PhysDeriv(
1282  const int dir,
1283  const Array<OneD, const NekDouble> &inarray,
1284  Array<OneD, NekDouble> &out_d);
1285 
1286  virtual void v_PhysDeriv(
1287  Direction edir,
1288  const Array<OneD, const NekDouble> &inarray,
1289  Array<OneD, NekDouble> &out_d);
1290 
1291  virtual void v_CurlCurl(
1292  Array<OneD, Array<OneD, NekDouble> > &Vel,
1293  Array<OneD, Array<OneD, NekDouble> > &Q);
1294 
1295  virtual void v_HomogeneousFwdTrans(
1296  const Array<OneD, const NekDouble> &inarray,
1297  Array<OneD, NekDouble> &outarray,
1298  CoeffState coeffstate = eLocal,
1299  bool Shuff = true,
1300  bool UnShuff = true);
1301 
1302  virtual void v_HomogeneousBwdTrans(
1303  const Array<OneD, const NekDouble> &inarray,
1304  Array<OneD, NekDouble> &outarray,
1305  CoeffState coeffstate = eLocal,
1306  bool Shuff = true,
1307  bool UnShuff = true);
1308 
1309  virtual void v_DealiasedProd(
1310  const Array<OneD, NekDouble> &inarray1,
1311  const Array<OneD, NekDouble> &inarray2,
1312  Array<OneD, NekDouble> &outarray,
1313  CoeffState coeffstate = eLocal);
1314 
1315  virtual void v_DealiasedDotProd(
1316  const Array<OneD, Array<OneD, NekDouble> > &inarray1,
1317  const Array<OneD, Array<OneD, NekDouble> > &inarray2,
1318  Array<OneD, Array<OneD, NekDouble> > &outarray,
1319  CoeffState coeffstate = eLocal);
1320 
1321  virtual void v_GetBCValues(
1322  Array<OneD, NekDouble> &BndVals,
1323  const Array<OneD, NekDouble> &TotField,
1324  int BndID);
1325 
1326  virtual void v_NormVectorIProductWRTBase(
1329  Array<OneD, NekDouble> &outarray,
1330  int BndID);
1331 
1332  virtual void v_NormVectorIProductWRTBase(
1333  Array<OneD, Array<OneD, NekDouble> > &V,
1334  Array<OneD, NekDouble> &outarray);
1335 
1336  virtual void v_SetUpPhysNormals();
1337 
1338  virtual void v_GetBoundaryToElmtMap(Array<OneD, int> &ElmtID,
1339  Array<OneD,int> &EdgeID);
1340 
1341  virtual void v_GetBndElmtExpansion(int i,
1342  boost::shared_ptr<ExpList> &result,
1343  const bool DeclareCoeffPhysArrays);
1344 
1345  virtual void v_ExtractElmtToBndPhys(int i,
1346  Array<OneD, NekDouble> &elmt,
1347  Array<OneD, NekDouble> &boundary);
1348 
1349  virtual void v_ExtractPhysToBndElmt(int i,
1350  const Array<OneD, const NekDouble> &phys,
1351  Array<OneD, NekDouble> &bndElmt);
1352 
1353  virtual void v_ExtractPhysToBnd(int i,
1354  const Array<OneD, const NekDouble> &phys,
1355  Array<OneD, NekDouble> &bnd);
1356 
1357  virtual void v_GetBoundaryNormals(int i,
1358  Array<OneD, Array<OneD, NekDouble> > &normals);
1359 
1360  virtual void v_ReadGlobalOptimizationParameters();
1361 
1362  virtual std::vector<LibUtilities::FieldDefinitionsSharedPtr>
1363  v_GetFieldDefinitions(void);
1364 
1365  virtual void v_GetFieldDefinitions(
1366  std::vector<LibUtilities::FieldDefinitionsSharedPtr> &fielddef);
1367 
1368  virtual void v_AppendFieldData(
1370  std::vector<NekDouble> &fielddata);
1371 
1372  virtual void v_AppendFieldData(
1374  std::vector<NekDouble> &fielddata,
1375  Array<OneD, NekDouble> &coeffs);
1376 
1377  virtual void v_ExtractDataToCoeffs(
1379  std::vector<NekDouble> &fielddata, std::string &field,
1380  Array<OneD, NekDouble> &coeffs);
1381 
1382  virtual void v_ExtractCoeffsToCoeffs(const boost::shared_ptr<ExpList> &fromExpList, const Array<OneD, const NekDouble> &fromCoeffs, Array<OneD, NekDouble> &toCoeffs);
1383 
1384  virtual void v_WriteTecplotHeader(std::ostream &outfile,
1385  std::string var = "");
1386  virtual void v_WriteTecplotZone(std::ostream &outfile,
1387  int expansion);
1388  virtual void v_WriteTecplotField(std::ostream &outfile,
1389  int expansion);
1390  virtual void v_WriteTecplotConnectivity(std::ostream &outfile,
1391  int expansion);
1392  virtual void v_WriteVtkPieceHeader(
1393  std::ostream &outfile,
1394  int expansion,
1395  int istrip);
1396 
1397  virtual void v_WriteVtkPieceData(
1398  std::ostream &outfile,
1399  int expansion,
1400  std::string var);
1401 
1402  virtual NekDouble v_L2(
1403  const Array<OneD, const NekDouble> &phys,
1405 
1406  virtual NekDouble v_Integral (
1407  const Array<OneD, const NekDouble> &inarray);
1408 
1411  virtual NekDouble v_GetHomoLen(void);
1414 
1415  // 1D Scaling and projection
1416  virtual void v_PhysInterp1DScaled(
1417  const NekDouble scale, const Array<OneD, NekDouble> &inarray,
1418  Array<OneD, NekDouble> &outarray);
1419 
1420  virtual void v_PhysGalerkinProjection1DScaled(
1421  const NekDouble scale,
1422  const Array<OneD, NekDouble> &inarray,
1423  Array<OneD, NekDouble> &outarray);
1424 
1425  virtual void v_ClearGlobalLinSysManager(void);
1426 
1427  void ExtractFileBCs(const std::string &fileName,
1429  const std::string &varName,
1430  const boost::shared_ptr<ExpList> locExpList);
1431 
1432  // Utility function for a common case of retrieving a
1433  // BoundaryCondition from a boundary condition collection.
1436  GetBoundaryCondition(const SpatialDomains::
1437  BoundaryConditionCollection& collection,
1438  unsigned int index, const std::string& variable);
1439 
1440 
1441  private:
1443 
1446 
1447  virtual void v_EvaluateBoundaryConditions(
1448  const NekDouble time = 0.0,
1449  const std::string varName = "",
1452 
1453  virtual std::map<int, RobinBCInfoSharedPtr> v_GetRobinBCInfo(void);
1454 
1455 
1456  virtual void v_GetPeriodicEntities(
1457  PeriodicMap &periodicVerts,
1458  PeriodicMap &periodicEdges,
1459  PeriodicMap &periodicFaces);
1460 
1461  // Homogeneous direction wrapper functions.
1463  {
1464  ASSERTL0(false,
1465  "This method is not defined or valid for this class type");
1467  }
1468 
1469  // wrapper function to set viscosity for Homo1D expansion
1471  {
1472  ASSERTL0(false,
1473  "This method is not defined or valid for this class type");
1474  }
1475 
1476 
1477  virtual boost::shared_ptr<ExpList> &v_GetPlane(int n);
1478  };
1479 
1480 
1481  /// Shared pointer to an ExpList object.
1482  typedef boost::shared_ptr<ExpList> ExpListSharedPtr;
1483  /// An empty ExpList object.
1485  static ExpListSharedPtr NullExpListSharedPtr;
1486 
1487  // Inline routines follow.
1488 
1489  /**
1490  * Returns the total number of local degrees of freedom
1491  * \f$N_{\mathrm{eof}}=\sum_{e=1}^{{N_{\mathrm{el}}}}N^{e}_m\f$.
1492  */
1493  inline int ExpList::GetNcoeffs() const
1494  {
1495  return m_ncoeffs;
1496  }
1497 
1498  inline int ExpList::GetNcoeffs(const int eid) const
1499  {
1500  return (*m_exp)[eid]->GetNcoeffs();
1501  }
1502 
1503  /**
1504  * Evaulates the maximum number of modes in the elemental basis
1505  * order over all elements
1506  */
1508  {
1509  unsigned int i;
1510  int returnval = 0;
1511 
1512  for(i= 0; i < (*m_exp).size(); ++i)
1513  {
1514  returnval = (std::max)(returnval,
1515  (*m_exp)[i]->EvalBasisNumModesMax());
1516  }
1517 
1518  return returnval;
1519  }
1520 
1521  /**
1522  *
1523  */
1525  const
1526  {
1527  unsigned int i;
1528  Array<OneD,int> returnval((*m_exp).size(),0);
1529 
1530  for(i= 0; i < (*m_exp).size(); ++i)
1531  {
1532  returnval[i]
1533  = (std::max)(returnval[i],(*m_exp)[i]->EvalBasisNumModesMax());
1534  }
1535 
1536  return returnval;
1537  }
1538 
1539 
1540  /**
1541  *
1542  */
1543  inline int ExpList::GetTotPoints() const
1544  {
1545  return m_npoints;
1546  }
1547 
1548  inline int ExpList::GetTotPoints(const int eid) const
1549  {
1550  return (*m_exp)[eid]->GetTotPoints();
1551  }
1552 
1553 
1554  inline int ExpList::Get1DScaledTotPoints(const NekDouble scale) const
1555  {
1556  int returnval = 0;
1557  int cnt;
1558  int nbase = (*m_exp)[0]->GetNumBases();
1559 
1560  for(int i = 0; i < (*m_exp).size(); ++i)
1561  {
1562  cnt = 1;
1563  for(int j = 0; j < nbase; ++j)
1564  {
1565  cnt *= (int)(scale*((*m_exp)[i]->GetNumPoints(j)));
1566  }
1567  returnval += cnt;
1568  }
1569  return returnval;
1570  }
1571 
1572  /**
1573  *
1574  */
1575  inline int ExpList::GetNpoints() const
1576  {
1577  return m_npoints;
1578  }
1579 
1580 
1581  /**
1582  *
1583  */
1584  inline void ExpList::SetWaveSpace(const bool wavespace)
1585  {
1586  m_WaveSpace = wavespace;
1587  }
1588 
1589  /**
1590  *
1591  */
1592  inline bool ExpList::GetWaveSpace() const
1593  {
1594  return m_WaveSpace;
1595  }
1596 
1597  /// Set the \a i th value of\a m_phys to value \a val
1598  inline void ExpList::SetPhys(int i, NekDouble val)
1599  {
1600  m_phys[i] = val;
1601  }
1602 
1603 
1604  /**
1605  * This function fills the array \f$\boldsymbol{u}_l\f$, the evaluation
1606  * of the expansion at the quadrature points (implemented as #m_phys),
1607  * with the values of the array \a inarray.
1608  *
1609  * @param inarray The array containing the values where
1610  * #m_phys should be filled with.
1611  */
1612  inline void ExpList::SetPhys(
1613  const Array<OneD, const NekDouble> &inarray)
1614  {
1615  ASSERTL0(inarray.num_elements() == m_npoints,
1616  "Input array does not have correct number of elements.");
1617 
1618  Vmath::Vcopy(m_npoints,&inarray[0],1,&m_phys[0],1);
1619  m_physState = true;
1620  }
1621 
1622 
1624  {
1625  m_phys = inarray;
1626  }
1627 
1628 
1629  /**
1630  * @param physState \a true (=filled) or \a false (=not filled).
1631  */
1632  inline void ExpList::SetPhysState(const bool physState)
1633  {
1634  m_physState = physState;
1635  }
1636 
1637 
1638  /**
1639  * @return physState \a true (=filled) or \a false (=not filled).
1640  */
1641  inline bool ExpList::GetPhysState() const
1642  {
1643  return m_physState;
1644  }
1645 
1646  /**
1647  *
1648  */
1650  const Array<OneD, const NekDouble> &inarray,
1651  Array<OneD, NekDouble> &outarray,
1652  CoeffState coeffstate)
1653  {
1654  v_IProductWRTBase(inarray,outarray, coeffstate);
1655  }
1656 
1657  /**
1658  *
1659  */
1661  const Array<OneD, const NekDouble> &inarray,
1662  Array<OneD, NekDouble> &outarray)
1663  {
1664  v_IProductWRTBase_IterPerExp(inarray,outarray);
1665  }
1666 
1667  /**
1668  *
1669  */
1670  inline void ExpList::FwdTrans(
1671  const Array<OneD, const NekDouble> &inarray,
1672  Array<OneD, NekDouble> &outarray,
1673  CoeffState coeffstate)
1674  {
1675  v_FwdTrans(inarray,outarray,coeffstate);
1676  }
1677 
1678  /**
1679  *
1680  */
1682  const Array<OneD, const NekDouble> &inarray,
1683  Array<OneD,NekDouble> &outarray)
1684  {
1685  v_FwdTrans_IterPerExp(inarray,outarray);
1686  }
1687 
1688  /**
1689  *
1690  */
1692  {
1693  v_SmoothField(field);
1694  }
1695 
1696  /**
1697  *
1698  */
1699  inline void ExpList::BwdTrans (
1700  const Array<OneD, const NekDouble> &inarray,
1701  Array<OneD, NekDouble> &outarray,
1702  CoeffState coeffstate)
1703  {
1704  v_BwdTrans(inarray,outarray,coeffstate);
1705  }
1706 
1707  /**
1708  *
1709  */
1711  const Array<OneD, const NekDouble> &inarray,
1712  Array<OneD, NekDouble> &outarray)
1713  {
1714  v_BwdTrans_IterPerExp(inarray,outarray);
1715  }
1716 
1717 
1718  /**
1719  *
1720  */
1722  const Array<OneD,const NekDouble> &inarray,
1723  Array<OneD, NekDouble> &outarray,
1724  CoeffState coeffstate)
1725  {
1726  v_MultiplyByInvMassMatrix(inarray,outarray,coeffstate);
1727  }
1728 
1729  /**
1730  *
1731  */
1732  inline void ExpList::HelmSolve(
1733  const Array<OneD, const NekDouble> &inarray,
1734  Array<OneD, NekDouble> &outarray,
1735  const FlagList &flags,
1736  const StdRegions::ConstFactorMap &factors,
1737  const StdRegions::VarCoeffMap &varcoeff,
1738  const Array<OneD, const NekDouble> &dirForcing,
1739  const bool PhysSpaceForcing)
1740 
1741  {
1742  v_HelmSolve(inarray, outarray, flags, factors, varcoeff,
1743  dirForcing, PhysSpaceForcing);
1744  }
1745 
1746 
1747  /**
1748  *
1749  */
1751  const Array<OneD, Array<OneD, NekDouble> > &velocity,
1752  const Array<OneD, const NekDouble> &inarray,
1753  Array<OneD, NekDouble> &outarray,
1754  const NekDouble lambda,
1755  CoeffState coeffstate,
1756  const Array<OneD, const NekDouble>& dirForcing)
1757  {
1758  v_LinearAdvectionDiffusionReactionSolve(velocity,inarray, outarray,
1759  lambda, coeffstate,dirForcing);
1760  }
1761 
1763  const Array<OneD, Array<OneD, NekDouble> > &velocity,
1764  const Array<OneD, const NekDouble> &inarray,
1765  Array<OneD, NekDouble> &outarray,
1766  const NekDouble lambda,
1767  CoeffState coeffstate,
1768  const Array<OneD, const NekDouble>& dirForcing)
1769  {
1770  v_LinearAdvectionReactionSolve(velocity,inarray, outarray,
1771  lambda, coeffstate,dirForcing);
1772  }
1773 
1774  /**
1775  *
1776  */
1778  Array<OneD, NekDouble> &coord_1,
1779  Array<OneD, NekDouble> &coord_2)
1780 
1781  {
1782  v_GetCoords(coord_0,coord_1,coord_2);
1783  }
1784 
1785  /**
1786  *
1787  */
1789  Array<OneD, NekDouble> &out_d0,
1790  Array<OneD, NekDouble> &out_d1,
1791  Array<OneD, NekDouble> &out_d2)
1792  {
1793  v_PhysDeriv(inarray,out_d0,out_d1,out_d2);
1794  }
1795 
1796  /**
1797  *
1798  */
1799  inline void ExpList::PhysDeriv(
1800  const int dir,
1801  const Array<OneD, const NekDouble> &inarray,
1802  Array<OneD, NekDouble> &out_d)
1803  {
1804  v_PhysDeriv(dir,inarray,out_d);
1805  }
1806 
1807  inline void ExpList::PhysDeriv(
1808  Direction edir,
1809  const Array<OneD, const NekDouble> &inarray,
1810  Array<OneD, NekDouble> &out_d)
1811  {
1812  v_PhysDeriv(edir, inarray,out_d);
1813  }
1814 
1815  inline void ExpList::CurlCurl(
1818  {
1819  v_CurlCurl(Vel, Q);
1820  }
1821 
1822  /**
1823  *
1824  */
1826  const Array<OneD, const NekDouble> &inarray,
1827  Array<OneD, NekDouble> &outarray,
1828  CoeffState coeffstate,
1829  bool Shuff,
1830  bool UnShuff)
1831  {
1832  v_HomogeneousFwdTrans(inarray,outarray,coeffstate,Shuff,UnShuff);
1833  }
1834 
1835  /**
1836  *
1837  */
1839  const Array<OneD, const NekDouble> &inarray,
1840  Array<OneD, NekDouble> &outarray,
1841  CoeffState coeffstate,
1842  bool Shuff,
1843  bool UnShuff)
1844  {
1845  v_HomogeneousBwdTrans(inarray,outarray,coeffstate,Shuff,UnShuff);
1846  }
1847 
1848  /**
1849  *
1850  */
1852  const Array<OneD, NekDouble> &inarray1,
1853  const Array<OneD, NekDouble> &inarray2,
1854  Array<OneD, NekDouble> &outarray,
1855  CoeffState coeffstate)
1856  {
1857  v_DealiasedProd(inarray1,inarray2,outarray,coeffstate);
1858  }
1859 
1860  /**
1861  *
1862  */
1864  const Array<OneD, Array<OneD, NekDouble> > &inarray1,
1865  const Array<OneD, Array<OneD, NekDouble> > &inarray2,
1866  Array<OneD, Array<OneD, NekDouble> > &outarray,
1867  CoeffState coeffstate)
1868  {
1869  v_DealiasedDotProd(inarray1,inarray2,outarray,coeffstate);
1870  }
1871 
1872  /**
1873  *
1874  */
1876  Array<OneD, NekDouble> &BndVals,
1877  const Array<OneD, NekDouble> &TotField,
1878  int BndID)
1879  {
1880  v_GetBCValues(BndVals,TotField,BndID);
1881  }
1882 
1883  /**
1884  *
1885  */
1889  Array<OneD, NekDouble> &outarray,
1890  int BndID)
1891  {
1892  v_NormVectorIProductWRTBase(V1,V2,outarray,BndID);
1893  }
1894 
1897  Array<OneD, NekDouble> &outarray)
1898  {
1899  v_NormVectorIProductWRTBase(V, outarray);
1900  }
1901 
1902  /**
1903  * @param eid The index of the element to be checked.
1904  * @return The dimension of the coordinates of the specific element.
1905  */
1906  inline int ExpList::GetCoordim(int eid)
1907  {
1908  ASSERTL2(eid <= (*m_exp).size(),
1909  "eid is larger than number of elements");
1910  return (*m_exp)[eid]->GetCoordim();
1911  }
1912 
1913  /**
1914  * @param i The index of m_coeffs to be set
1915  * @param val The value which m_coeffs[i] is to be set to.
1916  */
1917  inline void ExpList::SetCoeff(int i, NekDouble val)
1918  {
1919  m_coeffs[i] = val;
1920  }
1921 
1922 
1923  /**
1924  * @param i The index of #m_coeffs to be set.
1925  * @param val The value which #m_coeffs[i] is to be set to.
1926  */
1927  inline void ExpList::SetCoeffs(int i, NekDouble val)
1928  {
1929  m_coeffs[i] = val;
1930  }
1931 
1932 
1934  {
1935  m_coeffs = inarray;
1936  }
1937 
1938  /**
1939  * As the function returns a constant reference to a
1940  * <em>const Array</em>, it is not possible to modify the
1941  * underlying data of the array #m_coeffs. In order to do
1942  * so, use the function #UpdateCoeffs instead.
1943  *
1944  * @return (A constant reference to) the array #m_coeffs.
1945  */
1947  {
1948  return m_coeffs;
1949  }
1950 
1952  Array<OneD,NekDouble>& outarray)
1953  {
1954  v_ImposeDirichletConditions(outarray);
1955  }
1956 
1958  {
1960  }
1961 
1962  inline void ExpList::FillBndCondFromField(const int nreg)
1963  {
1964  v_FillBndCondFromField(nreg);
1965  }
1966 
1967  inline void ExpList::LocalToGlobal(bool useComm)
1968  {
1969  v_LocalToGlobal(useComm);
1970  }
1971 
1973  const Array<OneD, const NekDouble> &inarray,
1974  Array<OneD,NekDouble> &outarray,
1975  bool useComm)
1976  {
1977  v_LocalToGlobal(inarray, outarray,useComm);
1978  }
1979 
1980  inline void ExpList::GlobalToLocal(void)
1981  {
1982  v_GlobalToLocal();
1983  }
1984 
1985  /**
1986  * This operation is evaluated as:
1987  * \f{tabbing}
1988  * \hspace{1cm} \= Do \= $e=$ $1, N_{\mathrm{el}}$ \\
1989  * > > Do \= $i=$ $0,N_m^e-1$ \\
1990  * > > > $\boldsymbol{\hat{u}}^{e}[i] = \mbox{sign}[e][i] \cdot
1991  * \boldsymbol{\hat{u}}_g[\mbox{map}[e][i]]$ \\
1992  * > > continue \\
1993  * > continue
1994  * \f}
1995  * where \a map\f$[e][i]\f$ is the mapping array and \a
1996  * sign\f$[e][i]\f$ is an array of similar dimensions ensuring the
1997  * correct modal connectivity between the different elements (both
1998  * these arrays are contained in the data member #m_locToGloMap). This
1999  * operation is equivalent to the scatter operation
2000  * \f$\boldsymbol{\hat{u}}_l=\mathcal{A}\boldsymbol{\hat{u}}_g\f$, where
2001  * \f$\mathcal{A}\f$ is the
2002  * \f$N_{\mathrm{eof}}\times N_{\mathrm{dof}}\f$ permutation matrix.
2003  *
2004  * @param inarray An array of size \f$N_\mathrm{dof}\f$
2005  * containing the global degrees of freedom
2006  * \f$\boldsymbol{x}_g\f$.
2007  * @param outarray The resulting local degrees of freedom
2008  * \f$\boldsymbol{x}_l\f$ will be stored in this
2009  * array of size \f$N_\mathrm{eof}\f$.
2010  */
2012  const Array<OneD, const NekDouble> &inarray,
2013  Array<OneD,NekDouble> &outarray)
2014  {
2015  v_GlobalToLocal(inarray, outarray);
2016  }
2017 
2018 
2019  /**
2020  * @param i The index of #m_coeffs to be returned
2021  * @return The NekDouble held in #m_coeffs[i].
2022  */
2024  {
2025  return m_coeffs[i];
2026  }
2027 
2028  /**
2029  * @param i The index of #m_coeffs to be returned
2030  * @return The NekDouble held in #m_coeffs[i].
2031  */
2033  {
2034  return m_coeffs[i];
2035  }
2036 
2037  /**
2038  * As the function returns a constant reference to a
2039  * <em>const Array</em> it is not possible to modify the
2040  * underlying data of the array #m_phys. In order to do
2041  * so, use the function #UpdatePhys instead.
2042  *
2043  * @return (A constant reference to) the array #m_phys.
2044  */
2046  {
2047  return m_phys;
2048  }
2049 
2050  /**
2051  * @return \f$N_{\mathrm{el}}\f$, the number of elements in the
2052  * expansion.
2053  */
2054  inline int ExpList::GetExpSize(void)
2055  {
2056  return (*m_exp).size();
2057  }
2058 
2059 
2060  /**
2061  * @param n The index of the element concerned.
2062  *
2063  * @return (A shared pointer to) the local expansion of the
2064  * \f$n^{\mathrm{th}}\f$ element.
2065  */
2067  {
2068  return (*m_exp)[n];
2069  }
2070 
2071  /**
2072  * @return (A const shared pointer to) the local expansion vector #m_exp
2073  */
2074  inline const boost::shared_ptr<LocalRegions::ExpansionVector>
2075  ExpList::GetExp(void) const
2076  {
2077  return m_exp;
2078  }
2079 
2080 
2081  /**
2082  *
2083  */
2084  inline int ExpList::GetCoeff_Offset(int n) const
2085  {
2086  return m_coeff_offset[n];
2087  }
2088 
2089  /**
2090  *
2091  */
2092  inline int ExpList::GetPhys_Offset(int n) const
2093  {
2094  return m_phys_offset[n];
2095  }
2096 
2097  /**
2098  *
2099  */
2100  inline int ExpList::GetOffset_Elmt_Id(int n) const
2101  {
2102  return m_offset_elmt_id[n];
2103  }
2104 
2105  /**
2106  * If one wants to get hold of the underlying data without modifying
2107  * them, rather use the function #GetCoeffs instead.
2108  *
2109  * @return (A reference to) the array #m_coeffs.
2110  */
2112  {
2113  return m_coeffs;
2114  }
2115 
2116  /**
2117  * If one wants to get hold of the underlying data without modifying
2118  * them, rather use the function #GetPhys instead.
2119  *
2120  * @return (A reference to) the array #m_phys.
2121  */
2123  {
2124  m_physState = true;
2125  return m_phys;
2126  }
2127 
2128 
2129  // functions associated with DisContField
2132  {
2133  return v_GetBndCondExpansions();
2134  }
2135 
2136  inline boost::shared_ptr<ExpList> &ExpList::UpdateBndCondExpansion(int i)
2137  {
2138  return v_UpdateBndCondExpansion(i);
2139  }
2140 
2141  inline void ExpList::Upwind(
2142  const Array<OneD, const Array<OneD, NekDouble> > &Vec,
2143  const Array<OneD, const NekDouble> &Fwd,
2144  const Array<OneD, const NekDouble> &Bwd,
2145  Array<OneD, NekDouble> &Upwind)
2146  {
2147  v_Upwind(Vec, Fwd, Bwd, Upwind);
2148  }
2149 
2150  inline void ExpList::Upwind(
2151  const Array<OneD, const NekDouble> &Vn,
2152  const Array<OneD, const NekDouble> &Fwd,
2153  const Array<OneD, const NekDouble> &Bwd,
2154  Array<OneD, NekDouble> &Upwind)
2155  {
2156  v_Upwind(Vn, Fwd, Bwd, Upwind);
2157  }
2158 
2159  inline boost::shared_ptr<ExpList> &ExpList::GetTrace()
2160  {
2161  return v_GetTrace();
2162  }
2163 
2164  inline boost::shared_ptr<AssemblyMapDG> &ExpList::GetTraceMap()
2165  {
2166  return v_GetTraceMap();
2167  }
2168 
2170  {
2171  return v_GetTraceBndMap();
2172  }
2173 
2174  inline void ExpList::GetNormals(
2175  Array<OneD, Array<OneD, NekDouble> > &normals)
2176  {
2177  v_GetNormals(normals);
2178  }
2179 
2181  const Array<OneD, const NekDouble> &Fx,
2182  const Array<OneD, const NekDouble> &Fy,
2183  Array<OneD, NekDouble> &outarray)
2184  {
2185  v_AddTraceIntegral(Fx,Fy,outarray);
2186  }
2187 
2189  const Array<OneD, const NekDouble> &Fn,
2190  Array<OneD, NekDouble> &outarray)
2191  {
2192  v_AddTraceIntegral(Fn,outarray);
2193  }
2194 
2196  const Array<OneD, const NekDouble> &Fwd,
2197  const Array<OneD, const NekDouble> &Bwd,
2198  Array<OneD, NekDouble> &outarray)
2199  {
2200  v_AddFwdBwdTraceIntegral(Fwd,Bwd,outarray);
2201  }
2202 
2204  Array<OneD,NekDouble> &Fwd,
2205  Array<OneD,NekDouble> &Bwd)
2206  {
2207  v_GetFwdBwdTracePhys(Fwd,Bwd);
2208  }
2209 
2211  const Array<OneD,const NekDouble> &field,
2212  Array<OneD,NekDouble> &Fwd,
2213  Array<OneD,NekDouble> &Bwd)
2214  {
2215  v_GetFwdBwdTracePhys(field,Fwd,Bwd);
2216  }
2217 
2218  inline const std::vector<bool> &ExpList::GetLeftAdjacentFaces(void) const
2219  {
2220  return v_GetLeftAdjacentFaces();
2221  }
2222 
2224  {
2225  v_ExtractTracePhys(outarray);
2226  }
2227 
2228 
2230  const Array<OneD, const NekDouble> &inarray,
2231  Array<OneD,NekDouble> &outarray)
2232  {
2233  v_ExtractTracePhys(inarray,outarray);
2234  }
2235 
2238  {
2239  return v_GetBndConditions();
2240  }
2241 
2242 
2245  {
2246  return v_UpdateBndConditions();
2247  }
2248 
2250  const NekDouble time,
2251  const std::string varName,
2252  const NekDouble x2_in,
2253  const NekDouble x3_in)
2254  {
2255  v_EvaluateBoundaryConditions(time, varName, x2_in, x3_in);
2256  }
2257 
2258  // Routines for continous matrix solution
2259  /**
2260  * This operation is equivalent to the evaluation of
2261  * \f$\underline{\boldsymbol{M}}^e\boldsymbol{\hat{u}}_l\f$, that is,
2262  * \f[ \left[
2263  * \begin{array}{cccc}
2264  * \boldsymbol{M}^1 & 0 & \hspace{3mm}0 \hspace{3mm}& 0 \\
2265  * 0 & \boldsymbol{M}^2 & 0 & 0 \\
2266  * 0 & 0 & \ddots & 0 \\
2267  * 0 & 0 & 0 & \boldsymbol{M}^{N_{\mathrm{el}}} \end{array} \right]
2268  *\left [ \begin{array}{c}
2269  * \boldsymbol{\hat{u}}^{1} \\
2270  * \boldsymbol{\hat{u}}^{2} \\
2271  * \vdots \\
2272  * \boldsymbol{\hat{u}}^{{{N_{\mathrm{el}}}}} \end{array} \right ]\f]
2273  * where \f$\boldsymbol{M}^e\f$ are the local matrices of type
2274  * specified by the key \a mkey. The decoupling of the local matrices
2275  * allows for a local evaluation of the operation. However, rather than
2276  * a local matrix-vector multiplication, the local operations are
2277  * evaluated as implemented in the function
2278  * StdRegions#StdExpansion#GeneralMatrixOp.
2279  *
2280  * @param mkey This key uniquely defines the type matrix
2281  * required for the operation.
2282  * @param inarray The vector \f$\boldsymbol{\hat{u}}_l\f$ of
2283  * size \f$N_{\mathrm{eof}}\f$.
2284  * @param outarray The resulting vector of size
2285  * \f$N_{\mathrm{eof}}\f$.
2286  */
2288  const GlobalMatrixKey &gkey,
2289  const Array<OneD,const NekDouble> &inarray,
2290  Array<OneD, NekDouble> &outarray,
2291  CoeffState coeffstate)
2292  {
2293  v_GeneralMatrixOp(gkey,inarray,outarray,coeffstate);
2294  }
2295 
2296 
2298  {
2300  }
2301 
2303  Array<OneD,int> &EdgeID)
2304  {
2305  v_GetBoundaryToElmtMap(ElmtID,EdgeID);
2306  }
2307 
2308  inline void ExpList::GetBndElmtExpansion(int i,
2309  boost::shared_ptr<ExpList> &result,
2310  const bool DeclareCoeffPhysArrays)
2311  {
2312  v_GetBndElmtExpansion(i, result, DeclareCoeffPhysArrays);
2313  }
2314 
2316  Array<OneD, NekDouble> &elmt,
2317  Array<OneD, NekDouble> &boundary)
2318  {
2319  v_ExtractElmtToBndPhys(i, elmt, boundary);
2320  }
2321 
2323  const Array<OneD, const NekDouble> &phys,
2324  Array<OneD, NekDouble> &bndElmt)
2325  {
2326  v_ExtractPhysToBndElmt(i, phys, bndElmt);
2327  }
2328 
2329  inline void ExpList::ExtractPhysToBnd(int i,
2330  const Array<OneD, const NekDouble> &phys,
2332  {
2333  v_ExtractPhysToBnd(i, phys, bnd);
2334  }
2335 
2336  inline void ExpList::GetBoundaryNormals(int i,
2337  Array<OneD, Array<OneD, NekDouble> > &normals)
2338  {
2339  v_GetBoundaryNormals(i, normals);
2340  }
2341 
2343 
2344  } //end of namespace
2345 } //end of namespace
2346 
2347 #endif // EXPLIST_H
2348 
virtual NekDouble v_L2(const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &soln=NullNekDouble1DArray)
Definition: ExpList.cpp:1885
virtual void v_GlobalToLocal(void)
Definition: ExpList.cpp:2669
virtual void v_AddFwdBwdTraceIntegral(const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:2436
virtual void v_SmoothField(Array< OneD, NekDouble > &field)
Definition: ExpList.cpp:712
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:490
void AddTraceIntegral(const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:2180
virtual void v_GeneralMatrixOp(const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
Definition: ExpList.cpp:2713
const DNekScalBlkMatSharedPtr & GetBlockMatrix(const GlobalMatrixKey &gkey)
Definition: ExpList.cpp:890
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:2238
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:1946
void SmoothField(Array< OneD, NekDouble > &field)
Smooth a field across elements.
Definition: ExpList.h:1691
void GetBndElmtExpansion(int i, boost::shared_ptr< ExpList > &result, const bool DeclareCoeffPhysArrays=true)
Definition: ExpList.h:2308
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
int GetExpIndex(const Array< OneD, const NekDouble > &gloCoord, NekDouble tol=0.0, bool returnNearestElmt=false)
Definition: ExpList.cpp:1340
virtual void v_Reset()
Reset geometry information, metrics, matrix managers and geometry information.
Definition: ExpList.cpp:1513
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:973
virtual void v_SetHomo1DSpecVanVisc(Array< OneD, NekDouble > visc)
Definition: ExpList.h:1470
virtual void v_GetBoundaryNormals(int i, Array< OneD, Array< OneD, NekDouble > > &normals)
Definition: ExpList.cpp:2922
void LocalToGlobal(bool useComm=true)
Gathers the global coefficients from the local coefficients .
Definition: ExpList.h:1967
void GeneralMatrixOp_IterPerExp(const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:905
virtual const Array< OneD, const boost::shared_ptr< ExpList > > & v_GetBndCondExpansions(void)
Definition: ExpList.cpp:2355
void ExtractFileBCs(const std::string &fileName, LibUtilities::CommSharedPtr comm, const std::string &varName, const boost::shared_ptr< ExpList > locExpList)
Definition: ExpList.cpp:1987
boost::shared_ptr< Transposition > TranspositionSharedPtr
void ExtractPhysToBnd(int i, const Array< OneD, const NekDouble > &phys, Array< OneD, NekDouble > &bnd)
Definition: ExpList.h:2329
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:2084
ExpList()
The default constructor.
Definition: ExpList.cpp:95
Local coefficients.
static ExpListSharedPtr NullExpListSharedPtr
Definition: ExpList.h:1485
static Array< OneD, NekDouble > NullNekDouble1DArray
boost::unordered_map< int, int > m_elmtToExpId
Mapping from geometry ID of element to index inside m_exp.
Definition: ExpList.h:1070
void EvaluateBoundaryConditions(const NekDouble time=0.0, const std::string varName="", const NekDouble=NekConstants::kNekUnsetDouble, const NekDouble=NekConstants::kNekUnsetDouble)
Definition: ExpList.h:2249
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:1554
NekOptimize::GlobalOptParamSharedPtr m_globalOptParam
Definition: ExpList.h:1060
boost::shared_ptr< AssemblyMap > AssemblyMapSharedPtr
Definition: AssemblyMap.h:53
boost::shared_ptr< ExpList > & GetPlane(int n)
Definition: ExpList.h:946
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:571
void MultiplyByInvMassMatrix(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
Definition: ExpList.h:1721
virtual boost::shared_ptr< AssemblyMapDG > & v_GetTraceMap()
Definition: ExpList.cpp:2399
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:2092
virtual void v_ExtractPhysToBnd(int i, const Array< OneD, const NekDouble > &phys, Array< OneD, NekDouble > &bnd)
Definition: ExpList.cpp:2884
Array< OneD, NekDouble > & UpdateCoeffs()
This function returns (a reference to) the array (implemented as m_coeffs) containing all local expa...
Definition: ExpList.h:2111
virtual void v_FillBndCondFromField()
Definition: ExpList.cpp:2639
boost::shared_ptr< DNekMat > GenGlobalMatrixFull(const GlobalLinSysKey &mkey, const boost::shared_ptr< AssemblyMapCG > &locToGloMap)
Definition: ExpList.cpp:1110
NekDouble PhysIntegral(void)
This function integrates a function over the domain consisting of all the elements of the expansion...
Definition: ExpList.cpp:278
virtual void v_WriteTecplotZone(std::ostream &outfile, int expansion)
Definition: ExpList.cpp:1579
ExpansionType GetExpType(void)
Returns the type of the expansion.
Definition: ExpList.cpp:245
const Array< OneD, const int > & GetTraceBndMap(void)
Definition: ExpList.h:2169
boost::shared_ptr< FieldDefinitions > FieldDefinitionsSharedPtr
Definition: FieldIO.h:181
virtual void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
Definition: ExpList.cpp:2698
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:1524
const boost::shared_ptr< LocalRegions::ExpansionVector > GetExp() const
This function returns the vector of elements in the expansion.
Definition: ExpList.h:2075
std::vector< LibUtilities::FieldDefinitionsSharedPtr > GetFieldDefinitions()
Definition: ExpList.h:853
virtual void v_BwdTrans_IterPerExp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:1307
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:2066
static BasisSharedPtr NullBasisSharedPtr
Definition: Basis.h:358
#define MULTI_REGIONS_EXPORT
boost::shared_ptr< ExpList > & UpdateBndCondExpansion(int i)
Definition: ExpList.h:2136
void BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
Definition: ExpList.h:1699
virtual void v_GetNormals(Array< OneD, Array< OneD, NekDouble > > &normals)
Definition: ExpList.cpp:2412
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:564
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:252
void ImposeDirichletConditions(Array< OneD, NekDouble > &outarray)
Impose Dirichlet Boundary Conditions onto Array.
Definition: ExpList.h:1951
void GetBoundaryToElmtMap(Array< OneD, int > &ElmtID, Array< OneD, int > &EdgeID)
Definition: ExpList.h:2302
BlockMatrixMapShPtr m_blockMat
Definition: ExpList.h:1062
Array< OneD, NekDouble > m_phys
The global expansion evaluated at the quadrature points.
Definition: ExpList.h:1015
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:1750
virtual const std::vector< bool > & v_GetLeftAdjacentFaces(void) const
Definition: ExpList.cpp:2461
virtual void v_GetBndElmtExpansion(int i, boost::shared_ptr< ExpList > &result, const bool DeclareCoeffPhysArrays)
Definition: ExpList.cpp:2794
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:1762
void ApplyGeomInfo()
Apply geometry information to each expansion.
Definition: ExpList.cpp:1500
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Definition: ExpList.h:998
virtual void v_WriteVtkPieceHeader(std::ostream &outfile, int expansion, int istrip)
Definition: ExpList.cpp:1798
NekDouble Integral(const Array< OneD, const NekDouble > &inarray)
Definition: ExpList.h:541
virtual void v_NormVectorIProductWRTBase(Array< OneD, const NekDouble > &V1, Array< OneD, const NekDouble > &V2, Array< OneD, NekDouble > &outarray, int BndID)
Definition: ExpList.cpp:2576
virtual LibUtilities::TranspositionSharedPtr v_GetTransposition(void)
Definition: ExpList.cpp:1937
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:2054
void PhysDeriv(Direction edir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d)
Definition: ExpList.h:1807
boost::shared_ptr< SessionReader > SessionReaderSharedPtr
Definition: MeshPartition.h:51
virtual boost::shared_ptr< ExpList > & v_UpdateBndCondExpansion(int i)
Definition: ExpList.cpp:2363
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:1270
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:2804
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:608
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:641
virtual void v_WriteTecplotHeader(std::ostream &outfile, std::string var="")
Definition: ExpList.cpp:1540
virtual NekDouble v_GetHomoLen(void)
Definition: ExpList.cpp:1945
void WriteVtkFooter(std::ostream &outfile)
Definition: ExpList.cpp:1792
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:1927
virtual void v_GetCoords(Array< OneD, NekDouble > &coord_0, Array< OneD, NekDouble > &coord_1, Array< OneD, NekDouble > &coord_2=NullNekDouble1DArray)
Definition: ExpList.cpp:2733
const Array< OneD, const boost::shared_ptr< ExpList > > & GetBndCondExpansions()
Definition: ExpList.h:2131
void WriteTecplotField(std::ostream &outfile, int expansion=-1)
Definition: ExpList.h:396
boost::shared_ptr< ExpList > GetSharedThisPtr()
Returns a shared pointer to the current object.
Definition: ExpList.h:918
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:2141
virtual Array< OneD, const NekDouble > v_HomogeneousEnergy(void)
Definition: ExpList.cpp:1929
virtual const Array< OneD, const int > & v_GetTraceBndMap()
Definition: ExpList.cpp:2407
virtual void v_WriteVtkPieceData(std::ostream &outfile, int expansion, std::string var)
Definition: ExpList.cpp:1809
static std::vector< NekDouble > NullNekDoubleVector
void ExtractPhysToBndElmt(int i, const Array< OneD, const NekDouble > &phys, Array< OneD, NekDouble > &bndElmt)
Definition: ExpList.h:2322
virtual void v_ExtractPhysToBndElmt(int i, const Array< OneD, const NekDouble > &phys, Array< OneD, NekDouble > &bndElmt)
Definition: ExpList.cpp:2843
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:2287
void NormVectorIProductWRTBase(Array< OneD, const NekDouble > &V1, Array< OneD, const NekDouble > &V2, Array< OneD, NekDouble > &outarray, int BndID)
Definition: ExpList.h:1886
virtual void v_GetBCValues(Array< OneD, NekDouble > &BndVals, const Array< OneD, NekDouble > &TotField, int BndID)
Definition: ExpList.cpp:2568
void GetNormals(Array< OneD, Array< OneD, NekDouble > > &normals)
Definition: ExpList.h:2174
virtual void v_LocalToGlobal(bool UseComm)
Definition: ExpList.cpp:2653
void Reset()
Reset geometry information and reset matrices.
Definition: ExpList.h:378
boost::shared_ptr< LibUtilities::Comm > GetComm()
Returns the comm object.
Definition: ExpList.h:930
Array< OneD, int > m_coeff_offset
Offset of elemental data into the array m_coeffs.
Definition: ExpList.h:1047
LibUtilities::BasisSharedPtr GetHomogeneousBasis(void)
Definition: ExpList.h:941
const NekOptimize::GlobalOptParamSharedPtr & GetGlobalOptParam(void)
Definition: ExpList.h:834
int GetNumElmts(void)
This function returns the number of elements in the expansion which may be different for a homogeoeno...
Definition: ExpList.h:623
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:55
boost::shared_ptr< ExpList > & GetTrace()
Definition: ExpList.h:2159
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:1777
virtual void v_ImposeDirichletConditions(Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:2631
virtual LibUtilities::BasisSharedPtr v_GetHomogeneousBasis(void)
Definition: ExpList.h:1462
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:2371
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:1846
boost::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:1036
int GetTotPoints(void) const
Returns the total number of quadrature points m_npoints .
Definition: ExpList.h:1543
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:1710
virtual void v_FwdTrans_IterPerExp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:680
virtual ~ExpList()
The default destructor.
Definition: ExpList.cpp:258
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:383
std::map< StdRegions::VarCoeffType, Array< OneD, NekDouble > > VarCoeffMap
Definition: StdRegions.hpp:227
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:1251
void CurlCurl(Array< OneD, Array< OneD, NekDouble > > &Vel, Array< OneD, Array< OneD, NekDouble > > &Q)
Definition: ExpList.h:1815
void WriteVtkPieceData(std::ostream &outfile, int expansion, std::string var="v")
Definition: ExpList.h:421
void AddFwdBwdTraceIntegral(const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:2195
Array< OneD, int > m_phys_offset
Offset of elemental data into the array m_phys.
Definition: ExpList.h:1050
void FillBndCondFromField(void)
Fill Bnd Condition expansion from the values stored in expansion.
Definition: ExpList.h:1957
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:1660
NekDouble GetHomoLen(void)
This function returns the Width of homogeneous direction associaed with the homogeneous expansion...
Definition: ExpList.h:578
virtual NekDouble v_Integral(const Array< OneD, const NekDouble > &inarray)
Definition: ExpList.cpp:1915
virtual std::vector< LibUtilities::FieldDefinitionsSharedPtr > v_GetFieldDefinitions(void)
Definition: ExpList.cpp:2174
bool m_physState
The state of the array m_phys.
Definition: ExpList.h:1024
virtual void v_WriteTecplotConnectivity(std::ostream &outfile, int expansion)
Definition: ExpList.cpp:1673
static PeriodicMap NullPeriodicMap
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
Definition: ExpList.h:976
void WriteTecplotConnectivity(std::ostream &outfile, int expansion=-1)
Definition: ExpList.h:402
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:1058
boost::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
Definition: NekTypeDefs.hpp:74
void ExtractTracePhys(Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:2223
boost::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
static std::vector< unsigned int > NullUnsignedIntVector
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:2529
void SetWaveSpace(const bool wavespace)
Sets the wave space to the one of the possible configuration true or false.
Definition: ExpList.h:1584
void SetHomo1DSpecVanVisc(Array< OneD, NekDouble > visc)
This function sets the Spectral Vanishing Viscosity in homogeneous1D expansion.
Definition: ExpList.h:555
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, const bool PhysSpaceForcing)
Definition: ExpList.cpp:2493
SpatialDomains::MeshGraphSharedPtr m_graph
Mesh associated with this expansion list.
Definition: ExpList.h:972
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:3006
LibUtilities::SessionReaderSharedPtr m_session
Session.
Definition: ExpList.h:969
boost::shared_ptr< LibUtilities::SessionReader > GetSession() const
Returns the session object.
Definition: ExpList.h:924
SpatialDomains::MeshGraphSharedPtr GetGraph()
Definition: ExpList.h:935
double NekDouble
void WriteTecplotHeader(std::ostream &outfile, std::string var="")
Definition: ExpList.h:383
const DNekScalBlkMatSharedPtr GenBlockMatrix(const GlobalMatrixKey &gkey)
This function assembles the block diagonal matrix of local matrices of the type mtype.
Definition: ExpList.cpp:754
boost::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:68
Defines a list of flags.
std::map< int, BoundaryConditionMapShPtr > BoundaryConditionCollection
Definition: Conditions.h:226
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:1681
void GetPeriodicEntities(PeriodicMap &periodicVerts, PeriodicMap &periodicEdges, PeriodicMap &periodicFaces=NullPeriodicMap)
Definition: ExpList.h:844
virtual Array< OneD, SpatialDomains::BoundaryConditionShPtr > & v_UpdateBndConditions()
Definition: ExpList.cpp:2996
std::map< int, std::vector< PeriodicEntity > > PeriodicMap
virtual void v_ClearGlobalLinSysManager(void)
Definition: ExpList.cpp:1981
static const NekDouble kNekUnsetDouble
Describe a linear system.
void FwdTrans_BndConstrained(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:690
virtual boost::shared_ptr< ExpList > & v_GetTrace()
Definition: ExpList.cpp:2391
void HomogeneousFwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
Definition: ExpList.h:1825
virtual void v_DealiasedProd(const Array< OneD, NekDouble > &inarray1, const Array< OneD, NekDouble > &inarray2, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
Definition: ExpList.cpp:2549
virtual void v_ExtractCoeffsToCoeffs(const boost::shared_ptr< ExpList > &fromExpList, const Array< OneD, const NekDouble > &fromCoeffs, Array< OneD, NekDouble > &toCoeffs)
Definition: ExpList.cpp:2329
virtual boost::shared_ptr< ExpList > & v_GetPlane(int n)
Definition: ExpList.cpp:3051
Describes a matrix with ordering defined by a local to global map.
const std::vector< bool > & GetLeftAdjacentFaces(void) const
Definition: ExpList.h:2218
static ExpList NullExpList
An empty ExpList object.
Definition: ExpList.h:1484
std::map< int, RobinBCInfoSharedPtr > GetRobinBCInfo()
Definition: ExpList.h:839
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:2045
int GetNpoints(void) const
Returns the total number of quadrature points m_npoints .
Definition: ExpList.h:1575
Array< OneD, SpatialDomains::BoundaryConditionShPtr > & UpdateBndConditions()
Definition: ExpList.h:2244
virtual void v_GetPeriodicEntities(PeriodicMap &periodicVerts, PeriodicMap &periodicEdges, PeriodicMap &periodicFaces)
Definition: ExpList.cpp:3028
virtual void v_PhysInterp1DScaled(const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:1970
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:2122
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:86
LibUtilities::CommSharedPtr m_comm
Communicator.
Definition: ExpList.h:966
virtual void v_AddTraceIntegral(const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:2419
void GetFieldDefinitions(std::vector< LibUtilities::FieldDefinitionsSharedPtr > &fielddef)
Definition: ExpList.h:859
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:2225
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:548
virtual void v_CurlCurl(Array< OneD, Array< OneD, NekDouble > > &Vel, Array< OneD, Array< OneD, NekDouble > > &Q)
Definition: ExpList.cpp:562
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:2517
virtual Array< OneD, const unsigned int > v_GetZIDs(void)
Definition: ExpList.cpp:1953
void ExtractElmtToBndPhys(int i, Array< OneD, NekDouble > &elmt, Array< OneD, NekDouble > &boundary)
Definition: ExpList.h:2315
std::map< GlobalMatrixKey, DNekScalBlkMatSharedPtr > BlockMatrixMap
A map between global matrix keys and their associated block matrices.
Definition: ExpList.h:95
void DealiasedProd(const Array< OneD, NekDouble > &inarray1, const Array< OneD, NekDouble > &inarray2, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
Definition: ExpList.h:1851
bool GetPhysState(void) const
This function indicates whether the array of physical values (implemented as m_phys) is filled or no...
Definition: ExpList.h:1641
static const Array< OneD, ExpListSharedPtr > NullExpListSharedPtrArray
Definition: ExpList.h:2342
void FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
Definition: ExpList.h:1670
virtual std::map< int, RobinBCInfoSharedPtr > v_GetRobinBCInfo(void)
Definition: ExpList.cpp:3018
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:250
virtual void v_DealiasedDotProd(const Array< OneD, Array< OneD, NekDouble > > &inarray1, const Array< OneD, Array< OneD, NekDouble > > &inarray2, Array< OneD, Array< OneD, NekDouble > > &outarray, CoeffState coeffstate=eLocal)
Definition: ExpList.cpp:2558
NekDouble GetCoeff(int i)
Get the i th value (coefficient) of m_coeffs.
Definition: ExpList.h:2023
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:2047
void SetCoeff(int i, NekDouble val)
Set the i th coefficiient in m_coeffs to value val.
Definition: ExpList.h:1917
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:2100
virtual void v_IProductWRTBase_IterPerExp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:356
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:587
void IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
Definition: ExpList.h:1649
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:1632
void WriteVtkPieceFooter(std::ostream &outfile, int expansion)
Definition: ExpList.cpp:1803
void AppendFieldData(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata)
Append the element data listed in elements fielddef->m_ElementIDs onto fielddata. ...
Definition: ExpList.h:868
int GetNcoeffs(void) const
Returns the total number of local degrees of freedom .
Definition: ExpList.h:1493
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:2505
std::vector< int > m_coll_phys_offset
Offset of elemental data into the array m_phys.
Definition: ExpList.h:1044
virtual void v_AppendFieldData(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata)
Definition: ExpList.cpp:2188
static SpatialDomains::BoundaryConditionShPtr GetBoundaryCondition(const SpatialDomains::BoundaryConditionCollection &collection, unsigned int index, const std::string &variable)
Definition: ExpList.cpp:3037
boost::shared_ptr< AssemblyMapDG > & GetTraceMap(void)
Definition: ExpList.h:2164
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:3129
const Array< OneD, const SpatialDomains::BoundaryConditionShPtr > & GetBndConditions()
Definition: ExpList.h:2237
void WriteVtkPieceHeader(std::ostream &outfile, int expansion, int istrip=0)
Definition: ExpList.h:411
void MultiplyByBlockMatrix(const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:326
void SetCoeffsArray(Array< OneD, NekDouble > &inarray)
Set the m_coeffs array to inarray.
Definition: ExpList.h:1933
virtual void v_MultiplyByInvMassMatrix(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
Definition: ExpList.cpp:2484
void DealiasedDotProd(const Array< OneD, Array< OneD, NekDouble > > &inarray1, const Array< OneD, Array< OneD, NekDouble > > &inarray2, Array< OneD, Array< OneD, NekDouble > > &outarray, CoeffState coeffstate=eLocal)
Definition: ExpList.h:1863
int EvalBasisNumModesMax(void) const
Evaulates the maximum number of modes in the elemental basis order over all elements.
Definition: ExpList.h:1507
void GetBoundaryNormals(int i, Array< OneD, Array< OneD, NekDouble > > &normals)
Definition: ExpList.h:2336
bool GetWaveSpace(void) const
This function returns the third direction expansion condition, which can be in wave space (coefficien...
Definition: ExpList.h:1592
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:596
virtual void v_GetBoundaryToElmtMap(Array< OneD, int > &ElmtID, Array< OneD, int > &EdgeID)
Definition: ExpList.cpp:2967
virtual Array< OneD, const unsigned int > v_GetYIDs(void)
Definition: ExpList.cpp:1961
virtual void v_SetUpPhysNormals()
Definition: ExpList.cpp:2786
void GetBCValues(Array< OneD, NekDouble > &BndVals, const Array< OneD, NekDouble > &TotField, int BndID)
Definition: ExpList.h:1875
int GetCoordim(int eid)
This function returns the dimension of the coordinates of the element eid.
Definition: ExpList.h:1906
boost::shared_ptr< BoundaryConditionBase > BoundaryConditionShPtr
Definition: Conditions.h:219
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:878
virtual const Array< OneD, const SpatialDomains::BoundaryConditionShPtr > & v_GetBndConditions()
Definition: ExpList.cpp:2985
void SetExpType(ExpansionType Type)
Returns the type of the expansion.
Definition: ExpList.cpp:253
void SetPhysArray(Array< OneD, NekDouble > &inarray)
Sets the array m_phys.
Definition: ExpList.h:1623
void WriteTecplotZone(std::ostream &outfile, int expansion=-1)
Definition: ExpList.h:389
std::vector< int > m_coll_coeff_offset
Offset of elemental data into the array m_coeffs.
Definition: ExpList.h:1041
void WriteVtkHeader(std::ostream &outfile)
Definition: ExpList.cpp:1784
virtual void v_ExtractTracePhys(Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:2470
virtual void v_PhysGalerkinProjection1DScaled(const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:1976
boost::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:442
void GlobalToLocal(void)
Scatters from the global coefficients to the local coefficients .
Definition: ExpList.h:1980
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:528
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1061
virtual void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
Definition: ExpList.cpp:2684
void GetFwdBwdTracePhys(Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
Definition: ExpList.h:2203
virtual void v_GetFwdBwdTracePhys(Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
Definition: ExpList.cpp:2445
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:2216
void HomogeneousBwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
Definition: ExpList.h:1838
void SetPhys(int i, NekDouble val)
Set the i th value of m_phys to value val.
Definition: ExpList.h:1598
virtual void v_WriteTecplotField(std::ostream &outfile, int expansion)
Definition: ExpList.cpp:1750
virtual int v_GetNumElmts(void)
Definition: ExpList.h:1118
Collections::CollectionVector m_collections
Definition: ExpList.h:1038
virtual void v_ReadGlobalOptimizationParameters()
Definition: ExpList.cpp:2976
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, const bool PhysSpaceForcing=true)
Solve helmholtz problem.
Definition: ExpList.h:1732
static VarCoeffMap NullVarCoeffMap
Definition: StdRegions.hpp:228
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:2539
virtual void v_FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
Definition: ExpList.cpp:2691