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 
45 namespace Nektar
46 {
47 namespace MultiRegions
48 {
49 // Forward declarations
50 class AssemblyMap;
51 class ExpList;
52 typedef std::shared_ptr<AssemblyMap> AssemblyMapSharedPtr;
54 
55 /// Base class for constructing local to global mapping of degrees of
56 /// freedom.
58 {
59 public:
60  /// Default constructor.
62  /// Constructor with a communicator
65  const LibUtilities::CommSharedPtr &comm,
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 hash of this map
79  MULTI_REGIONS_EXPORT size_t GetHash() const;
80 
81  MULTI_REGIONS_EXPORT int GetLocalToGlobalMap(const int i) const;
82 
83  MULTI_REGIONS_EXPORT int GetGlobalToUniversalMap(const int i) const;
84 
86 
88 
91 
94 
96 
98  const;
99 
102  bool useComm = true) const;
103 
105  NekVector<NekDouble> &global,
106  bool useComm = true) const;
107 
109  const Array<OneD, const NekDouble> &global,
110  Array<OneD, NekDouble> &loc) const;
111 
113  NekVector<NekDouble> &loc) const;
114 
116  Array<OneD, NekDouble> &global) const;
117 
119  NekVector<NekDouble> &global) const;
120 
122  Array<OneD, NekDouble> &pGlobal) const;
123 
125  NekVector<NekDouble> &pGlobal) const;
126 
128  int offset) const;
129 
131  Array<OneD, NekDouble> &bndvals);
132 
133  /// Retrieve the global index of a given local boundary mode.
134  MULTI_REGIONS_EXPORT int GetLocalToGlobalBndMap(const int i) const;
135  /// Retrieve the global indices of the local boundary modes.
137 
140 
143 
144  /// Returns true if using a modal expansion requiring a change of
145  /// sign of some modes.
147 
148  /// Retrieve the sign change of a given local boundary mode.
150  /// Retrieve the sign change for all local boundary modes.
152  const;
153  /// Retrieves the local indices corresponding to the
154  /// boundary expansion modes.
157  /// Returns the modal sign associated with a given
158  /// boundary expansion mode.
161 
162  /// Retrieves the local indices corresponding to the
163  /// boundary expansion modes to global trace
166 
167  /// Returns the global index of the boundary trace giving the
168  /// index on the boundary expansion
172 
173  /// Returns the number of global Dirichlet boundary coefficients.
175  /// Returns the number of local Dirichlet boundary coefficients.
177  /// Returns the total number of global boundary coefficients.
179  /// Returns the total number of local boundary coefficients.
181  /// Returns the total number of local coefficients.
183  /// Returns the total number of global coefficients.
185  /// Retrieves if the system is singular (true) or not (false)
187 
188  ///
191  int offset) const;
192 
194  const NekVector<NekDouble> &global, NekVector<NekDouble> &loc) const;
195 
198  int offset) const;
199 
201  const Array<OneD, const NekDouble> &global,
202  Array<OneD, NekDouble> &loc) const;
203 
206  int offset, bool UseComm = true) const;
207 
210  bool UseComm = true) const;
211 
213  const Array<OneD, const NekDouble> &local,
214  Array<OneD, NekDouble> &locbnd) const;
215 
217  const Array<OneD, const NekDouble> &local,
218  Array<OneD, NekDouble> &locint) const;
219 
221  const Array<OneD, const NekDouble> &locbnd,
222  Array<OneD, NekDouble> &local) const;
223 
225  const Array<OneD, const NekDouble> &locbnd,
226  Array<OneD, NekDouble> &local) const;
227 
229  NekVector<NekDouble> &global,
230  int offset) const;
231 
233  NekVector<NekDouble> &global) const;
234 
237  int offset) const;
238 
241  Array<OneD, NekDouble> &global) const;
242 
244  Array<OneD, NekDouble> &pGlobal) const;
245 
247  NekVector<NekDouble> &pGlobal) const;
248 
250  Array<OneD, NekDouble> &pGlobal, int offset) const;
251 
253 
255 
257 
259 
261 
263 
265 
267 
268  MULTI_REGIONS_EXPORT void PrintStats(std::ostream &out,
269  std::string variable,
270  bool printHeader = true) const;
271 
273 
274  MULTI_REGIONS_EXPORT std::shared_ptr<AssemblyMap> LinearSpaceMap(
275  const ExpList &locexp, GlobalSysSolnType solnType);
276 
277  /// Returns the bandwidth of the boundary system.
279  /// Returns the level of static condensation for this map.
281  /// Returns the number of patches in this static condensation level.
283  /// Returns the number of local boundary coefficients in each patch.
286  /// Returns the number of local interior coefficients in each patch.
289  /// Returns the local to global mapping for the next level in the
290  /// multi-level static condensation.
293 
295  AssemblyMapSharedPtr pNextLevelLocalToGlobalMap);
296 
297  /// Returns the patch map from the previous level
298  /// of the multi-level static condensation.
300  void) const;
301 
302  /// Returns true if this is the last level in the multi-level
303  /// static condensation.
304  MULTI_REGIONS_EXPORT bool AtLastLevel() const;
305  /// Returns the method of solving global systems.
311  MULTI_REGIONS_EXPORT std::string GetLinSysIterSolver() const;
312 
314  {
316  }
317 
320  Array<OneD, NekDouble> &global) const;
321 
323  const Array<OneD, const NekDouble> &global,
324  Array<OneD, NekDouble> &loc) const;
325 
328  Array<OneD, NekDouble> &global) const;
329 
330 protected:
331  /// Session object
333 
334  /// Communicator
336 
337  /// Hash for map
338  size_t m_hash;
339 
340  /// Number of local boundary coefficients
342  /// Total number of global boundary coefficients
344  /// Number of Local Dirichlet Boundary Coefficients
346  /// Number of Global Dirichlet Boundary Coefficients
348  /// Flag indicating if the system is singular or not
350 
351  /// Total number of local coefficients
352  /** This corresponds to the number of total number of coefficients
353  * - For CG this corresponds to the total of bnd + int DOFs
354  * - For DG this corresponds to the number of bnd DOFs.
355  * This means that #m_numLocalCoeffs = #m_numLocalBndCoeffs
356  * This way, we can consider the trace-system solve as a
357  * statically condensed solve without interior DOFs. This allows
358  * us to use the same global system solver for both cases.
359  */
361 
362  /// Total number of global coefficients
363  /** This corresponds to the number of total number of coefficients
364  * - For CG this corresponds to the total of bnd + int DOFs.
365  * - For DG this corresponds to the number of bnd DOFs.
366  * This means that #m_numGlobalCoeffs = #m_numGlobalBndCoeffs
367  * This way, we can consider the trace-system solve as a
368  * statically condensed solve without interior DOFs. This allows
369  * us to use the same global system solver for both cases.
370  */
372 
373  /// Flag indicating if modes require sign reversal.
375 
376  /// Integer map of local coeffs to global Boundary Dofs
378  /// Integer sign of local boundary coeffs to global space
380  /// Integer map of local boundary coeffs to local boundary system numbering
382  /// Integer map of local boundary coeffs to local interior system numbering
384  /// Integer map of bnd cond coeffs to local coefficients
386  /// Integer map of sign of bnd cond coeffs to local coefficients
388  /// Integer map of bnd cond coeff to local trace coeff
390  /// Integer map of bnd cond trace number to global trace number
392  /// Integer map of process coeffs to universal space
394  /// Integer map of unique process coeffs to universal space (signed)
396 
397  /// The solution type of the global system
399  /// The bandwith of the global bnd system
401 
402  /// Type type of preconditioner to use in iterative solver.
404 
405  /// Maximum iterations for iterative solver
407 
408  /// Tolerance for iterative solver
410 
411  /// sucessive RHS for iterative solver
413 
414  /// Iterative solver: Conjugate Gradient, GMRES
415  std::string m_linSysIterSolver;
416 
419  /// gs gather communication to impose Dirhichlet BCs.
421 
422  /// The level of recursion in the case of multi-level static
423  /// condensation.
425  /// The number of patches (~elements) in the current level
427  /// The number of bnd dofs per patch
429  /// The number of int dofs per patch
431  /// Map from the patches of the previous level to the patches of
432  /// the current level
433 
434  /// The local to global mapping of the next level of recursion
436  /// Lowest static condensation level.
438 
439  /// Calculates the bandwidth of the boundary system.
441 
444 
445  virtual int v_GetLocalToGlobalMap(const int i) const;
446 
447  virtual int v_GetGlobalToUniversalMap(const int i) const;
448 
449  virtual int v_GetGlobalToUniversalMapUnique(const int i) const;
450 
452 
454 
456 
457  virtual NekDouble v_GetLocalToGlobalSign(const int i) const;
458 
459  virtual const Array<OneD, NekDouble> &v_GetLocalToGlobalSign() const;
460 
462  Array<OneD, NekDouble> &global,
463  bool useComm) const;
464 
465  virtual void v_LocalToGlobal(const NekVector<NekDouble> &loc,
466  NekVector<NekDouble> &global,
467  bool useComm) const;
468 
469  virtual void v_GlobalToLocal(const Array<OneD, const NekDouble> &global,
470  Array<OneD, NekDouble> &loc) const;
471 
472  virtual void v_GlobalToLocal(const NekVector<NekDouble> &global,
473  NekVector<NekDouble> &loc) const;
474 
475  virtual void v_Assemble(const Array<OneD, const NekDouble> &loc,
476  Array<OneD, NekDouble> &global) const;
477 
478  virtual void v_Assemble(const NekVector<NekDouble> &loc,
479  NekVector<NekDouble> &global) const;
480 
481  virtual void v_UniversalAssemble(Array<OneD, NekDouble> &pGlobal) const;
482 
483  virtual void v_UniversalAssemble(NekVector<NekDouble> &pGlobal) const;
484 
485  virtual void v_UniversalAssemble(Array<OneD, NekDouble> &pGlobal,
486  int offset) const;
487 
488  virtual int v_GetFullSystemBandWidth() const;
489 
490  virtual int v_GetNumNonDirVertexModes() const;
491 
492  virtual int v_GetNumNonDirEdgeModes() const;
493 
494  virtual int v_GetNumNonDirFaceModes() const;
495 
496  virtual int v_GetNumDirEdges() const;
497 
498  virtual int v_GetNumDirFaces() const;
499 
500  virtual int v_GetNumNonDirEdges() const;
501 
502  virtual int v_GetNumNonDirFaces() const;
503 
505 
506  /// Generate a linear space mapping from existing mapping
507  virtual std::shared_ptr<AssemblyMap> v_LinearSpaceMap(
508  const ExpList &locexp, GlobalSysSolnType solnType);
509 
510 private:
511  /// Mapping information for previous level in MultiLevel Solver
513 };
514 
515 } // namespace MultiRegions
516 } // namespace Nektar
517 
518 #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:437
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:398
int GetNumPatches() const
Returns the number of patches in this static condensation level.
PreconditionerType m_preconType
Type type of preconditioner to use in iterative solver.
Definition: AssemblyMap.h:403
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:360
void SetNextLevelLocalToGlobalMap(AssemblyMapSharedPtr pNextLevelLocalToGlobalMap)
Array< OneD, int > m_bndCondCoeffsToLocalTraceMap
Integer map of bnd cond coeff to local trace coeff.
Definition: AssemblyMap.h:389
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:385
virtual const Array< OneD, const int > & v_GetExtraDirEdges()
int GetNumLocalCoeffs() const
Returns the total number of local coefficients.
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:374
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
PreconditionerType GetPreconType() const
PatchMapSharedPtr m_patchMapFromPrevLevel
Mapping information for previous level in MultiLevel Solver.
Definition: AssemblyMap.h:512
int m_maxIterations
Maximum iterations for iterative solver.
Definition: AssemblyMap.h:406
Array< OneD, int > m_localToLocalIntMap
Integer map of local boundary coeffs to local interior system numbering.
Definition: AssemblyMap.h:383
virtual const Array< OneD, const int > & v_GetGlobalToUniversalMapUnique()
NekDouble m_iterativeTolerance
Tolerance for iterative solver.
Definition: AssemblyMap.h:409
NekDouble GetIterativeTolerance() const
int m_numGlobalCoeffs
Total number of global coefficients.
Definition: AssemblyMap.h:371
std::string m_linSysIterSolver
Iterative solver: Conjugate Gradient, GMRES.
Definition: AssemblyMap.h:415
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:393
void PatchAssemble(const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global) const
int m_successiveRHS
sucessive RHS for iterative solver
Definition: AssemblyMap.h:412
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:379
void UniversalAbsMaxBnd(Array< OneD, NekDouble > &bndvals)
Array< OneD, unsigned int > m_numLocalBndCoeffsPerPatch
The number of bnd dofs per patch.
Definition: AssemblyMap.h:428
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:332
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:341
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:435
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:424
int m_bndSystemBandWidth
The bandwith of the global bnd system.
Definition: AssemblyMap.h:400
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:345
bool GetSingularSystem() const
Retrieves if the system is singular (true) or not (false)
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:347
LibUtilities::CommSharedPtr GetComm()
Retrieves the communicator.
Array< OneD, unsigned int > m_numLocalIntCoeffsPerPatch
The number of int dofs per patch.
Definition: AssemblyMap.h:430
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:349
Array< OneD, int > m_localToGlobalBndMap
Integer map of local coeffs to global Boundary Dofs.
Definition: AssemblyMap.h:377
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:420
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:81
Array< OneD, int > m_globalToUniversalBndMapUnique
Integer map of unique process coeffs to universal space (signed)
Definition: AssemblyMap.h:395
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:335
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:381
Array< OneD, NekDouble > m_bndCondCoeffsToLocalCoeffsSign
Integer map of sign of bnd cond coeffs to local coefficients.
Definition: AssemblyMap.h:387
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:426
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
int GetNumLocalBndCoeffs() const
Returns the total number of local boundary coefficients.
int m_numGlobalBndCoeffs
Total number of global boundary coefficients.
Definition: AssemblyMap.h:343
Array< OneD, int > m_bndCondIDToGlobalTraceID
Integer map of bnd cond trace number to global trace number.
Definition: AssemblyMap.h:391
Base class for all multi-elemental spectral/hp expansions.
Definition: ExpList.h:103
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:54
static AssemblyMapSharedPtr NullAssemblyMapSharedPtr
Definition: AssemblyMap.h:53
std::shared_ptr< BottomUpSubStructuredGraph > BottomUpSubStructuredGraphSharedPtr
std::shared_ptr< AssemblyMap > AssemblyMapSharedPtr
Definition: AssemblyMap.h:51
std::shared_ptr< PatchMap > PatchMapSharedPtr
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
double NekDouble