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 
47 
48 namespace Nektar
49 {
50 namespace MultiRegions
51 {
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.
56 class ContField : public DisContField
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  const Collections::ImplementationType ImpType =
72 
73  /// Construct a global continuous field with solution type based on
74  /// another field but using a separate input mesh and boundary
75  /// conditions.
77  const ContField &In, const SpatialDomains::MeshGraphSharedPtr &graph2D,
78  const std::string &variable, const bool DeclareCoeffPhysArrays = true,
79  const bool CheckIfSingularSystem = false);
80 
81  /// The copy constructor.
83  bool DeclareCoeffPhysArrays = true);
84 
85  /// Copy constructor.
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(const Array<OneD, const NekDouble> &inarray,
100  Array<OneD, NekDouble> &outarray) const;
101 
102  /// Returns the map from local to global level.
103  inline const AssemblyMapCGSharedPtr &GetLocalToGlobalMap() const;
104 
105  /// Performs the global forward transformation of a function
106  /// \f$f(\boldsymbol{x})\f$, subject to the boundary conditions
107  /// specified.
109  const Array<OneD, const NekDouble> &inarray,
110  Array<OneD, NekDouble> &outarray);
111 
112  /// Multiply a solution by the inverse mass matrix.
114  const Array<OneD, const NekDouble> &inarray,
115  Array<OneD, NekDouble> &outarray);
116 
117  /// Solves the two-dimensional Laplace equation, subject to the
118  /// boundary conditions specified.
120  const Array<OneD, const NekDouble> &inarray,
121  Array<OneD, NekDouble> &outarray,
123  const Array<OneD, Array<OneD, NekDouble>> &variablecoeffs =
125  NekDouble time = 0.0);
126 
127  /// Compute the eigenvalues of the linear advection operator.
129  const NekDouble ax, const NekDouble ay, Array<OneD, NekDouble> &Real,
132 
133  /// Returns the boundary conditions expansion.
136 
137  /// Returns the boundary conditions.
139  &GetBndConditions();
140 
141  inline int GetGlobalMatrixNnz(const GlobalMatrixKey &gkey);
142 
143  /// Solves the linear system specified by the key \a key.
145  const GlobalLinSysKey &key, const Array<OneD, const NekDouble> &rhs,
146  Array<OneD, NekDouble> &inout,
148 
150  {
151  // initialize if required
152  if (!m_GJPData)
153  {
155  GetSharedThisPtr());
156  }
157 
158  return m_GJPData;
159  }
160 
162  const GJPStabilisationSharedPtr &GJPData)
163  {
164  m_GJPData = GJPData;
165  }
166 
167 protected:
168  // private:
169  /// (A shared pointer to) the object which contains all the
170  /// required information for the transformation from local to
171  /// global degrees of freedom.
173 
174  /// (A shared pointer to) a list which collects all the global
175  /// matrices being assembled, such that they should be constructed
176  /// only once.
178 
179  /// A manager which collects all the global
180  /// linear systems being assembled, such that they should be
181  /// constructed only once.
184 
185  /// Data for Gradient Jump Penalisation (GJP) stabilisaiton
187 
188  /// Returns the global matrix specified by \a mkey.
190  GetGlobalMatrix(const GlobalMatrixKey &mkey);
191 
192  /// Returns the linear system specified by the key \a mkey.
194  GetGlobalLinSys(const GlobalLinSysKey &mkey);
195 
197  GenGlobalLinSys(const GlobalLinSysKey &mkey);
198 
199  /// Impose the Dirichlet Boundary Conditions on outarray
201  Array<OneD, NekDouble> &outarray) override;
202 
204  const Array<OneD, NekDouble> coeffs) override;
205 
207  const int nreg, const Array<OneD, NekDouble> coeffs) override;
208 
209  /// Gathers the global coefficients \f$\boldsymbol{\hat{u}}_g\f$
210  /// from the local coefficients \f$\boldsymbol{\hat{u}}_l\f$.
212  const Array<OneD, const NekDouble> &inarray,
213  Array<OneD, NekDouble> &outarray, bool useComm) override;
214 
215  MULTI_REGIONS_EXPORT virtual void v_LocalToGlobal(bool useComm) override;
216 
217  /// Scatters from the global coefficients
218  /// \f$\boldsymbol{\hat{u}}_g\f$ to the local coefficients
219  /// \f$\boldsymbol{\hat{u}}_l\f$.
221  const Array<OneD, const NekDouble> &inarray,
222  Array<OneD, NekDouble> &outarray) override;
223 
224  MULTI_REGIONS_EXPORT virtual void v_GlobalToLocal(void) override;
225 
226  // /// Template method virtual forwarder for FwdTrans().
227  // MULTI_REGIONS_EXPORT virtual void v_BwdTrans(
228  // const Array<OneD, const NekDouble> &inarray,
229  // Array<OneD, NekDouble> &outarray);
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) override;
235 
236  /// Template method virtual forwarded for SmoothField().
237  MULTI_REGIONS_EXPORT virtual void v_SmoothField(
238  Array<OneD, NekDouble> &field) override;
239 
240  /// Template method virtual forwarder for MultiplyByInvMassMatrix().
242  const Array<OneD, const NekDouble> &inarray,
243  Array<OneD, NekDouble> &outarray) override;
244 
245  /// Solves the two-dimensional Helmholtz equation, subject to the
246  /// boundary conditions specified.
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) override;
255 
256  // Solve the linear advection problem assuming that m_coeffs
257  // vector contains an intial estimate for solution
260  const Array<OneD, Array<OneD, NekDouble>> &velocity,
261  const Array<OneD, const NekDouble> &inarray,
262  Array<OneD, NekDouble> &outarray, const NekDouble lambda,
263  const Array<OneD, const NekDouble> &dirForcing =
264  NullNekDouble1DArray) override;
265 
266  // Solve the linear advection problem assuming that m_coeff
267  // vector contains an intial estimate for solution
269  const Array<OneD, Array<OneD, NekDouble>> &velocity,
270  const Array<OneD, const NekDouble> &inarray,
271  Array<OneD, NekDouble> &outarray, const NekDouble lambda,
272  const Array<OneD, const NekDouble> &dirForcing =
273  NullNekDouble1DArray) override;
274 
275  /// Template method virtual forwarder for GetBndConditions().
276  MULTI_REGIONS_EXPORT virtual const Array<
278  &v_GetBndConditions() override;
279  MULTI_REGIONS_EXPORT virtual void v_ClearGlobalLinSysManager(void) override;
280 
281  // Get manager pool count; intended for unit tests
282  MULTI_REGIONS_EXPORT int v_GetPoolCount(std::string) override;
283 
284  // Remove GlobalLinSys, StaticCond Blocks and LocalMatrix Blocks
286  bool) override;
287 };
288 
289 typedef std::shared_ptr<ContField> ContFieldSharedPtr;
290 
291 /**
292  * This operation is evaluated as:
293  * \f{tabbing}
294  * \hspace{1cm} \= Do \= $e=$ $1, N_{\mathrm{el}}$ \\
295  * > > Do \= $i=$ $0,N_m^e-1$ \\
296  * > > > $\boldsymbol{\hat{u}}_g[\mbox{map}[e][i]] =
297  * \boldsymbol{\hat{u}}_g[\mbox{map}[e][i]]+\mbox{sign}[e][i] \cdot
298  * \boldsymbol{\hat{u}}^{e}[i]$\\
299  * > > continue\\
300  * > continue
301  * \f}
302  * where \a map\f$[e][i]\f$ is the mapping array and \a
303  * sign\f$[e][i]\f$ is an array of similar dimensions ensuring the
304  * correct modal connectivity between the different elements (both
305  * these arrays are contained in the data member #m_locToGloMap). This
306  * operation is equivalent to the gather operation
307  * \f$\boldsymbol{\hat{u}}_g=\mathcal{A}^{T}\boldsymbol{\hat{u}}_l\f$,
308  * where \f$\mathcal{A}\f$ is the
309  * \f$N_{\mathrm{eof}}\times N_{\mathrm{dof}}\f$ permutation matrix.
310  *
311  * @note The array #m_coeffs should be filled with the local
312  * coefficients \f$\boldsymbol{\hat{u}}_l\f$ and that the
313  * resulting global coefficients \f$\boldsymbol{\hat{u}}_g\f$
314  * will be stored in #m_coeffs.
315  */
316 inline void ContField::Assemble()
317 {
318  m_locToGloMap->Assemble(m_coeffs, m_coeffs);
319 }
320 
321 /**
322  * This operation is evaluated as:
323  * \f{tabbing}
324  * \hspace{1cm} \= Do \= $e=$ $1, N_{\mathrm{el}}$ \\
325  * > > Do \= $i=$ $0,N_m^e-1$ \\
326  * > > > $\boldsymbol{\hat{u}}_g[\mbox{map}[e][i]] =
327  * \boldsymbol{\hat{u}}_g[\mbox{map}[e][i]]+\mbox{sign}[e][i] \cdot
328  * \boldsymbol{\hat{u}}^{e}[i]$\\
329  * > > continue\\
330  * > continue
331  * \f}
332  * where \a map\f$[e][i]\f$ is the mapping array and \a
333  * sign\f$[e][i]\f$ is an array of similar dimensions ensuring the
334  * correct modal connectivity between the different elements (both
335  * these arrays are contained in the data member #m_locToGloMap). This
336  * operation is equivalent to the gather operation
337  * \f$\boldsymbol{\hat{u}}_g=\mathcal{A}^{T}\boldsymbol{\hat{u}}_l\f$,
338  * where \f$\mathcal{A}\f$ is the
339  * \f$N_{\mathrm{eof}}\times N_{\mathrm{dof}}\f$ permutation matrix.
340  *
341  * @param inarray An array of size \f$N_\mathrm{eof}\f$
342  * containing the local degrees of freedom
343  * \f$\boldsymbol{x}_l\f$.
344  * @param outarray The resulting global degrees of freedom
345  * \f$\boldsymbol{x}_g\f$ will be stored in this
346  * array of size \f$N_\mathrm{dof}\f$.
347  */
349  Array<OneD, NekDouble> &outarray) const
350 {
351  m_locToGloMap->Assemble(inarray, outarray);
352 }
353 
355 {
356  return m_locToGloMap;
357 }
358 
361 {
362  return m_bndCondExpansions;
363 }
364 
367 {
368  return m_bndConditions;
369 }
370 
372 {
374  "To use method must have a AssemblyMap "
375  "attached to key");
376 
377  auto matrixIter = m_globalMat->find(gkey);
378 
379  if (matrixIter == m_globalMat->end())
380  {
381  return 0;
382  }
383  else
384  {
385  return matrixIter->second->GetNumNonZeroEntries();
386  }
387 
388  return 0;
389 }
390 
391 } // namespace MultiRegions
392 } // namespace Nektar
393 
394 #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:249
#define MULTI_REGIONS_EXPORT
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
This class is the abstraction of a global continuous two- dimensional spectral/hp element expansion w...
Definition: ContField.h:57
virtual 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) override
Definition: ContField.cpp:985
virtual void v_FillBndCondFromField(const Array< OneD, NekDouble > coeffs) override
Definition: ContField.cpp:628
virtual ~ContField()
The default destructor.
Definition: ContField.cpp:227
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:183
GlobalLinSysSharedPtr GenGlobalLinSys(const GlobalLinSysKey &mkey)
Definition: ContField.cpp:551
virtual GlobalLinSysKey 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) override
Solves the two-dimensional Helmholtz equation, subject to the boundary conditions specified.
Definition: ContField.cpp:794
const GJPStabilisationSharedPtr GetGJPForcing()
Definition: ContField.h:149
AssemblyMapCGSharedPtr m_locToGloMap
(A shared pointer to) the object which contains all the required information for the transformation f...
Definition: ContField.h:172
virtual void v_FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Template method virtual forwarder for FwdTrans().
Definition: ContField.cpp:562
int v_GetPoolCount(std::string) override
Definition: ContField.cpp:1027
virtual void v_SmoothField(Array< OneD, NekDouble > &field) override
Template method virtual forwarded for SmoothField().
Definition: ContField.cpp:261
virtual const Array< OneD, const SpatialDomains ::BoundaryConditionShPtr > & v_GetBndConditions() override
Template method virtual forwarder for GetBndConditions().
Definition: ContField.cpp:1011
const Array< OneD, const SpatialDomains::BoundaryConditionShPtr > & GetBndConditions()
Returns the boundary conditions.
Definition: ContField.h:366
virtual void v_ClearGlobalLinSysManager(void) override
Definition: ContField.cpp:1019
void MultiplyByInvMassMatrix(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Multiply a solution by the inverse mass matrix.
Definition: ContField.cpp:279
virtual GlobalLinSysKey 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) override
Definition: ContField.cpp:905
virtual void v_GlobalToLocal(void) override
Definition: ContField.cpp:725
const AssemblyMapCGSharedPtr & GetLocalToGlobalMap() const
Returns the map from local to global level.
Definition: ContField.h:354
void v_UnsetGlobalLinSys(GlobalLinSysKey, bool) override
Definition: ContField.cpp:1037
int GetGlobalMatrixNnz(const GlobalMatrixKey &gkey)
Definition: ContField.h:371
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:339
virtual void v_LocalToGlobal(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool useComm) override
Gathers the global coefficients from the local coefficients .
Definition: ContField.cpp:751
void SetGJPForcing(const GJPStabilisationSharedPtr &GJPData)
Definition: ContField.h:161
GlobalMatrixMapShPtr m_globalMat
(A shared pointer to) a list which collects all the global matrices being assembled,...
Definition: ContField.h:177
GlobalLinSysSharedPtr GetGlobalLinSys(const GlobalLinSysKey &mkey)
Returns the linear system specified by the key mkey.
Definition: ContField.cpp:546
virtual void v_ImposeDirichletConditions(Array< OneD, NekDouble > &outarray) override
Impose the Dirichlet Boundary Conditions on outarray.
Definition: ContField.cpp:568
virtual void v_MultiplyByInvMassMatrix(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Template method virtual forwarder for MultiplyByInvMassMatrix().
Definition: ContField.cpp:766
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:244
GlobalMatrixSharedPtr GetGlobalMatrix(const GlobalMatrixKey &mkey)
Returns the global matrix specified by mkey.
Definition: ContField.cpp:516
ContField()
The default constructor.
Definition: ContField.cpp:89
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:417
const Array< OneD, const MultiRegions::ExpListSharedPtr > & GetBndCondExpansions()
Returns the boundary conditions expansion.
Definition: ContField.h:360
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:489
void Assemble()
Assembles the global coefficients from the local coefficients .
Definition: ContField.h:316
GJPStabilisationSharedPtr m_GJPData
Data for Gradient Jump Penalisation (GJP) stabilisaiton.
Definition: ContField.h:186
This class is the abstractio n of a global discontinuous two- dimensional spectral/hp element expansi...
Definition: DisContField.h:58
Array< OneD, SpatialDomains::BoundaryConditionShPtr > m_bndConditions
An array which contains the information about the boundary condition structure definition on the diff...
Definition: DisContField.h:137
Array< OneD, MultiRegions::ExpListSharedPtr > m_bndCondExpansions
An object which contains the discretised boundary conditions.
Definition: DisContField.h:150
Base class for all multi-elemental spectral/hp expansions.
Definition: ExpList.h:103
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Definition: ExpList.h:1076
std::shared_ptr< ExpList > GetSharedThisPtr()
Returns a shared pointer to the current object.
Definition: ExpList.h:950
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:86
std::shared_ptr< GlobalMatrixMap > GlobalMatrixMapShPtr
Shared pointer to a global matrix map.
Definition: GlobalMatrix.h:90
std::map< StdRegions::ConstFactorType, Array< OneD, NekDouble > > VarFactorsMap
std::shared_ptr< GJPStabilisation > GJPStabilisationSharedPtr
std::shared_ptr< ContField > ContFieldSharedPtr
Definition: ContField.h:289
std::shared_ptr< BoundaryConditionBase > BoundaryConditionShPtr
Definition: Conditions.h:212
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:172
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:399
std::map< StdRegions::VarCoeffType, VarCoeffEntry > VarCoeffMap
Definition: StdRegions.hpp:343
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
static Array< OneD, Array< OneD, NekDouble > > NullNekDoubleArrayOfArray
static Array< OneD, NekDouble > NullNekDouble1DArray
double NekDouble