Nektar++
AssemblyMap.h
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: AssemblyMap.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: Assembly (e.g. local to global) base mapping routines
32//
33///////////////////////////////////////////////////////////////////////////////
34
35#ifndef MULTIREGIONS_ASSEMBLY_MAP_H
36#define MULTIREGIONS_ASSEMBLY_MAP_H
37
43#include <vector>
44
45namespace Nektar
46{
47namespace MultiRegions
48{
49// Forward declarations
50class AssemblyMap;
51class ExpList;
52typedef std::shared_ptr<AssemblyMap> AssemblyMapSharedPtr;
54
55/// Base class for constructing local to global mapping of degrees of
56/// freedom.
58{
59public:
60 /// Default constructor.
62 /// Constructor with a communicator
66 const std::string variable = "DefaultVar");
67
68 /// Constructor for next level in multi-level static condensation.
70 AssemblyMap *oldLevelMap,
71 const BottomUpSubStructuredGraphSharedPtr &multiLevelGraph);
72 /// Destructor.
74
75 /// Retrieves the communicator
77
78 /// Retrieves the variable string
80
81 /// Retrieves the hash of this map
82 MULTI_REGIONS_EXPORT size_t GetHash() const;
83
84 MULTI_REGIONS_EXPORT int GetLocalToGlobalMap(const int i) const;
85
86 MULTI_REGIONS_EXPORT int GetGlobalToUniversalMap(const int i) const;
87
89
91
94
97
99
101 const;
102
105 bool useComm = true) const;
106
108 NekVector<NekDouble> &global,
109 bool useComm = true) const;
110
112 const Array<OneD, const NekDouble> &global,
114
117
119 Array<OneD, NekDouble> &global) const;
120
122 NekVector<NekDouble> &global) const;
123
125 Array<OneD, NekDouble> &pGlobal) const;
126
128 NekVector<NekDouble> &pGlobal) const;
129
131 int offset) const;
132
134 Array<OneD, NekDouble> &bndvals);
135
136 /// Retrieve the global index of a given local boundary mode.
137 MULTI_REGIONS_EXPORT int GetLocalToGlobalBndMap(const int i) const;
138 /// Retrieve the global indices of the local boundary modes.
140
143
146
147 /// Returns true if using a modal expansion requiring a change of
148 /// sign of some modes.
150
151 /// Retrieve the sign change of a given local boundary mode.
153 /// Retrieve the sign change for all local boundary modes.
155 const;
156 /// Retrieves the local indices corresponding to the
157 /// boundary expansion modes.
160 /// Returns the modal sign associated with a given
161 /// boundary expansion mode.
164
165 /// Retrieves the local indices corresponding to the
166 /// boundary expansion modes to global trace
169
170 /// Returns the global index of the boundary trace giving the
171 /// index on the boundary expansion
175
176 /// Returns the number of global Dirichlet boundary coefficients.
178 /// Returns the number of local Dirichlet boundary coefficients.
180 /// Returns the total number of global boundary coefficients.
182 /// Returns the total number of local boundary coefficients.
184 /// Returns the total number of local coefficients.
186 /// Returns the total number of global coefficients.
188 /// Retrieves if the system is singular (true) or not (false)
190
191 ///
194 int offset) const;
195
197 const NekVector<NekDouble> &global, NekVector<NekDouble> &loc) const;
198
201 int offset) const;
202
204 const Array<OneD, const NekDouble> &global,
206
209 int offset, bool UseComm = true) const;
210
213 bool UseComm = true) const;
214
216 const Array<OneD, const NekDouble> &local,
217 Array<OneD, NekDouble> &locbnd) const;
218
220 const Array<OneD, const NekDouble> &local,
221 Array<OneD, NekDouble> &locint) const;
222
224 const Array<OneD, const NekDouble> &locbnd,
225 Array<OneD, NekDouble> &local) const;
226
228 const Array<OneD, const NekDouble> &locbnd,
229 Array<OneD, NekDouble> &local) const;
230
232 NekVector<NekDouble> &global,
233 int offset) const;
234
236 NekVector<NekDouble> &global) const;
237
240 int offset) const;
241
244 Array<OneD, NekDouble> &global) const;
245
247 Array<OneD, NekDouble> &pGlobal) const;
248
250 NekVector<NekDouble> &pGlobal) const;
251
253 Array<OneD, NekDouble> &pGlobal, int offset) const;
254
256
258
260
262
264
266
268
270
271 MULTI_REGIONS_EXPORT void PrintStats(std::ostream &out,
272 std::string variable,
273 bool printHeader = true) const;
274
276
277 MULTI_REGIONS_EXPORT std::shared_ptr<AssemblyMap> LinearSpaceMap(
278 const ExpList &locexp, GlobalSysSolnType solnType);
279
280 /// Returns the bandwidth of the boundary system.
282 /// Returns the level of static condensation for this map.
284 /// Returns the number of patches in this static condensation level.
286 /// Returns the number of local boundary coefficients in each patch.
289 /// Returns the number of local interior coefficients in each patch.
292 /// Returns the local to global mapping for the next level in the
293 /// multi-level static condensation.
296
298 AssemblyMapSharedPtr pNextLevelLocalToGlobalMap);
299
300 /// Returns the patch map from the previous level
301 /// of the multi-level static condensation.
303 void) const;
304
305 /// Returns true if this is the last level in the multi-level
306 /// static condensation.
308 /// Returns the method of solving global systems.
310 MULTI_REGIONS_EXPORT std::string GetPreconType() const;
315 MULTI_REGIONS_EXPORT std::string GetLinSysIterSolver() const;
316
318 {
320 }
321
324 Array<OneD, NekDouble> &global) const;
325
327 const Array<OneD, const NekDouble> &global,
329
332 Array<OneD, NekDouble> &global) const;
333
334protected:
335 /// Session object
337
338 /// Communicator
340
341 /// Variable string identifier
342 std::string m_variable;
343
344 /// Hash for map
345 size_t m_hash;
346
347 /// Number of local boundary coefficients
349 /// Total number of global boundary coefficients
351 /// Number of Local Dirichlet Boundary Coefficients
353 /// Number of Global Dirichlet Boundary Coefficients
355 /// Flag indicating if the system is singular or not
357
358 /// Total number of local coefficients
359 /** This corresponds to the number of total number of coefficients
360 * - For CG this corresponds to the total of bnd + int DOFs
361 * - For DG this corresponds to the number of bnd DOFs.
362 * This means that #m_numLocalCoeffs = #m_numLocalBndCoeffs
363 * This way, we can consider the trace-system solve as a
364 * statically condensed solve without interior DOFs. This allows
365 * us to use the same global system solver for both cases.
366 */
368
369 /// Total number of global coefficients
370 /** This corresponds to the number of total number of coefficients
371 * - For CG this corresponds to the total of bnd + int DOFs.
372 * - For DG this corresponds to the number of bnd DOFs.
373 * This means that #m_numGlobalCoeffs = #m_numGlobalBndCoeffs
374 * This way, we can consider the trace-system solve as a
375 * statically condensed solve without interior DOFs. This allows
376 * us to use the same global system solver for both cases.
377 */
379
380 /// Flag indicating if modes require sign reversal.
382
383 /// Integer map of local coeffs to global Boundary Dofs
385 /// Integer sign of local boundary coeffs to global space
387 /// Integer map of local boundary coeffs to local boundary system numbering
389 /// Integer map of local boundary coeffs to local interior system numbering
391 /// Integer map of bnd cond coeffs to local coefficients
393 /// Integer map of sign of bnd cond coeffs to local coefficients
395 /// Integer map of bnd cond coeff to local trace coeff
397 /// Integer map of bnd cond trace number to global trace number
399 /// Integer map of process coeffs to universal space
401 /// Integer map of unique process coeffs to universal space (signed)
403
404 /// The solution type of the global system
406 /// The bandwith of the global bnd system
408
409 /// Type type of preconditioner to use in iterative solver.
410 std::string m_preconType;
411
412 /// Maximum iterations for iterative solver
414
415 /// Tolerance for iterative solver
418
419 /// sucessive RHS for iterative solver
421
422 /// Iterative solver: Conjugate Gradient, GMRES
424
427 /// gs gather communication to impose Dirhichlet BCs.
429
430 /// The level of recursion in the case of multi-level static
431 /// condensation.
433 /// The number of patches (~elements) in the current level
435 /// The number of bnd dofs per patch
437 /// The number of int dofs per patch
439 /// Map from the patches of the previous level to the patches of
440 /// the current level
441
442 /// The local to global mapping of the next level of recursion
444 /// Lowest static condensation level.
446
447 /// Calculates the bandwidth of the boundary system.
449
452
453 virtual int v_GetLocalToGlobalMap(const int i) const;
454
455 virtual int v_GetGlobalToUniversalMap(const int i) const;
456
457 virtual int v_GetGlobalToUniversalMapUnique(const int i) const;
458
460
462
464
465 virtual NekDouble v_GetLocalToGlobalSign(const int i) const;
466
467 virtual const Array<OneD, NekDouble> &v_GetLocalToGlobalSign() const;
468
471 bool useComm) const;
472
473 virtual void v_LocalToGlobal(const NekVector<NekDouble> &loc,
474 NekVector<NekDouble> &global,
475 bool useComm) const;
476
477 virtual void v_GlobalToLocal(const Array<OneD, const NekDouble> &global,
479
480 virtual void v_GlobalToLocal(const NekVector<NekDouble> &global,
482
483 virtual void v_Assemble(const Array<OneD, const NekDouble> &loc,
484 Array<OneD, NekDouble> &global) const;
485
486 virtual void v_Assemble(const NekVector<NekDouble> &loc,
487 NekVector<NekDouble> &global) const;
488
489 virtual void v_UniversalAssemble(Array<OneD, NekDouble> &pGlobal) const;
490
491 virtual void v_UniversalAssemble(NekVector<NekDouble> &pGlobal) const;
492
493 virtual void v_UniversalAssemble(Array<OneD, NekDouble> &pGlobal,
494 int offset) const;
495
496 virtual int v_GetFullSystemBandWidth() const;
497
498 virtual int v_GetNumNonDirVertexModes() const;
499
500 virtual int v_GetNumNonDirEdgeModes() const;
501
502 virtual int v_GetNumNonDirFaceModes() const;
503
504 virtual int v_GetNumDirEdges() const;
505
506 virtual int v_GetNumDirFaces() const;
507
508 virtual int v_GetNumNonDirEdges() const;
509
510 virtual int v_GetNumNonDirFaces() const;
511
513
514 /// Generate a linear space mapping from existing mapping
515 virtual std::shared_ptr<AssemblyMap> v_LinearSpaceMap(
516 const ExpList &locexp, GlobalSysSolnType solnType);
517
518private:
519 /// Mapping information for previous level in MultiLevel Solver
521};
522
523} // namespace MultiRegions
524} // namespace Nektar
525
526#endif // MULTIREGIONS_ASSEMBLY_MAP_H
#define MULTI_REGIONS_EXPORT
Base class for constructing local to global mapping of degrees of freedom.
Definition: AssemblyMap.h:58
int m_lowestStaticCondLevel
Lowest static condensation level.
Definition: AssemblyMap.h:445
void PatchLocalToGlobal(const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global) const
const Array< OneD, const int > & GetExtraDirEdges()
GlobalSysSolnType m_solnType
The solution type of the global system.
Definition: AssemblyMap.h:405
int GetNumPatches() const
Returns the number of patches in this static condensation level.
const Array< OneD, const int > & GetGlobalToUniversalMapUnique()
int GetNumGlobalCoeffs() const
Returns the total number of global coefficients.
int m_numLocalCoeffs
Total number of local coefficients.
Definition: AssemblyMap.h:367
void SetNextLevelLocalToGlobalMap(AssemblyMapSharedPtr pNextLevelLocalToGlobalMap)
Array< OneD, int > m_bndCondCoeffsToLocalTraceMap
Integer map of bnd cond coeff to local trace coeff.
Definition: AssemblyMap.h:396
const Array< OneD, const unsigned int > & GetNumLocalBndCoeffsPerPatch()
Returns the number of local boundary coefficients in each patch.
Array< OneD, int > m_bndCondCoeffsToLocalCoeffsMap
Integer map of bnd cond coeffs to local coefficients.
Definition: AssemblyMap.h:392
virtual const Array< OneD, const int > & v_GetExtraDirEdges()
int GetNumLocalCoeffs() const
Returns the total number of local coefficients.
std::string GetVariable()
Retrieves the variable string.
const Array< OneD, const unsigned int > & GetNumLocalIntCoeffsPerPatch()
Returns the number of local interior coefficients in each patch.
virtual void v_Assemble(const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global) const
bool m_signChange
Flag indicating if modes require sign reversal.
Definition: AssemblyMap.h:381
void Assemble(const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global) const
virtual std::shared_ptr< AssemblyMap > v_LinearSpaceMap(const ExpList &locexp, GlobalSysSolnType solnType)
Generate a linear space mapping from existing mapping.
void AssembleBnd(const NekVector< NekDouble > &loc, NekVector< NekDouble > &global, int offset) const
virtual void v_LocalToGlobal(const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global, bool useComm) const
PatchMapSharedPtr m_patchMapFromPrevLevel
Mapping information for previous level in MultiLevel Solver.
Definition: AssemblyMap.h:520
int m_maxIterations
Maximum iterations for iterative solver.
Definition: AssemblyMap.h:413
Array< OneD, int > m_localToLocalIntMap
Integer map of local boundary coeffs to local interior system numbering.
Definition: AssemblyMap.h:390
virtual const Array< OneD, const int > & v_GetGlobalToUniversalMapUnique()
NekDouble m_iterativeTolerance
Tolerance for iterative solver.
Definition: AssemblyMap.h:416
NekDouble GetIterativeTolerance() const
int m_numGlobalCoeffs
Total number of global coefficients.
Definition: AssemblyMap.h:378
std::string m_linSysIterSolver
Iterative solver: Conjugate Gradient, GMRES.
Definition: AssemblyMap.h:423
bool AtLastLevel() const
Returns true if this is the last level in the multi-level static condensation.
void UniversalAssembleBnd(Array< OneD, NekDouble > &pGlobal) const
const Array< OneD, const int > & GetLocalToGlobalMap()
const Array< OneD, const int > & GetBndCondCoeffsToLocalTraceMap()
Retrieves the local indices corresponding to the boundary expansion modes to global trace.
virtual const Array< OneD, const int > & v_GetGlobalToUniversalMap()
Array< OneD, int > m_globalToUniversalBndMap
Integer map of process coeffs to universal space.
Definition: AssemblyMap.h:400
void PatchAssemble(const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global) const
int m_successiveRHS
sucessive RHS for iterative solver
Definition: AssemblyMap.h:420
std::shared_ptr< AssemblyMap > LinearSpaceMap(const ExpList &locexp, GlobalSysSolnType solnType)
virtual int v_GetNumNonDirVertexModes() const
void CalculateBndSystemBandWidth()
Calculates the bandwidth of the boundary system.
void LocalToGlobal(const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global, bool useComm=true) const
Array< OneD, NekDouble > m_localToGlobalBndSign
Integer sign of local boundary coeffs to global space.
Definition: AssemblyMap.h:386
void UniversalAbsMaxBnd(Array< OneD, NekDouble > &bndvals)
Array< OneD, unsigned int > m_numLocalBndCoeffsPerPatch
The number of bnd dofs per patch.
Definition: AssemblyMap.h:436
virtual const Array< OneD, const int > & v_GetLocalToGlobalMap()
const Array< OneD, const int > & GetGlobalToUniversalMap()
std::string GetLinSysIterSolver() const
virtual int v_GetNumNonDirEdgeModes() const
virtual int v_GetNumDirFaces() const
LibUtilities::SessionReaderSharedPtr m_session
Session object.
Definition: AssemblyMap.h:336
virtual const Array< OneD, NekDouble > & v_GetLocalToGlobalSign() const
int GetNumGlobalDirBndCoeffs() const
Returns the number of global Dirichlet boundary coefficients.
const Array< OneD, NekDouble > & GetBndCondCoeffsToLocalCoeffsSign()
Returns the modal sign associated with a given boundary expansion mode.
int m_numLocalBndCoeffs
Number of local boundary coefficients.
Definition: AssemblyMap.h:348
void LocalToLocalBnd(const Array< OneD, const NekDouble > &local, Array< OneD, NekDouble > &locbnd) const
const Array< OneD, const int > & GetBndCondCoeffsToLocalCoeffsMap()
Retrieves the local indices corresponding to the boundary expansion modes.
void PrintStats(std::ostream &out, std::string variable, bool printHeader=true) const
AssemblyMapSharedPtr m_nextLevelLocalToGlobalMap
Map from the patches of the previous level to the patches of the current level.
Definition: AssemblyMap.h:443
const Array< OneD, const int > & GetLocalToGlobalBndMap()
Retrieve the global indices of the local boundary modes.
int m_staticCondLevel
The level of recursion in the case of multi-level static condensation.
Definition: AssemblyMap.h:432
int m_bndSystemBandWidth
The bandwith of the global bnd system.
Definition: AssemblyMap.h:407
const Array< OneD, const int > & GetGlobalToUniversalBndMapUnique()
void GlobalToLocal(const Array< OneD, const NekDouble > &global, Array< OneD, NekDouble > &loc) const
GlobalSysSolnType GetGlobalSysSolnType() const
Returns the method of solving global systems.
int m_numLocalDirBndCoeffs
Number of Local Dirichlet Boundary Coefficients.
Definition: AssemblyMap.h:352
bool GetSingularSystem() const
Retrieves if the system is singular (true) or not (false)
std::string m_variable
Variable string identifier.
Definition: AssemblyMap.h:342
virtual int v_GetNumNonDirFaceModes() const
void PatchGlobalToLocal(const Array< OneD, const NekDouble > &global, Array< OneD, NekDouble > &loc) const
virtual int v_GetFullSystemBandWidth() const
int m_numGlobalDirBndCoeffs
Number of Global Dirichlet Boundary Coefficients.
Definition: AssemblyMap.h:354
LibUtilities::CommSharedPtr GetComm()
Retrieves the communicator.
Array< OneD, unsigned int > m_numLocalIntCoeffsPerPatch
The number of int dofs per patch.
Definition: AssemblyMap.h:438
const AssemblyMapSharedPtr GetNextLevelLocalToGlobalMap() const
Returns the local to global mapping for the next level in the multi-level static condensation.
int GetBndSystemBandWidth() const
Returns the bandwidth of the boundary system.
virtual void v_UniversalAssemble(Array< OneD, NekDouble > &pGlobal) const
int GetStaticCondLevel() const
Returns the level of static condensation for this map.
void GlobalToLocalBndWithoutSign(const Array< OneD, const NekDouble > &global, Array< OneD, NekDouble > &loc)
bool m_systemSingular
Flag indicating if the system is singular or not.
Definition: AssemblyMap.h:356
Array< OneD, int > m_localToGlobalBndMap
Integer map of local coeffs to global Boundary Dofs.
Definition: AssemblyMap.h:384
size_t GetHash() const
Retrieves the hash of this map.
virtual ~AssemblyMap()
Destructor.
Gs::gs_data * m_dirBndGsh
gs gather communication to impose Dirhichlet BCs.
Definition: AssemblyMap.h:428
virtual void v_GlobalToLocal(const Array< OneD, const NekDouble > &global, Array< OneD, NekDouble > &loc) const
const Array< OneD, NekDouble > & GetLocalToGlobalSign() const
int GetNumLocalDirBndCoeffs() const
Returns the number of local Dirichlet boundary coefficients.
Array< OneD, const NekDouble > GetLocalToGlobalBndSign() const
Retrieve the sign change for all local boundary modes.
void LocalBndToLocal(const Array< OneD, const NekDouble > &locbnd, Array< OneD, NekDouble > &local) const
AssemblyMap()
Default constructor.
Definition: AssemblyMap.cpp:82
Array< OneD, int > m_globalToUniversalBndMapUnique
Integer map of unique process coeffs to universal space (signed)
Definition: AssemblyMap.h:402
virtual int v_GetNumNonDirEdges() const
virtual int v_GetNumDirEdges() const
const PatchMapSharedPtr & GetPatchMapFromPrevLevel(void) const
Returns the patch map from the previous level of the multi-level static condensation.
LibUtilities::CommSharedPtr m_comm
Communicator.
Definition: AssemblyMap.h:339
int GetNumGlobalBndCoeffs() const
Returns the total number of global boundary coefficients.
const Array< OneD, const int > & GetGlobalToUniversalBndMap()
Array< OneD, int > m_localToLocalBndMap
Integer map of local boundary coeffs to local boundary system numbering.
Definition: AssemblyMap.h:388
Array< OneD, NekDouble > m_bndCondCoeffsToLocalCoeffsSign
Integer map of sign of bnd cond coeffs to local coefficients.
Definition: AssemblyMap.h:394
const Array< OneD, const int > & GetBndCondIDToGlobalTraceID()
void LocalBndToGlobal(const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global, int offset, bool UseComm=true) const
void UniversalAssemble(Array< OneD, NekDouble > &pGlobal) const
int m_numPatches
The number of patches (~elements) in the current level.
Definition: AssemblyMap.h:434
void LocalToLocalInt(const Array< OneD, const NekDouble > &local, Array< OneD, NekDouble > &locint) const
virtual int v_GetNumNonDirFaces() const
bool GetSignChange()
Returns true if using a modal expansion requiring a change of sign of some modes.
void LocalIntToLocal(const Array< OneD, const NekDouble > &locbnd, Array< OneD, NekDouble > &local) const
void GlobalToLocalBnd(const NekVector< NekDouble > &global, NekVector< NekDouble > &loc, int offset) const
std::string m_preconType
Type type of preconditioner to use in iterative solver.
Definition: AssemblyMap.h:410
int GetNumLocalBndCoeffs() const
Returns the total number of local boundary coefficients.
int m_numGlobalBndCoeffs
Total number of global boundary coefficients.
Definition: AssemblyMap.h:350
Array< OneD, int > m_bndCondIDToGlobalTraceID
Integer map of bnd cond trace number to global trace number.
Definition: AssemblyMap.h:398
Base class for all multi-elemental spectral/hp expansions.
Definition: ExpList.h:102
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:57
static AssemblyMapSharedPtr NullAssemblyMapSharedPtr
Definition: AssemblyMap.h:53
std::shared_ptr< BottomUpSubStructuredGraph > BottomUpSubStructuredGraphSharedPtr
std::shared_ptr< AssemblyMap > AssemblyMapSharedPtr
Definition: AssemblyMap.h:52
std::shared_ptr< PatchMap > PatchMapSharedPtr
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
double NekDouble