Nektar++
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  /// Scatters from the global coefficients
87  /// \f$\boldsymbol{\hat{u}}_g\f$ to the local coefficients
88  /// \f$\boldsymbol{\hat{u}}_l\f$.
89  inline void GlobalToLocal(
90  Array<OneD,NekDouble> &outarray) const;
91 
92  /// Scatters from the global coefficients
93  /// \f$\boldsymbol{\hat{u}}_g\f$ to the local coefficients
94  /// \f$\boldsymbol{\hat{u}}_l\f$.
95  inline void GlobalToLocal(
96  const Array<OneD, const NekDouble> &inarray,
97  Array<OneD, NekDouble> &outarray) const;
98 
99  inline void LocalToGlobal(
100  const Array<OneD, const NekDouble> &inarray,
101  Array<OneD, NekDouble> &outarray) const;
102 
103  /// Assembles the global coefficients \f$\boldsymbol{\hat{u}}_g\f$
104  /// from the local coefficients \f$\boldsymbol{\hat{u}}_l\f$.
105  inline void Assemble();
106 
107  /// Assembles the global coefficients \f$\boldsymbol{\hat{u}}_g\f$
108  /// from the local coefficients \f$\boldsymbol{\hat{u}}_l\f$.
109  inline void Assemble(
110  const Array<OneD, const NekDouble> &inarray,
111  Array<OneD,NekDouble> &outarray) const;
112 
113  /// Returns the map from local to global level.
115  const;
116 
117 
118  /// Calculates the inner product of a function
119  /// \f$f(\boldsymbol{x})\f$ with respect to all <em>global</em>
120  /// expansion modes \f$\phi_n^e(\boldsymbol{x})\f$.
121  inline void IProductWRTBase(
122  const Array<OneD, const NekDouble> &inarray,
123  Array<OneD, NekDouble> &outarray,
124  CoeffState coeffstate = eLocal);
125 
126  /// Performs the global forward transformation of a function
127  /// \f$f(\boldsymbol{x})\f$, subject to the boundary conditions
128  /// specified.
130  Array<OneD, NekDouble> &outarray,
131  CoeffState coeffstate = eLocal);
132 
133  /// Performs the backward transformation of the spectral/hp
134  /// element expansion.
135  inline void BwdTrans(
136  const Array<OneD, const NekDouble> &inarray,
137  Array<OneD, NekDouble> &outarray,
138  CoeffState coeffstate = eLocal);
139 
140  /// Multiply a solution by the inverse mass matrix.
142  const Array<OneD, const NekDouble> &inarray,
143  Array<OneD, NekDouble> &outarray,
144  CoeffState coeffstate = eLocal);
145 
146  /// Solves the two-dimensional Laplace equation, subject to the
147  /// boundary conditions specified.
149  Array<OneD, NekDouble> &outarray,
150  const Array<OneD, const NekDouble> &dirForcing
153  variablecoeffs = NullNekDoubleArrayofArray,
154  NekDouble time = 0.0,
155  CoeffState coeffstate = eLocal);
156 
157 
158  /// Compute the eigenvalues of the linear advection operator.
160  const NekDouble ay,
165 
166  /// Returns the boundary conditions expansion.
169 
170  /// Returns the boundary conditions.
173 
174  inline int GetGlobalMatrixNnz(const GlobalMatrixKey &gkey);
175 
176  protected:
177 
178  private:
179  /// (A shared pointer to) the object which contains all the
180  /// required information for the transformation from local to
181  /// global degrees of freedom.
183 
184  /// (A shared pointer to) a list which collects all the global
185  /// matrices being assembled, such that they should be constructed
186  /// only once.
188 
189  /// A manager which collects all the global
190  /// linear systems being assembled, such that they should be
191  /// constructed only once.
193 
194  /// Solves the linear system specified by the key \a key.
196  const Array<OneD, const NekDouble> &rhs,
197  Array<OneD, NekDouble> &inout,
198  const Array<OneD, const NekDouble> &dirForcing
200 
201  /// Returns the global matrix specified by \a mkey.
203 
204  /// Returns the linear system specified by the key \a mkey.
206 
208 
209  /// Impose the Dirichlet Boundary Conditions on outarray
211 
213 
214  /// Gathers the global coefficients \f$\boldsymbol{\hat{u}}_g\f$
215  /// from the local coefficients \f$\boldsymbol{\hat{u}}_l\f$.
216  MULTI_REGIONS_EXPORT virtual void v_LocalToGlobal(void);
217 
218 
219  /// Scatters from the global coefficients
220  /// \f$\boldsymbol{\hat{u}}_g\f$ to the local coefficients
221  /// \f$\boldsymbol{\hat{u}}_l\f$.
222  MULTI_REGIONS_EXPORT virtual void v_GlobalToLocal(void);
223 
224  /// Template method virtual forwarder for FwdTrans().
225  MULTI_REGIONS_EXPORT virtual void v_BwdTrans(
226  const Array<OneD, const NekDouble> &inarray,
227  Array<OneD, NekDouble> &outarray,
228  CoeffState coeffstate);
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  CoeffState coeffstate);
236 
237  /// Template method virtual forwarded for SmoothField().
238  MULTI_REGIONS_EXPORT virtual void v_SmoothField(
239  Array<OneD,NekDouble> &field);
240 
241  /// Template method virtual forwarder for MultiplyByInvMassMatrix().
243  const Array<OneD, const NekDouble> &inarray,
244  Array<OneD, NekDouble> &outarray,
245  CoeffState coeffstate);
246 
247  /// Solves the two-dimensional Helmholtz equation, subject to the
248  /// boundary conditions specified.
249  MULTI_REGIONS_EXPORT virtual void v_HelmSolve(
250  const Array<OneD, const NekDouble> &inarray,
251  Array<OneD, NekDouble> &outarray,
252  const FlagList &flags,
253  const StdRegions::ConstFactorMap &factors,
254  const StdRegions::VarCoeffMap &varcoeff,
255  const Array<OneD, const NekDouble> &dirForcing);
256 
257  /// Calculates the result of the multiplication of a global
258  /// matrix of type specified by \a mkey with a vector given by \a
259  /// inarray.
260  virtual void v_GeneralMatrixOp(
261  const GlobalMatrixKey &gkey,
262  const Array<OneD,const NekDouble> &inarray,
263  Array<OneD, NekDouble> &outarray,
264  CoeffState coeffstate);
265 
266  // Solve the linear advection problem assuming that m_coeffs
267  // vector contains an intial estimate for solution
269  const Array<OneD, const NekDouble> &inarray,
270  Array<OneD, NekDouble> &outarray,
271  const NekDouble lambda,
272  CoeffState coeffstate = eLocal,
274  dirForcing = NullNekDouble1DArray);
275 
276  // Solve the linear advection problem assuming that m_coeff
277  // vector contains an intial estimate for solution
279  const Array<OneD, const NekDouble> &inarray,
280  Array<OneD, NekDouble> &outarray,
281  const NekDouble lambda,
282  CoeffState coeffstate = eLocal,
284  dirForcing = NullNekDouble1DArray);
285 
286  /// Template method virtual forwarder for GetBndConditions().
287  MULTI_REGIONS_EXPORT virtual const Array<OneD,const SpatialDomains
289  };
290 
291  typedef boost::shared_ptr<ContField2D> ContField2DSharedPtr;
292 
293 
294 
295  /**
296  * This operation is evaluated as:
297  * \f{tabbing}
298  * \hspace{1cm} \= Do \= $e=$ $1, N_{\mathrm{el}}$ \\
299  * > > Do \= $i=$ $0,N_m^e-1$ \\
300  * > > > $\boldsymbol{\hat{u}}^{e}[i] = \mbox{sign}[e][i] \cdot
301  * \boldsymbol{\hat{u}}_g[\mbox{map}[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 scatter operation
310  * \f$\boldsymbol{\hat{u}}_l=\mathcal{A}\boldsymbol{\hat{u}}_g\f$,
311  * where \f$\mathcal{A}\f$ is the
312  * \f$N_{\mathrm{eof}}\times N_{\mathrm{dof}}\f$ permutation matrix.
313  *
314  * @param outarray The resulting local degrees of freedom
315  * \f$\boldsymbol{x}_l\f$ will be stored in this
316  * array of size \f$N_\mathrm{eof}\f$.
317  */
319  Array<OneD,NekDouble> &outarray) const
320  {
321  m_locToGloMap->GlobalToLocal(m_coeffs,outarray);
322  }
323 
324 
325  /**
326  * This operation is evaluated as:
327  * \f{tabbing}
328  * \hspace{1cm} \= Do \= $e=$ $1, N_{\mathrm{el}}$ \\
329  * > > Do \= $i=$ $0,N_m^e-1$ \\
330  * > > > $\boldsymbol{\hat{u}}^{e}[i] = \mbox{sign}[e][i] \cdot
331  * \boldsymbol{\hat{u}}_g[\mbox{map}[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 scatter operation
340  * \f$\boldsymbol{\hat{u}}_l=\mathcal{A}\boldsymbol{\hat{u}}_g\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{dof}\f$
345  * containing the global degrees of freedom
346  * \f$\boldsymbol{x}_g\f$.
347  * @param outarray The resulting local degrees of freedom
348  * \f$\boldsymbol{x}_l\f$ will be stored in this
349  * array of size \f$N_\mathrm{eof}\f$.
350  */
352  const Array<OneD, const NekDouble> &inarray,
353  Array<OneD, NekDouble> &outarray) const
354  {
355  m_locToGloMap->GlobalToLocal(inarray,outarray);
356  }
357 
359  const Array<OneD, const NekDouble> &inarray,
360  Array<OneD,NekDouble> &outarray) const
361  {
362  m_locToGloMap->LocalToGlobal(inarray, outarray);
363  }
364 
365  /**
366  * This operation is evaluated as:
367  * \f{tabbing}
368  * \hspace{1cm} \= Do \= $e=$ $1, N_{\mathrm{el}}$ \\
369  * > > Do \= $i=$ $0,N_m^e-1$ \\
370  * > > > $\boldsymbol{\hat{u}}_g[\mbox{map}[e][i]] =
371  * \boldsymbol{\hat{u}}_g[\mbox{map}[e][i]]+\mbox{sign}[e][i] \cdot
372  * \boldsymbol{\hat{u}}^{e}[i]$\\
373  * > > continue\\
374  * > continue
375  * \f}
376  * where \a map\f$[e][i]\f$ is the mapping array and \a
377  * sign\f$[e][i]\f$ is an array of similar dimensions ensuring the
378  * correct modal connectivity between the different elements (both
379  * these arrays are contained in the data member #m_locToGloMap). This
380  * operation is equivalent to the gather operation
381  * \f$\boldsymbol{\hat{u}}_g=\mathcal{A}^{T}\boldsymbol{\hat{u}}_l\f$,
382  * where \f$\mathcal{A}\f$ is the
383  * \f$N_{\mathrm{eof}}\times N_{\mathrm{dof}}\f$ permutation matrix.
384  *
385  * @note The array #m_coeffs should be filled with the local
386  * coefficients \f$\boldsymbol{\hat{u}}_l\f$ and that the
387  * resulting global coefficients \f$\boldsymbol{\hat{u}}_g\f$
388  * will be stored in #m_coeffs.
389  */
390  inline void ContField2D::Assemble()
391  {
392  m_locToGloMap->Assemble(m_coeffs,m_coeffs);
393  }
394 
395  /**
396  * This operation is evaluated as:
397  * \f{tabbing}
398  * \hspace{1cm} \= Do \= $e=$ $1, N_{\mathrm{el}}$ \\
399  * > > Do \= $i=$ $0,N_m^e-1$ \\
400  * > > > $\boldsymbol{\hat{u}}_g[\mbox{map}[e][i]] =
401  * \boldsymbol{\hat{u}}_g[\mbox{map}[e][i]]+\mbox{sign}[e][i] \cdot
402  * \boldsymbol{\hat{u}}^{e}[i]$\\
403  * > > continue\\
404  * > continue
405  * \f}
406  * where \a map\f$[e][i]\f$ is the mapping array and \a
407  * sign\f$[e][i]\f$ is an array of similar dimensions ensuring the
408  * correct modal connectivity between the different elements (both
409  * these arrays are contained in the data member #m_locToGloMap). This
410  * operation is equivalent to the gather operation
411  * \f$\boldsymbol{\hat{u}}_g=\mathcal{A}^{T}\boldsymbol{\hat{u}}_l\f$,
412  * where \f$\mathcal{A}\f$ is the
413  * \f$N_{\mathrm{eof}}\times N_{\mathrm{dof}}\f$ permutation matrix.
414  *
415  * @param inarray An array of size \f$N_\mathrm{eof}\f$
416  * containing the local degrees of freedom
417  * \f$\boldsymbol{x}_l\f$.
418  * @param outarray The resulting global degrees of freedom
419  * \f$\boldsymbol{x}_g\f$ will be stored in this
420  * array of size \f$N_\mathrm{dof}\f$.
421  */
423  const Array<OneD, const NekDouble> &inarray,
424  Array<OneD,NekDouble> &outarray) const
425  {
426  m_locToGloMap->Assemble(inarray,outarray);
427  }
428 
429 
430  inline const AssemblyMapCGSharedPtr&
432  {
433  return m_locToGloMap;
434  }
435 
436 
437  /**
438  * The operation is evaluated locally (i.e. with respect to all local
439  * expansion modes) by the function ExpList#IProductWRTBase. The inner
440  * product with respect to the global expansion modes is than obtained
441  * by a global assembly operation.
442  *
443  * The values of the function \f$f(\boldsymbol{x})\f$ evaluated at the
444  * quadrature points \f$\boldsymbol{x}_i\f$ should be contained in the
445  * variable #m_phys of the ExpList object \a in. The result is stored
446  * in the array #m_coeffs.
447  *
448  * @param In An ExpList, containing the discrete evaluation
449  * of \f$f(\boldsymbol{x})\f$ at the quadrature
450  * points in its array #m_phys.
451  */
453  const Array<OneD, const NekDouble> &inarray,
454  Array<OneD, NekDouble> &outarray,
455  CoeffState coeffstate)
456 
457  {
458  if(coeffstate == eGlobal)
459  {
460  bool doGlobalOp = m_globalOptParam->DoGlobalMatOp(
462 
463  if(doGlobalOp)
464  {
466  m_locToGloMap);
468  mat->Multiply(inarray,outarray);
469  m_locToGloMap->UniversalAssemble(outarray);
470  }
471  else
472  {
474  IProductWRTBase_IterPerExp(inarray,wsp);
475  Assemble(wsp,outarray);
476  }
477  }
478  else
479  {
480  IProductWRTBase_IterPerExp(inarray,outarray);
481  }
482  }
483 
484  /**
485  * Given the coefficients of an expansion, this function evaluates the
486  * spectral/hp expansion \f$u^{\delta}(\boldsymbol{x})\f$ at the
487  * quadrature points \f$\boldsymbol{x}_i\f$. This operation is
488  * evaluated locally by the function ExpList#BwdTrans.
489  *
490  * The coefficients of the expansion should be contained in the variable
491  * #m_coeffs of the ExpList object \a In. The resulting physical values
492  * at the quadrature points \f$u^{\delta}(\boldsymbol{x}_i)\f$ are
493  * stored in the array #m_phys.
494  *
495  * @param In An ExpList, containing the local coefficients
496  * \f$\hat{u}_n^e\f$ in its array #m_coeffs.
497  */
499  const Array<OneD, const NekDouble> &inarray,
500  Array<OneD, NekDouble> &outarray,
501  CoeffState coeffstate)
502  {
503  if(coeffstate == eGlobal)
504  {
505  bool doGlobalOp = m_globalOptParam->DoGlobalMatOp(
507 
508  if(doGlobalOp)
509  {
512  mat->Multiply(inarray,outarray);
513  }
514  else
515  {
517  GlobalToLocal(inarray,wsp);
518  BwdTrans_IterPerExp(wsp,outarray);
519  }
520  }
521  else
522  {
523  BwdTrans_IterPerExp(inarray,outarray);
524  }
525  }
526 
529  {
530  return m_bndCondExpansions;
531  }
532 
535  {
536  return m_bndConditions;
537  }
538 
540  {
542  "To use method must have a AssemblyMap "
543  "attached to key");
544 
545  GlobalMatrixMap::iterator matrixIter = m_globalMat->find(gkey);
546 
547  if(matrixIter == m_globalMat->end())
548  {
549  return 0;
550  }
551  else
552  {
553  return matrixIter->second->GetNumNonZeroEntries();
554  }
555 
556  return 0;
557  }
558 
559  } //end of namespace
560 } //end of namespace
561 
562 #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:182
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:971
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:534
virtual void v_LocalToGlobal(void)
Gathers the global coefficients from the local coefficients .
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:498
#define MULTI_REGIONS_EXPORT
boost::shared_ptr< ContField2D > ContField2DSharedPtr
Definition: ContField2D.h:291
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:248
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Definition: ExpList.h:909
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:192
boost::shared_ptr< SessionReader > SessionReaderSharedPtr
Definition: MeshPartition.h:50
Global 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:1568
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:225
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:1518
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:887
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:187
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:86
virtual void v_GlobalToLocal(void)
Scatters from the global coefficients to the local coefficients .
void Assemble()
Assembles the global coefficients from the local coefficients .
Definition: ContField2D.h:390
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:452
Describes a matrix with ordering defined by a local to global map.
void LocalToGlobal(void)
Put the coefficients into global ordering using m_coeffs.
Definition: ExpList.h:1791
const AssemblyMapCGSharedPtr & GetLocalToGlobalMap() const
Returns the map from local to global level.
Definition: ContField2D.h:431
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs typedef NekMatrix< LhsDataType, StandardMatrixTag >::iterator iterator
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.
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)
Solves the two-dimensional Helmholtz equation, subject to the boundary conditions specified...
boost::shared_ptr< GlobalLinSys > GlobalLinSysSharedPtr
Pointer to a GlobalLinSys object.
Definition: GlobalLinSys.h:51
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.
static Array< OneD, Array< OneD, NekDouble > > NullNekDoubleArrayofArray
boost::shared_ptr< BoundaryConditionBase > BoundaryConditionShPtr
Definition: Conditions.h:206
int GetGlobalMatrixNnz(const GlobalMatrixKey &gkey)
Definition: ContField2D.h:539
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:165
boost::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:432
void GlobalToLocal(void)
Put the coefficients into local ordering and place in m_coeffs.
Definition: ExpList.h:1796
const Array< OneD, const MultiRegions::ExpListSharedPtr > & GetBndCondExpansions()
Returns the boundary conditions expansion.
Definition: ContField2D.h:528
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)