Nektar++
ContField.h
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File ContField.h
4 //
5 // For more information, please see: http://www.nektar.info
6 //
7 // The MIT License
8 //
9 // Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10 // Department of Aeronautics, Imperial College London (UK), and Scientific
11 // Computing and Imaging Institute, University of Utah (USA).
12 //
13 // Permission is hereby granted, free of charge, to any person obtaining a
14 // copy of this software and associated documentation files (the "Software"),
15 // to deal in the Software without restriction, including without limitation
16 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
17 // and/or sell copies of the Software, and to permit persons to whom the
18 // Software is furnished to do so, subject to the following conditions:
19 //
20 // The above copyright notice and this permission notice shall be included
21 // in all copies or substantial portions of the Software.
22 //
23 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29 // DEALINGS IN THE SOFTWARE.
30 //
31 // Description: Field definition in tow-dimensions
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
35 #ifndef NEKTAR_LIBS_MULTIREGIONS_CONTFIELD2D_H
36 #define NEKTAR_LIBS_MULTIREGIONS_CONTFIELD2D_H
37 
46 
47 
48 namespace Nektar
49 {
50  namespace MultiRegions
51  {
52  /// This class is the abstraction of a global continuous two-
53  /// dimensional spectral/hp element expansion which approximates the
54  /// solution of a set of partial differential equations.
55  class ContField: public DisContField
56  {
57  public:
58  /// The default constructor.
60 
61  /// This constructor sets up global continuous field based on an
62  /// input mesh and boundary conditions.
66  const std::string &variable = "DefaultVar",
67  const bool DeclareCoeffPhysArrays = true,
68  const bool CheckIfSingularSystem = false,
71 
72  /// Construct a global continuous field with solution type based on
73  /// another field but using a separate input mesh and boundary
74  /// conditions.
77  const std::string &variable,
78  const bool DeclareCoeffPhysArrays = true,
79  const bool CheckIfSingularSystem = false);
80 
81  /// The copy constructor.
83  (const ContField &In, bool DeclareCoeffPhysArrays = true);
84 
85  /// Copy constructor.
87  (const LibUtilities::SessionReaderSharedPtr &pSession,
88  const ExpList & In);
89 
90  /// The default destructor.
92 
93  /// Assembles the global coefficients \f$\boldsymbol{\hat{u}}_g\f$
94  /// from the local coefficients \f$\boldsymbol{\hat{u}}_l\f$.
95  inline void Assemble();
96 
97  /// Assembles the global coefficients \f$\boldsymbol{\hat{u}}_g\f$
98  /// from the local coefficients \f$\boldsymbol{\hat{u}}_l\f$.
99  inline void Assemble(
100  const Array<OneD, const NekDouble> &inarray,
101  Array<OneD,NekDouble> &outarray) const;
102 
103  /// Returns the map from local to global level.
105  const;
106 
107  /// Calculates the inner product of a function
108  /// \f$f(\boldsymbol{x})\f$ with respect to all <em>global</em>
109  /// expansion modes \f$\phi_n^e(\boldsymbol{x})\f$.
110  inline void IProductWRTBase(
111  const Array<OneD, const NekDouble> &inarray,
112  Array<OneD, NekDouble> &outarray);
113 
114  /// Performs the global forward transformation of a function
115  /// \f$f(\boldsymbol{x})\f$, subject to the boundary conditions
116  /// specified.
118  const NekDouble> &inarray,
119  Array<OneD,NekDouble> &outarray);
120 
121  /// Performs the backward transformation of the spectral/hp
122  /// element expansion.
123  inline void BwdTrans(
124  const Array<OneD, const NekDouble> &inarray,
125  Array<OneD, NekDouble> &outarray);
126 
127  /// Multiply a solution by the inverse mass matrix.
129  const Array<OneD, const NekDouble> &inarray,
130  Array<OneD, NekDouble> &outarray);
131 
132  /// Solves the two-dimensional Laplace equation, subject to the
133  /// boundary conditions specified.
135  const Array<OneD, const NekDouble> &inarray,
136  Array<OneD, NekDouble> &outarray,
137  const Array<OneD, const NekDouble> &dirForcing
140  variablecoeffs = NullNekDoubleArrayOfArray,
141  NekDouble time = 0.0);
142 
143  /// Compute the eigenvalues of the linear advection operator.
145  const NekDouble ay,
150 
151  /// Returns the boundary conditions expansion.
154 
155  /// Returns the boundary conditions.
158 
159  inline int GetGlobalMatrixNnz(const GlobalMatrixKey &gkey);
160 
161  protected:
162 
163  //private:
164  /// (A shared pointer to) the object which contains all the
165  /// required information for the transformation from local to
166  /// global degrees of freedom.
168 
169  /// (A shared pointer to) a list which collects all the global
170  /// matrices being assembled, such that they should be constructed
171  /// only once.
173 
174  /// A manager which collects all the global
175  /// linear systems being assembled, such that they should be
176  /// constructed only once.
178 
179  /// Solves the linear system specified by the key \a key.
181  const Array<OneD, const NekDouble> &rhs,
182  Array<OneD, NekDouble> &inout,
183  const Array<OneD, const NekDouble> &dirForcing
185 
186  /// Returns the global matrix specified by \a mkey.
188  GetGlobalMatrix(const GlobalMatrixKey &mkey);
189 
190  /// Returns the linear system specified by the key \a mkey.
192  GetGlobalLinSys(const GlobalLinSysKey &mkey);
193 
195  GenGlobalLinSys(const GlobalLinSysKey &mkey);
196 
197  /// Impose the Dirichlet Boundary Conditions on outarray
198  MULTI_REGIONS_EXPORT virtual void
200 
201  MULTI_REGIONS_EXPORT virtual void
203 
204  MULTI_REGIONS_EXPORT virtual void
205  v_FillBndCondFromField(const int nreg);
206 
207  /// Gathers the global coefficients \f$\boldsymbol{\hat{u}}_g\f$
208  /// from the local coefficients \f$\boldsymbol{\hat{u}}_l\f$.
210  const Array<OneD, const NekDouble> &inarray,
211  Array<OneD,NekDouble> &outarray,
212  bool useComm);
213 
214  MULTI_REGIONS_EXPORT virtual void v_LocalToGlobal(bool useComm);
215 
216  /// Scatters from the global coefficients
217  /// \f$\boldsymbol{\hat{u}}_g\f$ to the local coefficients
218  /// \f$\boldsymbol{\hat{u}}_l\f$.
220  const Array<OneD, const NekDouble> &inarray,
221  Array<OneD,NekDouble> &outarray);
222 
223  MULTI_REGIONS_EXPORT virtual void v_GlobalToLocal(void);
224 
225  // /// Template method virtual forwarder for FwdTrans().
226  // MULTI_REGIONS_EXPORT virtual void v_BwdTrans(
227  // const Array<OneD, const NekDouble> &inarray,
228  // Array<OneD, NekDouble> &outarray);
229 
230 
231  /// Template method virtual forwarder for FwdTrans().
232  MULTI_REGIONS_EXPORT virtual void v_FwdTrans(
233  const Array<OneD, const NekDouble> &inarray,
234  Array<OneD, NekDouble> &outarray);
235 
236  /// Template method virtual forwarded for SmoothField().
237  MULTI_REGIONS_EXPORT virtual void v_SmoothField(
238  Array<OneD,NekDouble> &field);
239 
240  /// Template method virtual forwarder for MultiplyByInvMassMatrix().
242  const Array<OneD, const NekDouble> &inarray,
243  Array<OneD, NekDouble> &outarray);
244 
245  /// Solves the two-dimensional Helmholtz equation, subject to the
246  /// boundary conditions specified.
247  MULTI_REGIONS_EXPORT virtual void v_HelmSolve(
248  const Array<OneD, const NekDouble> &inarray,
249  Array<OneD, NekDouble> &outarray,
250  const StdRegions::ConstFactorMap &factors,
251  const StdRegions::VarCoeffMap &varcoeff,
252  const MultiRegions::VarFactorsMap &varfactors,
253  const Array<OneD, const NekDouble> &dirForcing,
254  const bool PhysSpaceForcing);
255 
256  /// Calculates the result of the multiplication of a global
257  /// matrix of type specified by \a mkey with a vector given by \a
258  /// inarray.
259  virtual void v_GeneralMatrixOp(
260  const GlobalMatrixKey &gkey,
261  const Array<OneD,const NekDouble> &inarray,
262  Array<OneD, NekDouble> &outarray);
263 
264  // Solve the linear advection problem assuming that m_coeffs
265  // vector contains an intial estimate for solution
266  MULTI_REGIONS_EXPORT virtual void
268  (const Array<OneD, Array<OneD, NekDouble> > &velocity,
269  const Array<OneD, const NekDouble> &inarray,
270  Array<OneD, NekDouble> &outarray,
271  const NekDouble lambda,
273  dirForcing = NullNekDouble1DArray);
274 
275  // Solve the linear advection problem assuming that m_coeff
276  // vector contains an intial estimate for solution
279  (const Array<OneD, Array<OneD, NekDouble> > &velocity,
280  const Array<OneD, const NekDouble> &inarray,
281  Array<OneD, NekDouble> &outarray,
282  const NekDouble lambda,
284  dirForcing = NullNekDouble1DArray);
285 
286  /// Template method virtual forwarder for GetBndConditions().
287  MULTI_REGIONS_EXPORT virtual const Array<OneD,const SpatialDomains
290  };
291 
292  typedef std::shared_ptr<ContField> ContFieldSharedPtr;
293 
294  /**
295  * This operation is evaluated as:
296  * \f{tabbing}
297  * \hspace{1cm} \= Do \= $e=$ $1, N_{\mathrm{el}}$ \\
298  * > > Do \= $i=$ $0,N_m^e-1$ \\
299  * > > > $\boldsymbol{\hat{u}}_g[\mbox{map}[e][i]] =
300  * \boldsymbol{\hat{u}}_g[\mbox{map}[e][i]]+\mbox{sign}[e][i] \cdot
301  * \boldsymbol{\hat{u}}^{e}[i]$\\
302  * > > continue\\
303  * > continue
304  * \f}
305  * where \a map\f$[e][i]\f$ is the mapping array and \a
306  * sign\f$[e][i]\f$ is an array of similar dimensions ensuring the
307  * correct modal connectivity between the different elements (both
308  * these arrays are contained in the data member #m_locToGloMap). This
309  * operation is equivalent to the gather operation
310  * \f$\boldsymbol{\hat{u}}_g=\mathcal{A}^{T}\boldsymbol{\hat{u}}_l\f$,
311  * where \f$\mathcal{A}\f$ is the
312  * \f$N_{\mathrm{eof}}\times N_{\mathrm{dof}}\f$ permutation matrix.
313  *
314  * @note The array #m_coeffs should be filled with the local
315  * coefficients \f$\boldsymbol{\hat{u}}_l\f$ and that the
316  * resulting global coefficients \f$\boldsymbol{\hat{u}}_g\f$
317  * will be stored in #m_coeffs.
318  */
319  inline void ContField::Assemble()
320  {
321  m_locToGloMap->Assemble(m_coeffs,m_coeffs);
322  }
323 
324  /**
325  * This operation is evaluated as:
326  * \f{tabbing}
327  * \hspace{1cm} \= Do \= $e=$ $1, N_{\mathrm{el}}$ \\
328  * > > Do \= $i=$ $0,N_m^e-1$ \\
329  * > > > $\boldsymbol{\hat{u}}_g[\mbox{map}[e][i]] =
330  * \boldsymbol{\hat{u}}_g[\mbox{map}[e][i]]+\mbox{sign}[e][i] \cdot
331  * \boldsymbol{\hat{u}}^{e}[i]$\\
332  * > > continue\\
333  * > continue
334  * \f}
335  * where \a map\f$[e][i]\f$ is the mapping array and \a
336  * sign\f$[e][i]\f$ is an array of similar dimensions ensuring the
337  * correct modal connectivity between the different elements (both
338  * these arrays are contained in the data member #m_locToGloMap). This
339  * operation is equivalent to the gather operation
340  * \f$\boldsymbol{\hat{u}}_g=\mathcal{A}^{T}\boldsymbol{\hat{u}}_l\f$,
341  * where \f$\mathcal{A}\f$ is the
342  * \f$N_{\mathrm{eof}}\times N_{\mathrm{dof}}\f$ permutation matrix.
343  *
344  * @param inarray An array of size \f$N_\mathrm{eof}\f$
345  * containing the local degrees of freedom
346  * \f$\boldsymbol{x}_l\f$.
347  * @param outarray The resulting global degrees of freedom
348  * \f$\boldsymbol{x}_g\f$ will be stored in this
349  * array of size \f$N_\mathrm{dof}\f$.
350  */
351  inline void ContField::Assemble(
352  const Array<OneD, const NekDouble> &inarray,
353  Array<OneD,NekDouble> &outarray) const
354  {
355  m_locToGloMap->Assemble(inarray,outarray);
356  }
357 
358 
359  inline const AssemblyMapCGSharedPtr&
361  {
362  return m_locToGloMap;
363  }
364 
365 
366  /**
367  * The operation is evaluated locally (i.e. with respect to all local
368  * expansion modes) by the function ExpList#IProductWRTBase. The inner
369  * product with respect to the global expansion modes is than obtained
370  * by a global assembly operation.
371  *
372  * The values of the function \f$f(\boldsymbol{x})\f$ evaluated at the
373  * quadrature points \f$\boldsymbol{x}_i\f$ should be contained in the
374  * variable #m_phys of the ExpList object \a in. The result is stored
375  * in the array #m_coeffs.
376  *
377  * @param In An ExpList, containing the discrete evaluation
378  * of \f$f(\boldsymbol{x})\f$ at the quadrature
379  * points in its array #m_phys.
380  */
382  const Array<OneD, const NekDouble> &inarray,
383  Array<OneD, NekDouble> &outarray)
384 
385  {
386  IProductWRTBase_IterPerExp(inarray,outarray);
387  }
388 
389  /**
390  * Given the coefficients of an expansion, this function evaluates the
391  * spectral/hp expansion \f$u^{\delta}(\boldsymbol{x})\f$ at the
392  * quadrature points \f$\boldsymbol{x}_i\f$. This operation is
393  * evaluated locally by the function ExpList#BwdTrans.
394  *
395  * The coefficients of the expansion should be contained in the variable
396  * #m_coeffs of the ExpList object \a In. The resulting physical values
397  * at the quadrature points \f$u^{\delta}(\boldsymbol{x}_i)\f$ are
398  * stored in the array #m_phys.
399  *
400  * @param In An ExpList, containing the local coefficients
401  * \f$\hat{u}_n^e\f$ in its array #m_coeffs.
402  */
403  inline void ContField::BwdTrans(
404  const Array<OneD, const NekDouble> &inarray,
405  Array<OneD, NekDouble> &outarray)
406  {
407  BwdTrans_IterPerExp(inarray,outarray);
408  }
409 
412  {
413  return m_bndCondExpansions;
414  }
415 
418  {
419  return m_bndConditions;
420  }
421 
423  {
425  "To use method must have a AssemblyMap "
426  "attached to key");
427 
428  auto matrixIter = m_globalMat->find(gkey);
429 
430  if(matrixIter == m_globalMat->end())
431  {
432  return 0;
433  }
434  else
435  {
436  return matrixIter->second->GetNumNonZeroEntries();
437  }
438 
439  return 0;
440  }
441 
442  } //end of namespace
443 } //end of namespace
444 
445 #endif // MULTIERGIONS_CONTFIELD2D_H
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:250
#define MULTI_REGIONS_EXPORT
This class is the abstraction of a global continuous two- dimensional spectral/hp element expansion w...
Definition: ContField.h:56
virtual void v_MultiplyByInvMassMatrix(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Template method virtual forwarder for MultiplyByInvMassMatrix().
Definition: ContField.cpp:833
virtual ~ContField()
The default destructor.
Definition: ContField.cpp:248
LibUtilities::NekManager< GlobalLinSysKey, GlobalLinSys > m_globalLinSysManager
A manager which collects all the global linear systems being assembled, such that they should be cons...
Definition: ContField.h:177
void LaplaceSolve(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray, const Array< OneD, Array< OneD, NekDouble > > &variablecoeffs=NullNekDoubleArrayOfArray, NekDouble time=0.0)
Solves the two-dimensional Laplace equation, subject to the boundary conditions specified.
Definition: ContField.cpp:364
GlobalLinSysSharedPtr GenGlobalLinSys(const GlobalLinSysKey &mkey)
Definition: ContField.cpp:585
AssemblyMapCGSharedPtr m_locToGloMap
(A shared pointer to) the object which contains all the required information for the transformation f...
Definition: ContField.h:167
virtual void v_SmoothField(Array< OneD, NekDouble > &field)
Template method virtual forwarded for SmoothField().
Definition: ContField.cpp:287
virtual void v_LocalToGlobal(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool useComm)
Gathers the global coefficients from the local coefficients .
Definition: ContField.cpp:815
const Array< OneD, const SpatialDomains::BoundaryConditionShPtr > & GetBndConditions()
Returns the boundary conditions.
Definition: ContField.h:417
virtual void v_FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Template method virtual forwarder for FwdTrans().
Definition: ContField.cpp:609
void MultiplyByInvMassMatrix(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Multiply a solution by the inverse mass matrix.
Definition: ContField.cpp:306
const AssemblyMapCGSharedPtr & GetLocalToGlobalMap() const
Returns the map from local to global level.
Definition: ContField.h:360
int GetGlobalMatrixNnz(const GlobalMatrixKey &gkey)
Definition: ContField.h:422
virtual void v_GeneralMatrixOp(const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Calculates the result of the multiplication of a global matrix of type specified by mkey with a vecto...
Definition: ContField.cpp:953
virtual void v_GlobalToLocal(void)
Definition: ContField.cpp:781
void IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Calculates the inner product of a function with respect to all global expansion modes .
Definition: ContField.h:381
GlobalMatrixMapShPtr m_globalMat
(A shared pointer to) a list which collects all the global matrices being assembled,...
Definition: ContField.h:172
void BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Performs the backward transformation of the spectral/hp element expansion.
Definition: ContField.h:403
virtual void v_LinearAdvectionDiffusionReactionSolve(const Array< OneD, Array< OneD, NekDouble > > &velocity, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const NekDouble lambda, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
Definition: ContField.cpp:973
virtual void v_HelmSolve(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::ConstFactorMap &factors, const StdRegions::VarCoeffMap &varcoeff, const MultiRegions::VarFactorsMap &varfactors, const Array< OneD, const NekDouble > &dirForcing, const bool PhysSpaceForcing)
Solves the two-dimensional Helmholtz equation, subject to the boundary conditions specified.
Definition: ContField.cpp:867
virtual void v_ImposeDirichletConditions(Array< OneD, NekDouble > &outarray)
Impose the Dirichlet Boundary Conditions on outarray.
Definition: ContField.cpp:617
virtual void v_ClearGlobalLinSysManager(void)
Definition: ContField.cpp:1085
GlobalLinSysSharedPtr GetGlobalLinSys(const GlobalLinSysKey &mkey)
Returns the linear system specified by the key mkey.
Definition: ContField.cpp:579
void FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Performs the global forward transformation of a function , subject to the boundary conditions specifi...
Definition: ContField.cpp:270
GlobalMatrixSharedPtr GetGlobalMatrix(const GlobalMatrixKey &mkey)
Returns the global matrix specified by mkey.
Definition: ContField.cpp:547
virtual void v_FillBndCondFromField()
Definition: ContField.cpp:678
ContField()
The default constructor.
Definition: ContField.cpp:88
virtual const Array< OneD, const SpatialDomains ::BoundaryConditionShPtr > & v_GetBndConditions()
Template method virtual forwarder for GetBndConditions().
Definition: ContField.cpp:1076
void LinearAdvectionEigs(const NekDouble ax, const NekDouble ay, Array< OneD, NekDouble > &Real, Array< OneD, NekDouble > &Imag, Array< OneD, NekDouble > &Evecs=NullNekDouble1DArray)
Compute the eigenvalues of the linear advection operator.
Definition: ContField.cpp:443
const Array< OneD, const MultiRegions::ExpListSharedPtr > & GetBndCondExpansions()
Returns the boundary conditions expansion.
Definition: ContField.h:411
void GlobalSolve(const GlobalLinSysKey &key, const Array< OneD, const NekDouble > &rhs, Array< OneD, NekDouble > &inout, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
Solves the linear system specified by the key key.
Definition: ContField.cpp:517
void Assemble()
Assembles the global coefficients from the local coefficients .
Definition: ContField.h:319
void v_LinearAdvectionReactionSolve(const Array< OneD, Array< OneD, NekDouble > > &velocity, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const NekDouble lambda, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
Definition: ContField.cpp:1049
This class is the abstractio n of a global discontinuous two- dimensional spectral/hp element expansi...
Definition: DisContField.h:56
Array< OneD, SpatialDomains::BoundaryConditionShPtr > m_bndConditions
An array which contains the information about the boundary condition on the different boundary region...
Definition: DisContField.h:148
Array< OneD, MultiRegions::ExpListSharedPtr > m_bndCondExpansions
An object which contains the discretised boundary conditions.
Definition: DisContField.h:142
Base class for all multi-elemental spectral/hp expansions.
Definition: ExpList.h:107
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Definition: ExpList.h:1252
void IProductWRTBase_IterPerExp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function calculates the inner product of a function with respect to all local expansion modes .
Definition: ExpList.h:1966
void BwdTrans_IterPerExp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function elementally evaluates the backward transformation of the global spectral/hp element exp...
Definition: ExpList.h:2025
Describes a matrix with ordering defined by a local to global map.
bool LocToGloMapIsDefined() const
Returns true if a local to global map is defined.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< AssemblyMapCG > AssemblyMapCGSharedPtr
Definition: AssemblyMapCG.h:51
std::shared_ptr< GlobalLinSys > GlobalLinSysSharedPtr
Pointer to a GlobalLinSys object.
Definition: GlobalLinSys.h:50
std::shared_ptr< GlobalMatrix > GlobalMatrixSharedPtr
Shared pointer to a GlobalMatrix object.
Definition: GlobalMatrix.h:88
std::shared_ptr< GlobalMatrixMap > GlobalMatrixMapShPtr
Shared pointer to a global matrix map.
Definition: GlobalMatrix.h:92
std::map< StdRegions::ConstFactorType, Array< OneD, NekDouble > > VarFactorsMap
std::shared_ptr< ContField > ContFieldSharedPtr
Definition: ContField.h:292
std::shared_ptr< BoundaryConditionBase > BoundaryConditionShPtr
Definition: Conditions.h:219
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174
std::map< StdRegions::VarCoeffType, Array< OneD, NekDouble > > VarCoeffMap
Definition: StdRegions.hpp:272
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:314
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
static Array< OneD, Array< OneD, NekDouble > > NullNekDoubleArrayOfArray
static Array< OneD, NekDouble > NullNekDouble1DArray
double NekDouble