Nektar++
Loading...
Searching...
No Matches
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
46{
47// Forward declarations
48class AssemblyMap;
49class ExpList;
50typedef std::shared_ptr<AssemblyMap> AssemblyMapSharedPtr;
52
53/// Base class for constructing local to global mapping of degrees of
54/// freedom.
56{
57public:
58 /// Default constructor.
60 /// Constructor with a communicator
64 const std::string variable = "DefaultVar");
65
66 /// Constructor for next level in multi-level static condensation.
68 AssemblyMap *oldLevelMap,
69 const BottomUpSubStructuredGraphSharedPtr &multiLevelGraph);
70 /// Destructor.
72
73 /// Retrieves the communicator
75
76 /// Retrieves the variable string
78
79 /// Retrieves the hash of this map
80 MULTI_REGIONS_EXPORT size_t GetHash() const;
81
82 MULTI_REGIONS_EXPORT int GetLocalToGlobalMap(const int i) const;
83
84 MULTI_REGIONS_EXPORT int GetGlobalToUniversalMap(const int i) const;
85
87
89
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
139
142
143 /// Returns true if using a modal expansion requiring a change of
144 /// sign of some modes.
146
147 /// Retrieve the sign change of a given local boundary mode.
149 /// Retrieve the sign change for all local boundary modes.
151 const;
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 global index of the rotational periodic boundary trace
173 /// giving the index on the rotational periodic boundary condition
177
178 /// Returns the number of global Dirichlet boundary coefficients.
180 /// Returns the number of local Dirichlet boundary coefficients.
182 /// Returns the total number of global boundary coefficients.
184 /// Returns the total number of local boundary coefficients.
186 /// Returns the total number of local coefficients.
188 /// Returns the total number of global coefficients.
190 /// Retrieves if the system is singular (true) or not (false)
192
193 ///
195 const NekVector<NekDouble> &global, NekVector<NekDouble> &loc,
196 int offset) const;
197
199 const NekVector<NekDouble> &global, NekVector<NekDouble> &loc) const;
200
203 int offset) const;
204
206 const Array<OneD, const NekDouble> &global,
207 Array<OneD, NekDouble> &loc) const;
208
211 int offset, bool UseComm = true) const;
212
215 bool UseComm = true) const;
216
218 const Array<OneD, const NekDouble> &local,
219 Array<OneD, NekDouble> &locbnd) const;
220
222 const Array<OneD, const NekDouble> &local,
223 Array<OneD, NekDouble> &locint) const;
224
226 const Array<OneD, const NekDouble> &locbnd,
227 Array<OneD, NekDouble> &local) const;
228
230 const Array<OneD, const NekDouble> &locbnd,
231 Array<OneD, NekDouble> &local) const;
232
234 NekVector<NekDouble> &global,
235 int offset) const;
236
238 NekVector<NekDouble> &global) const;
239
242 int offset) const;
243
246 Array<OneD, NekDouble> &global) const;
247
249 Array<OneD, NekDouble> &pGlobal) const;
250
252 NekVector<NekDouble> &pGlobal) const;
253
255 Array<OneD, NekDouble> &pGlobal, int offset) const;
256
258
260
262
264
266
268
270
272
273 MULTI_REGIONS_EXPORT void PrintStats(std::ostream &out,
274 std::string variable,
275 bool printHeader = true) const;
276
278
279 MULTI_REGIONS_EXPORT std::shared_ptr<AssemblyMap> LinearSpaceMap(
280 const ExpList &locexp, GlobalSysSolnType solnType);
281
282 /// Returns the bandwidth of the boundary system.
284 /// Returns the level of static condensation for this map.
286 /// Returns the number of patches in this static condensation level.
288 /// Returns the number of local boundary coefficients in each patch.
291 /// Returns the number of local interior coefficients in each patch.
294 /// Returns the local to global mapping for the next level in the
295 /// multi-level static condensation.
298
300 AssemblyMapSharedPtr pNextLevelLocalToGlobalMap);
301
302 /// Returns the patch map from the previous level
303 /// of the multi-level static condensation.
305 void) const;
306
307 /// Returns true if this is the last level in the multi-level
308 /// static condensation.
310 /// Returns the method of solving global systems.
312 MULTI_REGIONS_EXPORT std::string GetPreconType() const;
315 MULTI_REGIONS_EXPORT std::string GetLinSysIterSolver() const;
316
321
324 Array<OneD, NekDouble> &global) const;
325
327 const Array<OneD, const NekDouble> &global,
328 Array<OneD, NekDouble> &loc) const;
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 rotational periodic bnd cond trace number to global trace
400 /// number
402 /// Integer map of process coeffs to universal space
404 /// Integer map of unique process coeffs to universal space (signed)
406
407 /// The solution type of the global system
409 /// The bandwith of the global bnd system
411
412 /// Type type of preconditioner to use in iterative solver.
413 std::string m_preconType;
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
469 virtual void v_LocalToGlobal(const Array<OneD, const NekDouble> &loc,
471 bool useComm) const;
472
473 virtual void v_GlobalToLocal(const Array<OneD, const NekDouble> &global,
474 Array<OneD, NekDouble> &loc) const;
475
476 virtual void v_GlobalToLocal(const NekVector<NekDouble> &global,
477 NekVector<NekDouble> &loc) const;
478
479 virtual void v_Assemble(const Array<OneD, const NekDouble> &loc,
480 Array<OneD, NekDouble> &global) const;
481
482 virtual void v_Assemble(const NekVector<NekDouble> &loc,
483 NekVector<NekDouble> &global) const;
484
485 virtual void v_UniversalAssemble(Array<OneD, NekDouble> &pGlobal) const;
486
487 virtual void v_UniversalAssemble(Array<OneD, NekDouble> &pGlobal,
488 int offset) const;
489
490 virtual int v_GetFullSystemBandWidth() const;
491
492 virtual int v_GetNumNonDirVertexModes() const;
493
494 virtual int v_GetNumNonDirEdgeModes() const;
495
496 virtual int v_GetNumNonDirFaceModes() const;
497
498 virtual int v_GetNumDirEdges() const;
499
500 virtual int v_GetNumDirFaces() const;
501
502 virtual int v_GetNumNonDirEdges() const;
503
504 virtual int v_GetNumNonDirFaces() const;
505
507
508 /// Generate a linear space mapping from existing mapping
509 virtual std::shared_ptr<AssemblyMap> v_LinearSpaceMap(
510 const ExpList &locexp, GlobalSysSolnType solnType);
511
512private:
513 /// Mapping information for previous level in MultiLevel Solver
515};
516
517} // namespace Nektar::MultiRegions
518
519#endif // MULTIREGIONS_ASSEMBLY_MAP_H
#define MULTI_REGIONS_EXPORT
Base class for constructing local to global mapping of degrees of freedom.
Definition AssemblyMap.h:56
int m_lowestStaticCondLevel
Lowest static condensation level.
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.
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.
void SetNextLevelLocalToGlobalMap(AssemblyMapSharedPtr pNextLevelLocalToGlobalMap)
Array< OneD, int > m_bndCondCoeffsToLocalTraceMap
Integer map of bnd cond coeff to local trace coeff.
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.
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.
Array< OneD, int > m_perbndCondIDToGlobalTraceID
Integer map of rotational periodic bnd cond trace number to global trace number.
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.
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.
Array< OneD, int > m_localToLocalIntMap
Integer map of local boundary coeffs to local interior system numbering.
virtual const Array< OneD, const int > & v_GetGlobalToUniversalMapUnique()
NekDouble m_iterativeTolerance
Tolerance for iterative solver.
int m_numGlobalCoeffs
Total number of global coefficients.
std::string m_linSysIterSolver
Iterative solver: Conjugate Gradient, GMRES.
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.
void PatchAssemble(const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global) const
int m_successiveRHS
sucessive RHS for iterative solver
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.
void UniversalAbsMaxBnd(Array< OneD, NekDouble > &bndvals)
Array< OneD, unsigned int > m_numLocalBndCoeffsPerPatch
The number of bnd dofs per patch.
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.
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.
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.
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.
int m_bndSystemBandWidth
The bandwith of the global bnd system.
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.
bool GetSingularSystem() const
Retrieves if the system is singular (true) or not (false)
std::string m_variable
Variable string identifier.
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.
LibUtilities::CommSharedPtr GetComm()
Retrieves the communicator.
Array< OneD, unsigned int > m_numLocalIntCoeffsPerPatch
The number of int dofs per patch.
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.
Array< OneD, int > m_localToGlobalBndMap
Integer map of local coeffs to global Boundary Dofs.
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.
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.
Array< OneD, int > m_globalToUniversalBndMapUnique
Integer map of unique process coeffs to universal space (signed)
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.
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.
Array< OneD, NekDouble > m_bndCondCoeffsToLocalCoeffsSign
Integer map of sign of bnd cond coeffs to local coefficients.
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.
const Array< OneD, const int > & GetPerBndCondIDToGlobalTraceID()
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.
int GetNumLocalBndCoeffs() const
Returns the total number of local boundary coefficients.
int m_numGlobalBndCoeffs
Total number of global boundary coefficients.
Array< OneD, int > m_bndCondIDToGlobalTraceID
Integer map of bnd cond trace number to global trace number.
Base class for all multi-elemental spectral/hp expansions.
Definition ExpList.h:99
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition Comm.h:55
static AssemblyMapSharedPtr NullAssemblyMapSharedPtr
Definition AssemblyMap.h:51
std::shared_ptr< BottomUpSubStructuredGraph > BottomUpSubStructuredGraphSharedPtr
std::shared_ptr< AssemblyMap > AssemblyMapSharedPtr
Definition AssemblyMap.h:50
std::shared_ptr< PatchMap > PatchMapSharedPtr