Nektar++
Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | 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 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...
 
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...
 
PreconditionerType GetPreconType () const
 
NekDouble GetIterativeTolerance () const
 
int GetMaxIterations () 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)
 

Protected Attributes

LibUtilities::SessionReaderSharedPtr m_session
 Session object. More...
 
LibUtilities::CommSharedPtr m_comm
 Communicator. 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...
 
PreconditionerType m_preconType
 Type type of preconditioner to use in iterative solver. More...
 
int m_maxIterations
 Maximum iterations for iterative solver. More...
 
NekDouble m_iterativeTolerance
 Tolerance for iterative solver. More...
 
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 Member Functions

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...
 

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 58 of file AssemblyMap.h.

Constructor & Destructor Documentation

◆ AssemblyMap() [1/3]

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

Default constructor.

Initialises an empty mapping.

Definition at line 80 of file AssemblyMap.cpp.

80  :
81  m_session(),
82  m_comm(),
83  m_hash(0),
90  m_successiveRHS(0),
91  m_linSysIterSolver("ConjugateGradient"),
92  m_gsh(0),
93  m_bndGsh(0)
94  {
95  }
GlobalSysSolnType m_solnType
The solution type of the global system.
Definition: AssemblyMap.h:395
std::string m_linSysIterSolver
Iterative solver: Conjugate Gradient, GMRES.
Definition: AssemblyMap.h:412
int m_successiveRHS
sucessive RHS for iterative solver
Definition: AssemblyMap.h:409
LibUtilities::SessionReaderSharedPtr m_session
Session object.
Definition: AssemblyMap.h:329
int m_numLocalBndCoeffs
Number of local boundary coefficients.
Definition: AssemblyMap.h:338
int m_bndSystemBandWidth
The bandwith of the global bnd system.
Definition: AssemblyMap.h:397
int m_numLocalDirBndCoeffs
Number of Local Dirichlet Boundary Coefficients.
Definition: AssemblyMap.h:342
int m_numGlobalDirBndCoeffs
Number of Global Dirichlet Boundary Coefficients.
Definition: AssemblyMap.h:344
LibUtilities::CommSharedPtr m_comm
Communicator.
Definition: AssemblyMap.h:332
int m_numGlobalBndCoeffs
Total number of global boundary coefficients.
Definition: AssemblyMap.h:340
@ eNoSolnType
No Solution type specified.

◆ AssemblyMap() [2/3]

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

Constructor with a communicator.

Definition at line 97 of file AssemblyMap.cpp.

99  :
100  m_session(pSession),
101  m_comm(pSession->GetComm()),
102  m_hash(0),
108  m_successiveRHS(0),
109  m_linSysIterSolver("ConjugateGradient"),
110  m_gsh(0),
111  m_bndGsh(0)
112  {
113  // Default value from Solver Info
114  m_solnType = pSession->GetSolverInfoAsEnum<GlobalSysSolnType>(
115  "GlobalSysSoln");
116  m_preconType = pSession->GetSolverInfoAsEnum<PreconditionerType>(
117  "Preconditioner");
118 
119  // Override values with data from GlobalSysSolnInfo section
120  if(pSession->DefinesGlobalSysSolnInfo(variable, "GlobalSysSoln"))
121  {
122  std::string sysSoln = pSession->GetGlobalSysSolnInfo(variable,
123  "GlobalSysSoln");
124  m_solnType = pSession->GetValueAsEnum<GlobalSysSolnType>(
125  "GlobalSysSoln", sysSoln);
126  }
127 
128  if(pSession->DefinesGlobalSysSolnInfo(variable, "Preconditioner"))
129  {
130  std::string precon = pSession->GetGlobalSysSolnInfo(variable,
131  "Preconditioner");
132  m_preconType = pSession->GetValueAsEnum<PreconditionerType>(
133  "Preconditioner", precon);
134  }
135 
136  if(pSession->DefinesGlobalSysSolnInfo(variable,
137  "IterativeSolverTolerance"))
138  {
139  m_iterativeTolerance = boost::lexical_cast<NekDouble>(
140  pSession->GetGlobalSysSolnInfo(variable,
141  "IterativeSolverTolerance").c_str());
142  }
143  else
144  {
145  pSession->LoadParameter("IterativeSolverTolerance",
148  }
149 
150 
151  if(pSession->DefinesGlobalSysSolnInfo(variable,
152  "MaxIterations"))
153  {
154  m_maxIterations = boost::lexical_cast<int>(
155  pSession->GetGlobalSysSolnInfo(variable,
156  "MaxIterations").c_str());
157  }
158  else
159  {
160  pSession->LoadParameter("MaxIterations",
162  5000);
163  }
164 
165 
166  if(pSession->DefinesGlobalSysSolnInfo(variable,"SuccessiveRHS"))
167  {
168  m_successiveRHS = boost::lexical_cast<int>(
169  pSession->GetGlobalSysSolnInfo(variable,
170  "SuccessiveRHS").c_str());
171  }
172  else
173  {
174  pSession->LoadParameter("SuccessiveRHS",
175  m_successiveRHS,0);
176  }
177 
178  if (pSession->DefinesGlobalSysSolnInfo(variable, "LinSysIterSolver"))
179  {
180  m_linSysIterSolver = pSession->GetGlobalSysSolnInfo(
181  variable,"LinSysIterSolver");
182  }
183  else if (pSession->DefinesSolverInfo("LinSysIterSolver"))
184  {
185  m_linSysIterSolver = pSession->GetSolverInfo("LinSysIterSolver");
186  }
187  else
188  {
189  m_linSysIterSolver = "ConjugateGradient";
190  }
191  }
PreconditionerType m_preconType
Type type of preconditioner to use in iterative solver.
Definition: AssemblyMap.h:400
int m_maxIterations
Maximum iterations for iterative solver.
Definition: AssemblyMap.h:403
NekDouble m_iterativeTolerance
Tolerance for iterative solver.
Definition: AssemblyMap.h:406
static const NekDouble kNekIterativeTol

References Nektar::NekConstants::kNekIterativeTol, m_iterativeTolerance, m_linSysIterSolver, m_maxIterations, 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 197 of file AssemblyMap.cpp.

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

492  {
493  }

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 848 of file AssemblyMap.cpp.

851  {
852  v_Assemble(loc,global);
853  }
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 855 of file AssemblyMap.cpp.

858  {
859  v_Assemble(loc,global);
860  }

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 1301 of file AssemblyMap.cpp.

1304  {
1305  ASSERTL1(loc.size() >= m_numLocalBndCoeffs,"Local vector is not of correct dimension");
1306  ASSERTL1(global.size() >= m_numGlobalBndCoeffs,"Global vector is not of correct dimension");
1307 
1308  Vmath::Zero(m_numGlobalBndCoeffs, global.get(), 1);
1309 
1310  if(m_signChange)
1311  {
1313  loc.get(), m_localToGlobalBndMap.get(), global.get());
1314  }
1315  else
1316  {
1317  Vmath::Assmb(m_numLocalBndCoeffs,loc.get(), m_localToGlobalBndMap.get(), global.get());
1318  }
1319  UniversalAssembleBnd(global);
1320  }
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.cpp:813
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:436

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 1280 of file AssemblyMap.cpp.

1283  {
1284  ASSERTL1(loc.size() >= m_numLocalBndCoeffs,"Local array is not of correct dimension");
1285  ASSERTL1(global.size() >= m_numGlobalBndCoeffs-offset,"Global array is not of correct dimension");
1286  Array<OneD,NekDouble> tmp(m_numGlobalBndCoeffs,0.0);
1287 
1288  if(m_signChange)
1289  {
1291  }
1292  else
1293  {
1294  Vmath::Assmb(m_numLocalBndCoeffs,loc.get(), m_localToGlobalBndMap.get(), tmp.get());
1295  }
1296  UniversalAssembleBnd(tmp);
1297  Vmath::Vcopy(m_numGlobalBndCoeffs-offset, tmp.get() + offset, 1, global.get(), 1);
1298  }
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1199

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 1272 of file AssemblyMap.cpp.

1275  {
1276  AssembleBnd(loc.GetPtr(), global.GetPtr());
1277  }
void AssembleBnd(const NekVector< NekDouble > &loc, NekVector< NekDouble > &global, int offset) const

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 1264 of file AssemblyMap.cpp.

1267  {
1268  AssembleBnd(loc.GetPtr(), global.GetPtr(), offset);
1269  }

References 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 1392 of file AssemblyMap.cpp.

1393  {
1394  return !( m_nextLevelLocalToGlobalMap.get() );
1395  }

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 510 of file AssemblyMap.cpp.

511  {
512  int i,j;
513  int cnt = 0;
514  int locSize;
515  int maxId;
516  int minId;
517  int bwidth = -1;
518  for(i = 0; i < m_numPatches; ++i)
519  {
520  locSize = m_numLocalBndCoeffsPerPatch[i];
521  maxId = -1;
522  minId = m_numLocalBndCoeffs+1;
523  for(j = 0; j < locSize; j++)
524  {
526  {
527  if(m_localToGlobalBndMap[cnt+j] > maxId)
528  {
529  maxId = m_localToGlobalBndMap[cnt+j];
530  }
531 
532  if(m_localToGlobalBndMap[cnt+j] < minId)
533  {
534  minId = m_localToGlobalBndMap[cnt+j];
535  }
536  }
537  }
538  bwidth = (bwidth>(maxId-minId))?bwidth:(maxId-minId);
539 
540  cnt+=locSize;
541  }
542 
543  m_bndSystemBandWidth = bwidth;
544  }

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 1049 of file AssemblyMap.cpp.

1050  {
1052  }
Array< OneD, int > m_bndCondCoeffsToLocalCoeffsMap
Integer map of bnd cond coeffs to local coefficients.
Definition: AssemblyMap.h:382

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 1054 of file AssemblyMap.cpp.

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

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 1060 of file AssemblyMap.cpp.

1061  {
1063  }
Array< OneD, int > m_bndCondCoeffsToLocalTraceMap
Integer map of bnd cond coeff to local trace coeff.
Definition: AssemblyMap.h:386

References m_bndCondCoeffsToLocalTraceMap.

◆ GetBndCondIDToGlobalTraceID() [1/2]

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

Definition at line 1074 of file AssemblyMap.cpp.

1075  {
1077  }
Array< OneD, int > m_bndCondIDToGlobalTraceID
Integer map of bnd cond trace number to global trace number.
Definition: AssemblyMap.h:388

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 1066 of file AssemblyMap.cpp.

1068  {
1070  "Index out of range.");
1071  return m_bndCondIDToGlobalTraceID[i];
1072  }

References ASSERTL1, and m_bndCondIDToGlobalTraceID.

◆ GetBndSystemBandWidth()

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

Returns the bandwidth of the boundary system.

Definition at line 1351 of file AssemblyMap.cpp.

1352  {
1353  return m_bndSystemBandWidth;
1354  }

References m_bndSystemBandWidth.

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

◆ GetComm()

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

Retrieves the communicator.

Definition at line 767 of file AssemblyMap.cpp.

768  {
769  return m_comm;
770  }

References m_comm.

◆ GetExtraDirEdges()

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

Definition at line 992 of file AssemblyMap.cpp.

993  {
994  return v_GetExtraDirEdges();
995  }
virtual const Array< OneD, const int > & v_GetExtraDirEdges()

References v_GetExtraDirEdges().

◆ GetFullSystemBandWidth()

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

Definition at line 952 of file AssemblyMap.cpp.

953  {
954  return v_GetFullSystemBandWidth();
955  }
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 1398 of file AssemblyMap.cpp.

1399  {
1400  return m_solnType;
1401  }

References m_solnType.

Referenced by AssemblyMap().

◆ GetGlobalToUniversalBndMap()

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

Definition at line 1025 of file AssemblyMap.cpp.

1026  {
1028  }

References m_globalToUniversalBndMap.

Referenced by AssemblyMap().

◆ GetGlobalToUniversalBndMapUnique()

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

Definition at line 1030 of file AssemblyMap.cpp.

1031  {
1033  }

References m_globalToUniversalBndMapUnique.

Referenced by AssemblyMap().

◆ GetGlobalToUniversalMap() [1/2]

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

Definition at line 797 of file AssemblyMap.cpp.

798  {
799  return v_GetGlobalToUniversalMap();
800  }
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 782 of file AssemblyMap.cpp.

783  {
784  return v_GetGlobalToUniversalMap(i);
785  }

References v_GetGlobalToUniversalMap().

◆ GetGlobalToUniversalMapUnique() [1/2]

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

Definition at line 802 of file AssemblyMap.cpp.

803  {
805  }
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 787 of file AssemblyMap.cpp.

788  {
790  }

References v_GetGlobalToUniversalMapUnique().

◆ GetHash()

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

Retrieves the hash of this map.

Definition at line 772 of file AssemblyMap.cpp.

773  {
774  return m_hash;
775  }

References m_hash.

◆ GetIterativeTolerance()

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

Definition at line 1408 of file AssemblyMap.cpp.

1409  {
1410  return m_iterativeTolerance;
1411  }

References m_iterativeTolerance.

◆ GetLinSysIterSolver()

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

Definition at line 1423 of file AssemblyMap.cpp.

1424  {
1425  return m_linSysIterSolver;
1426  }

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 1008 of file AssemblyMap.cpp.

1009  {
1010  return m_localToGlobalBndMap;
1011  }

References m_localToGlobalBndMap.

Referenced by AssemblyMap().

◆ 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 1002 of file AssemblyMap.cpp.

1003  {
1004  return m_localToGlobalBndMap[i];
1005  }

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 1020 of file AssemblyMap.cpp.

1021  {
1022  return m_localToGlobalBndSign;
1023  }

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 1035 of file AssemblyMap.cpp.

1036  {
1037  if(m_signChange)
1038  {
1039  return m_localToGlobalBndSign[i];
1040  }
1041  else
1042  {
1043  return 1.0;
1044  }
1045  }

References m_localToGlobalBndSign, and m_signChange.

Referenced by AssemblyMap().

◆ GetLocalToGlobalMap() [1/2]

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

Definition at line 792 of file AssemblyMap.cpp.

793  {
794  return v_GetLocalToGlobalMap();
795  }
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 777 of file AssemblyMap.cpp.

778  {
779  return v_GetLocalToGlobalMap(i);
780  }

References v_GetLocalToGlobalMap().

◆ GetLocalToGlobalSign() [1/2]

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

Definition at line 812 of file AssemblyMap.cpp.

813  {
814  return v_GetLocalToGlobalSign();
815  }
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 807 of file AssemblyMap.cpp.

808  {
809  return v_GetLocalToGlobalSign(i);
810  }

References v_GetLocalToGlobalSign().

◆ GetLowestStaticCondLevel()

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

Definition at line 311 of file AssemblyMap.h.

312  {
314  }

References m_lowestStaticCondLevel.

◆ GetMaxIterations()

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

Definition at line 1413 of file AssemblyMap.cpp.

1414  {
1415  return m_maxIterations;
1416  }

References m_maxIterations.

◆ 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 1380 of file AssemblyMap.cpp.

1381  {
1383  }

References m_nextLevelLocalToGlobalMap.

◆ GetNumDirEdges()

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

Definition at line 972 of file AssemblyMap.cpp.

973  {
974  return v_GetNumDirEdges();
975  }
virtual int v_GetNumDirEdges() const

References v_GetNumDirEdges().

◆ GetNumDirFaces()

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

Definition at line 977 of file AssemblyMap.cpp.

978  {
979  return v_GetNumDirFaces();
980  }
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 1096 of file AssemblyMap.cpp.

1097  {
1098  return m_numGlobalBndCoeffs;
1099  }

References m_numGlobalBndCoeffs.

Referenced by AssemblyMap().

◆ GetNumGlobalCoeffs()

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

Returns the total number of global coefficients.

Definition at line 1106 of file AssemblyMap.cpp.

1107  {
1108  return m_numGlobalCoeffs;
1109  }

References m_numGlobalCoeffs.

◆ GetNumGlobalDirBndCoeffs()

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

Returns the number of global Dirichlet boundary coefficients.

Definition at line 1080 of file AssemblyMap.cpp.

1081  {
1082  return m_numGlobalDirBndCoeffs;
1083  }

References m_numGlobalDirBndCoeffs.

Referenced by AssemblyMap(), and Nektar::MultiRegions::PreconditionerBlock::BlockPreconditionerHDG().

◆ GetNumLocalBndCoeffs()

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

Returns the total number of local boundary coefficients.

Definition at line 1091 of file AssemblyMap.cpp.

1092  {
1093  return m_numLocalBndCoeffs;
1094  }

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 1367 of file AssemblyMap.cpp.

1368  {
1370  }

References m_numLocalBndCoeffsPerPatch.

Referenced by AssemblyMap().

◆ GetNumLocalCoeffs()

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

Returns the total number of local coefficients.

Definition at line 1101 of file AssemblyMap.cpp.

1102  {
1103  return m_numLocalCoeffs;
1104  }

References m_numLocalCoeffs.

◆ GetNumLocalDirBndCoeffs()

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

Returns the number of local Dirichlet boundary coefficients.

Definition at line 1086 of file AssemblyMap.cpp.

1087  {
1088  return m_numLocalDirBndCoeffs;
1089  }

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 1374 of file AssemblyMap.cpp.

1375  {
1377  }

References m_numLocalIntCoeffsPerPatch.

◆ GetNumNonDirEdgeModes()

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

Definition at line 962 of file AssemblyMap.cpp.

963  {
964  return v_GetNumNonDirEdgeModes();
965  }
virtual int v_GetNumNonDirEdgeModes() const

References v_GetNumNonDirEdgeModes().

◆ GetNumNonDirEdges()

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

Definition at line 982 of file AssemblyMap.cpp.

983  {
984  return v_GetNumNonDirEdges();
985  }
virtual int v_GetNumNonDirEdges() const

References v_GetNumNonDirEdges().

◆ GetNumNonDirFaceModes()

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

Definition at line 967 of file AssemblyMap.cpp.

968  {
969  return v_GetNumNonDirFaceModes();
970  }
virtual int v_GetNumNonDirFaceModes() const

References v_GetNumNonDirFaceModes().

◆ GetNumNonDirFaces()

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

Definition at line 987 of file AssemblyMap.cpp.

988  {
989  return v_GetNumNonDirFaces();
990  }
virtual int v_GetNumNonDirFaces() const

References v_GetNumNonDirFaces().

◆ GetNumNonDirVertexModes()

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

Definition at line 957 of file AssemblyMap.cpp.

958  {
959  return v_GetNumNonDirVertexModes();
960  }
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 1361 of file AssemblyMap.cpp.

1362  {
1363  return m_numPatches;
1364  }

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 1386 of file AssemblyMap.cpp.

1388  {
1389  return m_patchMapFromPrevLevel;
1390  }

References m_patchMapFromPrevLevel.

◆ GetPreconType()

PreconditionerType Nektar::MultiRegions::AssemblyMap::GetPreconType ( ) const

Definition at line 1403 of file AssemblyMap.cpp.

1404  {
1405  return m_preconType;
1406  }

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 1013 of file AssemblyMap.cpp.

1014  {
1015  return m_signChange;
1016  }

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 1111 of file AssemblyMap.cpp.

1112  {
1113  return m_systemSingular;
1114  }
bool m_systemSingular
Flag indicating if the system is singular or not.
Definition: AssemblyMap.h:346

References m_systemSingular.

◆ GetStaticCondLevel()

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

Returns the level of static condensation for this map.

Definition at line 1356 of file AssemblyMap.cpp.

1357  {
1358  return m_staticCondLevel;
1359  }

References m_staticCondLevel.

Referenced by AssemblyMap().

◆ GetSuccessiveRHS()

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

Definition at line 1418 of file AssemblyMap.cpp.

1419  {
1420  return m_successiveRHS;
1421  }

References m_successiveRHS.

◆ GlobalToLocal() [1/2]

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

Definition at line 834 of file AssemblyMap.cpp.

837  {
838  v_GlobalToLocal(global,loc);
839  }
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 841 of file AssemblyMap.cpp.

844  {
845  v_GlobalToLocal(global,loc);
846  }

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 1155 of file AssemblyMap.cpp.

1158  {
1159  ASSERTL1(loc.size() >= m_numLocalBndCoeffs,"Local vector is not of correct dimension");
1160  ASSERTL1(global.size() >= m_numGlobalBndCoeffs,"Global vector is not of correct dimension");
1161 
1162  if(m_signChange)
1163  {
1165  }
1166  else
1167  {
1168  Vmath::Gathr(m_numLocalBndCoeffs, global.get(), m_localToGlobalBndMap.get(), loc.get());
1169  }
1170  }
void Gathr(int n, const T *sign, const T *x, const int *y, T *z)
Gather vector z[i] = sign[i]*x[y[i]].
Definition: Vmath.cpp:756

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 1133 of file AssemblyMap.cpp.

1136  {
1137  ASSERTL1(loc.size() >= m_numLocalBndCoeffs,"Local vector is not of correct dimension");
1138  ASSERTL1(global.size() >= m_numGlobalBndCoeffs-offset,"Global vector is not of correct dimension");
1139 
1140  // offset input data by length "offset" for Dirichlet boundary conditions.
1141  Array<OneD,NekDouble> tmp(m_numGlobalBndCoeffs,0.0);
1142  Vmath::Vcopy(m_numGlobalBndCoeffs-offset, global.get(), 1, tmp.get() + offset, 1);
1143 
1144  if(m_signChange)
1145  {
1147  }
1148  else
1149  {
1150  Vmath::Gathr(m_numLocalBndCoeffs, tmp.get(), m_localToGlobalBndMap.get(), loc.get());
1151  }
1152  }

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 1125 of file AssemblyMap.cpp.

1128  {
1129  GlobalToLocalBnd(global.GetPtr(), loc.GetPtr());
1130  }
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

Definition at line 1116 of file AssemblyMap.cpp.

1120  {
1121  GlobalToLocalBnd(global.GetPtr(), loc.GetPtr(), offset);
1122  }

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

Referenced by GlobalToLocalBnd(), and Nektar::MultiRegions::AssemblyMapDG::v_GlobalToLocal().

◆ GlobalToLocalBndWithoutSign()

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

Definition at line 1428 of file AssemblyMap.cpp.

1431  {
1432  ASSERTL1(loc.size() >= m_numLocalBndCoeffs,"Local vector is not of correct dimension");
1433  ASSERTL1(global.size() >= m_numGlobalBndCoeffs,"Global vector is not of correct dimension");
1434 
1435  Vmath::Gathr(m_numLocalBndCoeffs, global.get(), m_localToGlobalBndMap.get(), loc.get());
1436  }

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

Referenced by AssemblyMap().

◆ LinearSpaceMap()

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

Definition at line 997 of file AssemblyMap.cpp.

998  {
999  return v_LinearSpaceMap(locexp, solnType);
1000  }
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 1200 of file AssemblyMap.cpp.

1203  {
1204  ASSERTL1(loc.size() >= m_numLocalBndCoeffs,"Local vector is not of correct dimension");
1205  ASSERTL1(global.size() >= m_numGlobalBndCoeffs,"Global vector is not of correct dimension");
1206 
1207  if(m_signChange)
1208  {
1210  }
1211  else
1212  {
1213  Vmath::Scatr(m_numLocalBndCoeffs, loc.get(), m_localToGlobalBndMap.get(), global.get());
1214  }
1215  if(UseComm)
1216  {
1217  Gs::Gather(global, Gs::gs_max, m_bndGsh);
1218  }
1219  }
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:245
@ gs_max
Definition: GsLib.hpp:53
void Scatr(int n, const T *x, const int *y, T *z)
Scatter vector z[y[i]] = x[i].
Definition: Vmath.cpp:772

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 1172 of file AssemblyMap.cpp.

1176  {
1177  ASSERTL1(loc.size() >= m_numLocalBndCoeffs,"Local vector is not of correct dimension");
1178  ASSERTL1(global.size() >= m_numGlobalBndCoeffs-offset,"Global vector is not of correct dimension");
1179 
1180  // offset input data by length "offset" for Dirichlet boundary conditions.
1181  Array<OneD,NekDouble> tmp(m_numGlobalBndCoeffs,0.0);
1182 
1183  if(m_signChange)
1184  {
1186  }
1187  else
1188  {
1189  Vmath::Scatr(m_numLocalBndCoeffs, loc.get(), m_localToGlobalBndMap.get(), tmp.get());
1190  }
1191 
1192  // Ensure each processor has unique value with a max gather.
1193  if(UseComm)
1194  {
1196  }
1197  Vmath::Vcopy(m_numGlobalBndCoeffs-offset, tmp.get()+offset, 1, global.get(), 1);
1198  }

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 1243 of file AssemblyMap.cpp.

1246  {
1247  ASSERTL1(locbnd.size() >= m_numLocalBndCoeffs,"LocBnd vector is not of correct dimension");
1248  ASSERTL1(local.size() >= m_numLocalCoeffs,"Local vector is not of correct dimension");
1249 
1250  Vmath::Scatr(m_numLocalBndCoeffs, locbnd.get(), m_localToLocalBndMap.get(), local.get());
1251  }

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 1253 of file AssemblyMap.cpp.

1256  {
1257  ASSERTL1(locint.size() >= m_numLocalCoeffs-m_numLocalBndCoeffs,"LocBnd vector is not of correct dimension");
1258  ASSERTL1(local.size() >= m_numLocalCoeffs,"Local vector is not of correct dimension");
1259 
1260  Vmath::Scatr(m_numLocalCoeffs-m_numLocalBndCoeffs, locint.get(), m_localToLocalIntMap.get(), local.get());
1261  }

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 818 of file AssemblyMap.cpp.

822  {
823  v_LocalToGlobal(loc,global,useComm);
824  }
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 826 of file AssemblyMap.cpp.

830  {
831  v_LocalToGlobal(loc,global,useComm);
832  }

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 1222 of file AssemblyMap.cpp.

1225  {
1226  ASSERTL1(locbnd.size() >= m_numLocalBndCoeffs,"LocBnd vector is not of correct dimension");
1227  ASSERTL1(local.size() >= m_numLocalCoeffs,"Local vector is not of correct dimension");
1228 
1229  Vmath::Gathr(m_numLocalBndCoeffs, local.get(), m_localToLocalBndMap.get(), locbnd.get());
1230  }

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 1232 of file AssemblyMap.cpp.

1235  {
1236  ASSERTL1(locint.size() >= m_numLocalCoeffs-m_numLocalBndCoeffs,"Locint vector is not of correct dimension");
1237  ASSERTL1(local.size() >= m_numLocalCoeffs,"Local vector is not of correct dimension");
1238 
1239  Vmath::Gathr(m_numLocalCoeffs-m_numLocalBndCoeffs, local.get(), m_localToLocalIntMap.get(), locint.get());
1240  }

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 926 of file AssemblyMap.cpp.

929  {
930  Array<OneD, const NekDouble> local;
931  Array<OneD, const int> map = m_patchMapFromPrevLevel->GetNewLevelMap();
932  Array<OneD, const NekDouble> sign = m_patchMapFromPrevLevel->GetSign();
933 
934  if(global.data() == loc.data())
935  {
936  local = Array<OneD, NekDouble>(map.size(),loc.data());
937  }
938  else
939  {
940  local = loc; // create reference
941  }
942 
943  // since we are calling mapping from level down from array
944  // the m_numLocaBndCoeffs represents the size of the
945  // boundary elements we need to assemble into
946  Vmath::Zero(m_numLocalCoeffs,global.get(), 1);
947 
948  Vmath::Assmb(map.size(), sign.get(), local.get(),
949  map.get(), global.get());
950  }

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 905 of file AssemblyMap.cpp.

908  {
909  Array<OneD, const NekDouble> glo;
910 
911  Array<OneD, const int> map = m_patchMapFromPrevLevel->GetNewLevelMap();
912  Array<OneD, const NekDouble> sign = m_patchMapFromPrevLevel->GetSign();
913 
914  if(global.data() == loc.data())
915  {
916  glo = Array<OneD, NekDouble>(global.size(),global.data());
917  }
918  else
919  {
920  glo = global; // create reference
921  }
922 
923  Vmath::Gathr(map.size(),sign.get(),glo.get(),map.get(),loc.get());
924  }

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 881 of file AssemblyMap.cpp.

884  {
885  Array<OneD, const NekDouble> local;
886 
887  Array<OneD, const int> map = m_patchMapFromPrevLevel->GetNewLevelMap();
888  Array<OneD, const NekDouble> sign = m_patchMapFromPrevLevel->GetSign();
889 
890  if(global.data() == loc.data())
891  {
892  local = Array<OneD, NekDouble>(map.size(),loc.data());
893  }
894  else
895  {
896  local = loc; // create reference
897  }
898 
899 
900  Vmath::Scatr(map.size(), sign.get(),
901  local.get(), map.get(), global.get());
902  }

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 1438 of file AssemblyMap.cpp.

1440  {
1442  = m_session->GetComm()->GetRowComm();
1443  bool isRoot = vRowComm->GetRank() == 0;
1444  int n = vRowComm->GetSize();
1445  int i;
1446 
1447  // Determine number of global degrees of freedom.
1448  int globBndCnt = 0, globDirCnt = 0;
1449 
1450  for (i = 0; i < m_numGlobalBndCoeffs; ++i)
1451  {
1453  {
1454  globBndCnt++;
1455 
1456  if (i < m_numGlobalDirBndCoeffs)
1457  {
1458  globDirCnt++;
1459  }
1460  }
1461  }
1462 
1463  int globCnt = m_numGlobalCoeffs - m_numGlobalBndCoeffs + globBndCnt;
1464 
1465  // Calculate maximum valency
1466  Array<OneD, NekDouble> tmpLoc (m_numLocalBndCoeffs, 1.0);
1467  Array<OneD, NekDouble> tmpGlob(m_numGlobalBndCoeffs, 0.0);
1468 
1469  Vmath::Assmb(
1470  m_numLocalBndCoeffs, tmpLoc.get(),
1471  m_localToGlobalBndMap.get(), tmpGlob.get());
1472  UniversalAssembleBnd(tmpGlob);
1473 
1474  int totGlobDof = globCnt;
1475  int totGlobBndDof = globBndCnt;
1476  int totGlobDirDof = globDirCnt;
1477  int totLocalDof = m_numLocalCoeffs;
1478  int totLocalBndDof = m_numLocalBndCoeffs;
1479  int totLocalDirDof = m_numLocalDirBndCoeffs;
1480 
1481  int meanValence = 0;
1482  int maxValence = 0;
1483  int minValence = 10000000;
1484  for (int i = 0; i < m_numGlobalBndCoeffs; ++i)
1485  {
1487  {
1488  continue;
1489  }
1490 
1491  if (tmpGlob[i] > maxValence)
1492  {
1493  maxValence = tmpGlob[i];
1494  }
1495  if (tmpGlob[i] < minValence)
1496  {
1497  minValence = tmpGlob[i];
1498  }
1499  meanValence += tmpGlob[i];
1500  }
1501 
1502  vRowComm->AllReduce(maxValence, LibUtilities::ReduceMax);
1503  vRowComm->AllReduce(minValence, LibUtilities::ReduceMin);
1504  vRowComm->AllReduce(meanValence, LibUtilities::ReduceSum);
1505  vRowComm->AllReduce(totGlobDof, LibUtilities::ReduceSum);
1506  vRowComm->AllReduce(totGlobBndDof, LibUtilities::ReduceSum);
1507  vRowComm->AllReduce(totGlobDirDof, LibUtilities::ReduceSum);
1508  vRowComm->AllReduce(totLocalDof, LibUtilities::ReduceSum);
1509  vRowComm->AllReduce(totLocalBndDof, LibUtilities::ReduceSum);
1510  vRowComm->AllReduce(totLocalDirDof, LibUtilities::ReduceSum);
1511 
1512  meanValence /= totGlobBndDof;
1513 
1514  if (isRoot)
1515  {
1516  if (printHeader)
1517  {
1518  out << "Assembly map statistics for field " << variable
1519  << ":" << endl;
1520  }
1521 
1522  out << " - Number of local/global dof : "
1523  << totLocalDof << " " << totGlobDof << endl;
1524  out << " - Number of local/global boundary dof : "
1525  << totLocalBndDof << " " << totGlobBndDof << endl;
1526  out << " - Number of local/global Dirichlet dof : "
1527  << totLocalDirDof << " " << totGlobDirDof << endl;
1528  out << " - dof valency (min/max/mean) : "
1529  << minValence << " " << maxValence << " " << meanValence
1530  << endl;
1531 
1532  if (n > 1)
1533  {
1534  NekDouble mean = m_numLocalCoeffs, mean2 = mean * mean;
1535  NekDouble minval = mean, maxval = mean;
1536  Array<OneD, NekDouble> tmp(1);
1537 
1538  for (i = 1; i < n; ++i)
1539  {
1540  vRowComm->Recv(i, tmp);
1541  mean += tmp[0];
1542  mean2 += tmp[0]*tmp[0];
1543 
1544  if (tmp[0] > maxval)
1545  {
1546  maxval = tmp[0];
1547  }
1548  if (tmp[0] < minval)
1549  {
1550  minval = tmp[0];
1551  }
1552  }
1553 
1554  if (maxval > 0.1)
1555  {
1556  out << " - Local dof dist. (min/max/mean/dev) : "
1557  << minval << " " << maxval << " " << (mean / n)
1558  << " " << sqrt(mean2/n - mean*mean/n/n) << endl;
1559  }
1560 
1561  vRowComm->Block();
1562 
1563  mean = minval = maxval = m_numLocalBndCoeffs;
1564  mean2 = mean * mean;
1565 
1566  for (i = 1; i < n; ++i)
1567  {
1568  vRowComm->Recv(i, tmp);
1569  mean += tmp[0];
1570  mean2 += tmp[0]*tmp[0];
1571 
1572  if (tmp[0] > maxval)
1573  {
1574  maxval = tmp[0];
1575  }
1576  if (tmp[0] < minval)
1577  {
1578  minval = tmp[0];
1579  }
1580  }
1581 
1582  out << " - Local bnd dof dist. (min/max/mean/dev) : "
1583  << minval << " " << maxval << " " << (mean / n) << " "
1584  << sqrt(mean2/n - mean*mean/n/n) << endl;
1585  }
1586  }
1587  else
1588  {
1589  Array<OneD, NekDouble> tmp(1);
1590  tmp[0] = m_numLocalCoeffs;
1591  vRowComm->Send(0, tmp);
1592  vRowComm->Block();
1593  tmp[0] = m_numLocalBndCoeffs;
1594  vRowComm->Send(0, tmp);
1595  }
1596 
1597  // Either we have no more levels in the static condensation, or we
1598  // are not multi-level.
1600  {
1601  return;
1602  }
1603 
1604  int level = 2;
1606  while (tmp->m_nextLevelLocalToGlobalMap)
1607  {
1608  tmp = tmp->m_nextLevelLocalToGlobalMap;
1609  ++level;
1610  }
1611 
1612  // Print out multi-level static condensation information.
1613  if (n > 1)
1614  {
1615  if (isRoot)
1616  {
1617  NekDouble mean = level, mean2 = mean * mean;
1618  int minval = level, maxval = level;
1619 
1620  Array<OneD, NekDouble> tmpRecv(1);
1621  for (i = 1; i < n; ++i)
1622  {
1623  vRowComm->Recv(i, tmpRecv);
1624  mean += tmpRecv[0];
1625  mean2 += tmpRecv[0]*tmpRecv[0];
1626 
1627  if (tmpRecv[0] > maxval)
1628  {
1629  maxval = (int)(tmpRecv[0] + 0.5);
1630  }
1631  if (tmpRecv[0] < minval)
1632  {
1633  minval = (int)(tmpRecv[0] + 0.5);
1634  }
1635  }
1636 
1637  out << " - M-level sc. dist. (min/max/mean/dev) : "
1638  << minval << " " << maxval << " " << (mean / n) << " "
1639  << sqrt(mean2/n - mean*mean/n/n) << endl;
1640  }
1641  else
1642  {
1643  Array<OneD, NekDouble> tmpSend(1);
1644  tmpSend[0] = level;
1645  vRowComm->Send(0, tmpSend);
1646  }
1647  }
1648  else
1649  {
1650  out << " - Number of static cond. levels : "
1651  << level << endl;
1652  }
1653 
1654  if (isRoot)
1655  {
1656  out << "Stats at lowest static cond. level:" << endl;
1657  }
1658  tmp->PrintStats(out, variable, false);
1659  }
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:54
std::shared_ptr< AssemblyMap > AssemblyMapSharedPtr
Definition: AssemblyMap.h:52
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:267

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 1346 of file AssemblyMap.cpp.

1347  {
1348  Gs::Gather(bndvals, Gs::gs_amax, m_dirBndGsh);
1349  }
Gs::gs_data * m_dirBndGsh
gs gather communication to impose Dirhichlet BCs.
Definition: AssemblyMap.h:417
@ gs_amax
Definition: GsLib.hpp:53

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 874 of file AssemblyMap.cpp.

877  {
878  v_UniversalAssemble(pGlobal, offset);
879  }

References v_UniversalAssemble().

◆ UniversalAssemble() [3/3]

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

Definition at line 868 of file AssemblyMap.cpp.

870  {
871  v_UniversalAssemble(pGlobal);
872  }

References v_UniversalAssemble().

◆ UniversalAssembleBnd() [1/3]

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

Definition at line 1322 of file AssemblyMap.cpp.

1324  {
1325  ASSERTL1(pGlobal.size() >= m_numGlobalBndCoeffs,
1326  "Wrong size.");
1327  Gs::Gather(pGlobal, Gs::gs_add, m_bndGsh);
1328  }
@ gs_add
Definition: GsLib.hpp:53

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 1336 of file AssemblyMap.cpp.

1339  {
1340  Array<OneD, NekDouble> tmp(offset);
1341  if (offset > 0) Vmath::Vcopy(offset, pGlobal, 1, tmp, 1);
1342  UniversalAssembleBnd(pGlobal);
1343  if (offset > 0) Vmath::Vcopy(offset, tmp, 1, pGlobal, 1);
1344  }

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

◆ UniversalAssembleBnd() [3/3]

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

Definition at line 1330 of file AssemblyMap.cpp.

1332  {
1333  UniversalAssembleBnd(pGlobal.GetPtr());
1334  }

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
privatevirtual

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

Definition at line 650 of file AssemblyMap.cpp.

653  {
654  boost::ignore_unused(loc, global);
656  "Not defined for this type of mapping.");
657  }
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:209

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

Referenced by Assemble().

◆ v_Assemble() [2/2]

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

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

Definition at line 659 of file AssemblyMap.cpp.

662  {
663  boost::ignore_unused(loc, global);
665  "Not defined for this type of mapping.");
666  }

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

◆ v_GetExtraDirEdges()

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

Reimplemented in Nektar::MultiRegions::AssemblyMapCG.

Definition at line 749 of file AssemblyMap.cpp.

750  {
752  "Not defined for this type of mapping.");
753  static Array<OneD, const int> result;
754  return result;
755  }

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

Referenced by GetExtraDirEdges().

◆ v_GetFullSystemBandWidth()

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

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

Definition at line 693 of file AssemblyMap.cpp.

694  {
696  "Not defined for this type of mapping.");
697  return 0;
698  }

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

Referenced by GetFullSystemBandWidth().

◆ v_GetGlobalToUniversalMap() [1/2]

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

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

Definition at line 579 of file AssemblyMap.cpp.

580  {
582  "Not defined for this type of mapping.");
583  static Array<OneD, const int> result;
584  return result;
585  }

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

Referenced by GetGlobalToUniversalMap().

◆ v_GetGlobalToUniversalMap() [2/2]

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

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

Definition at line 555 of file AssemblyMap.cpp.

556  {
557  boost::ignore_unused(i);
559  "Not defined for this type of mapping.");
560  return 0;
561  }

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

◆ v_GetGlobalToUniversalMapUnique() [1/2]

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

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

Definition at line 587 of file AssemblyMap.cpp.

588  {
590  "Not defined for this type of mapping.");
591  static Array<OneD, const int> result;
592  return result;
593  }

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

Referenced by GetGlobalToUniversalMapUnique().

◆ v_GetGlobalToUniversalMapUnique() [2/2]

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

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

Definition at line 563 of file AssemblyMap.cpp.

564  {
565  boost::ignore_unused(i);
567  "Not defined for this type of mapping.");
568  return 0;
569  }

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

◆ v_GetLocalToGlobalMap() [1/2]

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

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

Definition at line 571 of file AssemblyMap.cpp.

572  {
574  "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 GetLocalToGlobalMap().

◆ v_GetLocalToGlobalMap() [2/2]

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

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

Definition at line 547 of file AssemblyMap.cpp.

548  {
549  boost::ignore_unused(i);
551  "Not defined for this type of mapping.");
552  return 0;
553  }

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

◆ v_GetLocalToGlobalSign() [1/2]

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

Reimplemented in Nektar::MultiRegions::AssemblyMapCG.

Definition at line 603 of file AssemblyMap.cpp.

604  {
606  "Not defined for this type of mapping.");
607  static Array<OneD, NekDouble> result;
608  return result;
609  }

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

Referenced by GetLocalToGlobalSign().

◆ v_GetLocalToGlobalSign() [2/2]

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

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

Definition at line 595 of file AssemblyMap.cpp.

596  {
597  boost::ignore_unused(i);
599  "Not defined for this type of mapping.");
600  return 0.0;
601  }

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

◆ v_GetNumDirEdges()

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

Reimplemented in Nektar::MultiRegions::AssemblyMapCG.

Definition at line 721 of file AssemblyMap.cpp.

722  {
724  "Not defined for this type of mapping.");
725  return 0;
726  }

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

Referenced by GetNumDirEdges().

◆ v_GetNumDirFaces()

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

Reimplemented in Nektar::MultiRegions::AssemblyMapCG.

Definition at line 728 of file AssemblyMap.cpp.

729  {
731  "Not defined for this type of mapping.");
732  return 0;
733  }

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

Referenced by GetNumDirFaces().

◆ v_GetNumNonDirEdgeModes()

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

Reimplemented in Nektar::MultiRegions::AssemblyMapCG.

Definition at line 707 of file AssemblyMap.cpp.

708  {
710  "Not defined for this type of mapping.");
711  return 0;
712  }

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

Referenced by GetNumNonDirEdgeModes().

◆ v_GetNumNonDirEdges()

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

Reimplemented in Nektar::MultiRegions::AssemblyMapCG.

Definition at line 735 of file AssemblyMap.cpp.

736  {
738  "Not defined for this type of mapping.");
739  return 0;
740  }

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

Referenced by GetNumNonDirEdges().

◆ v_GetNumNonDirFaceModes()

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

Reimplemented in Nektar::MultiRegions::AssemblyMapCG.

Definition at line 714 of file AssemblyMap.cpp.

715  {
717  "Not defined for this type of mapping.");
718  return 0;
719  }

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

Referenced by GetNumNonDirFaceModes().

◆ v_GetNumNonDirFaces()

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

Reimplemented in Nektar::MultiRegions::AssemblyMapCG.

Definition at line 742 of file AssemblyMap.cpp.

743  {
745  "Not defined for this type of mapping.");
746  return 0;
747  }

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

Referenced by GetNumNonDirFaces().

◆ v_GetNumNonDirVertexModes()

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

Reimplemented in Nektar::MultiRegions::AssemblyMapCG.

Definition at line 700 of file AssemblyMap.cpp.

701  {
703  "Not defined for this type of mapping.");
704  return 0;
705  }

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
privatevirtual

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

Definition at line 632 of file AssemblyMap.cpp.

635  {
636  boost::ignore_unused(loc, global);
638  "Not defined for this type of mapping.");
639  }

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

Referenced by GlobalToLocal().

◆ v_GlobalToLocal() [2/2]

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

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

Definition at line 641 of file AssemblyMap.cpp.

644  {
645  boost::ignore_unused(loc, global);
647  "Not defined for this type of mapping.");
648  }

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

◆ v_LinearSpaceMap()

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

Generate a linear space mapping from existing mapping.

Reimplemented in Nektar::MultiRegions::AssemblyMapCG.

Definition at line 757 of file AssemblyMap.cpp.

759  {
760  boost::ignore_unused(locexp, solnType);
762  "Not defined for this type of mapping.");
763  static std::shared_ptr<AssemblyMap> result;
764  return result;
765  }

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
privatevirtual

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

Definition at line 611 of file AssemblyMap.cpp.

615  {
616  boost::ignore_unused(loc, global, useComm);
618  "Not defined for this type of mapping.");
619  }

References Nektar::ErrorUtil::efatal, CG_Iterations::loc, 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
privatevirtual

Reimplemented in Nektar::MultiRegions::AssemblyMapCG.

Definition at line 621 of file AssemblyMap.cpp.

625  {
626  boost::ignore_unused(loc, global, useComm);
628  "Not defined for this type of mapping.");
629  }

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

◆ v_UniversalAssemble() [1/3]

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

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

Definition at line 668 of file AssemblyMap.cpp.

670  {
671  boost::ignore_unused(pGlobal);
672  // Do nothing here since multi-level static condensation uses a
673  // AssemblyMap and thus will call this routine in serial.
674  }

Referenced by UniversalAssemble().

◆ v_UniversalAssemble() [2/3]

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

Reimplemented in Nektar::MultiRegions::AssemblyMapCG.

Definition at line 684 of file AssemblyMap.cpp.

687  {
688  boost::ignore_unused(pGlobal, offset);
689  // Do nothing here since multi-level static condensation uses a
690  // AssemblyMap and thus will call this routine in serial.
691  }

◆ v_UniversalAssemble() [3/3]

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

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

Definition at line 676 of file AssemblyMap.cpp.

678  {
679  boost::ignore_unused(pGlobal);
680  // Do nothing here since multi-level static condensation uses a
681  // AssemblyMap and thus will call this routine in serial.
682  }

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 382 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 384 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 386 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 388 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 397 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 417 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_iterativeTolerance

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

Tolerance for iterative solver.

Definition at line 406 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 412 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_maxIterations

int Nektar::MultiRegions::AssemblyMap::m_maxIterations
protected

Maximum iterations for iterative solver.

Definition at line 403 of file AssemblyMap.h.

Referenced by AssemblyMap(), and GetMaxIterations().

◆ 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 432 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 368 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 357 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 445 of file AssemblyMap.h.

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

◆ m_preconType

PreconditionerType Nektar::MultiRegions::AssemblyMap::m_preconType
protected

Type type of preconditioner to use in iterative solver.

Definition at line 400 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 409 of file AssemblyMap.h.

Referenced by AssemblyMap(), and GetSuccessiveRHS().

◆ m_systemSingular

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