Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ContField2D.h
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File ContField2D.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: Field definition in tow-dimensions
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
36 #ifndef NEKTAR_LIBS_MULTIREGIONS_CONTFIELD2D_H
37 #define NEKTAR_LIBS_MULTIREGIONS_CONTFIELD2D_H
38 
47 
48 
49 namespace Nektar
50 {
51  namespace MultiRegions
52  {
53  /// This class is the abstraction of a global continuous two-
54  /// dimensional spectral/hp element expansion which approximates the
55  /// solution of a set of partial differential equations.
57  {
58  public:
59  /// The default constructor.
61 
62  /// This constructor sets up global continuous field based on an
63  /// input mesh and boundary conditions.
67  const std::string &variable = "DefaultVar",
68  const bool DeclareCoeffPhysArrays = true,
69  const bool CheckIfSingularSystem = false);
70 
71  /// Construct a global continuous field with solution type based on
72  /// another field but using a separate input mesh and boundary
73  /// conditions.
74  MULTI_REGIONS_EXPORT ContField2D(const ContField2D &In,
76  const std::string &variable,
77  const bool DeclareCoeffPhysArrays = true,
78  const bool CheckIfSingularSystem = false);
79 
80  /// The copy constructor.
81  MULTI_REGIONS_EXPORT ContField2D(const ContField2D &In, bool DeclareCoeffPhysArrays = true);
82 
83  /// The default destructor.
85 
86  /// Assembles the global coefficients \f$\boldsymbol{\hat{u}}_g\f$
87  /// from the local coefficients \f$\boldsymbol{\hat{u}}_l\f$.
88  inline void Assemble();
89 
90  /// Assembles the global coefficients \f$\boldsymbol{\hat{u}}_g\f$
91  /// from the local coefficients \f$\boldsymbol{\hat{u}}_l\f$.
92  inline void Assemble(
93  const Array<OneD, const NekDouble> &inarray,
94  Array<OneD,NekDouble> &outarray) const;
95 
96  /// Returns the map from local to global level.
98  const;
99 
100 
101  /// Calculates the inner product of a function
102  /// \f$f(\boldsymbol{x})\f$ with respect to all <em>global</em>
103  /// expansion modes \f$\phi_n^e(\boldsymbol{x})\f$.
104  inline void IProductWRTBase(
105  const Array<OneD, const NekDouble> &inarray,
106  Array<OneD, NekDouble> &outarray,
107  CoeffState coeffstate = eLocal);
108 
109  /// Performs the global forward transformation of a function
110  /// \f$f(\boldsymbol{x})\f$, subject to the boundary conditions
111  /// specified.
113  Array<OneD, NekDouble> &outarray,
114  CoeffState coeffstate = eLocal);
115 
116  /// Performs the backward transformation of the spectral/hp
117  /// element expansion.
118  inline void BwdTrans(
119  const Array<OneD, const NekDouble> &inarray,
120  Array<OneD, NekDouble> &outarray,
121  CoeffState coeffstate = eLocal);
122 
123  /// Multiply a solution by the inverse mass matrix.
125  const Array<OneD, const NekDouble> &inarray,
126  Array<OneD, NekDouble> &outarray,
127  CoeffState coeffstate = eLocal);
128 
129  /// Solves the two-dimensional Laplace equation, subject to the
130  /// boundary conditions specified.
132  Array<OneD, NekDouble> &outarray,
133  const Array<OneD, const NekDouble> &dirForcing
136  variablecoeffs = NullNekDoubleArrayofArray,
137  NekDouble time = 0.0,
138  CoeffState coeffstate = eLocal);
139 
140 
141  /// Compute the eigenvalues of the linear advection operator.
143  const NekDouble ay,
148 
149  /// Returns the boundary conditions expansion.
152 
153  /// Returns the boundary conditions.
156 
157  inline int GetGlobalMatrixNnz(const GlobalMatrixKey &gkey);
158 
159  protected:
160 
161  private:
162  /// (A shared pointer to) the object which contains all the
163  /// required information for the transformation from local to
164  /// global degrees of freedom.
166 
167  /// (A shared pointer to) a list which collects all the global
168  /// matrices being assembled, such that they should be constructed
169  /// only once.
171 
172  /// A manager which collects all the global
173  /// linear systems being assembled, such that they should be
174  /// constructed only once.
176 
177  /// Solves the linear system specified by the key \a key.
179  const Array<OneD, const NekDouble> &rhs,
180  Array<OneD, NekDouble> &inout,
181  const Array<OneD, const NekDouble> &dirForcing
183 
184  /// Returns the global matrix specified by \a mkey.
186 
187  /// Returns the linear system specified by the key \a mkey.
189 
191 
192  /// Impose the Dirichlet Boundary Conditions on outarray
194 
196 
197  MULTI_REGIONS_EXPORT virtual void v_FillBndCondFromField(const int nreg);
198 
199  /// Gathers the global coefficients \f$\boldsymbol{\hat{u}}_g\f$
200  /// from the local coefficients \f$\boldsymbol{\hat{u}}_l\f$.
202  const Array<OneD, const NekDouble> &inarray,
203  Array<OneD,NekDouble> &outarray,
204  bool useComm);
205 
206  MULTI_REGIONS_EXPORT virtual void v_LocalToGlobal(bool useComm);
207 
208  /// Scatters from the global coefficients
209  /// \f$\boldsymbol{\hat{u}}_g\f$ to the local coefficients
210  /// \f$\boldsymbol{\hat{u}}_l\f$.
212  const Array<OneD, const NekDouble> &inarray,
213  Array<OneD,NekDouble> &outarray);
214 
215  MULTI_REGIONS_EXPORT virtual void v_GlobalToLocal(void);
216 
217  /// Template method virtual forwarder for FwdTrans().
218  MULTI_REGIONS_EXPORT virtual void v_BwdTrans(
219  const Array<OneD, const NekDouble> &inarray,
220  Array<OneD, NekDouble> &outarray,
221  CoeffState coeffstate);
222 
223 
224  /// Template method virtual forwarder for FwdTrans().
225  MULTI_REGIONS_EXPORT virtual void v_FwdTrans(
226  const Array<OneD, const NekDouble> &inarray,
227  Array<OneD, NekDouble> &outarray,
228  CoeffState coeffstate);
229 
230  /// Template method virtual forwarded for SmoothField().
231  MULTI_REGIONS_EXPORT virtual void v_SmoothField(
232  Array<OneD,NekDouble> &field);
233 
234  /// Template method virtual forwarder for MultiplyByInvMassMatrix().
236  const Array<OneD, const NekDouble> &inarray,
237  Array<OneD, NekDouble> &outarray,
238  CoeffState coeffstate);
239 
240  /// Solves the two-dimensional Helmholtz equation, subject to the
241  /// boundary conditions specified.
242  MULTI_REGIONS_EXPORT virtual void v_HelmSolve(
243  const Array<OneD, const NekDouble> &inarray,
244  Array<OneD, NekDouble> &outarray,
245  const FlagList &flags,
246  const StdRegions::ConstFactorMap &factors,
247  const StdRegions::VarCoeffMap &varcoeff,
248  const Array<OneD, const NekDouble> &dirForcing,
249  const bool PhysSpaceForcing);
250 
251  /// Calculates the result of the multiplication of a global
252  /// matrix of type specified by \a mkey with a vector given by \a
253  /// inarray.
254  virtual void v_GeneralMatrixOp(
255  const GlobalMatrixKey &gkey,
256  const Array<OneD,const NekDouble> &inarray,
257  Array<OneD, NekDouble> &outarray,
258  CoeffState coeffstate);
259 
260  // Solve the linear advection problem assuming that m_coeffs
261  // vector contains an intial estimate for solution
263  const Array<OneD, const NekDouble> &inarray,
264  Array<OneD, NekDouble> &outarray,
265  const NekDouble lambda,
266  CoeffState coeffstate = eLocal,
268  dirForcing = NullNekDouble1DArray);
269 
270  // Solve the linear advection problem assuming that m_coeff
271  // vector contains an intial estimate for solution
273  const Array<OneD, const NekDouble> &inarray,
274  Array<OneD, NekDouble> &outarray,
275  const NekDouble lambda,
276  CoeffState coeffstate = eLocal,
278  dirForcing = NullNekDouble1DArray);
279 
280  /// Template method virtual forwarder for GetBndConditions().
281  MULTI_REGIONS_EXPORT virtual const Array<OneD,const SpatialDomains
284 
285  };
286 
287  typedef boost::shared_ptr<ContField2D> ContField2DSharedPtr;
288 
289  /**
290  * This operation is evaluated as:
291  * \f{tabbing}
292  * \hspace{1cm} \= Do \= $e=$ $1, N_{\mathrm{el}}$ \\
293  * > > Do \= $i=$ $0,N_m^e-1$ \\
294  * > > > $\boldsymbol{\hat{u}}_g[\mbox{map}[e][i]] =
295  * \boldsymbol{\hat{u}}_g[\mbox{map}[e][i]]+\mbox{sign}[e][i] \cdot
296  * \boldsymbol{\hat{u}}^{e}[i]$\\
297  * > > continue\\
298  * > continue
299  * \f}
300  * where \a map\f$[e][i]\f$ is the mapping array and \a
301  * sign\f$[e][i]\f$ is an array of similar dimensions ensuring the
302  * correct modal connectivity between the different elements (both
303  * these arrays are contained in the data member #m_locToGloMap). This
304  * operation is equivalent to the gather operation
305  * \f$\boldsymbol{\hat{u}}_g=\mathcal{A}^{T}\boldsymbol{\hat{u}}_l\f$,
306  * where \f$\mathcal{A}\f$ is the
307  * \f$N_{\mathrm{eof}}\times N_{\mathrm{dof}}\f$ permutation matrix.
308  *
309  * @note The array #m_coeffs should be filled with the local
310  * coefficients \f$\boldsymbol{\hat{u}}_l\f$ and that the
311  * resulting global coefficients \f$\boldsymbol{\hat{u}}_g\f$
312  * will be stored in #m_coeffs.
313  */
314  inline void ContField2D::Assemble()
315  {
316  m_locToGloMap->Assemble(m_coeffs,m_coeffs);
317  }
318 
319  /**
320  * This operation is evaluated as:
321  * \f{tabbing}
322  * \hspace{1cm} \= Do \= $e=$ $1, N_{\mathrm{el}}$ \\
323  * > > Do \= $i=$ $0,N_m^e-1$ \\
324  * > > > $\boldsymbol{\hat{u}}_g[\mbox{map}[e][i]] =
325  * \boldsymbol{\hat{u}}_g[\mbox{map}[e][i]]+\mbox{sign}[e][i] \cdot
326  * \boldsymbol{\hat{u}}^{e}[i]$\\
327  * > > continue\\
328  * > continue
329  * \f}
330  * where \a map\f$[e][i]\f$ is the mapping array and \a
331  * sign\f$[e][i]\f$ is an array of similar dimensions ensuring the
332  * correct modal connectivity between the different elements (both
333  * these arrays are contained in the data member #m_locToGloMap). This
334  * operation is equivalent to the gather operation
335  * \f$\boldsymbol{\hat{u}}_g=\mathcal{A}^{T}\boldsymbol{\hat{u}}_l\f$,
336  * where \f$\mathcal{A}\f$ is the
337  * \f$N_{\mathrm{eof}}\times N_{\mathrm{dof}}\f$ permutation matrix.
338  *
339  * @param inarray An array of size \f$N_\mathrm{eof}\f$
340  * containing the local degrees of freedom
341  * \f$\boldsymbol{x}_l\f$.
342  * @param outarray The resulting global degrees of freedom
343  * \f$\boldsymbol{x}_g\f$ will be stored in this
344  * array of size \f$N_\mathrm{dof}\f$.
345  */
347  const Array<OneD, const NekDouble> &inarray,
348  Array<OneD,NekDouble> &outarray) const
349  {
350  m_locToGloMap->Assemble(inarray,outarray);
351  }
352 
353 
354  inline const AssemblyMapCGSharedPtr&
356  {
357  return m_locToGloMap;
358  }
359 
360 
361  /**
362  * The operation is evaluated locally (i.e. with respect to all local
363  * expansion modes) by the function ExpList#IProductWRTBase. The inner
364  * product with respect to the global expansion modes is than obtained
365  * by a global assembly operation.
366  *
367  * The values of the function \f$f(\boldsymbol{x})\f$ evaluated at the
368  * quadrature points \f$\boldsymbol{x}_i\f$ should be contained in the
369  * variable #m_phys of the ExpList object \a in. The result is stored
370  * in the array #m_coeffs.
371  *
372  * @param In An ExpList, containing the discrete evaluation
373  * of \f$f(\boldsymbol{x})\f$ at the quadrature
374  * points in its array #m_phys.
375  */
377  const Array<OneD, const NekDouble> &inarray,
378  Array<OneD, NekDouble> &outarray,
379  CoeffState coeffstate)
380 
381  {
382  if(coeffstate == eGlobal)
383  {
384  bool doGlobalOp = m_globalOptParam->DoGlobalMatOp(
386 
387  if(doGlobalOp)
388  {
390  m_locToGloMap);
392  mat->Multiply(inarray,outarray);
393  m_locToGloMap->UniversalAssemble(outarray);
394  }
395  else
396  {
398  IProductWRTBase_IterPerExp(inarray,wsp);
399  Assemble(wsp,outarray);
400  }
401  }
402  else
403  {
404  IProductWRTBase_IterPerExp(inarray,outarray);
405  }
406  }
407 
408  /**
409  * Given the coefficients of an expansion, this function evaluates the
410  * spectral/hp expansion \f$u^{\delta}(\boldsymbol{x})\f$ at the
411  * quadrature points \f$\boldsymbol{x}_i\f$. This operation is
412  * evaluated locally by the function ExpList#BwdTrans.
413  *
414  * The coefficients of the expansion should be contained in the variable
415  * #m_coeffs of the ExpList object \a In. The resulting physical values
416  * at the quadrature points \f$u^{\delta}(\boldsymbol{x}_i)\f$ are
417  * stored in the array #m_phys.
418  *
419  * @param In An ExpList, containing the local coefficients
420  * \f$\hat{u}_n^e\f$ in its array #m_coeffs.
421  */
423  const Array<OneD, const NekDouble> &inarray,
424  Array<OneD, NekDouble> &outarray,
425  CoeffState coeffstate)
426  {
427  if(coeffstate == eGlobal)
428  {
429  bool doGlobalOp = m_globalOptParam->DoGlobalMatOp(
431 
432  if(doGlobalOp)
433  {
436  mat->Multiply(inarray,outarray);
437  }
438  else
439  {
441  GlobalToLocal(inarray,wsp);
442  BwdTrans_IterPerExp(wsp,outarray);
443  }
444  }
445  else
446  {
447  BwdTrans_IterPerExp(inarray,outarray);
448  }
449  }
450 
453  {
454  return m_bndCondExpansions;
455  }
456 
459  {
460  return m_bndConditions;
461  }
462 
464  {
466  "To use method must have a AssemblyMap "
467  "attached to key");
468 
469  GlobalMatrixMap::iterator matrixIter = m_globalMat->find(gkey);
470 
471  if(matrixIter == m_globalMat->end())
472  {
473  return 0;
474  }
475  else
476  {
477  return matrixIter->second->GetNumNonZeroEntries();
478  }
479 
480  return 0;
481  }
482 
483  } //end of namespace
484 } //end of namespace
485 
486 #endif // MULTIERGIONS_CONTFIELD2D_H
virtual void v_FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
Template method virtual forwarder for FwdTrans().
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, CoeffState coeffstate=eLocal)
Solves the two-dimensional Laplace equation, subject to the boundary conditions specified.
AssemblyMapCGSharedPtr m_locToGloMap
(A shared pointer to) the object which contains all the required information for the transformation f...
Definition: ContField2D.h:165
Local coefficients.
static Array< OneD, NekDouble > NullNekDouble1DArray
void FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
Performs the global forward transformation of a function , subject to the boundary conditions specifi...
NekOptimize::GlobalOptParamSharedPtr m_globalOptParam
Definition: ExpList.h:1060
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)
const Array< OneD, const SpatialDomains::BoundaryConditionShPtr > & GetBndConditions()
Returns the boundary conditions.
Definition: ContField2D.h:458
virtual void v_ImposeDirichletConditions(Array< OneD, NekDouble > &outarray)
Impose the Dirichlet Boundary Conditions on outarray.
GlobalLinSysSharedPtr GenGlobalLinSys(const GlobalLinSysKey &mkey)
void BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
Performs the backward transformation of the spectral/hp element expansion.
Definition: ContField2D.h:422
#define MULTI_REGIONS_EXPORT
boost::shared_ptr< ContField2D > ContField2DSharedPtr
Definition: ContField2D.h:287
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:252
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Definition: ExpList.h:998
boost::shared_ptr< GlobalMatrix > GlobalMatrixSharedPtr
Shared pointer to a GlobalMatrix object.
Definition: GlobalMatrix.h:89
LibUtilities::NekManager< GlobalLinSysKey, GlobalLinSys > m_globalLinSysManager
A manager which collects all the global linear systems being assembled, such that they should be cons...
Definition: ContField2D.h:175
boost::shared_ptr< SessionReader > SessionReaderSharedPtr
Definition: MeshPartition.h:51
Global coefficients.
virtual void v_LocalToGlobal(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool useComm)
Gathers the global coefficients from the local coefficients .
Array< OneD, MultiRegions::ExpListSharedPtr > m_bndCondExpansions
An object which contains the discretised boundary conditions.
This class is the abstraction of a global continuous two- dimensional spectral/hp element expansion w...
Definition: ContField2D.h:56
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
boost::shared_ptr< GlobalMatrixMap > GlobalMatrixMapShPtr
Shared pointer to a global matrix map.
Definition: GlobalMatrix.h:93
std::map< StdRegions::VarCoeffType, Array< OneD, NekDouble > > VarCoeffMap
Definition: StdRegions.hpp:227
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
virtual void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
Template method virtual forwarder for FwdTrans().
virtual void v_SmoothField(Array< OneD, NekDouble > &field)
Template method virtual forwarded for SmoothField().
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
Definition: ExpList.h:976
GlobalMatrixMapShPtr m_globalMat
(A shared pointer to) a list which collects all the global matrices being assembled, such that they should be constructed only once.
Definition: ContField2D.h:170
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.
GlobalMatrixSharedPtr GetGlobalMatrix(const GlobalMatrixKey &mkey)
Returns the global matrix specified by mkey.
virtual ~ContField2D()
The default destructor.
virtual void v_GeneralMatrixOp(const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
Calculates the result of the multiplication of a global matrix of type specified by mkey with a vecto...
ContField2D()
The default constructor.
Definition: ContField2D.cpp:88
virtual void v_GlobalToLocal(void)
void Assemble()
Assembles the global coefficients from the local coefficients .
Definition: ContField2D.h:314
double NekDouble
Defines a list of flags.
Array< OneD, SpatialDomains::BoundaryConditionShPtr > m_bndConditions
An array which contains the information about the boundary condition on the different boundary region...
virtual const Array< OneD, const SpatialDomains::BoundaryConditionShPtr > & v_GetBndConditions()
Template method virtual forwarder for GetBndConditions().
Describe a linear system.
void IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
Calculates the inner product of a function with respect to all global expansion modes ...
Definition: ContField2D.h:376
Describes a matrix with ordering defined by a local to global map.
const AssemblyMapCGSharedPtr & GetLocalToGlobalMap() const
Returns the map from local to global level.
Definition: ContField2D.h:355
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs typedef NekMatrix< LhsDataType, StandardMatrixTag >::iterator iterator
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)
Solves the two-dimensional Helmholtz equation, subject to the boundary conditions specified...
bool LocToGloMapIsDefined() const
Returns true if a local to global map is defined.
virtual void v_MultiplyByInvMassMatrix(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
Template method virtual forwarder for MultiplyByInvMassMatrix().
void MultiplyByInvMassMatrix(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
Multiply a solution by the inverse mass matrix.
boost::shared_ptr< GlobalLinSys > GlobalLinSysSharedPtr
Pointer to a GlobalLinSys object.
Definition: GlobalLinSys.h:52
GlobalLinSysSharedPtr GetGlobalLinSys(const GlobalLinSysKey &mkey)
Returns the linear system specified by the key mkey.
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.
virtual void v_ClearGlobalLinSysManager(void)
static Array< OneD, Array< OneD, NekDouble > > NullNekDoubleArrayofArray
boost::shared_ptr< BoundaryConditionBase > BoundaryConditionShPtr
Definition: Conditions.h:219
int GetGlobalMatrixNnz(const GlobalMatrixKey &gkey)
Definition: ContField2D.h:463
boost::shared_ptr< AssemblyMapCG > AssemblyMapCGSharedPtr
Definition: AssemblyMapCG.h:52
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:228
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
const Array< OneD, const MultiRegions::ExpListSharedPtr > & GetBndCondExpansions()
Returns the boundary conditions expansion.
Definition: ContField2D.h:452
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)