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 
41 #include <vector>
44 
45 
46 namespace Nektar
47 {
48  namespace MultiRegions
49  {
50  // Forward declarations
51  class AssemblyMap;
52  class ExpList;
53  typedef std::shared_ptr<AssemblyMap> AssemblyMapSharedPtr;
55 
56  /// Base class for constructing local to global mapping of degrees of
57  /// freedom.
59  {
60  public:
61  /// Default constructor.
63  /// Constructor with a communicator
66  const std::string variable = "DefaultVar");
67 
68  /// Constructor for next level in multi-level static condensation.
70  const BottomUpSubStructuredGraphSharedPtr& multiLevelGraph);
71  /// Destructor.
73 
74  /// Retrieves the communicator
76 
77  /// Retrieves the hash of this map
78  MULTI_REGIONS_EXPORT size_t GetHash() const;
79 
80  MULTI_REGIONS_EXPORT int GetLocalToGlobalMap(const int i) const;
81 
82  MULTI_REGIONS_EXPORT int GetGlobalToUniversalMap(const int i) const;
83 
85 
87 
89 
91 
93 
95 
98  Array<OneD, NekDouble>& global,
99  bool useComm = true) const;
100 
102  const NekVector<NekDouble>& loc,
103  NekVector< NekDouble>& global,
104  bool useComm = true) const;
105 
107  const Array<OneD, const NekDouble>& global,
108  Array<OneD, NekDouble>& loc) const;
109 
111  const NekVector<NekDouble>& global,
112  NekVector< NekDouble>& loc) const;
113 
116  Array<OneD, NekDouble> &global) const;
117 
119  const NekVector<NekDouble>& loc,
120  NekVector< NekDouble>& global) const;
121 
123  Array<OneD, NekDouble>& pGlobal) const;
124 
126  NekVector< NekDouble>& pGlobal) const;
127 
129  Array<OneD, NekDouble>& pGlobal,
130  int offset) const;
131 
133  Array<OneD, NekDouble> &bndvals);
134 
135  /// Retrieve the global index of a given local boundary mode.
136  MULTI_REGIONS_EXPORT int GetLocalToGlobalBndMap(const int i) const;
137  /// Retrieve the global indices of the local boundary modes.
139 
141 
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  /// Retrieves the local indices corresponding to the
153  /// boundary expansion modes.
156  /// Returns the modal sign associated with a given
157  /// boundary expansion mode.
160 
161  /// Retrieves the local indices corresponding to the
162  /// boundary expansion modes to global trace
165 
166  /// Returns the global index of the boundary trace giving the
167  /// index on the boundary expansion
171 
172  /// Returns the number of global Dirichlet boundary coefficients.
174  /// Returns the number of local Dirichlet boundary coefficients.
176  /// Returns the total number of global boundary coefficients.
178  /// Returns the total number of local boundary coefficients.
180  /// Returns the total number of local coefficients.
182  /// Returns the total number of global coefficients.
184  /// Retrieves if the system is singular (true) or not (false)
186 
187  ///
189  const NekVector<NekDouble>& global,
191  int offset) const;
192 
194  const NekVector<NekDouble>& global,
195  NekVector<NekDouble>& loc) const;
196 
198  const Array<OneD, const NekDouble>& global,
200  int offset) const;
201 
203  const Array<OneD, const NekDouble>& global,
204  Array<OneD,NekDouble>& loc) const;
205 
208  Array<OneD,NekDouble>& global,
209  int offset, bool UseComm = true) const;
210 
213  Array<OneD,NekDouble>& global, 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, int offset) const;
233 
235  NekVector<NekDouble>& global) const;
236 
238  Array<OneD, NekDouble>& global, int offset) const;
239 
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,
251  int offset) const;
252 
254 
256 
258 
260 
262 
264 
266 
268 
269  MULTI_REGIONS_EXPORT void PrintStats(std::ostream &out, std::string variable, bool printHeader = true) const;
270 
273 
274  MULTI_REGIONS_EXPORT std::shared_ptr<AssemblyMap> LinearSpaceMap(const ExpList &locexp, GlobalSysSolnType solnType);
275 
276  /// Returns the bandwidth of the boundary system.
278  /// Returns the level of static condensation for this map.
280  /// Returns the number of patches in this static condensation level.
282  /// Returns the number of local boundary coefficients in each patch.
285  /// Returns the number of local interior coefficients in each patch.
288  /// Returns the local to global mapping for the next level in the
289  /// multi-level static condensation.
292 
294 
295  /// Returns the patch map from the previous level
296  /// of the multi-level static condensation.
298  GetPatchMapFromPrevLevel(void) const;
299 
300  /// Returns true if this is the last level in the multi-level
301  /// static condensation.
302  MULTI_REGIONS_EXPORT bool AtLastLevel() const;
303  /// Returns the method of solving global systems.
309  MULTI_REGIONS_EXPORT std::string GetLinSysIterSolver() const;
310 
312  {
314  }
315 
318  Array<OneD, NekDouble>& global) const;
319 
321  const Array<OneD, const NekDouble>& global,
322  Array<OneD, NekDouble>& loc) const;
323 
326  Array<OneD, NekDouble> &global) const;
327  protected:
328  /// Session object
330 
331  /// Communicator
333 
334  /// Hash for map
335  size_t m_hash;
336 
337  /// Number of local boundary coefficients
339  /// Total number of global boundary coefficients
341  /// Number of Local Dirichlet Boundary Coefficients
343  /// Number of Global Dirichlet Boundary Coefficients
345  /// Flag indicating if the system is singular or not
347 
348  /// Total number of local coefficients
349  /** This corresponds to the number of total number of coefficients
350  * - For CG this corresponds to the total of bnd + int DOFs
351  * - For DG this corresponds to the number of bnd DOFs.
352  * This means that #m_numLocalCoeffs = #m_numLocalBndCoeffs
353  * This way, we can consider the trace-system solve as a
354  * statically condensed solve without interior DOFs. This allows
355  * us to use the same global system solver for both cases.
356  */
358 
359  /// Total number of global coefficients
360  /** This corresponds to the number of total number of coefficients
361  * - For CG this corresponds to the total of bnd + int DOFs.
362  * - For DG this corresponds to the number of bnd DOFs.
363  * This means that #m_numGlobalCoeffs = #m_numGlobalBndCoeffs
364  * This way, we can consider the trace-system solve as a
365  * statically condensed solve without interior DOFs. This allows
366  * us to use the same global system solver for both cases.
367  */
369 
370  /// Flag indicating if modes require sign reversal.
372 
373  /// Integer map of local coeffs to global Boundary Dofs
375  /// Integer sign of local boundary coeffs to global space
377  /// Integer map of local boundary coeffs to local boundary system numbering
379  /// Integer map of local boundary coeffs to local interior system numbering
381  /// Integer map of bnd cond coeffs to local coefficients
383  /// Integer map of sign of bnd cond coeffs to local coefficients
385  /// Integer map of bnd cond coeff to local trace coeff
387  /// Integer map of bnd cond trace number to global trace number
389  /// Integer map of process coeffs to universal space
391  /// Integer map of unique process coeffs to universal space (signed)
393 
394  /// The solution type of the global system
396  /// The bandwith of the global bnd system
398 
399  /// Type type of preconditioner to use in iterative solver.
401 
402  /// Maximum iterations for iterative solver
404 
405  /// Tolerance for iterative solver
407 
408  /// sucessive RHS for iterative solver
410 
411  /// Iterative solver: Conjugate Gradient, GMRES
412  std::string m_linSysIterSolver;
413 
416  /// gs gather communication to impose Dirhichlet BCs.
418 
419  /// The level of recursion in the case of multi-level static
420  /// condensation.
422  /// The number of patches (~elements) in the current level
424  /// The number of bnd dofs per patch
426  /// The number of int dofs per patch
428  /// Map from the patches of the previous level to the patches of
429  /// the current level
430 
431  /// The local to global mapping of the next level of recursion
433  /// Lowest static condensation level.
435 
436  /// Calculates the bandwidth of the boundary system.
438 
440  const Array<OneD, const NekDouble>& global,
442 
443  private:
444  /// Mapping information for previous level in MultiLevel Solver
446 
447  virtual int v_GetLocalToGlobalMap(const int i) const;
448 
449  virtual int v_GetGlobalToUniversalMap(const int i) const;
450 
451  virtual int v_GetGlobalToUniversalMapUnique(const int i) const;
452 
454 
456 
458 
459  virtual NekDouble v_GetLocalToGlobalSign(const int i) const;
460 
461  virtual const Array<OneD, NekDouble>& v_GetLocalToGlobalSign() const;
462 
463  virtual void v_LocalToGlobal(
465  Array<OneD, NekDouble>& global,
466  bool useComm) const;
467 
468  virtual void v_LocalToGlobal(
469  const NekVector<NekDouble>& loc,
470  NekVector< NekDouble>& global,
471  bool useComm) const;
472 
473  virtual void v_GlobalToLocal(
474  const Array<OneD, const NekDouble>& global,
475  Array<OneD, NekDouble>& loc) const;
476 
477  virtual void v_GlobalToLocal(
478  const NekVector<NekDouble>& global,
479  NekVector< NekDouble>& loc) const;
480 
481  virtual void v_Assemble(
483  Array<OneD, NekDouble> &global) const;
484 
485  virtual void v_Assemble(
486  const NekVector<NekDouble>& loc,
487  NekVector< NekDouble>& global) const;
488 
489  virtual void v_UniversalAssemble(
490  Array<OneD, NekDouble>& pGlobal) const;
491 
492  virtual void v_UniversalAssemble(
493  NekVector< NekDouble>& pGlobal) const;
494 
495  virtual void v_UniversalAssemble(
496  Array<OneD, NekDouble>& pGlobal,
497  int offset) const;
498 
499  virtual int v_GetFullSystemBandWidth() const;
500 
501  virtual int v_GetNumNonDirVertexModes() const;
502 
503  virtual int v_GetNumNonDirEdgeModes() const;
504 
505  virtual int v_GetNumNonDirFaceModes() const;
506 
507  virtual int v_GetNumDirEdges() const;
508 
509  virtual int v_GetNumDirFaces() const;
510 
511  virtual int v_GetNumNonDirEdges() const;
512 
513  virtual int v_GetNumNonDirFaces() const;
514 
515  virtual const Array<OneD, const int>&
517 
518  /// Generate a linear space mapping from existing mapping
519  virtual std::shared_ptr<AssemblyMap> v_LinearSpaceMap(
520  const ExpList &locexp, GlobalSysSolnType solnType);
521  };
522 
523 
524  } // end of namespace
525 } // end of namespace
526 
527 #endif //MULTIREGIONS_ASSEMBLY_MAP_H
528 
529 
#define MULTI_REGIONS_EXPORT
Base class for constructing local to global mapping of degrees of freedom.
Definition: AssemblyMap.h:59
int m_lowestStaticCondLevel
Lowest static condensation level.
Definition: AssemblyMap.h:434
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:395
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:400
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:357
void SetNextLevelLocalToGlobalMap(AssemblyMapSharedPtr pNextLevelLocalToGlobalMap)
const Array< OneD, const unsigned int > & GetNumLocalBndCoeffsPerPatch()
Returns the number of local boundary coefficients in each patch.
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
Array< OneD, int > m_globalToUniversalBndMap
Integer map of process coeffs to universal space.
Definition: AssemblyMap.h:390
bool m_signChange
Flag indicating if modes require sign reversal.
Definition: AssemblyMap.h:371
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
Array< OneD, int > m_bndCondCoeffsToLocalCoeffsMap
Integer map of bnd cond coeffs to local coefficients.
Definition: AssemblyMap.h:382
PatchMapSharedPtr m_patchMapFromPrevLevel
Mapping information for previous level in MultiLevel Solver.
Definition: AssemblyMap.h:445
int m_maxIterations
Maximum iterations for iterative solver.
Definition: AssemblyMap.h:403
virtual const Array< OneD, const int > & v_GetGlobalToUniversalMapUnique()
NekDouble m_iterativeTolerance
Tolerance for iterative solver.
Definition: AssemblyMap.h:406
NekDouble GetIterativeTolerance() const
int m_numGlobalCoeffs
Total number of global coefficients.
Definition: AssemblyMap.h:368
std::string m_linSysIterSolver
Iterative solver: Conjugate Gradient, GMRES.
Definition: AssemblyMap.h:412
bool AtLastLevel() const
Returns true if this is the last level in the multi-level static condensation.
Array< OneD, int > m_localToGlobalBndMap
Integer map of local coeffs to global Boundary Dofs.
Definition: AssemblyMap.h:374
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()
void PatchAssemble(const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global) const
int m_successiveRHS
sucessive RHS for iterative solver
Definition: AssemblyMap.h:409
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
void UniversalAbsMaxBnd(Array< OneD, NekDouble > &bndvals)
Array< OneD, unsigned int > m_numLocalBndCoeffsPerPatch
The number of bnd dofs per patch.
Definition: AssemblyMap.h:425
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:329
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:338
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:432
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:421
Array< OneD, NekDouble > m_localToGlobalBndSign
Integer sign of local boundary coeffs to global space.
Definition: AssemblyMap.h:376
int m_bndSystemBandWidth
The bandwith of the global bnd system.
Definition: AssemblyMap.h:397
Array< OneD, int > m_localToLocalIntMap
Integer map of local boundary coeffs to local interior system numbering.
Definition: AssemblyMap.h:380
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:342
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:344
LibUtilities::CommSharedPtr GetComm()
Retrieves the communicator.
Array< OneD, unsigned int > m_numLocalIntCoeffsPerPatch
The number of int dofs per patch.
Definition: AssemblyMap.h:427
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:346
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:417
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:80
Array< OneD, NekDouble > m_bndCondCoeffsToLocalCoeffsSign
Integer map of sign of bnd cond coeffs to local coefficients.
Definition: AssemblyMap.h:384
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:332
int GetNumGlobalBndCoeffs() const
Returns the total number of global boundary coefficients.
Array< OneD, int > m_localToLocalBndMap
Integer map of local boundary coeffs to local boundary system numbering.
Definition: AssemblyMap.h:378
Array< OneD, int > m_globalToUniversalBndMapUnique
Integer map of unique process coeffs to universal space (signed)
Definition: AssemblyMap.h:392
const Array< OneD, const int > & GetGlobalToUniversalBndMap()
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:423
void LocalToLocalInt(const Array< OneD, const NekDouble > &local, Array< OneD, NekDouble > &locint) const
Array< OneD, int > m_bndCondIDToGlobalTraceID
Integer map of bnd cond trace number to global trace number.
Definition: AssemblyMap.h:388
virtual int v_GetNumNonDirFaces() const
bool GetSignChange()
Returns true if using a modal expansion requiring a change of sign of some modes.
Array< OneD, int > m_bndCondCoeffsToLocalTraceMap
Integer map of bnd cond coeff to local trace coeff.
Definition: AssemblyMap.h:386
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:340
Base class for all multi-elemental spectral/hp expansions.
Definition: ExpList.h:107
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:54
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:1
double NekDouble