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
92
95
97
99 const;
100
103 bool useComm = true) const;
104
106 NekVector<NekDouble> &global,
107 bool useComm = true) const;
108
110 const Array<OneD, const NekDouble> &global,
112
115
117 Array<OneD, NekDouble> &global) const;
118
120 NekVector<NekDouble> &global) const;
121
123 Array<OneD, NekDouble> &pGlobal) const;
124
126 NekVector<NekDouble> &pGlobal) const;
127
129 int offset) const;
130
132 Array<OneD, NekDouble> &bndvals);
133
134 /// Retrieve the global index of a given local boundary mode.
135 MULTI_REGIONS_EXPORT int GetLocalToGlobalBndMap(const int i) const;
136 /// Retrieve the global indices of the local boundary modes.
138
141
144
145 /// Returns true if using a modal expansion requiring a change of
146 /// sign of some modes.
148
149 /// Retrieve the sign change of a given local boundary mode.
151 /// Retrieve the sign change for all local boundary modes.
153 const;
154 /// Retrieves the local indices corresponding to the
155 /// boundary expansion modes.
158 /// Returns the modal sign associated with a given
159 /// boundary expansion mode.
162
163 /// Retrieves the local indices corresponding to the
164 /// boundary expansion modes to global trace
167
168 /// Returns the global index of the boundary trace giving the
169 /// index on the boundary expansion
173
174 /// Returns the number of global Dirichlet boundary coefficients.
176 /// Returns the number of local Dirichlet boundary coefficients.
178 /// Returns the total number of global boundary coefficients.
180 /// Returns the total number of local boundary coefficients.
182 /// Returns the total number of local coefficients.
184 /// Returns the total number of global coefficients.
186 /// Retrieves if the system is singular (true) or not (false)
188
189 ///
192 int offset) const;
193
195 const NekVector<NekDouble> &global, NekVector<NekDouble> &loc) const;
196
199 int offset) const;
200
202 const Array<OneD, const NekDouble> &global,
204
207 int offset, bool UseComm = true) const;
208
211 bool UseComm = true) const;
212
214 const Array<OneD, const NekDouble> &local,
215 Array<OneD, NekDouble> &locbnd) const;
216
218 const Array<OneD, const NekDouble> &local,
219 Array<OneD, NekDouble> &locint) const;
220
222 const Array<OneD, const NekDouble> &locbnd,
223 Array<OneD, NekDouble> &local) const;
224
226 const Array<OneD, const NekDouble> &locbnd,
227 Array<OneD, NekDouble> &local) const;
228
230 NekVector<NekDouble> &global,
231 int offset) const;
232
234 NekVector<NekDouble> &global) const;
235
238 int offset) const;
239
242 Array<OneD, NekDouble> &global) const;
243
245 Array<OneD, NekDouble> &pGlobal) const;
246
248 NekVector<NekDouble> &pGlobal) const;
249
251 Array<OneD, NekDouble> &pGlobal, int offset) const;
252
254
256
258
260
262
264
266
268
269 MULTI_REGIONS_EXPORT void PrintStats(std::ostream &out,
270 std::string variable,
271 bool printHeader = true) const;
272
274
275 MULTI_REGIONS_EXPORT std::shared_ptr<AssemblyMap> LinearSpaceMap(
276 const ExpList &locexp, GlobalSysSolnType solnType);
277
278 /// Returns the bandwidth of the boundary system.
280 /// Returns the level of static condensation for this map.
282 /// Returns the number of patches in this static condensation level.
284 /// Returns the number of local boundary coefficients in each patch.
287 /// Returns the number of local interior coefficients in each patch.
290 /// Returns the local to global mapping for the next level in the
291 /// multi-level static condensation.
294
296 AssemblyMapSharedPtr pNextLevelLocalToGlobalMap);
297
298 /// Returns the patch map from the previous level
299 /// of the multi-level static condensation.
301 void) const;
302
303 /// Returns true if this is the last level in the multi-level
304 /// static condensation.
306 /// Returns the method of solving global systems.
308 MULTI_REGIONS_EXPORT std::string GetPreconType() const;
312 MULTI_REGIONS_EXPORT std::string GetLinSysIterSolver() const;
313
315 {
317 }
318
321 Array<OneD, NekDouble> &global) const;
322
324 const Array<OneD, const NekDouble> &global,
326
329 Array<OneD, NekDouble> &global) const;
330
331protected:
332 /// Session object
334
335 /// Communicator
337
338 /// Variable string identifier
339 std::string m_variable;
340
341 /// Hash for map
342 size_t m_hash;
343
344 /// Number of local boundary coefficients
346 /// Total number of global boundary coefficients
348 /// Number of Local Dirichlet Boundary Coefficients
350 /// Number of Global Dirichlet Boundary Coefficients
352 /// Flag indicating if the system is singular or not
354
355 /// Total number of local coefficients
356 /** This corresponds to the number of total number of coefficients
357 * - For CG this corresponds to the total of bnd + int DOFs
358 * - For DG this corresponds to the number of bnd DOFs.
359 * This means that #m_numLocalCoeffs = #m_numLocalBndCoeffs
360 * This way, we can consider the trace-system solve as a
361 * statically condensed solve without interior DOFs. This allows
362 * us to use the same global system solver for both cases.
363 */
365
366 /// Total number of global coefficients
367 /** This corresponds to the number of total number of coefficients
368 * - For CG this corresponds to the total of bnd + int DOFs.
369 * - For DG this corresponds to the number of bnd DOFs.
370 * This means that #m_numGlobalCoeffs = #m_numGlobalBndCoeffs
371 * This way, we can consider the trace-system solve as a
372 * statically condensed solve without interior DOFs. This allows
373 * us to use the same global system solver for both cases.
374 */
376
377 /// Flag indicating if modes require sign reversal.
379
380 /// Integer map of local coeffs to global Boundary Dofs
382 /// Integer sign of local boundary coeffs to global space
384 /// Integer map of local boundary coeffs to local boundary system numbering
386 /// Integer map of local boundary coeffs to local interior system numbering
388 /// Integer map of bnd cond coeffs to local coefficients
390 /// Integer map of sign of bnd cond coeffs to local coefficients
392 /// Integer map of bnd cond coeff to local trace coeff
394 /// Integer map of bnd cond trace number to global trace number
396 /// Integer map of process coeffs to universal space
398 /// Integer map of unique process coeffs to universal space (signed)
400
401 /// The solution type of the global system
403 /// The bandwith of the global bnd system
405
406 /// Type type of preconditioner to use in iterative solver.
407 std::string m_preconType;
408
409 /// Tolerance for iterative solver
412
413 /// sucessive RHS for iterative solver
415
416 /// Iterative solver: Conjugate Gradient, GMRES
418
421 /// gs gather communication to impose Dirhichlet BCs.
423
424 /// The level of recursion in the case of multi-level static
425 /// condensation.
427 /// The number of patches (~elements) in the current level
429 /// The number of bnd dofs per patch
431 /// The number of int dofs per patch
433 /// Map from the patches of the previous level to the patches of
434 /// the current level
435
436 /// The local to global mapping of the next level of recursion
438 /// Lowest static condensation level.
440
441 /// Calculates the bandwidth of the boundary system.
443
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
465 bool useComm) const;
466
467 virtual void v_LocalToGlobal(const NekVector<NekDouble> &loc,
468 NekVector<NekDouble> &global,
469 bool useComm) const;
470
471 virtual void v_GlobalToLocal(const Array<OneD, const NekDouble> &global,
473
474 virtual void v_GlobalToLocal(const NekVector<NekDouble> &global,
476
477 virtual void v_Assemble(const Array<OneD, const NekDouble> &loc,
478 Array<OneD, NekDouble> &global) const;
479
480 virtual void v_Assemble(const NekVector<NekDouble> &loc,
481 NekVector<NekDouble> &global) const;
482
483 virtual void v_UniversalAssemble(Array<OneD, NekDouble> &pGlobal) const;
484
485 virtual void v_UniversalAssemble(NekVector<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.
Definition: AssemblyMap.h:439
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:402
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:364
void SetNextLevelLocalToGlobalMap(AssemblyMapSharedPtr pNextLevelLocalToGlobalMap)
Array< OneD, int > m_bndCondCoeffsToLocalTraceMap
Integer map of bnd cond coeff to local trace coeff.
Definition: AssemblyMap.h:393
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:389
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:378
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:514
Array< OneD, int > m_localToLocalIntMap
Integer map of local boundary coeffs to local interior system numbering.
Definition: AssemblyMap.h:387
virtual const Array< OneD, const int > & v_GetGlobalToUniversalMapUnique()
NekDouble m_iterativeTolerance
Tolerance for iterative solver.
Definition: AssemblyMap.h:410
NekDouble GetIterativeTolerance() const
int m_numGlobalCoeffs
Total number of global coefficients.
Definition: AssemblyMap.h:375
std::string m_linSysIterSolver
Iterative solver: Conjugate Gradient, GMRES.
Definition: AssemblyMap.h:417
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:397
void PatchAssemble(const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global) const
int m_successiveRHS
sucessive RHS for iterative solver
Definition: AssemblyMap.h:414
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:383
void UniversalAbsMaxBnd(Array< OneD, NekDouble > &bndvals)
Array< OneD, unsigned int > m_numLocalBndCoeffsPerPatch
The number of bnd dofs per patch.
Definition: AssemblyMap.h:430
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:333
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:345
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:437
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:426
int m_bndSystemBandWidth
The bandwith of the global bnd system.
Definition: AssemblyMap.h:404
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:349
bool GetSingularSystem() const
Retrieves if the system is singular (true) or not (false)
std::string m_variable
Variable string identifier.
Definition: AssemblyMap.h:339
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:351
LibUtilities::CommSharedPtr GetComm()
Retrieves the communicator.
Array< OneD, unsigned int > m_numLocalIntCoeffsPerPatch
The number of int dofs per patch.
Definition: AssemblyMap.h:432
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:353
Array< OneD, int > m_localToGlobalBndMap
Integer map of local coeffs to global Boundary Dofs.
Definition: AssemblyMap.h:381
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:422
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:399
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:336
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:385
Array< OneD, NekDouble > m_bndCondCoeffsToLocalCoeffsSign
Integer map of sign of bnd cond coeffs to local coefficients.
Definition: AssemblyMap.h:391
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:428
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:407
int GetNumLocalBndCoeffs() const
Returns the total number of local boundary coefficients.
int m_numGlobalBndCoeffs
Total number of global boundary coefficients.
Definition: AssemblyMap.h:347
Array< OneD, int > m_bndCondIDToGlobalTraceID
Integer map of bnd cond trace number to global trace number.
Definition: AssemblyMap.h:395
Base class for all multi-elemental spectral/hp expansions.
Definition: ExpList.h:98
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