Nektar++
Public Member Functions | Protected Member Functions | Protected Attributes | Private Attributes | List of all members
Nektar::MultiRegions::AssemblyMap Class Reference

Base class for constructing local to global mapping of degrees of freedom. More...

#include <AssemblyMap.h>

Inheritance diagram for Nektar::MultiRegions::AssemblyMap:
[legend]

Public Member Functions

 AssemblyMap ()
 Default constructor. More...
 
 AssemblyMap (const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &comm, const std::string variable="DefaultVar")
 Constructor with a communicator. More...
 
 AssemblyMap (AssemblyMap *oldLevelMap, const BottomUpSubStructuredGraphSharedPtr &multiLevelGraph)
 Constructor for next level in multi-level static condensation. More...
 
virtual ~AssemblyMap ()
 Destructor. More...
 
LibUtilities::CommSharedPtr GetComm ()
 Retrieves the communicator. More...
 
std::string GetVariable ()
 Retrieves the variable string. More...
 
size_t GetHash () const
 Retrieves the hash of this map. More...
 
int GetLocalToGlobalMap (const int i) const
 
int GetGlobalToUniversalMap (const int i) const
 
int GetGlobalToUniversalMapUnique (const int i) const
 
const Array< OneD, const int > & GetLocalToGlobalMap ()
 
const Array< OneD, const int > & GetGlobalToUniversalMap ()
 
const Array< OneD, const int > & GetGlobalToUniversalMapUnique ()
 
NekDouble GetLocalToGlobalSign (const int i) const
 
const Array< OneD, NekDouble > & GetLocalToGlobalSign () const
 
void LocalToGlobal (const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global, bool useComm=true) const
 
void LocalToGlobal (const NekVector< NekDouble > &loc, NekVector< NekDouble > &global, bool useComm=true) const
 
void GlobalToLocal (const Array< OneD, const NekDouble > &global, Array< OneD, NekDouble > &loc) const
 
void GlobalToLocal (const NekVector< NekDouble > &global, NekVector< NekDouble > &loc) const
 
void Assemble (const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global) const
 
void Assemble (const NekVector< NekDouble > &loc, NekVector< NekDouble > &global) const
 
void UniversalAssemble (Array< OneD, NekDouble > &pGlobal) const
 
void UniversalAssemble (NekVector< NekDouble > &pGlobal) const
 
void UniversalAssemble (Array< OneD, NekDouble > &pGlobal, int offset) const
 
void UniversalAbsMaxBnd (Array< OneD, NekDouble > &bndvals)
 
int GetLocalToGlobalBndMap (const int i) const
 Retrieve the global index of a given local boundary mode. More...
 
const Array< OneD, const int > & GetLocalToGlobalBndMap ()
 Retrieve the global indices of the local boundary modes. More...
 
const Array< OneD, const int > & GetGlobalToUniversalBndMap ()
 
const Array< OneD, const int > & GetGlobalToUniversalBndMapUnique ()
 
bool GetSignChange ()
 Returns true if using a modal expansion requiring a change of sign of some modes. More...
 
NekDouble GetLocalToGlobalBndSign (const int i) const
 Retrieve the sign change of a given local boundary mode. More...
 
Array< OneD, const NekDoubleGetLocalToGlobalBndSign () const
 Retrieve the sign change for all local boundary modes. More...
 
const Array< OneD, const int > & GetBndCondCoeffsToLocalCoeffsMap ()
 Retrieves the local indices corresponding to the boundary expansion modes. More...
 
const Array< OneD, NekDouble > & GetBndCondCoeffsToLocalCoeffsSign ()
 Returns the modal sign associated with a given boundary expansion mode. More...
 
const Array< OneD, const int > & GetBndCondCoeffsToLocalTraceMap ()
 Retrieves the local indices corresponding to the boundary expansion modes to global trace. More...
 
int GetBndCondIDToGlobalTraceID (const int i)
 Returns the global index of the boundary trace giving the index on the boundary expansion. More...
 
const Array< OneD, const int > & GetBndCondIDToGlobalTraceID ()
 
int GetNumGlobalDirBndCoeffs () const
 Returns the number of global Dirichlet boundary coefficients. More...
 
int GetNumLocalDirBndCoeffs () const
 Returns the number of local Dirichlet boundary coefficients. More...
 
int GetNumGlobalBndCoeffs () const
 Returns the total number of global boundary coefficients. More...
 
int GetNumLocalBndCoeffs () const
 Returns the total number of local boundary coefficients. More...
 
int GetNumLocalCoeffs () const
 Returns the total number of local coefficients. More...
 
int GetNumGlobalCoeffs () const
 Returns the total number of global coefficients. More...
 
bool GetSingularSystem () const
 Retrieves if the system is singular (true) or not (false) More...
 
void GlobalToLocalBnd (const NekVector< NekDouble > &global, NekVector< NekDouble > &loc, int offset) const
 
void GlobalToLocalBnd (const NekVector< NekDouble > &global, NekVector< NekDouble > &loc) const
 
void GlobalToLocalBnd (const Array< OneD, const NekDouble > &global, Array< OneD, NekDouble > &loc, int offset) const
 
void GlobalToLocalBnd (const Array< OneD, const NekDouble > &global, Array< OneD, NekDouble > &loc) const
 
void LocalBndToGlobal (const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global, int offset, bool UseComm=true) const
 
void LocalBndToGlobal (const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global, bool UseComm=true) const
 
void LocalToLocalBnd (const Array< OneD, const NekDouble > &local, Array< OneD, NekDouble > &locbnd) const
 
void LocalToLocalInt (const Array< OneD, const NekDouble > &local, Array< OneD, NekDouble > &locint) const
 
void LocalBndToLocal (const Array< OneD, const NekDouble > &locbnd, Array< OneD, NekDouble > &local) const
 
void LocalIntToLocal (const Array< OneD, const NekDouble > &locbnd, Array< OneD, NekDouble > &local) const
 
void AssembleBnd (const NekVector< NekDouble > &loc, NekVector< NekDouble > &global, int offset) const
 
void AssembleBnd (const NekVector< NekDouble > &loc, NekVector< NekDouble > &global) const
 
void AssembleBnd (const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global, int offset) const
 
void AssembleBnd (const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global) const
 
void UniversalAssembleBnd (Array< OneD, NekDouble > &pGlobal) const
 
void UniversalAssembleBnd (NekVector< NekDouble > &pGlobal) const
 
void UniversalAssembleBnd (Array< OneD, NekDouble > &pGlobal, int offset) const
 
int GetFullSystemBandWidth () const
 
int GetNumNonDirVertexModes () const
 
int GetNumNonDirEdgeModes () const
 
int GetNumNonDirFaceModes () const
 
int GetNumDirEdges () const
 
int GetNumDirFaces () const
 
int GetNumNonDirEdges () const
 
int GetNumNonDirFaces () const
 
void PrintStats (std::ostream &out, std::string variable, bool printHeader=true) const
 
const Array< OneD, const int > & GetExtraDirEdges ()
 
std::shared_ptr< AssemblyMapLinearSpaceMap (const ExpList &locexp, GlobalSysSolnType solnType)
 
int GetBndSystemBandWidth () const
 Returns the bandwidth of the boundary system. More...
 
int GetStaticCondLevel () const
 Returns the level of static condensation for this map. More...
 
int GetNumPatches () const
 Returns the number of patches in this static condensation level. More...
 
const Array< OneD, const unsigned int > & GetNumLocalBndCoeffsPerPatch ()
 Returns the number of local boundary coefficients in each patch. More...
 
const Array< OneD, const unsigned int > & GetNumLocalIntCoeffsPerPatch ()
 Returns the number of local interior coefficients in each patch. More...
 
const AssemblyMapSharedPtr GetNextLevelLocalToGlobalMap () const
 Returns the local to global mapping for the next level in the multi-level static condensation. More...
 
void SetNextLevelLocalToGlobalMap (AssemblyMapSharedPtr pNextLevelLocalToGlobalMap)
 
const PatchMapSharedPtrGetPatchMapFromPrevLevel (void) const
 Returns the patch map from the previous level of the multi-level static condensation. More...
 
bool AtLastLevel () const
 Returns true if this is the last level in the multi-level static condensation. More...
 
GlobalSysSolnType GetGlobalSysSolnType () const
 Returns the method of solving global systems. More...
 
std::string GetPreconType () const
 
NekDouble GetIterativeTolerance () const
 
bool IsAbsoluteTolerance () const
 
int GetSuccessiveRHS () const
 
std::string GetLinSysIterSolver () const
 
int GetLowestStaticCondLevel () const
 
void PatchLocalToGlobal (const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global) const
 
void PatchGlobalToLocal (const Array< OneD, const NekDouble > &global, Array< OneD, NekDouble > &loc) const
 
void PatchAssemble (const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global) const
 

Protected Member Functions

void CalculateBndSystemBandWidth ()
 Calculates the bandwidth of the boundary system. More...
 
void GlobalToLocalBndWithoutSign (const Array< OneD, const NekDouble > &global, Array< OneD, NekDouble > &loc)
 
virtual int v_GetLocalToGlobalMap (const int i) const
 
virtual int v_GetGlobalToUniversalMap (const int i) const
 
virtual int v_GetGlobalToUniversalMapUnique (const int i) const
 
virtual const Array< OneD, const int > & v_GetLocalToGlobalMap ()
 
virtual const Array< OneD, const int > & v_GetGlobalToUniversalMap ()
 
virtual const Array< OneD, const int > & v_GetGlobalToUniversalMapUnique ()
 
virtual NekDouble v_GetLocalToGlobalSign (const int i) const
 
virtual const Array< OneD, NekDouble > & v_GetLocalToGlobalSign () const
 
virtual void v_LocalToGlobal (const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global, bool useComm) const
 
virtual void v_LocalToGlobal (const NekVector< NekDouble > &loc, NekVector< NekDouble > &global, bool useComm) const
 
virtual void v_GlobalToLocal (const Array< OneD, const NekDouble > &global, Array< OneD, NekDouble > &loc) const
 
virtual void v_GlobalToLocal (const NekVector< NekDouble > &global, NekVector< NekDouble > &loc) const
 
virtual void v_Assemble (const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global) const
 
virtual void v_Assemble (const NekVector< NekDouble > &loc, NekVector< NekDouble > &global) const
 
virtual void v_UniversalAssemble (Array< OneD, NekDouble > &pGlobal) const
 
virtual void v_UniversalAssemble (NekVector< NekDouble > &pGlobal) const
 
virtual void v_UniversalAssemble (Array< OneD, NekDouble > &pGlobal, int offset) const
 
virtual int v_GetFullSystemBandWidth () const
 
virtual int v_GetNumNonDirVertexModes () const
 
virtual int v_GetNumNonDirEdgeModes () const
 
virtual int v_GetNumNonDirFaceModes () const
 
virtual int v_GetNumDirEdges () const
 
virtual int v_GetNumDirFaces () const
 
virtual int v_GetNumNonDirEdges () const
 
virtual int v_GetNumNonDirFaces () const
 
virtual const Array< OneD, const int > & v_GetExtraDirEdges ()
 
virtual std::shared_ptr< AssemblyMapv_LinearSpaceMap (const ExpList &locexp, GlobalSysSolnType solnType)
 Generate a linear space mapping from existing mapping. More...
 

Protected Attributes

LibUtilities::SessionReaderSharedPtr m_session
 Session object. More...
 
LibUtilities::CommSharedPtr m_comm
 Communicator. More...
 
std::string m_variable
 Variable string identifier. More...
 
size_t m_hash
 Hash for map. More...
 
int m_numLocalBndCoeffs
 Number of local boundary coefficients. More...
 
int m_numGlobalBndCoeffs
 Total number of global boundary coefficients. More...
 
int m_numLocalDirBndCoeffs
 Number of Local Dirichlet Boundary Coefficients. More...
 
int m_numGlobalDirBndCoeffs
 Number of Global Dirichlet Boundary Coefficients. More...
 
bool m_systemSingular
 Flag indicating if the system is singular or not. More...
 
int m_numLocalCoeffs
 Total number of local coefficients. More...
 
int m_numGlobalCoeffs
 Total number of global coefficients. More...
 
bool m_signChange
 Flag indicating if modes require sign reversal. More...
 
Array< OneD, int > m_localToGlobalBndMap
 Integer map of local coeffs to global Boundary Dofs. More...
 
Array< OneD, NekDoublem_localToGlobalBndSign
 Integer sign of local boundary coeffs to global space. More...
 
Array< OneD, int > m_localToLocalBndMap
 Integer map of local boundary coeffs to local boundary system numbering. More...
 
Array< OneD, int > m_localToLocalIntMap
 Integer map of local boundary coeffs to local interior system numbering. More...
 
Array< OneD, int > m_bndCondCoeffsToLocalCoeffsMap
 Integer map of bnd cond coeffs to local coefficients. More...
 
Array< OneD, NekDoublem_bndCondCoeffsToLocalCoeffsSign
 Integer map of sign of bnd cond coeffs to local coefficients. More...
 
Array< OneD, int > m_bndCondCoeffsToLocalTraceMap
 Integer map of bnd cond coeff to local trace coeff. More...
 
Array< OneD, int > m_bndCondIDToGlobalTraceID
 Integer map of bnd cond trace number to global trace number. More...
 
Array< OneD, int > m_globalToUniversalBndMap
 Integer map of process coeffs to universal space. More...
 
Array< OneD, int > m_globalToUniversalBndMapUnique
 Integer map of unique process coeffs to universal space (signed) More...
 
GlobalSysSolnType m_solnType
 The solution type of the global system. More...
 
int m_bndSystemBandWidth
 The bandwith of the global bnd system. More...
 
std::string m_preconType
 Type type of preconditioner to use in iterative solver. More...
 
NekDouble m_iterativeTolerance
 Tolerance for iterative solver. More...
 
bool m_isAbsoluteTolerance
 
int m_successiveRHS
 sucessive RHS for iterative solver More...
 
std::string m_linSysIterSolver
 Iterative solver: Conjugate Gradient, GMRES. More...
 
Gs::gs_datam_gsh
 
Gs::gs_datam_bndGsh
 
Gs::gs_datam_dirBndGsh
 gs gather communication to impose Dirhichlet BCs. More...
 
int m_staticCondLevel
 The level of recursion in the case of multi-level static condensation. More...
 
int m_numPatches
 The number of patches (~elements) in the current level. More...
 
Array< OneD, unsigned int > m_numLocalBndCoeffsPerPatch
 The number of bnd dofs per patch. More...
 
Array< OneD, unsigned int > m_numLocalIntCoeffsPerPatch
 The number of int dofs per patch. More...
 
AssemblyMapSharedPtr m_nextLevelLocalToGlobalMap
 Map from the patches of the previous level to the patches of the current level. More...
 
int m_lowestStaticCondLevel
 Lowest static condensation level. More...
 

Private Attributes

PatchMapSharedPtr m_patchMapFromPrevLevel
 Mapping information for previous level in MultiLevel Solver. More...
 

Detailed Description

Base class for constructing local to global mapping of degrees of freedom.

This class acts as a base class for constructing mappings between local, global and boundary degrees of freedom. It holds the storage for the maps and provides the accessors needed to retrieve them.

There are two derived classes: AssemblyMapCG and AssemblyMapDG. These perform the actual construction of the maps within their specific contexts.

Definition at line 55 of file AssemblyMap.h.

Constructor & Destructor Documentation

◆ AssemblyMap() [1/3]

Nektar::MultiRegions::AssemblyMap::AssemblyMap ( )

Default constructor.

Initialises an empty mapping.

Definition at line 79 of file AssemblyMap.cpp.

84 m_linSysIterSolver("ConjugateGradient"), m_gsh(nullptr), m_bndGsh(nullptr)
85{
86}
GlobalSysSolnType m_solnType
The solution type of the global system.
Definition: AssemblyMap.h:402
std::string m_linSysIterSolver
Iterative solver: Conjugate Gradient, GMRES.
Definition: AssemblyMap.h:417
int m_successiveRHS
sucessive RHS for iterative solver
Definition: AssemblyMap.h:414
LibUtilities::SessionReaderSharedPtr m_session
Session object.
Definition: AssemblyMap.h:333
int m_numLocalBndCoeffs
Number of local boundary coefficients.
Definition: AssemblyMap.h:345
int m_bndSystemBandWidth
The bandwith of the global bnd system.
Definition: AssemblyMap.h:404
int m_numLocalDirBndCoeffs
Number of Local Dirichlet Boundary Coefficients.
Definition: AssemblyMap.h:349
int m_numGlobalDirBndCoeffs
Number of Global Dirichlet Boundary Coefficients.
Definition: AssemblyMap.h:351
LibUtilities::CommSharedPtr m_comm
Communicator.
Definition: AssemblyMap.h:336
int m_numGlobalBndCoeffs
Total number of global boundary coefficients.
Definition: AssemblyMap.h:347
@ eNoSolnType
No Solution type specified.

◆ AssemblyMap() [2/3]

Nektar::MultiRegions::AssemblyMap::AssemblyMap ( const LibUtilities::SessionReaderSharedPtr pSession,
const LibUtilities::CommSharedPtr comm,
const std::string  variable = "DefaultVar" 
)

Constructor with a communicator.

Definition at line 88 of file AssemblyMap.cpp.

91 : m_session(pSession), m_comm(comm), m_variable(variable), m_hash(0),
95 m_linSysIterSolver("ConjugateGradient"), m_gsh(nullptr), m_bndGsh(nullptr)
96{
97 // Default value from Solver Info
99 pSession->GetSolverInfoAsEnum<GlobalSysSolnType>("GlobalSysSoln");
100
101 // Override values with data from GlobalSysSolnInfo section
102 if (pSession->DefinesGlobalSysSolnInfo(variable, "GlobalSysSoln"))
103 {
104 std::string sysSoln =
105 pSession->GetGlobalSysSolnInfo(variable, "GlobalSysSoln");
106 m_solnType = pSession->GetValueAsEnum<GlobalSysSolnType>(
107 "GlobalSysSoln", sysSoln);
108 }
109
110 // Set up preconditioner with SysSoln as override then solverinfo otherwise
111 // set a default as diagonal
112 if (pSession->DefinesGlobalSysSolnInfo(variable, "Preconditioner"))
113 {
115 pSession->GetGlobalSysSolnInfo(variable, "Preconditioner");
116 }
117 else if (pSession->DefinesSolverInfo("Preconditioner"))
118 {
119 m_preconType = pSession->GetSolverInfo("Preconditioner");
120 }
121 else
122 { // Possibly preconditioner is default registered in
123 // GlobLinSysIterative.cpp
124 m_preconType = "Diagonal";
125 }
126
127 if (pSession->DefinesGlobalSysSolnInfo(variable,
128 "IterativeSolverTolerance"))
129 {
130 m_iterativeTolerance = boost::lexical_cast<NekDouble>(
131 pSession->GetGlobalSysSolnInfo(variable, "IterativeSolverTolerance")
132 .c_str());
133 }
134 else
135 {
136 pSession->LoadParameter("IterativeSolverTolerance",
139 }
140
141 if (pSession->DefinesGlobalSysSolnInfo(variable, "AbsoluteTolerance"))
142 {
143 std::string abstol =
144 pSession->GetGlobalSysSolnInfo(variable, "AbsoluteTolerance");
146 boost::iequals(boost::to_upper_copy(abstol), "TRUE");
147 }
148 else
149 {
150 m_isAbsoluteTolerance = false;
151 }
152
153 if (pSession->DefinesGlobalSysSolnInfo(variable, "SuccessiveRHS"))
154 {
155 m_successiveRHS = boost::lexical_cast<int>(
156 pSession->GetGlobalSysSolnInfo(variable, "SuccessiveRHS").c_str());
157 }
158 else
159 {
160 pSession->LoadParameter("SuccessiveRHS", m_successiveRHS, 0);
161 }
162
163 if (pSession->DefinesGlobalSysSolnInfo(variable, "LinSysIterSolver"))
164 {
166 pSession->GetGlobalSysSolnInfo(variable, "LinSysIterSolver");
167 }
168 else if (pSession->DefinesSolverInfo("LinSysIterSolver"))
169 {
170 m_linSysIterSolver = pSession->GetSolverInfo("LinSysIterSolver");
171 }
172 else
173 {
174 m_linSysIterSolver = "ConjugateGradient";
175 }
176}
NekDouble m_iterativeTolerance
Tolerance for iterative solver.
Definition: AssemblyMap.h:410
std::string m_variable
Variable string identifier.
Definition: AssemblyMap.h:339
std::string m_preconType
Type type of preconditioner to use in iterative solver.
Definition: AssemblyMap.h:407
static const NekDouble kNekIterativeTol

References Nektar::NekConstants::kNekIterativeTol, m_isAbsoluteTolerance, m_iterativeTolerance, m_linSysIterSolver, m_preconType, m_solnType, and m_successiveRHS.

◆ AssemblyMap() [3/3]

Nektar::MultiRegions::AssemblyMap::AssemblyMap ( AssemblyMap oldLevelMap,
const BottomUpSubStructuredGraphSharedPtr multiLevelGraph 
)

Constructor for next level in multi-level static condensation.

Create a new level of mapping using the information in multiLevelGraph and performing the following steps:

  • STEP 1: setup a mask array to determine to which patch of the new level every patch of the current level belongs. To do so we make four arrays, #gloPatchMask, #globHomPatchMask, #locPatchMask_NekDouble and #locPatchMask. These arrays are then used to check which local dofs of the old level belong to which patch of the new level
  • STEP 2: We calculate how many local bnd dofs of the old level belong to the boundaries of each patch at the new level. We need this information to set up the mapping between different levels.

Definition at line 182 of file AssemblyMap.cpp.

185 : m_session(oldLevelMap->m_session), m_comm(oldLevelMap->GetComm()),
186 m_hash(0), m_solnType(oldLevelMap->m_solnType),
187 m_preconType(oldLevelMap->m_preconType),
188 m_iterativeTolerance(oldLevelMap->m_iterativeTolerance),
189 m_successiveRHS(oldLevelMap->m_successiveRHS),
190 m_linSysIterSolver(oldLevelMap->m_linSysIterSolver),
191 m_gsh(oldLevelMap->m_gsh), m_bndGsh(oldLevelMap->m_bndGsh),
192 m_lowestStaticCondLevel(oldLevelMap->m_lowestStaticCondLevel)
193{
194 int i;
195 int j;
196 int cnt;
197
198 //--------------------------------------------------------------
199 // -- Extract information from the input argument
200 int numGlobalBndCoeffsOld = oldLevelMap->GetNumGlobalBndCoeffs();
201 int numGlobalDirBndCoeffsOld = oldLevelMap->GetNumGlobalDirBndCoeffs();
202 int numLocalBndCoeffsOld = oldLevelMap->GetNumLocalBndCoeffs();
203 int numLocalDirBndCoeffsOld = oldLevelMap->GetNumLocalDirBndCoeffs();
204 bool signChangeOld = oldLevelMap->GetSignChange();
205
206 int staticCondLevelOld = oldLevelMap->GetStaticCondLevel();
207 int numPatchesOld = oldLevelMap->GetNumPatches();
208 GlobalSysSolnType solnTypeOld = oldLevelMap->GetGlobalSysSolnType();
209 const Array<OneD, const unsigned int> &numLocalBndCoeffsPerPatchOld =
210 oldLevelMap->GetNumLocalBndCoeffsPerPatch();
211 //--------------------------------------------------------------
212
213 //--------------------------------------------------------------
214 int newLevel = staticCondLevelOld + 1;
215 /** - STEP 1: setup a mask array to determine to which patch
216 * of the new level every patch of the current
217 * level belongs. To do so we make four arrays,
218 * #gloPatchMask, #globHomPatchMask,
219 * #locPatchMask_NekDouble and #locPatchMask.
220 * These arrays are then used to check which local
221 * dofs of the old level belong to which patch of
222 * the new level
223 */
224 Array<OneD, NekDouble> globPatchMask(numGlobalBndCoeffsOld, -1.0);
225 Array<OneD, NekDouble> globHomPatchMask(globPatchMask +
226 numGlobalDirBndCoeffsOld);
227 Array<OneD, NekDouble> locPatchMask_NekDouble(numLocalBndCoeffsOld, -3.0);
228 Array<OneD, int> locPatchMask(numLocalBndCoeffsOld);
229
230 // Fill the array globPatchMask as follows:
231 // - The first part (i.e. the glob bnd dofs) is filled with the
232 // value -1
233 // - The second part (i.e. the glob interior dofs) is numbered
234 // according to the patch it belongs to (i.e. dofs in first block
235 // all are numbered 0, the second block numbered are 1, etc...)
236 multiLevelGraph->MaskPatches(newLevel, globHomPatchMask);
237
238 // Map from Global Dofs to Local Dofs
239 // As a result, we know for each local dof whether
240 // it is mapped to the boundary of the next level, or to which
241 // patch. Based upon this, we can than later associate every patch
242 // of the current level with a patch in the next level.
243 oldLevelMap->GlobalToLocalBndWithoutSign(globPatchMask,
244 locPatchMask_NekDouble);
245
246 // Convert the result to an array of integers rather than doubles
247 RoundNekDoubleToInt(locPatchMask_NekDouble, locPatchMask);
248
249 /** - STEP 2: We calculate how many local bnd dofs of the
250 * old level belong to the boundaries of each patch at
251 * the new level. We need this information to set up the
252 * mapping between different levels.
253 */
254
255 // Retrieve the number of patches at the next level
256 int numPatchesWithIntNew =
257 multiLevelGraph->GetNpatchesWithInterior(newLevel);
258 int numPatchesNew = numPatchesWithIntNew;
259
260 // Allocate memory to store the number of local dofs associated to
261 // each of elemental boundaries of these patches
262 std::map<int, int> numLocalBndCoeffsPerPatchNew;
263 for (int i = 0; i < numPatchesNew; i++)
264 {
265 numLocalBndCoeffsPerPatchNew[i] = 0;
266 }
267
268 int minval;
269 int maxval;
270 int curPatch;
271 for (i = cnt = 0; i < numPatchesOld; i++)
272 {
273 // For every patch at the current level, the mask array
274 // locPatchMask should be filled with
275 // - the same (positive) number for each entry (which will
276 // correspond to the patch at the next level it belongs to)
277 // - the same (positive) number for each entry, except some
278 // entries that are -1 (the enties correspond to -1, will be
279 // mapped to the local boundary of the next level patch given
280 // by the positive number)
281 // - -1 for all entries. In this case, we will make an
282 // additional patch only consisting of boundaries at the next
283 // level
284 minval =
285 *min_element(&locPatchMask[cnt],
286 &locPatchMask[cnt] + numLocalBndCoeffsPerPatchOld[i]);
287 maxval =
288 *max_element(&locPatchMask[cnt],
289 &locPatchMask[cnt] + numLocalBndCoeffsPerPatchOld[i]);
290 ASSERTL0((minval == maxval) || (minval == -1),
291 "These values should never be the same");
292
293 if (maxval == -1)
294 {
295 curPatch = numPatchesNew;
296 numLocalBndCoeffsPerPatchNew[curPatch] = 0;
297 numPatchesNew++;
298 }
299 else
300 {
301 curPatch = maxval;
302 }
303
304 for (j = 0; j < numLocalBndCoeffsPerPatchOld[i]; j++)
305 {
306 ASSERTL0((locPatchMask[cnt] == maxval) ||
307 (locPatchMask[cnt] == minval),
308 "These values should never be the same");
309 if (locPatchMask[cnt] == -1)
310 {
311 ++numLocalBndCoeffsPerPatchNew[curPatch];
312 }
313 cnt++;
314 }
315 }
316
317 // Count how many local dofs of the old level are mapped
318 // to the local boundary dofs of the new level
321 m_numPatches = numLocalBndCoeffsPerPatchNew.size();
322 m_numLocalBndCoeffsPerPatch = Array<OneD, unsigned int>(m_numPatches);
323 m_numLocalIntCoeffsPerPatch = Array<OneD, unsigned int>(m_numPatches, 0u);
324 multiLevelGraph->GetNintDofsPerPatch(newLevel, m_numLocalIntCoeffsPerPatch);
325
326 for (int i = 0; i < m_numPatches; i++)
327 {
329 (unsigned int)numLocalBndCoeffsPerPatchNew[i];
333 }
334
335 // Also initialise some more data members
336 m_solnType = solnTypeOld;
341 "This method should only be called for in "
342 "case of multi-level static condensation.");
343 m_staticCondLevel = newLevel;
344 m_signChange = signChangeOld;
345 m_numLocalDirBndCoeffs = numLocalDirBndCoeffsOld;
346 m_numGlobalDirBndCoeffs = numGlobalDirBndCoeffsOld;
348 multiLevelGraph->GetInteriorOffset(newLevel) + m_numGlobalDirBndCoeffs;
350 multiLevelGraph->GetNumGlobalDofs(newLevel) + m_numGlobalDirBndCoeffs;
351 m_localToGlobalBndMap = Array<OneD, int>(m_numLocalBndCoeffs);
352
353 m_localToLocalBndMap = Array<OneD, int>(m_numLocalBndCoeffs);
355 Array<OneD, int>(m_numLocalCoeffs - m_numLocalBndCoeffs);
356
357 // local to bnd map is just a copy
358 for (int i = 0; i < m_numLocalBndCoeffs; ++i)
359 {
361 }
362
363 // local to int map is just a copy plus offset
364 for (int i = m_numLocalBndCoeffs; i < m_numLocalCoeffs; ++i)
365 {
367 }
368
369 if (m_signChange)
370 {
371 m_localToGlobalBndSign = Array<OneD, NekDouble>(m_numLocalBndCoeffs);
372 }
373
375 MemoryManager<PatchMap>::AllocateSharedPtr(numLocalBndCoeffsOld);
376
377 m_globalToUniversalBndMap = Array<OneD, int>(
378 m_numGlobalBndCoeffs, oldLevelMap->GetGlobalToUniversalBndMap());
379 m_globalToUniversalBndMapUnique = Array<OneD, int>(
380 m_numGlobalBndCoeffs, oldLevelMap->GetGlobalToUniversalBndMapUnique());
381
382 // Set up an offset array that denotes the offset of the local
383 // boundary degrees of freedom of the next level
384 Array<OneD, int> numLocalBndCoeffsPerPatchOffset(m_numPatches + 1, 0);
385 for (int i = 1; i < m_numPatches + 1; i++)
386 {
387 numLocalBndCoeffsPerPatchOffset[i] +=
388 numLocalBndCoeffsPerPatchOffset[i - 1] +
389 numLocalBndCoeffsPerPatchNew[i - 1];
390 }
391
392 int additionalPatchCnt = numPatchesWithIntNew;
393 int newid;
394 int blockid;
395 bool isBndDof;
397 Array<OneD, int> bndDofPerPatchCnt(m_numPatches, 0);
398
399 for (i = cnt = 0; i < numPatchesOld; i++)
400 {
401 minval =
402 *min_element(&locPatchMask[cnt],
403 &locPatchMask[cnt] + numLocalBndCoeffsPerPatchOld[i]);
404 maxval =
405 *max_element(&locPatchMask[cnt],
406 &locPatchMask[cnt] + numLocalBndCoeffsPerPatchOld[i]);
407 ASSERTL0((minval == maxval) || (minval == -1),
408 "These values should never be the same");
409
410 if (maxval == -1)
411 {
412 curPatch = additionalPatchCnt;
413 additionalPatchCnt++;
414 }
415 else
416 {
417 curPatch = maxval;
418 }
419
420 for (j = 0; j < numLocalBndCoeffsPerPatchOld[i]; j++)
421 {
422 ASSERTL0((locPatchMask[cnt] == maxval) ||
423 (locPatchMask[cnt] == minval),
424 "These values should never be the same");
425
426 sign = oldLevelMap->GetLocalToGlobalBndSign(cnt);
427
428 if (locPatchMask[cnt] == -1)
429 {
430 newid = numLocalBndCoeffsPerPatchOffset[curPatch];
431
432 m_localToGlobalBndMap[newid] =
433 oldLevelMap->GetLocalToGlobalBndMap(cnt);
434
435 if (m_signChange)
436 {
438 }
439
440 blockid = bndDofPerPatchCnt[curPatch];
441 isBndDof = true;
442
443 numLocalBndCoeffsPerPatchOffset[curPatch]++;
444 bndDofPerPatchCnt[curPatch]++;
445 }
446 else
447 {
448 newid = oldLevelMap->GetLocalToGlobalBndMap(cnt) -
450
451 blockid =
452 oldLevelMap->GetLocalToGlobalBndMap(cnt) -
454 multiLevelGraph->GetInteriorOffset(newLevel, curPatch);
455 isBndDof = false;
456 }
457
458 sign = isBndDof ? 1.0 : sign;
459
460 m_patchMapFromPrevLevel->SetPatchMap(cnt, curPatch, blockid,
461 isBndDof, sign);
462 cnt++;
463 }
464 }
465
466 // set up local to local mapping from previous to new level
469
471
472 // Postprocess the computed information - Update the old
473 // level with the mapping to new level
474 // oldLevelMap->SetLocalBndToNextLevelMap(oldLocalBndToNewLevelMap,oldLocalBndToNewLevelSign);
475
476 // - Construct the next level mapping object
477 if (m_staticCondLevel < (multiLevelGraph->GetNlevels() - 1))
478 {
481 multiLevelGraph);
482 }
483}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:47
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
int m_lowestStaticCondLevel
Lowest static condensation level.
Definition: AssemblyMap.h:439
int m_numLocalCoeffs
Total number of local coefficients.
Definition: AssemblyMap.h:364
bool m_signChange
Flag indicating if modes require sign reversal.
Definition: AssemblyMap.h:378
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
int m_numGlobalCoeffs
Total number of global coefficients.
Definition: AssemblyMap.h:375
Array< OneD, int > m_globalToUniversalBndMap
Integer map of process coeffs to universal space.
Definition: AssemblyMap.h:397
void CalculateBndSystemBandWidth()
Calculates the bandwidth of the boundary system.
Array< OneD, NekDouble > m_localToGlobalBndSign
Integer sign of local boundary coeffs to global space.
Definition: AssemblyMap.h:383
Array< OneD, unsigned int > m_numLocalBndCoeffsPerPatch
The number of bnd dofs per patch.
Definition: AssemblyMap.h:430
AssemblyMapSharedPtr m_nextLevelLocalToGlobalMap
Map from the patches of the previous level to the patches of the current level.
Definition: AssemblyMap.h:437
int m_staticCondLevel
The level of recursion in the case of multi-level static condensation.
Definition: AssemblyMap.h:426
Array< OneD, unsigned int > m_numLocalIntCoeffsPerPatch
The number of int dofs per patch.
Definition: AssemblyMap.h:432
Array< OneD, int > m_localToGlobalBndMap
Integer map of local coeffs to global Boundary Dofs.
Definition: AssemblyMap.h:381
Array< OneD, int > m_globalToUniversalBndMapUnique
Integer map of unique process coeffs to universal space (signed)
Definition: AssemblyMap.h:399
Array< OneD, int > m_localToLocalBndMap
Integer map of local boundary coeffs to local boundary system numbering.
Definition: AssemblyMap.h:385
int m_numPatches
The number of patches (~elements) in the current level.
Definition: AssemblyMap.h:428
int RoundNekDoubleToInt(NekDouble x)
Rounds a double precision number to an integer.
Definition: AssemblyMap.cpp:56
double NekDouble

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, ASSERTL1, CalculateBndSystemBandWidth(), Nektar::MultiRegions::eDirectMultiLevelStaticCond, Nektar::MultiRegions::eIterativeMultiLevelStaticCond, Nektar::MultiRegions::ePETScMultiLevelStaticCond, Nektar::MultiRegions::eXxtMultiLevelStaticCond, GetGlobalSysSolnType(), GetGlobalToUniversalBndMap(), GetGlobalToUniversalBndMapUnique(), GetLocalToGlobalBndMap(), GetLocalToGlobalBndSign(), GetNumGlobalBndCoeffs(), GetNumGlobalDirBndCoeffs(), GetNumLocalBndCoeffs(), GetNumLocalBndCoeffsPerPatch(), GetNumLocalDirBndCoeffs(), GetNumPatches(), GetSignChange(), GetStaticCondLevel(), GlobalToLocalBndWithoutSign(), m_globalToUniversalBndMap, m_globalToUniversalBndMapUnique, m_localToGlobalBndMap, m_localToGlobalBndSign, m_localToLocalBndMap, m_localToLocalIntMap, m_nextLevelLocalToGlobalMap, m_numGlobalBndCoeffs, m_numGlobalCoeffs, m_numGlobalDirBndCoeffs, m_numLocalBndCoeffs, m_numLocalBndCoeffsPerPatch, m_numLocalCoeffs, m_numLocalDirBndCoeffs, m_numLocalIntCoeffsPerPatch, m_numPatches, m_patchMapFromPrevLevel, m_signChange, m_solnType, m_staticCondLevel, Nektar::MultiRegions::RoundNekDoubleToInt(), and sign.

◆ ~AssemblyMap()

Nektar::MultiRegions::AssemblyMap::~AssemblyMap ( void  )
virtual

Destructor.

Definition at line 485 of file AssemblyMap.cpp.

486{
487}

Member Function Documentation

◆ Assemble() [1/2]

void Nektar::MultiRegions::AssemblyMap::Assemble ( const Array< OneD, const NekDouble > &  loc,
Array< OneD, NekDouble > &  global 
) const

Definition at line 804 of file AssemblyMap.cpp.

806{
807 v_Assemble(loc, global);
808}
virtual void v_Assemble(const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global) const

References CG_Iterations::loc, and v_Assemble().

Referenced by Nektar::MultiRegions::AssemblyMapCG::v_Assemble().

◆ Assemble() [2/2]

void Nektar::MultiRegions::AssemblyMap::Assemble ( const NekVector< NekDouble > &  loc,
NekVector< NekDouble > &  global 
) const

Definition at line 810 of file AssemblyMap.cpp.

812{
813 v_Assemble(loc, global);
814}

References CG_Iterations::loc, and v_Assemble().

◆ AssembleBnd() [1/4]

void Nektar::MultiRegions::AssemblyMap::AssembleBnd ( const Array< OneD, const NekDouble > &  loc,
Array< OneD, NekDouble > &  global 
) const

Definition at line 1259 of file AssemblyMap.cpp.

1261{
1263 "Local vector is not of correct dimension");
1264 ASSERTL1(global.size() >= m_numGlobalBndCoeffs,
1265 "Global vector is not of correct dimension");
1266
1267 Vmath::Zero(m_numGlobalBndCoeffs, global.get(), 1);
1268
1269 if (m_signChange)
1270 {
1272 loc.get(), m_localToGlobalBndMap.get(), global.get());
1273 }
1274 else
1275 {
1277 m_localToGlobalBndMap.get(), global.get());
1278 }
1279 UniversalAssembleBnd(global);
1280}
void UniversalAssembleBnd(Array< OneD, NekDouble > &pGlobal) const
void Assmb(int n, const T *x, const int *y, T *z)
Assemble z[y[i]] += x[i]; z should be zero'd first.
Definition: Vmath.hpp:577
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273

References ASSERTL1, Vmath::Assmb(), CG_Iterations::loc, m_localToGlobalBndMap, m_localToGlobalBndSign, m_numGlobalBndCoeffs, m_numLocalBndCoeffs, m_signChange, UniversalAssembleBnd(), and Vmath::Zero().

◆ AssembleBnd() [2/4]

void Nektar::MultiRegions::AssemblyMap::AssembleBnd ( const Array< OneD, const NekDouble > &  loc,
Array< OneD, NekDouble > &  global,
int  offset 
) const

Definition at line 1235 of file AssemblyMap.cpp.

1237{
1239 "Local array is not of correct dimension");
1240 ASSERTL1(global.size() >= m_numGlobalBndCoeffs - offset,
1241 "Global array is not of correct dimension");
1242 Array<OneD, NekDouble> tmp(m_numGlobalBndCoeffs, 0.0);
1243
1244 if (m_signChange)
1245 {
1247 loc.get(), m_localToGlobalBndMap.get(), tmp.get());
1248 }
1249 else
1250 {
1252 m_localToGlobalBndMap.get(), tmp.get());
1253 }
1255 Vmath::Vcopy(m_numGlobalBndCoeffs - offset, tmp.get() + offset, 1,
1256 global.get(), 1);
1257}
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825

References ASSERTL1, Vmath::Assmb(), CG_Iterations::loc, m_localToGlobalBndMap, m_localToGlobalBndSign, m_numGlobalBndCoeffs, m_numLocalBndCoeffs, m_signChange, UniversalAssembleBnd(), and Vmath::Vcopy().

◆ AssembleBnd() [3/4]

void Nektar::MultiRegions::AssemblyMap::AssembleBnd ( const NekVector< NekDouble > &  loc,
NekVector< NekDouble > &  global 
) const

Definition at line 1229 of file AssemblyMap.cpp.

1231{
1232 AssembleBnd(loc.GetPtr(), global.GetPtr());
1233}
void AssembleBnd(const NekVector< NekDouble > &loc, NekVector< NekDouble > &global, int offset) const
Array< OneD, DataType > & GetPtr()
Definition: NekVector.cpp:217

References AssembleBnd(), Nektar::NekVector< DataType >::GetPtr(), and CG_Iterations::loc.

◆ AssembleBnd() [4/4]

void Nektar::MultiRegions::AssemblyMap::AssembleBnd ( const NekVector< NekDouble > &  loc,
NekVector< NekDouble > &  global,
int  offset 
) const

Definition at line 1223 of file AssemblyMap.cpp.

1225{
1226 AssembleBnd(loc.GetPtr(), global.GetPtr(), offset);
1227}

References AssembleBnd(), Nektar::NekVector< DataType >::GetPtr(), and CG_Iterations::loc.

Referenced by AssembleBnd(), and Nektar::MultiRegions::AssemblyMapDG::v_Assemble().

◆ AtLastLevel()

bool Nektar::MultiRegions::AssemblyMap::AtLastLevel ( ) const

Returns true if this is the last level in the multi-level static condensation.

Definition at line 1350 of file AssemblyMap.cpp.

1351{
1352 return !(m_nextLevelLocalToGlobalMap.get());
1353}

References m_nextLevelLocalToGlobalMap.

◆ CalculateBndSystemBandWidth()

void Nektar::MultiRegions::AssemblyMap::CalculateBndSystemBandWidth ( )
protected

Calculates the bandwidth of the boundary system.

The bandwidth calculated corresponds to what is referred to as half-bandwidth. If the elements of the matrix are designated as a_ij, it corresponds to the maximum value of |i-j| for non-zero a_ij. As a result, the value also corresponds to the number of sub- or super-diagonals.

The bandwith can be calculated elementally as it corresponds to the maximal elemental bandwith (i.e. the maximal difference in global DOF index for every element).

We here calculate the bandwith of the global boundary system (as used for static condensation).

Definition at line 503 of file AssemblyMap.cpp.

504{
505 int i, j;
506 int cnt = 0;
507 int locSize;
508 int maxId;
509 int minId;
510 int bwidth = -1;
511 for (i = 0; i < m_numPatches; ++i)
512 {
513 locSize = m_numLocalBndCoeffsPerPatch[i];
514 maxId = -1;
515 minId = m_numLocalBndCoeffs + 1;
516 for (j = 0; j < locSize; j++)
517 {
519 {
520 if (m_localToGlobalBndMap[cnt + j] > maxId)
521 {
522 maxId = m_localToGlobalBndMap[cnt + j];
523 }
524
525 if (m_localToGlobalBndMap[cnt + j] < minId)
526 {
527 minId = m_localToGlobalBndMap[cnt + j];
528 }
529 }
530 }
531 bwidth = (bwidth > (maxId - minId)) ? bwidth : (maxId - minId);
532
533 cnt += locSize;
534 }
535
536 m_bndSystemBandWidth = bwidth;
537}

References m_bndSystemBandWidth, m_localToGlobalBndMap, m_numGlobalDirBndCoeffs, m_numLocalBndCoeffs, m_numLocalBndCoeffsPerPatch, and m_numPatches.

Referenced by AssemblyMap(), Nektar::MultiRegions::AssemblyMapCG::AssemblyMapCG(), and Nektar::MultiRegions::AssemblyMapDG::AssemblyMapDG().

◆ GetBndCondCoeffsToLocalCoeffsMap()

const Array< OneD, const int > & Nektar::MultiRegions::AssemblyMap::GetBndCondCoeffsToLocalCoeffsMap ( )

Retrieves the local indices corresponding to the boundary expansion modes.

Definition at line 989 of file AssemblyMap.cpp.

990{
992}
Array< OneD, int > m_bndCondCoeffsToLocalCoeffsMap
Integer map of bnd cond coeffs to local coefficients.
Definition: AssemblyMap.h:389

References m_bndCondCoeffsToLocalCoeffsMap.

◆ GetBndCondCoeffsToLocalCoeffsSign()

const Array< OneD, NekDouble > & Nektar::MultiRegions::AssemblyMap::GetBndCondCoeffsToLocalCoeffsSign ( )

Returns the modal sign associated with a given boundary expansion mode.

Definition at line 994 of file AssemblyMap.cpp.

995{
997}
Array< OneD, NekDouble > m_bndCondCoeffsToLocalCoeffsSign
Integer map of sign of bnd cond coeffs to local coefficients.
Definition: AssemblyMap.h:391

References m_bndCondCoeffsToLocalCoeffsSign.

◆ GetBndCondCoeffsToLocalTraceMap()

const Array< OneD, const int > & Nektar::MultiRegions::AssemblyMap::GetBndCondCoeffsToLocalTraceMap ( )

Retrieves the local indices corresponding to the boundary expansion modes to global trace.

Definition at line 999 of file AssemblyMap.cpp.

1000{
1002}
Array< OneD, int > m_bndCondCoeffsToLocalTraceMap
Integer map of bnd cond coeff to local trace coeff.
Definition: AssemblyMap.h:393

References m_bndCondCoeffsToLocalTraceMap.

◆ GetBndCondIDToGlobalTraceID() [1/2]

const Array< OneD, const int > & Nektar::MultiRegions::AssemblyMap::GetBndCondIDToGlobalTraceID ( )

Definition at line 1010 of file AssemblyMap.cpp.

1011{
1013}
Array< OneD, int > m_bndCondIDToGlobalTraceID
Integer map of bnd cond trace number to global trace number.
Definition: AssemblyMap.h:395

References m_bndCondIDToGlobalTraceID.

◆ GetBndCondIDToGlobalTraceID() [2/2]

int Nektar::MultiRegions::AssemblyMap::GetBndCondIDToGlobalTraceID ( const int  i)

Returns the global index of the boundary trace giving the index on the boundary expansion.

Definition at line 1004 of file AssemblyMap.cpp.

1005{
1006 ASSERTL1(i < m_bndCondIDToGlobalTraceID.size(), "Index out of range.");
1008}

References ASSERTL1, and m_bndCondIDToGlobalTraceID.

◆ GetBndSystemBandWidth()

int Nektar::MultiRegions::AssemblyMap::GetBndSystemBandWidth ( ) const

Returns the bandwidth of the boundary system.

Definition at line 1313 of file AssemblyMap.cpp.

1314{
1315 return m_bndSystemBandWidth;
1316}

References m_bndSystemBandWidth.

Referenced by Nektar::MultiRegions::AssemblyMapDG::v_GetFullSystemBandWidth().

◆ GetComm()

LibUtilities::CommSharedPtr Nektar::MultiRegions::AssemblyMap::GetComm ( )

Retrieves the communicator.

Definition at line 723 of file AssemblyMap.cpp.

724{
725 return m_comm;
726}

References m_comm.

◆ GetExtraDirEdges()

const Array< OneD, const int > & Nektar::MultiRegions::AssemblyMap::GetExtraDirEdges ( )

Definition at line 936 of file AssemblyMap.cpp.

937{
938 return v_GetExtraDirEdges();
939}
virtual const Array< OneD, const int > & v_GetExtraDirEdges()

References v_GetExtraDirEdges().

◆ GetFullSystemBandWidth()

int Nektar::MultiRegions::AssemblyMap::GetFullSystemBandWidth ( ) const

Definition at line 896 of file AssemblyMap.cpp.

897{
899}
virtual int v_GetFullSystemBandWidth() const

References v_GetFullSystemBandWidth().

◆ GetGlobalSysSolnType()

GlobalSysSolnType Nektar::MultiRegions::AssemblyMap::GetGlobalSysSolnType ( ) const

Returns the method of solving global systems.

Definition at line 1355 of file AssemblyMap.cpp.

1356{
1357 return m_solnType;
1358}

References m_solnType.

Referenced by AssemblyMap().

◆ GetGlobalToUniversalBndMap()

const Array< OneD, const int > & Nektar::MultiRegions::AssemblyMap::GetGlobalToUniversalBndMap ( )

Definition at line 967 of file AssemblyMap.cpp.

968{
970}

References m_globalToUniversalBndMap.

Referenced by AssemblyMap().

◆ GetGlobalToUniversalBndMapUnique()

const Array< OneD, const int > & Nektar::MultiRegions::AssemblyMap::GetGlobalToUniversalBndMapUnique ( )

Definition at line 972 of file AssemblyMap.cpp.

973{
975}

References m_globalToUniversalBndMapUnique.

Referenced by AssemblyMap().

◆ GetGlobalToUniversalMap() [1/2]

const Array< OneD, const int > & Nektar::MultiRegions::AssemblyMap::GetGlobalToUniversalMap ( )

Definition at line 758 of file AssemblyMap.cpp.

759{
761}
virtual const Array< OneD, const int > & v_GetGlobalToUniversalMap()

References v_GetGlobalToUniversalMap().

◆ GetGlobalToUniversalMap() [2/2]

int Nektar::MultiRegions::AssemblyMap::GetGlobalToUniversalMap ( const int  i) const

Definition at line 743 of file AssemblyMap.cpp.

744{
746}

References v_GetGlobalToUniversalMap().

◆ GetGlobalToUniversalMapUnique() [1/2]

const Array< OneD, const int > & Nektar::MultiRegions::AssemblyMap::GetGlobalToUniversalMapUnique ( )

Definition at line 763 of file AssemblyMap.cpp.

764{
766}
virtual const Array< OneD, const int > & v_GetGlobalToUniversalMapUnique()

References v_GetGlobalToUniversalMapUnique().

◆ GetGlobalToUniversalMapUnique() [2/2]

int Nektar::MultiRegions::AssemblyMap::GetGlobalToUniversalMapUnique ( const int  i) const

Definition at line 748 of file AssemblyMap.cpp.

749{
751}

References v_GetGlobalToUniversalMapUnique().

◆ GetHash()

size_t Nektar::MultiRegions::AssemblyMap::GetHash ( ) const

Retrieves the hash of this map.

Definition at line 733 of file AssemblyMap.cpp.

734{
735 return m_hash;
736}

References m_hash.

◆ GetIterativeTolerance()

NekDouble Nektar::MultiRegions::AssemblyMap::GetIterativeTolerance ( ) const

Definition at line 1365 of file AssemblyMap.cpp.

1366{
1367 return m_iterativeTolerance;
1368}

References m_iterativeTolerance.

◆ GetLinSysIterSolver()

std::string Nektar::MultiRegions::AssemblyMap::GetLinSysIterSolver ( ) const

Definition at line 1380 of file AssemblyMap.cpp.

1381{
1382 return m_linSysIterSolver;
1383}

References m_linSysIterSolver.

◆ GetLocalToGlobalBndMap() [1/2]

const Array< OneD, const int > & Nektar::MultiRegions::AssemblyMap::GetLocalToGlobalBndMap ( void  )

Retrieve the global indices of the local boundary modes.

Definition at line 952 of file AssemblyMap.cpp.

953{
955}

References m_localToGlobalBndMap.

◆ GetLocalToGlobalBndMap() [2/2]

int Nektar::MultiRegions::AssemblyMap::GetLocalToGlobalBndMap ( const int  i) const

Retrieve the global index of a given local boundary mode.

Definition at line 947 of file AssemblyMap.cpp.

948{
949 return m_localToGlobalBndMap[i];
950}

References m_localToGlobalBndMap.

Referenced by AssemblyMap().

◆ GetLocalToGlobalBndSign() [1/2]

Array< OneD, const NekDouble > Nektar::MultiRegions::AssemblyMap::GetLocalToGlobalBndSign ( void  ) const

Retrieve the sign change for all local boundary modes.

Definition at line 962 of file AssemblyMap.cpp.

963{
965}

References m_localToGlobalBndSign.

Referenced by Nektar::MultiRegions::AssemblyMapDG::v_GetLocalToGlobalSign().

◆ GetLocalToGlobalBndSign() [2/2]

NekDouble Nektar::MultiRegions::AssemblyMap::GetLocalToGlobalBndSign ( const int  i) const

Retrieve the sign change of a given local boundary mode.

Definition at line 977 of file AssemblyMap.cpp.

978{
979 if (m_signChange)
980 {
981 return m_localToGlobalBndSign[i];
982 }
983 else
984 {
985 return 1.0;
986 }
987}

References m_localToGlobalBndSign, and m_signChange.

Referenced by AssemblyMap().

◆ GetLocalToGlobalMap() [1/2]

const Array< OneD, const int > & Nektar::MultiRegions::AssemblyMap::GetLocalToGlobalMap ( )

Definition at line 753 of file AssemblyMap.cpp.

754{
755 return v_GetLocalToGlobalMap();
756}
virtual const Array< OneD, const int > & v_GetLocalToGlobalMap()

References v_GetLocalToGlobalMap().

◆ GetLocalToGlobalMap() [2/2]

int Nektar::MultiRegions::AssemblyMap::GetLocalToGlobalMap ( const int  i) const

Definition at line 738 of file AssemblyMap.cpp.

739{
740 return v_GetLocalToGlobalMap(i);
741}

References v_GetLocalToGlobalMap().

◆ GetLocalToGlobalSign() [1/2]

const Array< OneD, NekDouble > & Nektar::MultiRegions::AssemblyMap::GetLocalToGlobalSign ( ) const

Definition at line 773 of file AssemblyMap.cpp.

774{
775 return v_GetLocalToGlobalSign();
776}
virtual const Array< OneD, NekDouble > & v_GetLocalToGlobalSign() const

References v_GetLocalToGlobalSign().

◆ GetLocalToGlobalSign() [2/2]

NekDouble Nektar::MultiRegions::AssemblyMap::GetLocalToGlobalSign ( const int  i) const

Definition at line 768 of file AssemblyMap.cpp.

769{
770 return v_GetLocalToGlobalSign(i);
771}

References v_GetLocalToGlobalSign().

◆ GetLowestStaticCondLevel()

int Nektar::MultiRegions::AssemblyMap::GetLowestStaticCondLevel ( ) const
inline

Definition at line 314 of file AssemblyMap.h.

315 {
317 }

References m_lowestStaticCondLevel.

◆ GetNextLevelLocalToGlobalMap()

const AssemblyMapSharedPtr Nektar::MultiRegions::AssemblyMap::GetNextLevelLocalToGlobalMap ( ) const

Returns the local to global mapping for the next level in the multi-level static condensation.

Definition at line 1340 of file AssemblyMap.cpp.

1341{
1343}

References m_nextLevelLocalToGlobalMap.

◆ GetNumDirEdges()

int Nektar::MultiRegions::AssemblyMap::GetNumDirEdges ( ) const

Definition at line 916 of file AssemblyMap.cpp.

917{
918 return v_GetNumDirEdges();
919}
virtual int v_GetNumDirEdges() const

References v_GetNumDirEdges().

◆ GetNumDirFaces()

int Nektar::MultiRegions::AssemblyMap::GetNumDirFaces ( ) const

Definition at line 921 of file AssemblyMap.cpp.

922{
923 return v_GetNumDirFaces();
924}
virtual int v_GetNumDirFaces() const

References v_GetNumDirFaces().

◆ GetNumGlobalBndCoeffs()

int Nektar::MultiRegions::AssemblyMap::GetNumGlobalBndCoeffs ( ) const

Returns the total number of global boundary coefficients.

Definition at line 1030 of file AssemblyMap.cpp.

1031{
1032 return m_numGlobalBndCoeffs;
1033}

References m_numGlobalBndCoeffs.

Referenced by AssemblyMap().

◆ GetNumGlobalCoeffs()

int Nektar::MultiRegions::AssemblyMap::GetNumGlobalCoeffs ( ) const

Returns the total number of global coefficients.

Definition at line 1040 of file AssemblyMap.cpp.

1041{
1042 return m_numGlobalCoeffs;
1043}

References m_numGlobalCoeffs.

◆ GetNumGlobalDirBndCoeffs()

int Nektar::MultiRegions::AssemblyMap::GetNumGlobalDirBndCoeffs ( ) const

Returns the number of global Dirichlet boundary coefficients.

Definition at line 1015 of file AssemblyMap.cpp.

1016{
1018}

References m_numGlobalDirBndCoeffs.

Referenced by AssemblyMap().

◆ GetNumLocalBndCoeffs()

int Nektar::MultiRegions::AssemblyMap::GetNumLocalBndCoeffs ( ) const

Returns the total number of local boundary coefficients.

Definition at line 1025 of file AssemblyMap.cpp.

1026{
1027 return m_numLocalBndCoeffs;
1028}

References m_numLocalBndCoeffs.

Referenced by AssemblyMap().

◆ GetNumLocalBndCoeffsPerPatch()

const Array< OneD, const unsigned int > & Nektar::MultiRegions::AssemblyMap::GetNumLocalBndCoeffsPerPatch ( )

Returns the number of local boundary coefficients in each patch.

Definition at line 1329 of file AssemblyMap.cpp.

1330{
1332}

References m_numLocalBndCoeffsPerPatch.

Referenced by AssemblyMap().

◆ GetNumLocalCoeffs()

int Nektar::MultiRegions::AssemblyMap::GetNumLocalCoeffs ( ) const

Returns the total number of local coefficients.

Definition at line 1035 of file AssemblyMap.cpp.

1036{
1037 return m_numLocalCoeffs;
1038}

References m_numLocalCoeffs.

◆ GetNumLocalDirBndCoeffs()

int Nektar::MultiRegions::AssemblyMap::GetNumLocalDirBndCoeffs ( ) const

Returns the number of local Dirichlet boundary coefficients.

Definition at line 1020 of file AssemblyMap.cpp.

1021{
1023}

References m_numLocalDirBndCoeffs.

Referenced by AssemblyMap().

◆ GetNumLocalIntCoeffsPerPatch()

const Array< OneD, const unsigned int > & Nektar::MultiRegions::AssemblyMap::GetNumLocalIntCoeffsPerPatch ( )

Returns the number of local interior coefficients in each patch.

Definition at line 1335 of file AssemblyMap.cpp.

1336{
1338}

References m_numLocalIntCoeffsPerPatch.

◆ GetNumNonDirEdgeModes()

int Nektar::MultiRegions::AssemblyMap::GetNumNonDirEdgeModes ( ) const

Definition at line 906 of file AssemblyMap.cpp.

907{
909}
virtual int v_GetNumNonDirEdgeModes() const

References v_GetNumNonDirEdgeModes().

◆ GetNumNonDirEdges()

int Nektar::MultiRegions::AssemblyMap::GetNumNonDirEdges ( ) const

Definition at line 926 of file AssemblyMap.cpp.

927{
928 return v_GetNumNonDirEdges();
929}
virtual int v_GetNumNonDirEdges() const

References v_GetNumNonDirEdges().

◆ GetNumNonDirFaceModes()

int Nektar::MultiRegions::AssemblyMap::GetNumNonDirFaceModes ( ) const

Definition at line 911 of file AssemblyMap.cpp.

912{
914}
virtual int v_GetNumNonDirFaceModes() const

References v_GetNumNonDirFaceModes().

◆ GetNumNonDirFaces()

int Nektar::MultiRegions::AssemblyMap::GetNumNonDirFaces ( ) const

Definition at line 931 of file AssemblyMap.cpp.

932{
933 return v_GetNumNonDirFaces();
934}
virtual int v_GetNumNonDirFaces() const

References v_GetNumNonDirFaces().

◆ GetNumNonDirVertexModes()

int Nektar::MultiRegions::AssemblyMap::GetNumNonDirVertexModes ( ) const

Definition at line 901 of file AssemblyMap.cpp.

902{
904}
virtual int v_GetNumNonDirVertexModes() const

References v_GetNumNonDirVertexModes().

◆ GetNumPatches()

int Nektar::MultiRegions::AssemblyMap::GetNumPatches ( ) const

Returns the number of patches in this static condensation level.

Definition at line 1323 of file AssemblyMap.cpp.

1324{
1325 return m_numPatches;
1326}

References m_numPatches.

Referenced by AssemblyMap().

◆ GetPatchMapFromPrevLevel()

const PatchMapSharedPtr & Nektar::MultiRegions::AssemblyMap::GetPatchMapFromPrevLevel ( void  ) const

Returns the patch map from the previous level of the multi-level static condensation.

Definition at line 1345 of file AssemblyMap.cpp.

1346{
1348}

References m_patchMapFromPrevLevel.

◆ GetPreconType()

std::string Nektar::MultiRegions::AssemblyMap::GetPreconType ( ) const

Definition at line 1360 of file AssemblyMap.cpp.

1361{
1362 return m_preconType;
1363}

References m_preconType.

◆ GetSignChange()

bool Nektar::MultiRegions::AssemblyMap::GetSignChange ( )

Returns true if using a modal expansion requiring a change of sign of some modes.

Definition at line 957 of file AssemblyMap.cpp.

958{
959 return m_signChange;
960}

References m_signChange.

Referenced by AssemblyMap().

◆ GetSingularSystem()

bool Nektar::MultiRegions::AssemblyMap::GetSingularSystem ( ) const

Retrieves if the system is singular (true) or not (false)

Definition at line 1045 of file AssemblyMap.cpp.

1046{
1047 return m_systemSingular;
1048}
bool m_systemSingular
Flag indicating if the system is singular or not.
Definition: AssemblyMap.h:353

References m_systemSingular.

◆ GetStaticCondLevel()

int Nektar::MultiRegions::AssemblyMap::GetStaticCondLevel ( ) const

Returns the level of static condensation for this map.

Definition at line 1318 of file AssemblyMap.cpp.

1319{
1320 return m_staticCondLevel;
1321}

References m_staticCondLevel.

Referenced by AssemblyMap().

◆ GetSuccessiveRHS()

int Nektar::MultiRegions::AssemblyMap::GetSuccessiveRHS ( ) const

Definition at line 1375 of file AssemblyMap.cpp.

1376{
1377 return m_successiveRHS;
1378}

References m_successiveRHS.

◆ GetVariable()

std::string Nektar::MultiRegions::AssemblyMap::GetVariable ( )

Retrieves the variable string.

Definition at line 728 of file AssemblyMap.cpp.

729{
730 return m_variable;
731}

References m_variable.

◆ GlobalToLocal() [1/2]

void Nektar::MultiRegions::AssemblyMap::GlobalToLocal ( const Array< OneD, const NekDouble > &  global,
Array< OneD, NekDouble > &  loc 
) const

Definition at line 792 of file AssemblyMap.cpp.

794{
795 v_GlobalToLocal(global, loc);
796}
virtual void v_GlobalToLocal(const Array< OneD, const NekDouble > &global, Array< OneD, NekDouble > &loc) const

References CG_Iterations::loc, and v_GlobalToLocal().

Referenced by Nektar::MultiRegions::AssemblyMapCG::v_GlobalToLocal().

◆ GlobalToLocal() [2/2]

void Nektar::MultiRegions::AssemblyMap::GlobalToLocal ( const NekVector< NekDouble > &  global,
NekVector< NekDouble > &  loc 
) const

Definition at line 798 of file AssemblyMap.cpp.

800{
801 v_GlobalToLocal(global, loc);
802}

References CG_Iterations::loc, and v_GlobalToLocal().

◆ GlobalToLocalBnd() [1/4]

void Nektar::MultiRegions::AssemblyMap::GlobalToLocalBnd ( const Array< OneD, const NekDouble > &  global,
Array< OneD, NekDouble > &  loc 
) const

Definition at line 1088 of file AssemblyMap.cpp.

1090{
1092 "Local vector is not of correct dimension");
1093 ASSERTL1(global.size() >= m_numGlobalBndCoeffs,
1094 "Global vector is not of correct dimension");
1095
1096 Array<OneD, const NekDouble> glo;
1097 if (global.data() == loc.data())
1098 {
1099 glo = Array<OneD, NekDouble>(m_numLocalBndCoeffs, global.data());
1100 }
1101 else
1102 {
1103 glo = global; // create reference
1104 }
1105
1106 if (m_signChange)
1107 {
1109 glo.get(), m_localToGlobalBndMap.get(), loc.get());
1110 }
1111 else
1112 {
1114 m_localToGlobalBndMap.get(), loc.get());
1115 }
1116}
void Gathr(I n, const T *x, const I *y, T *z)
Gather vector z[i] = x[y[i]].
Definition: Vmath.hpp:507

References ASSERTL1, Vmath::Gathr(), CG_Iterations::loc, m_localToGlobalBndMap, m_localToGlobalBndSign, m_numGlobalBndCoeffs, m_numLocalBndCoeffs, and m_signChange.

◆ GlobalToLocalBnd() [2/4]

void Nektar::MultiRegions::AssemblyMap::GlobalToLocalBnd ( const Array< OneD, const NekDouble > &  global,
Array< OneD, NekDouble > &  loc,
int  offset 
) const

Definition at line 1062 of file AssemblyMap.cpp.

1065{
1067 "Local vector is not of correct dimension");
1068 ASSERTL1(global.size() >= m_numGlobalBndCoeffs - offset,
1069 "Global vector is not of correct dimension");
1070
1071 // offset input data by length "offset" for Dirichlet boundary conditions.
1072 Array<OneD, NekDouble> tmp(m_numGlobalBndCoeffs, 0.0);
1073 Vmath::Vcopy(m_numGlobalBndCoeffs - offset, global.get(), 1,
1074 tmp.get() + offset, 1);
1075
1076 if (m_signChange)
1077 {
1079 tmp.get(), m_localToGlobalBndMap.get(), loc.get());
1080 }
1081 else
1082 {
1084 m_localToGlobalBndMap.get(), loc.get());
1085 }
1086}

References ASSERTL1, Vmath::Gathr(), CG_Iterations::loc, m_localToGlobalBndMap, m_localToGlobalBndSign, m_numGlobalBndCoeffs, m_numLocalBndCoeffs, m_signChange, and Vmath::Vcopy().

◆ GlobalToLocalBnd() [3/4]

void Nektar::MultiRegions::AssemblyMap::GlobalToLocalBnd ( const NekVector< NekDouble > &  global,
NekVector< NekDouble > &  loc 
) const

Definition at line 1056 of file AssemblyMap.cpp.

1058{
1059 GlobalToLocalBnd(global.GetPtr(), loc.GetPtr());
1060}
void GlobalToLocalBnd(const NekVector< NekDouble > &global, NekVector< NekDouble > &loc, int offset) const

References Nektar::NekVector< DataType >::GetPtr(), GlobalToLocalBnd(), and CG_Iterations::loc.

◆ GlobalToLocalBnd() [4/4]

void Nektar::MultiRegions::AssemblyMap::GlobalToLocalBnd ( const NekVector< NekDouble > &  global,
NekVector< NekDouble > &  loc,
int  offset 
) const

◆ GlobalToLocalBndWithoutSign()

void Nektar::MultiRegions::AssemblyMap::GlobalToLocalBndWithoutSign ( const Array< OneD, const NekDouble > &  global,
Array< OneD, NekDouble > &  loc 
)
protected

Definition at line 1385 of file AssemblyMap.cpp.

1387{
1389 "Local vector is not of correct dimension");
1390 ASSERTL1(global.size() >= m_numGlobalBndCoeffs,
1391 "Global vector is not of correct dimension");
1392
1394 loc.get());
1395}

References ASSERTL1, Vmath::Gathr(), CG_Iterations::loc, m_localToGlobalBndMap, m_numGlobalBndCoeffs, and m_numLocalBndCoeffs.

Referenced by AssemblyMap().

◆ IsAbsoluteTolerance()

bool Nektar::MultiRegions::AssemblyMap::IsAbsoluteTolerance ( ) const

Definition at line 1370 of file AssemblyMap.cpp.

1371{
1372 return m_isAbsoluteTolerance;
1373}

References m_isAbsoluteTolerance.

◆ LinearSpaceMap()

std::shared_ptr< AssemblyMap > Nektar::MultiRegions::AssemblyMap::LinearSpaceMap ( const ExpList locexp,
GlobalSysSolnType  solnType 
)

Definition at line 941 of file AssemblyMap.cpp.

943{
944 return v_LinearSpaceMap(locexp, solnType);
945}
virtual std::shared_ptr< AssemblyMap > v_LinearSpaceMap(const ExpList &locexp, GlobalSysSolnType solnType)
Generate a linear space mapping from existing mapping.

References v_LinearSpaceMap().

◆ LocalBndToGlobal() [1/2]

void Nektar::MultiRegions::AssemblyMap::LocalBndToGlobal ( const Array< OneD, const NekDouble > &  loc,
Array< OneD, NekDouble > &  global,
bool  UseComm = true 
) const

Definition at line 1150 of file AssemblyMap.cpp.

1153{
1155 "Local vector is not of correct dimension");
1156 ASSERTL1(global.size() >= m_numGlobalBndCoeffs,
1157 "Global vector is not of correct dimension");
1158
1159 if (m_signChange)
1160 {
1162 loc.get(), m_localToGlobalBndMap.get(), global.get());
1163 }
1164 else
1165 {
1167 m_localToGlobalBndMap.get(), global.get());
1168 }
1169 if (UseComm)
1170 {
1171 Gs::Gather(global, Gs::gs_max, m_bndGsh);
1172 }
1173}
static void Gather(Nektar::Array< OneD, NekDouble > pU, gs_op pOp, gs_data *pGsh, Nektar::Array< OneD, NekDouble > pBuffer=NullNekDouble1DArray)
Performs a gather-scatter operation of the provided values.
Definition: GsLib.hpp:265
@ gs_max
Definition: GsLib.hpp:63
void Scatr(int n, const T *x, const int *y, T *z)
Scatter vector z[y[i]] = x[i].
Definition: Vmath.hpp:539

References ASSERTL1, Gs::Gather(), Gs::gs_max, CG_Iterations::loc, m_bndGsh, m_localToGlobalBndMap, m_localToGlobalBndSign, m_numGlobalBndCoeffs, m_numLocalBndCoeffs, m_signChange, and Vmath::Scatr().

◆ LocalBndToGlobal() [2/2]

void Nektar::MultiRegions::AssemblyMap::LocalBndToGlobal ( const Array< OneD, const NekDouble > &  loc,
Array< OneD, NekDouble > &  global,
int  offset,
bool  UseComm = true 
) const

Definition at line 1118 of file AssemblyMap.cpp.

1121{
1123 "Local vector is not of correct dimension");
1124 ASSERTL1(global.size() >= m_numGlobalBndCoeffs - offset,
1125 "Global vector is not of correct dimension");
1126
1127 // offset input data by length "offset" for Dirichlet boundary conditions.
1128 Array<OneD, NekDouble> tmp(m_numGlobalBndCoeffs, 0.0);
1129
1130 if (m_signChange)
1131 {
1133 loc.get(), m_localToGlobalBndMap.get(), tmp.get());
1134 }
1135 else
1136 {
1138 m_localToGlobalBndMap.get(), tmp.get());
1139 }
1140
1141 // Ensure each processor has unique value with a max gather.
1142 if (UseComm)
1143 {
1145 }
1146 Vmath::Vcopy(m_numGlobalBndCoeffs - offset, tmp.get() + offset, 1,
1147 global.get(), 1);
1148}

References ASSERTL1, Gs::Gather(), Gs::gs_max, CG_Iterations::loc, m_bndGsh, m_localToGlobalBndMap, m_localToGlobalBndSign, m_numGlobalBndCoeffs, m_numLocalBndCoeffs, m_signChange, Vmath::Scatr(), and Vmath::Vcopy().

Referenced by Nektar::MultiRegions::AssemblyMapDG::v_LocalToGlobal().

◆ LocalBndToLocal()

void Nektar::MultiRegions::AssemblyMap::LocalBndToLocal ( const Array< OneD, const NekDouble > &  locbnd,
Array< OneD, NekDouble > &  local 
) const

Definition at line 1199 of file AssemblyMap.cpp.

1201{
1202 ASSERTL1(locbnd.size() >= m_numLocalBndCoeffs,
1203 "LocBnd vector is not of correct dimension");
1204 ASSERTL1(local.size() >= m_numLocalCoeffs,
1205 "Local vector is not of correct dimension");
1206
1208 local.get());
1209}

References ASSERTL1, m_localToLocalBndMap, m_numLocalBndCoeffs, m_numLocalCoeffs, and Vmath::Scatr().

◆ LocalIntToLocal()

void Nektar::MultiRegions::AssemblyMap::LocalIntToLocal ( const Array< OneD, const NekDouble > &  locbnd,
Array< OneD, NekDouble > &  local 
) const

Definition at line 1211 of file AssemblyMap.cpp.

1213{
1215 "LocBnd vector is not of correct dimension");
1216 ASSERTL1(local.size() >= m_numLocalCoeffs,
1217 "Local vector is not of correct dimension");
1218
1220 m_localToLocalIntMap.get(), local.get());
1221}

References ASSERTL1, m_localToLocalIntMap, m_numLocalBndCoeffs, m_numLocalCoeffs, and Vmath::Scatr().

◆ LocalToGlobal() [1/2]

void Nektar::MultiRegions::AssemblyMap::LocalToGlobal ( const Array< OneD, const NekDouble > &  loc,
Array< OneD, NekDouble > &  global,
bool  useComm = true 
) const

Definition at line 778 of file AssemblyMap.cpp.

781{
782 v_LocalToGlobal(loc, global, useComm);
783}
virtual void v_LocalToGlobal(const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global, bool useComm) const

References CG_Iterations::loc, and v_LocalToGlobal().

Referenced by Nektar::MultiRegions::AssemblyMapCG::v_LocalToGlobal().

◆ LocalToGlobal() [2/2]

void Nektar::MultiRegions::AssemblyMap::LocalToGlobal ( const NekVector< NekDouble > &  loc,
NekVector< NekDouble > &  global,
bool  useComm = true 
) const

Definition at line 785 of file AssemblyMap.cpp.

788{
789 v_LocalToGlobal(loc, global, useComm);
790}

References CG_Iterations::loc, and v_LocalToGlobal().

◆ LocalToLocalBnd()

void Nektar::MultiRegions::AssemblyMap::LocalToLocalBnd ( const Array< OneD, const NekDouble > &  local,
Array< OneD, NekDouble > &  locbnd 
) const

Definition at line 1175 of file AssemblyMap.cpp.

1177{
1178 ASSERTL1(locbnd.size() >= m_numLocalBndCoeffs,
1179 "LocBnd vector is not of correct dimension");
1180 ASSERTL1(local.size() >= m_numLocalCoeffs,
1181 "Local vector is not of correct dimension");
1182
1184 locbnd.get());
1185}

References ASSERTL1, Vmath::Gathr(), m_localToLocalBndMap, m_numLocalBndCoeffs, and m_numLocalCoeffs.

◆ LocalToLocalInt()

void Nektar::MultiRegions::AssemblyMap::LocalToLocalInt ( const Array< OneD, const NekDouble > &  local,
Array< OneD, NekDouble > &  locint 
) const

Definition at line 1187 of file AssemblyMap.cpp.

1189{
1191 "Locint vector is not of correct dimension");
1192 ASSERTL1(local.size() >= m_numLocalCoeffs,
1193 "Local vector is not of correct dimension");
1194
1196 m_localToLocalIntMap.get(), locint.get());
1197}

References ASSERTL1, Vmath::Gathr(), m_localToLocalIntMap, m_numLocalBndCoeffs, and m_numLocalCoeffs.

◆ PatchAssemble()

void Nektar::MultiRegions::AssemblyMap::PatchAssemble ( const Array< OneD, const NekDouble > &  loc,
Array< OneD, NekDouble > &  global 
) const

Definition at line 872 of file AssemblyMap.cpp.

874{
875 Array<OneD, const NekDouble> local;
876 Array<OneD, const int> map = m_patchMapFromPrevLevel->GetNewLevelMap();
877 Array<OneD, const NekDouble> sign = m_patchMapFromPrevLevel->GetSign();
878
879 if (global.data() == loc.data())
880 {
881 local = Array<OneD, NekDouble>(map.size(), loc.data());
882 }
883 else
884 {
885 local = loc; // create reference
886 }
887
888 // since we are calling mapping from level down from array
889 // the m_numLocaBndCoeffs represents the size of the
890 // boundary elements we need to assemble into
891 Vmath::Zero(m_numLocalCoeffs, global.get(), 1);
892
893 Vmath::Assmb(map.size(), sign.get(), local.get(), map.get(), global.get());
894}

References Vmath::Assmb(), CG_Iterations::loc, m_numLocalCoeffs, m_patchMapFromPrevLevel, sign, and Vmath::Zero().

◆ PatchGlobalToLocal()

void Nektar::MultiRegions::AssemblyMap::PatchGlobalToLocal ( const Array< OneD, const NekDouble > &  global,
Array< OneD, NekDouble > &  loc 
) const

Definition at line 852 of file AssemblyMap.cpp.

854{
855 Array<OneD, const NekDouble> glo;
856
857 Array<OneD, const int> map = m_patchMapFromPrevLevel->GetNewLevelMap();
858 Array<OneD, const NekDouble> sign = m_patchMapFromPrevLevel->GetSign();
859
860 if (global.data() == loc.data())
861 {
862 glo = Array<OneD, NekDouble>(global.size(), global.data());
863 }
864 else
865 {
866 glo = global; // create reference
867 }
868
869 Vmath::Gathr(map.size(), sign.get(), glo.get(), map.get(), loc.get());
870}

References Vmath::Gathr(), CG_Iterations::loc, m_patchMapFromPrevLevel, and sign.

◆ PatchLocalToGlobal()

void Nektar::MultiRegions::AssemblyMap::PatchLocalToGlobal ( const Array< OneD, const NekDouble > &  loc,
Array< OneD, NekDouble > &  global 
) const

Definition at line 832 of file AssemblyMap.cpp.

834{
835 Array<OneD, const NekDouble> local;
836
837 Array<OneD, const int> map = m_patchMapFromPrevLevel->GetNewLevelMap();
838 Array<OneD, const NekDouble> sign = m_patchMapFromPrevLevel->GetSign();
839
840 if (global.data() == loc.data())
841 {
842 local = Array<OneD, NekDouble>(map.size(), loc.data());
843 }
844 else
845 {
846 local = loc; // create reference
847 }
848
849 Vmath::Scatr(map.size(), sign.get(), local.get(), map.get(), global.get());
850}

References CG_Iterations::loc, m_patchMapFromPrevLevel, Vmath::Scatr(), and sign.

◆ PrintStats()

void Nektar::MultiRegions::AssemblyMap::PrintStats ( std::ostream &  out,
std::string  variable,
bool  printHeader = true 
) const

Definition at line 1397 of file AssemblyMap.cpp.

1399{
1400 LibUtilities::CommSharedPtr vRowComm = m_session->GetComm()->GetRowComm();
1401 bool isRoot = m_session->GetComm()->IsParallelInTime()
1402 ? m_session->GetComm()->GetRank() == 0
1403 : vRowComm->GetRank() == 0;
1404 int n = vRowComm->GetSize();
1405 int i;
1406
1407 // Determine number of global degrees of freedom.
1408 int globBndCnt = 0, globDirCnt = 0;
1409
1410 for (i = 0; i < m_numGlobalBndCoeffs; ++i)
1411 {
1413 {
1414 globBndCnt++;
1415
1417 {
1418 globDirCnt++;
1419 }
1420 }
1421 }
1422
1423 int globCnt = m_numGlobalCoeffs - m_numGlobalBndCoeffs + globBndCnt;
1424
1425 // Calculate maximum valency
1426 Array<OneD, NekDouble> tmpLoc(m_numLocalBndCoeffs, 1.0);
1427 Array<OneD, NekDouble> tmpGlob(m_numGlobalBndCoeffs, 0.0);
1428
1430 tmpGlob.get());
1431 UniversalAssembleBnd(tmpGlob);
1432
1433 int totGlobDof = globCnt;
1434 int totGlobBndDof = globBndCnt;
1435 int totGlobDirDof = globDirCnt;
1436 int totLocalDof = m_numLocalCoeffs;
1437 int totLocalBndDof = m_numLocalBndCoeffs;
1438 int totLocalDirDof = m_numLocalDirBndCoeffs;
1439
1440 int meanValence = 0;
1441 int maxValence = 0;
1442 int minValence = 10000000;
1443 for (int i = 0; i < m_numGlobalBndCoeffs; ++i)
1444 {
1446 {
1447 continue;
1448 }
1449
1450 if (tmpGlob[i] > maxValence)
1451 {
1452 maxValence = tmpGlob[i];
1453 }
1454 if (tmpGlob[i] < minValence)
1455 {
1456 minValence = tmpGlob[i];
1457 }
1458 meanValence += tmpGlob[i];
1459 }
1460
1461 vRowComm->AllReduce(maxValence, LibUtilities::ReduceMax);
1462 vRowComm->AllReduce(minValence, LibUtilities::ReduceMin);
1463 vRowComm->AllReduce(meanValence, LibUtilities::ReduceSum);
1464 vRowComm->AllReduce(totGlobDof, LibUtilities::ReduceSum);
1465 vRowComm->AllReduce(totGlobBndDof, LibUtilities::ReduceSum);
1466 vRowComm->AllReduce(totGlobDirDof, LibUtilities::ReduceSum);
1467 vRowComm->AllReduce(totLocalDof, LibUtilities::ReduceSum);
1468 vRowComm->AllReduce(totLocalBndDof, LibUtilities::ReduceSum);
1469 vRowComm->AllReduce(totLocalDirDof, LibUtilities::ReduceSum);
1470
1471 meanValence /= totGlobBndDof;
1472
1473 if (isRoot)
1474 {
1475 if (printHeader)
1476 {
1477 out << "Assembly map statistics for field " << variable << ":"
1478 << endl;
1479 }
1480
1481 out << " - Number of local/global dof : " << totLocalDof
1482 << " " << totGlobDof << endl;
1483 out << " - Number of local/global boundary dof : " << totLocalBndDof
1484 << " " << totGlobBndDof << endl;
1485 out << " - Number of local/global Dirichlet dof : " << totLocalDirDof
1486 << " " << totGlobDirDof << endl;
1487 out << " - dof valency (min/max/mean) : " << minValence
1488 << " " << maxValence << " " << meanValence << endl;
1489
1490 if (n > 1)
1491 {
1492 NekDouble mean = m_numLocalCoeffs, mean2 = mean * mean;
1493 NekDouble minval = mean, maxval = mean;
1494 Array<OneD, NekDouble> tmp(1);
1495
1496 for (i = 1; i < n; ++i)
1497 {
1498 vRowComm->Recv(i, tmp);
1499 mean += tmp[0];
1500 mean2 += tmp[0] * tmp[0];
1501
1502 if (tmp[0] > maxval)
1503 {
1504 maxval = tmp[0];
1505 }
1506 if (tmp[0] < minval)
1507 {
1508 minval = tmp[0];
1509 }
1510 }
1511
1512 if (maxval > 0.1)
1513 {
1514 out << " - Local dof dist. (min/max/mean/dev) : " << minval
1515 << " " << maxval << " " << (mean / n) << " "
1516 << sqrt(mean2 / n - mean * mean / n / n) << endl;
1517 }
1518
1519 vRowComm->Block();
1520
1521 mean = minval = maxval = m_numLocalBndCoeffs;
1522 mean2 = mean * mean;
1523
1524 for (i = 1; i < n; ++i)
1525 {
1526 vRowComm->Recv(i, tmp);
1527 mean += tmp[0];
1528 mean2 += tmp[0] * tmp[0];
1529
1530 if (tmp[0] > maxval)
1531 {
1532 maxval = tmp[0];
1533 }
1534 if (tmp[0] < minval)
1535 {
1536 minval = tmp[0];
1537 }
1538 }
1539
1540 out << " - Local bnd dof dist. (min/max/mean/dev) : " << minval
1541 << " " << maxval << " " << (mean / n) << " "
1542 << sqrt(mean2 / n - mean * mean / n / n) << endl;
1543 }
1544 }
1545 else
1546 {
1547 Array<OneD, NekDouble> tmp(1);
1548 tmp[0] = m_numLocalCoeffs;
1549 vRowComm->Send(0, tmp);
1550 vRowComm->Block();
1551 tmp[0] = m_numLocalBndCoeffs;
1552 vRowComm->Send(0, tmp);
1553 }
1554
1555 // Either we have no more levels in the static condensation, or we
1556 // are not multi-level.
1558 {
1559 return;
1560 }
1561
1562 int level = 2;
1564 while (tmp->m_nextLevelLocalToGlobalMap)
1565 {
1566 tmp = tmp->m_nextLevelLocalToGlobalMap;
1567 ++level;
1568 }
1569
1570 // Print out multi-level static condensation information.
1571 if (n > 1)
1572 {
1573 if (isRoot)
1574 {
1575 NekDouble mean = level, mean2 = mean * mean;
1576 int minval = level, maxval = level;
1577
1578 Array<OneD, NekDouble> tmpRecv(1);
1579 for (i = 1; i < n; ++i)
1580 {
1581 vRowComm->Recv(i, tmpRecv);
1582 mean += tmpRecv[0];
1583 mean2 += tmpRecv[0] * tmpRecv[0];
1584
1585 if (tmpRecv[0] > maxval)
1586 {
1587 maxval = (int)(tmpRecv[0] + 0.5);
1588 }
1589 if (tmpRecv[0] < minval)
1590 {
1591 minval = (int)(tmpRecv[0] + 0.5);
1592 }
1593 }
1594
1595 out << " - M-level sc. dist. (min/max/mean/dev) : " << minval
1596 << " " << maxval << " " << (mean / n) << " "
1597 << sqrt(mean2 / n - mean * mean / n / n) << endl;
1598 }
1599 else
1600 {
1601 Array<OneD, NekDouble> tmpSend(1);
1602 tmpSend[0] = level;
1603 vRowComm->Send(0, tmpSend);
1604 }
1605 }
1606 else
1607 {
1608 if (isRoot)
1609 {
1610 out << " - Number of static cond. levels : " << level
1611 << endl;
1612 }
1613 }
1614
1615 if (isRoot)
1616 {
1617 out << "Stats at lowest static cond. level:" << endl;
1618 }
1619 tmp->PrintStats(out, variable, false);
1620}
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:55
std::shared_ptr< AssemblyMap > AssemblyMapSharedPtr
Definition: AssemblyMap.h:50
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294

References Vmath::Assmb(), CellMLToNektar.pycml::level, m_globalToUniversalBndMapUnique, m_localToGlobalBndMap, m_nextLevelLocalToGlobalMap, m_numGlobalBndCoeffs, m_numGlobalCoeffs, m_numGlobalDirBndCoeffs, m_numLocalBndCoeffs, m_numLocalCoeffs, m_numLocalDirBndCoeffs, m_session, Nektar::LibUtilities::ReduceMax, Nektar::LibUtilities::ReduceMin, Nektar::LibUtilities::ReduceSum, tinysimd::sqrt(), and UniversalAssembleBnd().

◆ SetNextLevelLocalToGlobalMap()

void Nektar::MultiRegions::AssemblyMap::SetNextLevelLocalToGlobalMap ( AssemblyMapSharedPtr  pNextLevelLocalToGlobalMap)

◆ UniversalAbsMaxBnd()

void Nektar::MultiRegions::AssemblyMap::UniversalAbsMaxBnd ( Array< OneD, NekDouble > &  bndvals)

Definition at line 1308 of file AssemblyMap.cpp.

1309{
1311}
Gs::gs_data * m_dirBndGsh
gs gather communication to impose Dirhichlet BCs.
Definition: AssemblyMap.h:422
@ gs_amax
Definition: GsLib.hpp:64

References Gs::Gather(), Gs::gs_amax, and m_dirBndGsh.

◆ UniversalAssemble() [1/3]

void Nektar::MultiRegions::AssemblyMap::UniversalAssemble ( Array< OneD, NekDouble > &  pGlobal) const

◆ UniversalAssemble() [2/3]

void Nektar::MultiRegions::AssemblyMap::UniversalAssemble ( Array< OneD, NekDouble > &  pGlobal,
int  offset 
) const

Definition at line 826 of file AssemblyMap.cpp.

828{
829 v_UniversalAssemble(pGlobal, offset);
830}

References v_UniversalAssemble().

◆ UniversalAssemble() [3/3]

void Nektar::MultiRegions::AssemblyMap::UniversalAssemble ( NekVector< NekDouble > &  pGlobal) const

Definition at line 821 of file AssemblyMap.cpp.

822{
823 v_UniversalAssemble(pGlobal);
824}

References v_UniversalAssemble().

◆ UniversalAssembleBnd() [1/3]

void Nektar::MultiRegions::AssemblyMap::UniversalAssembleBnd ( Array< OneD, NekDouble > &  pGlobal) const

Definition at line 1282 of file AssemblyMap.cpp.

1283{
1284 ASSERTL1(pGlobal.size() >= m_numGlobalBndCoeffs, "Wrong size.");
1285 Gs::Gather(pGlobal, Gs::gs_add, m_bndGsh);
1286}
@ gs_add
Definition: GsLib.hpp:60

References ASSERTL1, Gs::Gather(), Gs::gs_add, m_bndGsh, and m_numGlobalBndCoeffs.

Referenced by AssembleBnd(), PrintStats(), and UniversalAssembleBnd().

◆ UniversalAssembleBnd() [2/3]

void Nektar::MultiRegions::AssemblyMap::UniversalAssembleBnd ( Array< OneD, NekDouble > &  pGlobal,
int  offset 
) const

Definition at line 1293 of file AssemblyMap.cpp.

1295{
1296 Array<OneD, NekDouble> tmp(offset);
1297 if (offset > 0)
1298 {
1299 Vmath::Vcopy(offset, pGlobal, 1, tmp, 1);
1300 }
1301 UniversalAssembleBnd(pGlobal);
1302 if (offset > 0)
1303 {
1304 Vmath::Vcopy(offset, tmp, 1, pGlobal, 1);
1305 }
1306}

References UniversalAssembleBnd(), and Vmath::Vcopy().

◆ UniversalAssembleBnd() [3/3]

void Nektar::MultiRegions::AssemblyMap::UniversalAssembleBnd ( NekVector< NekDouble > &  pGlobal) const

Definition at line 1288 of file AssemblyMap.cpp.

1289{
1290 UniversalAssembleBnd(pGlobal.GetPtr());
1291}

References Nektar::NekVector< DataType >::GetPtr(), and UniversalAssembleBnd().

◆ v_Assemble() [1/2]

void Nektar::MultiRegions::AssemblyMap::v_Assemble ( const Array< OneD, const NekDouble > &  loc,
Array< OneD, NekDouble > &  global 
) const
protectedvirtual

Reimplemented in Nektar::MultiRegions::AssemblyMapCG, and Nektar::MultiRegions::AssemblyMapDG.

Definition at line 623 of file AssemblyMap.cpp.

626{
627 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
628}
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:202

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by Assemble().

◆ v_Assemble() [2/2]

void Nektar::MultiRegions::AssemblyMap::v_Assemble ( const NekVector< NekDouble > &  loc,
NekVector< NekDouble > &  global 
) const
protectedvirtual

Reimplemented in Nektar::MultiRegions::AssemblyMapCG, and Nektar::MultiRegions::AssemblyMapDG.

Definition at line 630 of file AssemblyMap.cpp.

633{
634 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
635}

References Nektar::ErrorUtil::efatal, and NEKERROR.

◆ v_GetExtraDirEdges()

const Array< OneD, const int > & Nektar::MultiRegions::AssemblyMap::v_GetExtraDirEdges ( )
protectedvirtual

Reimplemented in Nektar::MultiRegions::AssemblyMapCG.

Definition at line 707 of file AssemblyMap.cpp.

708{
709 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
710 static Array<OneD, const int> result;
711 return result;
712}

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by GetExtraDirEdges().

◆ v_GetFullSystemBandWidth()

int Nektar::MultiRegions::AssemblyMap::v_GetFullSystemBandWidth ( ) const
protectedvirtual

Reimplemented in Nektar::MultiRegions::AssemblyMapCG, and Nektar::MultiRegions::AssemblyMapDG.

Definition at line 659 of file AssemblyMap.cpp.

660{
661 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
662 return 0;
663}

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by GetFullSystemBandWidth().

◆ v_GetGlobalToUniversalMap() [1/2]

const Array< OneD, const int > & Nektar::MultiRegions::AssemblyMap::v_GetGlobalToUniversalMap ( void  )
protectedvirtual

Reimplemented in Nektar::MultiRegions::AssemblyMapCG, and Nektar::MultiRegions::AssemblyMapDG.

Definition at line 565 of file AssemblyMap.cpp.

566{
567 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
568 static Array<OneD, const int> result;
569 return result;
570}

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by GetGlobalToUniversalMap().

◆ v_GetGlobalToUniversalMap() [2/2]

int Nektar::MultiRegions::AssemblyMap::v_GetGlobalToUniversalMap ( const int  i) const
protectedvirtual

Reimplemented in Nektar::MultiRegions::AssemblyMapCG, and Nektar::MultiRegions::AssemblyMapDG.

Definition at line 545 of file AssemblyMap.cpp.

546{
547 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
548 return 0;
549}

References Nektar::ErrorUtil::efatal, and NEKERROR.

◆ v_GetGlobalToUniversalMapUnique() [1/2]

const Array< OneD, const int > & Nektar::MultiRegions::AssemblyMap::v_GetGlobalToUniversalMapUnique ( void  )
protectedvirtual

Reimplemented in Nektar::MultiRegions::AssemblyMapCG, and Nektar::MultiRegions::AssemblyMapDG.

Definition at line 572 of file AssemblyMap.cpp.

573{
574 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
575 static Array<OneD, const int> result;
576 return result;
577}

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by GetGlobalToUniversalMapUnique().

◆ v_GetGlobalToUniversalMapUnique() [2/2]

int Nektar::MultiRegions::AssemblyMap::v_GetGlobalToUniversalMapUnique ( const int  i) const
protectedvirtual

Reimplemented in Nektar::MultiRegions::AssemblyMapCG, and Nektar::MultiRegions::AssemblyMapDG.

Definition at line 551 of file AssemblyMap.cpp.

553{
554 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
555 return 0;
556}

References Nektar::ErrorUtil::efatal, and NEKERROR.

◆ v_GetLocalToGlobalMap() [1/2]

const Array< OneD, const int > & Nektar::MultiRegions::AssemblyMap::v_GetLocalToGlobalMap ( void  )
protectedvirtual

Reimplemented in Nektar::MultiRegions::AssemblyMapCG, and Nektar::MultiRegions::AssemblyMapDG.

Definition at line 558 of file AssemblyMap.cpp.

559{
560 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
561 static Array<OneD, const int> result;
562 return result;
563}

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by GetLocalToGlobalMap().

◆ v_GetLocalToGlobalMap() [2/2]

int Nektar::MultiRegions::AssemblyMap::v_GetLocalToGlobalMap ( const int  i) const
protectedvirtual

Reimplemented in Nektar::MultiRegions::AssemblyMapCG, and Nektar::MultiRegions::AssemblyMapDG.

Definition at line 539 of file AssemblyMap.cpp.

540{
541 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
542 return 0;
543}

References Nektar::ErrorUtil::efatal, and NEKERROR.

◆ v_GetLocalToGlobalSign() [1/2]

const Array< OneD, NekDouble > & Nektar::MultiRegions::AssemblyMap::v_GetLocalToGlobalSign ( ) const
protectedvirtual

Reimplemented in Nektar::MultiRegions::AssemblyMapCG.

Definition at line 586 of file AssemblyMap.cpp.

587{
588 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
589 static Array<OneD, NekDouble> result;
590 return result;
591}

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by GetLocalToGlobalSign().

◆ v_GetLocalToGlobalSign() [2/2]

NekDouble Nektar::MultiRegions::AssemblyMap::v_GetLocalToGlobalSign ( const int  i) const
protectedvirtual

Reimplemented in Nektar::MultiRegions::AssemblyMapCG, and Nektar::MultiRegions::AssemblyMapDG.

Definition at line 579 of file AssemblyMap.cpp.

581{
582 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
583 return 0.0;
584}

References Nektar::ErrorUtil::efatal, and NEKERROR.

◆ v_GetNumDirEdges()

int Nektar::MultiRegions::AssemblyMap::v_GetNumDirEdges ( ) const
protectedvirtual

Reimplemented in Nektar::MultiRegions::AssemblyMapCG.

Definition at line 683 of file AssemblyMap.cpp.

684{
685 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
686 return 0;
687}

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by GetNumDirEdges().

◆ v_GetNumDirFaces()

int Nektar::MultiRegions::AssemblyMap::v_GetNumDirFaces ( ) const
protectedvirtual

Reimplemented in Nektar::MultiRegions::AssemblyMapCG.

Definition at line 689 of file AssemblyMap.cpp.

690{
691 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
692 return 0;
693}

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by GetNumDirFaces().

◆ v_GetNumNonDirEdgeModes()

int Nektar::MultiRegions::AssemblyMap::v_GetNumNonDirEdgeModes ( ) const
protectedvirtual

Reimplemented in Nektar::MultiRegions::AssemblyMapCG.

Definition at line 671 of file AssemblyMap.cpp.

672{
673 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
674 return 0;
675}

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by GetNumNonDirEdgeModes().

◆ v_GetNumNonDirEdges()

int Nektar::MultiRegions::AssemblyMap::v_GetNumNonDirEdges ( ) const
protectedvirtual

Reimplemented in Nektar::MultiRegions::AssemblyMapCG.

Definition at line 695 of file AssemblyMap.cpp.

696{
697 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
698 return 0;
699}

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by GetNumNonDirEdges().

◆ v_GetNumNonDirFaceModes()

int Nektar::MultiRegions::AssemblyMap::v_GetNumNonDirFaceModes ( ) const
protectedvirtual

Reimplemented in Nektar::MultiRegions::AssemblyMapCG.

Definition at line 677 of file AssemblyMap.cpp.

678{
679 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
680 return 0;
681}

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by GetNumNonDirFaceModes().

◆ v_GetNumNonDirFaces()

int Nektar::MultiRegions::AssemblyMap::v_GetNumNonDirFaces ( ) const
protectedvirtual

Reimplemented in Nektar::MultiRegions::AssemblyMapCG.

Definition at line 701 of file AssemblyMap.cpp.

702{
703 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
704 return 0;
705}

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by GetNumNonDirFaces().

◆ v_GetNumNonDirVertexModes()

int Nektar::MultiRegions::AssemblyMap::v_GetNumNonDirVertexModes ( ) const
protectedvirtual

Reimplemented in Nektar::MultiRegions::AssemblyMapCG.

Definition at line 665 of file AssemblyMap.cpp.

666{
667 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
668 return 0;
669}

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by GetNumNonDirVertexModes().

◆ v_GlobalToLocal() [1/2]

void Nektar::MultiRegions::AssemblyMap::v_GlobalToLocal ( const Array< OneD, const NekDouble > &  global,
Array< OneD, NekDouble > &  loc 
) const
protectedvirtual

Reimplemented in Nektar::MultiRegions::AssemblyMapCG, and Nektar::MultiRegions::AssemblyMapDG.

Definition at line 609 of file AssemblyMap.cpp.

612{
613 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
614}

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by GlobalToLocal().

◆ v_GlobalToLocal() [2/2]

void Nektar::MultiRegions::AssemblyMap::v_GlobalToLocal ( const NekVector< NekDouble > &  global,
NekVector< NekDouble > &  loc 
) const
protectedvirtual

Reimplemented in Nektar::MultiRegions::AssemblyMapCG, and Nektar::MultiRegions::AssemblyMapDG.

Definition at line 616 of file AssemblyMap.cpp.

619{
620 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
621}

References Nektar::ErrorUtil::efatal, and NEKERROR.

◆ v_LinearSpaceMap()

std::shared_ptr< AssemblyMap > Nektar::MultiRegions::AssemblyMap::v_LinearSpaceMap ( const ExpList locexp,
GlobalSysSolnType  solnType 
)
protectedvirtual

Generate a linear space mapping from existing mapping.

Reimplemented in Nektar::MultiRegions::AssemblyMapCG.

Definition at line 714 of file AssemblyMap.cpp.

717{
718 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
719 static std::shared_ptr<AssemblyMap> result;
720 return result;
721}

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by LinearSpaceMap().

◆ v_LocalToGlobal() [1/2]

void Nektar::MultiRegions::AssemblyMap::v_LocalToGlobal ( const Array< OneD, const NekDouble > &  loc,
Array< OneD, NekDouble > &  global,
bool  useComm 
) const
protectedvirtual

Reimplemented in Nektar::MultiRegions::AssemblyMapCG, and Nektar::MultiRegions::AssemblyMapDG.

Definition at line 593 of file AssemblyMap.cpp.

597{
598 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
599}

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by LocalToGlobal().

◆ v_LocalToGlobal() [2/2]

void Nektar::MultiRegions::AssemblyMap::v_LocalToGlobal ( const NekVector< NekDouble > &  loc,
NekVector< NekDouble > &  global,
bool  useComm 
) const
protectedvirtual

Reimplemented in Nektar::MultiRegions::AssemblyMapCG.

Definition at line 601 of file AssemblyMap.cpp.

605{
606 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
607}

References Nektar::ErrorUtil::efatal, and NEKERROR.

◆ v_UniversalAssemble() [1/3]

void Nektar::MultiRegions::AssemblyMap::v_UniversalAssemble ( Array< OneD, NekDouble > &  pGlobal) const
protectedvirtual

Reimplemented in Nektar::MultiRegions::AssemblyMapCG, and Nektar::MultiRegions::AssemblyMapDG.

Definition at line 637 of file AssemblyMap.cpp.

639{
640 // Do nothing here since multi-level static condensation uses a
641 // AssemblyMap and thus will call this routine in serial.
642}

Referenced by UniversalAssemble().

◆ v_UniversalAssemble() [2/3]

void Nektar::MultiRegions::AssemblyMap::v_UniversalAssemble ( Array< OneD, NekDouble > &  pGlobal,
int  offset 
) const
protectedvirtual

Reimplemented in Nektar::MultiRegions::AssemblyMapCG.

Definition at line 651 of file AssemblyMap.cpp.

654{
655 // Do nothing here since multi-level static condensation uses a
656 // AssemblyMap and thus will call this routine in serial.
657}

◆ v_UniversalAssemble() [3/3]

void Nektar::MultiRegions::AssemblyMap::v_UniversalAssemble ( NekVector< NekDouble > &  pGlobal) const
protectedvirtual

Reimplemented in Nektar::MultiRegions::AssemblyMapCG, and Nektar::MultiRegions::AssemblyMapDG.

Definition at line 644 of file AssemblyMap.cpp.

646{
647 // Do nothing here since multi-level static condensation uses a
648 // AssemblyMap and thus will call this routine in serial.
649}

Member Data Documentation

◆ m_bndCondCoeffsToLocalCoeffsMap

Array<OneD, int> Nektar::MultiRegions::AssemblyMap::m_bndCondCoeffsToLocalCoeffsMap
protected

Integer map of bnd cond coeffs to local coefficients.

Definition at line 389 of file AssemblyMap.h.

Referenced by Nektar::MultiRegions::AssemblyMapCG::AssemblyMapCG(), and GetBndCondCoeffsToLocalCoeffsMap().

◆ m_bndCondCoeffsToLocalCoeffsSign

Array<OneD, NekDouble> Nektar::MultiRegions::AssemblyMap::m_bndCondCoeffsToLocalCoeffsSign
protected

Integer map of sign of bnd cond coeffs to local coefficients.

Definition at line 391 of file AssemblyMap.h.

Referenced by Nektar::MultiRegions::AssemblyMapCG::AssemblyMapCG(), and GetBndCondCoeffsToLocalCoeffsSign().

◆ m_bndCondCoeffsToLocalTraceMap

Array<OneD, int> Nektar::MultiRegions::AssemblyMap::m_bndCondCoeffsToLocalTraceMap
protected

Integer map of bnd cond coeff to local trace coeff.

Definition at line 393 of file AssemblyMap.h.

Referenced by Nektar::MultiRegions::AssemblyMapDG::AssemblyMapDG(), and GetBndCondCoeffsToLocalTraceMap().

◆ m_bndCondIDToGlobalTraceID

Array<OneD, int> Nektar::MultiRegions::AssemblyMap::m_bndCondIDToGlobalTraceID
protected

Integer map of bnd cond trace number to global trace number.

Definition at line 395 of file AssemblyMap.h.

Referenced by Nektar::MultiRegions::AssemblyMapDG::AssemblyMapDG(), Nektar::CoupledLocalToGlobalC0ContMap::CoupledLocalToGlobalC0ContMap(), and GetBndCondIDToGlobalTraceID().

◆ m_bndGsh

Gs::gs_data* Nektar::MultiRegions::AssemblyMap::m_bndGsh
protected

◆ m_bndSystemBandWidth

int Nektar::MultiRegions::AssemblyMap::m_bndSystemBandWidth
protected

The bandwith of the global bnd system.

Definition at line 404 of file AssemblyMap.h.

Referenced by CalculateBndSystemBandWidth(), and GetBndSystemBandWidth().

◆ m_comm

LibUtilities::CommSharedPtr Nektar::MultiRegions::AssemblyMap::m_comm
protected

◆ m_dirBndGsh

Gs::gs_data* Nektar::MultiRegions::AssemblyMap::m_dirBndGsh
protected

gs gather communication to impose Dirhichlet BCs.

Definition at line 422 of file AssemblyMap.h.

Referenced by Nektar::MultiRegions::AssemblyMapCG::AssemblyMapCG(), and UniversalAbsMaxBnd().

◆ m_globalToUniversalBndMap

Array<OneD, int> Nektar::MultiRegions::AssemblyMap::m_globalToUniversalBndMap
protected

◆ m_globalToUniversalBndMapUnique

Array<OneD, int> Nektar::MultiRegions::AssemblyMap::m_globalToUniversalBndMapUnique
protected

◆ m_gsh

Gs::gs_data* Nektar::MultiRegions::AssemblyMap::m_gsh
protected

◆ m_hash

size_t Nektar::MultiRegions::AssemblyMap::m_hash
protected

◆ m_isAbsoluteTolerance

bool Nektar::MultiRegions::AssemblyMap::m_isAbsoluteTolerance
protected

Definition at line 411 of file AssemblyMap.h.

Referenced by AssemblyMap(), and IsAbsoluteTolerance().

◆ m_iterativeTolerance

NekDouble Nektar::MultiRegions::AssemblyMap::m_iterativeTolerance
protected

Tolerance for iterative solver.

Definition at line 410 of file AssemblyMap.h.

Referenced by AssemblyMap(), and GetIterativeTolerance().

◆ m_linSysIterSolver

std::string Nektar::MultiRegions::AssemblyMap::m_linSysIterSolver
protected

Iterative solver: Conjugate Gradient, GMRES.

Definition at line 417 of file AssemblyMap.h.

Referenced by AssemblyMap(), and GetLinSysIterSolver().

◆ m_localToGlobalBndMap

Array<OneD, int> Nektar::MultiRegions::AssemblyMap::m_localToGlobalBndMap
protected

◆ m_localToGlobalBndSign

Array<OneD, NekDouble> Nektar::MultiRegions::AssemblyMap::m_localToGlobalBndSign
protected

◆ m_localToLocalBndMap

Array<OneD, int> Nektar::MultiRegions::AssemblyMap::m_localToLocalBndMap
protected

◆ m_localToLocalIntMap

Array<OneD, int> Nektar::MultiRegions::AssemblyMap::m_localToLocalIntMap
protected

◆ m_lowestStaticCondLevel

int Nektar::MultiRegions::AssemblyMap::m_lowestStaticCondLevel
protected

◆ m_nextLevelLocalToGlobalMap

AssemblyMapSharedPtr Nektar::MultiRegions::AssemblyMap::m_nextLevelLocalToGlobalMap
protected

Map from the patches of the previous level to the patches of the current level.

The local to global mapping of the next level of recursion

Definition at line 437 of file AssemblyMap.h.

Referenced by AssemblyMap(), Nektar::MultiRegions::AssemblyMapCG::AssemblyMapCG(), Nektar::MultiRegions::AssemblyMapDG::AssemblyMapDG(), AtLastLevel(), Nektar::CoupledLocalToGlobalC0ContMap::CoupledLocalToGlobalC0ContMap(), GetNextLevelLocalToGlobalMap(), and PrintStats().

◆ m_numGlobalBndCoeffs

int Nektar::MultiRegions::AssemblyMap::m_numGlobalBndCoeffs
protected

◆ m_numGlobalCoeffs

int Nektar::MultiRegions::AssemblyMap::m_numGlobalCoeffs
protected

Total number of global coefficients.

This corresponds to the number of total number of coefficients

  • For CG this corresponds to the total of bnd + int DOFs.
  • For DG this corresponds to the number of bnd DOFs. This means that m_numGlobalCoeffs = m_numGlobalBndCoeffs This way, we can consider the trace-system solve as a statically condensed solve without interior DOFs. This allows us to use the same global system solver for both cases.

Definition at line 375 of file AssemblyMap.h.

Referenced by AssemblyMap(), Nektar::MultiRegions::AssemblyMapCG::AssemblyMapCG(), Nektar::MultiRegions::AssemblyMapDG::AssemblyMapDG(), Nektar::CoupledAssemblyMap::CoupledAssemblyMap(), Nektar::CoupledLocalToGlobalC0ContMap::CoupledLocalToGlobalC0ContMap(), GetNumGlobalCoeffs(), PrintStats(), Nektar::MultiRegions::AssemblyMapCG::SetUpUniversalC0ContMap(), Nektar::MultiRegions::AssemblyMapCG::v_Assemble(), Nektar::MultiRegions::AssemblyMapCG::v_GlobalToLocal(), and Nektar::MultiRegions::AssemblyMapCG::v_LinearSpaceMap().

◆ m_numGlobalDirBndCoeffs

int Nektar::MultiRegions::AssemblyMap::m_numGlobalDirBndCoeffs
protected

◆ m_numLocalBndCoeffs

int Nektar::MultiRegions::AssemblyMap::m_numLocalBndCoeffs
protected

◆ m_numLocalBndCoeffsPerPatch

Array<OneD, unsigned int> Nektar::MultiRegions::AssemblyMap::m_numLocalBndCoeffsPerPatch
protected

◆ m_numLocalCoeffs

int Nektar::MultiRegions::AssemblyMap::m_numLocalCoeffs
protected

Total number of local coefficients.

This corresponds to the number of total number of coefficients

  • For CG this corresponds to the total of bnd + int DOFs
  • For DG this corresponds to the number of bnd DOFs. This means that m_numLocalCoeffs = m_numLocalBndCoeffs This way, we can consider the trace-system solve as a statically condensed solve without interior DOFs. This allows us to use the same global system solver for both cases.

Definition at line 364 of file AssemblyMap.h.

Referenced by AssemblyMap(), Nektar::MultiRegions::AssemblyMapCG::AssemblyMapCG(), Nektar::MultiRegions::AssemblyMapDG::AssemblyMapDG(), Nektar::MultiRegions::AssemblyMapCG::CalculateFullSystemBandWidth(), Nektar::CoupledAssemblyMap::CoupledAssemblyMap(), Nektar::CoupledLocalToGlobalC0ContMap::CoupledLocalToGlobalC0ContMap(), GetNumLocalCoeffs(), LocalBndToLocal(), LocalIntToLocal(), LocalToLocalBnd(), LocalToLocalInt(), PatchAssemble(), PrintStats(), Nektar::MultiRegions::AssemblyMapCG::v_Assemble(), Nektar::MultiRegions::AssemblyMapCG::v_GlobalToLocal(), and Nektar::MultiRegions::AssemblyMapCG::v_LocalToGlobal().

◆ m_numLocalDirBndCoeffs

int Nektar::MultiRegions::AssemblyMap::m_numLocalDirBndCoeffs
protected

◆ m_numLocalIntCoeffsPerPatch

Array<OneD, unsigned int> Nektar::MultiRegions::AssemblyMap::m_numLocalIntCoeffsPerPatch
protected

◆ m_numPatches

int Nektar::MultiRegions::AssemblyMap::m_numPatches
protected

◆ m_patchMapFromPrevLevel

PatchMapSharedPtr Nektar::MultiRegions::AssemblyMap::m_patchMapFromPrevLevel
private

Mapping information for previous level in MultiLevel Solver.

Definition at line 514 of file AssemblyMap.h.

Referenced by AssemblyMap(), GetPatchMapFromPrevLevel(), PatchAssemble(), PatchGlobalToLocal(), and PatchLocalToGlobal().

◆ m_preconType

std::string Nektar::MultiRegions::AssemblyMap::m_preconType
protected

Type type of preconditioner to use in iterative solver.

Definition at line 407 of file AssemblyMap.h.

Referenced by AssemblyMap(), and GetPreconType().

◆ m_session

LibUtilities::SessionReaderSharedPtr Nektar::MultiRegions::AssemblyMap::m_session
protected

◆ m_signChange

bool Nektar::MultiRegions::AssemblyMap::m_signChange
protected

◆ m_solnType

GlobalSysSolnType Nektar::MultiRegions::AssemblyMap::m_solnType
protected

◆ m_staticCondLevel

int Nektar::MultiRegions::AssemblyMap::m_staticCondLevel
protected

◆ m_successiveRHS

int Nektar::MultiRegions::AssemblyMap::m_successiveRHS
protected

sucessive RHS for iterative solver

Definition at line 414 of file AssemblyMap.h.

Referenced by AssemblyMap(), and GetSuccessiveRHS().

◆ m_systemSingular

bool Nektar::MultiRegions::AssemblyMap::m_systemSingular
protected

◆ m_variable

std::string Nektar::MultiRegions::AssemblyMap::m_variable
protected

Variable string identifier.

Definition at line 339 of file AssemblyMap.h.

Referenced by GetVariable().