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);
202 
204 
205  MULTI_REGIONS_EXPORT virtual void 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, bool useComm);
212 
213  MULTI_REGIONS_EXPORT virtual void v_LocalToGlobal(bool useComm);
214 
215  /// Scatters from the global coefficients
216  /// \f$\boldsymbol{\hat{u}}_g\f$ to the local coefficients
217  /// \f$\boldsymbol{\hat{u}}_l\f$.
219  const Array<OneD, const NekDouble> &inarray,
220  Array<OneD, NekDouble> &outarray);
221 
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 
229  /// Template method virtual forwarder for FwdTrans().
230  MULTI_REGIONS_EXPORT virtual void v_FwdTrans(
231  const Array<OneD, const NekDouble> &inarray,
232  Array<OneD, NekDouble> &outarray);
233 
234  /// Template method virtual forwarded for SmoothField().
235  MULTI_REGIONS_EXPORT virtual void v_SmoothField(
236  Array<OneD, NekDouble> &field);
237 
238  /// Template method virtual forwarder for MultiplyByInvMassMatrix().
240  const Array<OneD, const NekDouble> &inarray,
241  Array<OneD, NekDouble> &outarray);
242 
243  /// Solves the two-dimensional Helmholtz equation, subject to the
244  /// boundary conditions specified.
245  MULTI_REGIONS_EXPORT virtual void v_HelmSolve(
246  const Array<OneD, const NekDouble> &inarray,
247  Array<OneD, NekDouble> &outarray,
248  const StdRegions::ConstFactorMap &factors,
249  const StdRegions::VarCoeffMap &varcoeff,
250  const MultiRegions::VarFactorsMap &varfactors,
251  const Array<OneD, const NekDouble> &dirForcing,
252  const bool PhysSpaceForcing);
253 
254  // Solve the linear advection problem assuming that m_coeffs
255  // vector contains an intial estimate for solution
257  const Array<OneD, Array<OneD, NekDouble>> &velocity,
258  const Array<OneD, const NekDouble> &inarray,
259  Array<OneD, NekDouble> &outarray, const NekDouble lambda,
261 
262  // Solve the linear advection problem assuming that m_coeff
263  // vector contains an intial estimate for solution
265  const Array<OneD, Array<OneD, NekDouble>> &velocity,
266  const Array<OneD, const NekDouble> &inarray,
267  Array<OneD, NekDouble> &outarray, const NekDouble lambda,
269 
270  /// Template method virtual forwarder for GetBndConditions().
271  MULTI_REGIONS_EXPORT virtual const Array<
275 };
276 
277 typedef std::shared_ptr<ContField> ContFieldSharedPtr;
278 
279 /**
280  * This operation is evaluated as:
281  * \f{tabbing}
282  * \hspace{1cm} \= Do \= $e=$ $1, N_{\mathrm{el}}$ \\
283  * > > Do \= $i=$ $0,N_m^e-1$ \\
284  * > > > $\boldsymbol{\hat{u}}_g[\mbox{map}[e][i]] =
285  * \boldsymbol{\hat{u}}_g[\mbox{map}[e][i]]+\mbox{sign}[e][i] \cdot
286  * \boldsymbol{\hat{u}}^{e}[i]$\\
287  * > > continue\\
288  * > continue
289  * \f}
290  * where \a map\f$[e][i]\f$ is the mapping array and \a
291  * sign\f$[e][i]\f$ is an array of similar dimensions ensuring the
292  * correct modal connectivity between the different elements (both
293  * these arrays are contained in the data member #m_locToGloMap). This
294  * operation is equivalent to the gather operation
295  * \f$\boldsymbol{\hat{u}}_g=\mathcal{A}^{T}\boldsymbol{\hat{u}}_l\f$,
296  * where \f$\mathcal{A}\f$ is the
297  * \f$N_{\mathrm{eof}}\times N_{\mathrm{dof}}\f$ permutation matrix.
298  *
299  * @note The array #m_coeffs should be filled with the local
300  * coefficients \f$\boldsymbol{\hat{u}}_l\f$ and that the
301  * resulting global coefficients \f$\boldsymbol{\hat{u}}_g\f$
302  * will be stored in #m_coeffs.
303  */
304 inline void ContField::Assemble()
305 {
306  m_locToGloMap->Assemble(m_coeffs, m_coeffs);
307 }
308 
309 /**
310  * This operation is evaluated as:
311  * \f{tabbing}
312  * \hspace{1cm} \= Do \= $e=$ $1, N_{\mathrm{el}}$ \\
313  * > > Do \= $i=$ $0,N_m^e-1$ \\
314  * > > > $\boldsymbol{\hat{u}}_g[\mbox{map}[e][i]] =
315  * \boldsymbol{\hat{u}}_g[\mbox{map}[e][i]]+\mbox{sign}[e][i] \cdot
316  * \boldsymbol{\hat{u}}^{e}[i]$\\
317  * > > continue\\
318  * > continue
319  * \f}
320  * where \a map\f$[e][i]\f$ is the mapping array and \a
321  * sign\f$[e][i]\f$ is an array of similar dimensions ensuring the
322  * correct modal connectivity between the different elements (both
323  * these arrays are contained in the data member #m_locToGloMap). This
324  * operation is equivalent to the gather operation
325  * \f$\boldsymbol{\hat{u}}_g=\mathcal{A}^{T}\boldsymbol{\hat{u}}_l\f$,
326  * where \f$\mathcal{A}\f$ is the
327  * \f$N_{\mathrm{eof}}\times N_{\mathrm{dof}}\f$ permutation matrix.
328  *
329  * @param inarray An array of size \f$N_\mathrm{eof}\f$
330  * containing the local degrees of freedom
331  * \f$\boldsymbol{x}_l\f$.
332  * @param outarray The resulting global degrees of freedom
333  * \f$\boldsymbol{x}_g\f$ will be stored in this
334  * array of size \f$N_\mathrm{dof}\f$.
335  */
337  Array<OneD, NekDouble> &outarray) const
338 {
339  m_locToGloMap->Assemble(inarray, outarray);
340 }
341 
343 {
344  return m_locToGloMap;
345 }
346 
349 {
350  return m_bndCondExpansions;
351 }
352 
355 {
356  return m_bndConditions;
357 }
358 
360 {
362  "To use method must have a AssemblyMap "
363  "attached to key");
364 
365  auto matrixIter = m_globalMat->find(gkey);
366 
367  if (matrixIter == m_globalMat->end())
368  {
369  return 0;
370  }
371  else
372  {
373  return matrixIter->second->GetNumNonZeroEntries();
374  }
375 
376  return 0;
377 }
378 
379 } // namespace MultiRegions
380 } // namespace Nektar
381 
382 #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_MultiplyByInvMassMatrix(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Template method virtual forwarder for MultiplyByInvMassMatrix().
Definition: ContField.cpp:774
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
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:994
GlobalLinSysSharedPtr GenGlobalLinSys(const GlobalLinSysKey &mkey)
Definition: ContField.cpp:551
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_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:917
virtual void v_SmoothField(Array< OneD, NekDouble > &field)
Template method virtual forwarded for SmoothField().
Definition: ContField.cpp:265
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:759
const Array< OneD, const SpatialDomains::BoundaryConditionShPtr > & GetBndConditions()
Returns the boundary conditions.
Definition: ContField.h:354
virtual void v_FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Template method virtual forwarder for FwdTrans().
Definition: ContField.cpp:562
void MultiplyByInvMassMatrix(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Multiply a solution by the inverse mass matrix.
Definition: ContField.cpp:283
const AssemblyMapCGSharedPtr & GetLocalToGlobalMap() const
Returns the map from local to global level.
Definition: ContField.h:342
int GetGlobalMatrixNnz(const GlobalMatrixKey &gkey)
Definition: ContField.h:359
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:340
virtual void v_GlobalToLocal(void)
Definition: ContField.cpp:726
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
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:807
virtual void v_ImposeDirichletConditions(Array< OneD, NekDouble > &outarray)
Impose the Dirichlet Boundary Conditions on outarray.
Definition: ContField.cpp:568
virtual void v_ClearGlobalLinSysManager(void)
Definition: ContField.cpp:1028
GlobalLinSysSharedPtr GetGlobalLinSys(const GlobalLinSysKey &mkey)
Returns the linear system specified by the key mkey.
Definition: ContField.cpp:546
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:248
GlobalMatrixSharedPtr GetGlobalMatrix(const GlobalMatrixKey &mkey)
Returns the global matrix specified by mkey.
Definition: ContField.cpp:516
virtual void v_FillBndCondFromField()
Definition: ContField.cpp:627
ContField()
The default constructor.
Definition: ContField.cpp:89
virtual const Array< OneD, const SpatialDomains ::BoundaryConditionShPtr > & v_GetBndConditions()
Template method virtual forwarder for GetBndConditions().
Definition: ContField.cpp:1020
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:348
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:304
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: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:149
Array< OneD, MultiRegions::ExpListSharedPtr > m_bndCondExpansions
An object which contains the discretised boundary conditions.
Definition: DisContField.h:143
Base class for all multi-elemental spectral/hp expansions.
Definition: ExpList.h:101
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Definition: ExpList.h:1158
std::shared_ptr< ExpList > GetSharedThisPtr()
Returns a shared pointer to the current object.
Definition: ExpList.h:1034
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:277
std::shared_ptr< BoundaryConditionBase > BoundaryConditionShPtr
Definition: Conditions.h:212
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:172
std::map< StdRegions::VarCoeffType, Array< OneD, NekDouble > > VarCoeffMap
Definition: StdRegions.hpp:240
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:282
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