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
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,
111
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 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 ///
190 int offset) const;
191
193 const NekVector<NekDouble> &global, NekVector<NekDouble> &loc) const;
194
197 int offset) const;
198
200 const Array<OneD, const NekDouble> &global,
202
205 int offset, bool UseComm = true) const;
206
209 bool UseComm = true) const;
210
212 const Array<OneD, const NekDouble> &local,
213 Array<OneD, NekDouble> &locbnd) const;
214
216 const Array<OneD, const NekDouble> &local,
217 Array<OneD, NekDouble> &locint) const;
218
220 const Array<OneD, const NekDouble> &locbnd,
221 Array<OneD, NekDouble> &local) const;
222
224 const Array<OneD, const NekDouble> &locbnd,
225 Array<OneD, NekDouble> &local) const;
226
228 NekVector<NekDouble> &global,
229 int offset) const;
230
232 NekVector<NekDouble> &global) const;
233
236 int offset) const;
237
240 Array<OneD, NekDouble> &global) const;
241
243 Array<OneD, NekDouble> &pGlobal) const;
244
246 NekVector<NekDouble> &pGlobal) const;
247
249 Array<OneD, NekDouble> &pGlobal, int offset) const;
250
252
254
256
258
260
262
264
266
267 MULTI_REGIONS_EXPORT void PrintStats(std::ostream &out,
268 std::string variable,
269 bool printHeader = true) const;
270
272
273 MULTI_REGIONS_EXPORT std::shared_ptr<AssemblyMap> LinearSpaceMap(
274 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 AssemblyMapSharedPtr pNextLevelLocalToGlobalMap);
295
296 /// Returns the patch map from the previous level
297 /// of the multi-level static condensation.
299 void) const;
300
301 /// Returns true if this is the last level in the multi-level
302 /// static condensation.
304 /// Returns the method of solving global systems.
306 MULTI_REGIONS_EXPORT std::string GetPreconType() const;
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,
323
326 Array<OneD, NekDouble> &global) const;
327
328protected:
329 /// Session object
331
332 /// Communicator
334
335 /// Variable string identifier
336 std::string m_variable;
337
338 /// Hash for map
339 size_t m_hash;
340
341 /// Number of local boundary coefficients
343 /// Total number of global boundary coefficients
345 /// Number of Local Dirichlet Boundary Coefficients
347 /// Number of Global Dirichlet Boundary Coefficients
349 /// Flag indicating if the system is singular or not
351
352 /// Total number of local coefficients
353 /** This corresponds to the number of total number of coefficients
354 * - For CG this corresponds to the total of bnd + int DOFs
355 * - For DG this corresponds to the number of bnd DOFs.
356 * This means that #m_numLocalCoeffs = #m_numLocalBndCoeffs
357 * This way, we can consider the trace-system solve as a
358 * statically condensed solve without interior DOFs. This allows
359 * us to use the same global system solver for both cases.
360 */
362
363 /// Total number of global coefficients
364 /** This corresponds to the number of total number of coefficients
365 * - For CG this corresponds to the total of bnd + int DOFs.
366 * - For DG this corresponds to the number of bnd DOFs.
367 * This means that #m_numGlobalCoeffs = #m_numGlobalBndCoeffs
368 * This way, we can consider the trace-system solve as a
369 * statically condensed solve without interior DOFs. This allows
370 * us to use the same global system solver for both cases.
371 */
373
374 /// Flag indicating if modes require sign reversal.
376
377 /// Integer map of local coeffs to global Boundary Dofs
379 /// Integer sign of local boundary coeffs to global space
381 /// Integer map of local boundary coeffs to local boundary system numbering
383 /// Integer map of local boundary coeffs to local interior system numbering
385 /// Integer map of bnd cond coeffs to local coefficients
387 /// Integer map of sign of bnd cond coeffs to local coefficients
389 /// Integer map of bnd cond coeff to local trace coeff
391 /// Integer map of bnd cond trace number to global trace number
393 /// Integer map of process coeffs to universal space
395 /// Integer map of unique process coeffs to universal space (signed)
397
398 /// The solution type of the global system
400 /// The bandwith of the global bnd system
402
403 /// Type type of preconditioner to use in iterative solver.
404 std::string m_preconType;
405
406 /// Tolerance for iterative solver
409
410 /// sucessive RHS for iterative solver
412
413 /// Iterative solver: Conjugate Gradient, GMRES
415
418 /// gs gather communication to impose Dirhichlet BCs.
420
421 /// The level of recursion in the case of multi-level static
422 /// condensation.
424 /// The number of patches (~elements) in the current level
426 /// The number of bnd dofs per patch
428 /// The number of int dofs per patch
430 /// Map from the patches of the previous level to the patches of
431 /// the current level
432
433 /// The local to global mapping of the next level of recursion
435 /// Lowest static condensation level.
437
438 /// Calculates the bandwidth of the boundary system.
440
443
444 virtual int v_GetLocalToGlobalMap(const int i) const;
445
446 virtual int v_GetGlobalToUniversalMap(const int i) const;
447
448 virtual int v_GetGlobalToUniversalMapUnique(const int i) const;
449
451
453
455
456 virtual NekDouble v_GetLocalToGlobalSign(const int i) const;
457
458 virtual const Array<OneD, NekDouble> &v_GetLocalToGlobalSign() const;
459
462 bool useComm) const;
463
464 virtual void v_LocalToGlobal(const NekVector<NekDouble> &loc,
465 NekVector<NekDouble> &global,
466 bool useComm) const;
467
468 virtual void v_GlobalToLocal(const Array<OneD, const NekDouble> &global,
470
471 virtual void v_GlobalToLocal(const NekVector<NekDouble> &global,
473
474 virtual void v_Assemble(const Array<OneD, const NekDouble> &loc,
475 Array<OneD, NekDouble> &global) const;
476
477 virtual void v_Assemble(const NekVector<NekDouble> &loc,
478 NekVector<NekDouble> &global) const;
479
480 virtual void v_UniversalAssemble(Array<OneD, NekDouble> &pGlobal) const;
481
482 virtual void v_UniversalAssemble(NekVector<NekDouble> &pGlobal) const;
483
484 virtual void v_UniversalAssemble(Array<OneD, NekDouble> &pGlobal,
485 int offset) const;
486
487 virtual int v_GetFullSystemBandWidth() const;
488
489 virtual int v_GetNumNonDirVertexModes() const;
490
491 virtual int v_GetNumNonDirEdgeModes() const;
492
493 virtual int v_GetNumNonDirFaceModes() const;
494
495 virtual int v_GetNumDirEdges() const;
496
497 virtual int v_GetNumDirFaces() const;
498
499 virtual int v_GetNumNonDirEdges() const;
500
501 virtual int v_GetNumNonDirFaces() const;
502
504
505 /// Generate a linear space mapping from existing mapping
506 virtual std::shared_ptr<AssemblyMap> v_LinearSpaceMap(
507 const ExpList &locexp, GlobalSysSolnType solnType);
508
509private:
510 /// Mapping information for previous level in MultiLevel Solver
512};
513
514} // namespace Nektar::MultiRegions
515
516#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.
Definition: AssemblyMap.h:436
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:399
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:361
void SetNextLevelLocalToGlobalMap(AssemblyMapSharedPtr pNextLevelLocalToGlobalMap)
Array< OneD, int > m_bndCondCoeffsToLocalTraceMap
Integer map of bnd cond coeff to local trace coeff.
Definition: AssemblyMap.h:390
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:386
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:375
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:511
Array< OneD, int > m_localToLocalIntMap
Integer map of local boundary coeffs to local interior system numbering.
Definition: AssemblyMap.h:384
virtual const Array< OneD, const int > & v_GetGlobalToUniversalMapUnique()
NekDouble m_iterativeTolerance
Tolerance for iterative solver.
Definition: AssemblyMap.h:407
int m_numGlobalCoeffs
Total number of global coefficients.
Definition: AssemblyMap.h:372
std::string m_linSysIterSolver
Iterative solver: Conjugate Gradient, GMRES.
Definition: AssemblyMap.h:414
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:394
void PatchAssemble(const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global) const
int m_successiveRHS
sucessive RHS for iterative solver
Definition: AssemblyMap.h:411
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:380
void UniversalAbsMaxBnd(Array< OneD, NekDouble > &bndvals)
Array< OneD, unsigned int > m_numLocalBndCoeffsPerPatch
The number of bnd dofs per patch.
Definition: AssemblyMap.h:427
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:330
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:342
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:434
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:423
int m_bndSystemBandWidth
The bandwith of the global bnd system.
Definition: AssemblyMap.h:401
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:346
bool GetSingularSystem() const
Retrieves if the system is singular (true) or not (false)
std::string m_variable
Variable string identifier.
Definition: AssemblyMap.h:336
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:348
LibUtilities::CommSharedPtr GetComm()
Retrieves the communicator.
Array< OneD, unsigned int > m_numLocalIntCoeffsPerPatch
The number of int dofs per patch.
Definition: AssemblyMap.h:429
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:350
Array< OneD, int > m_localToGlobalBndMap
Integer map of local coeffs to global Boundary Dofs.
Definition: AssemblyMap.h:378
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:419
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:79
Array< OneD, int > m_globalToUniversalBndMapUnique
Integer map of unique process coeffs to universal space (signed)
Definition: AssemblyMap.h:396
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:333
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:382
Array< OneD, NekDouble > m_bndCondCoeffsToLocalCoeffsSign
Integer map of sign of bnd cond coeffs to local coefficients.
Definition: AssemblyMap.h:388
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:425
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:404
int GetNumLocalBndCoeffs() const
Returns the total number of local boundary coefficients.
int m_numGlobalBndCoeffs
Total number of global boundary coefficients.
Definition: AssemblyMap.h:344
Array< OneD, int > m_bndCondIDToGlobalTraceID
Integer map of bnd cond trace number to global trace number.
Definition: AssemblyMap.h:392
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
double NekDouble