Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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:
Inheritance graph
[legend]
Collaboration diagram for Nektar::MultiRegions::AssemblyMap:
Collaboration graph
[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) const
 
void LocalToGlobal (const NekVector< NekDouble > &loc, NekVector< NekDouble > &global) 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
 
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...
 
int GetBndCondCoeffsToGlobalCoeffsMap (const int i)
 Retrieves the global index corresponding to a boundary expansion mode. More...
 
const Array< OneD, const int > & GetBndCondCoeffsToGlobalCoeffsMap ()
 Retrieves the global indices corresponding to the boundary expansion modes. More...
 
NekDouble GetBndCondCoeffsToGlobalCoeffsSign (const int i)
 Returns the modal sign associated with a given boundary expansion mode. More...
 
int GetBndCondTraceToGlobalTraceMap (const int i)
 Returns the global index of the boundary trace giving the index on the boundary expansion. More...
 
const Array< OneD, const int > & GetBndCondTraceToGlobalTraceMap ()
 
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 NekVector< NekDouble > &loc, NekVector< NekDouble > &global, int offset) const
 
void LocalBndToGlobal (const NekVector< NekDouble > &loc, NekVector< NekDouble > &global) const
 
void LocalBndToGlobal (const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global, int offset) const
 
void LocalBndToGlobal (const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global) 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) const
 
const Array< OneD, const int > & GetExtraDirEdges ()
 
boost::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
 
int GetLowestStaticCondLevel () 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 boundary coeffs to global space. More...
 
Array< OneD, NekDoublem_localToGlobalBndSign
 Integer sign of local boundary coeffs to global space. More...
 
Array< OneD, int > m_bndCondCoeffsToGlobalCoeffsMap
 Integer map of bnd cond coeffs to global coefficients. More...
 
Array< OneD, NekDoublem_bndCondCoeffsToGlobalCoeffsSign
 Integer map of bnd cond coeffs to global coefficients. More...
 
Array< OneD, int > m_bndCondTraceToGlobalTraceMap
 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...
 
Gs::gs_datam_gsh
 
Gs::gs_datam_bndGsh
 
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) const
 
virtual void v_LocalToGlobal (const NekVector< NekDouble > &loc, NekVector< NekDouble > &global) 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 boost::shared_ptr
< AssemblyMap
v_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 59 of file AssemblyMap.h.

Constructor & Destructor Documentation

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

Default constructor.

Initialises an empty mapping.

Definition at line 79 of file AssemblyMap.cpp.

79  :
80  m_session(),
81  m_comm(),
82  m_hash(0),
89  m_successiveRHS(0),
90  m_gsh(0),
91  m_bndGsh(0)
92  {
93  }
int m_numGlobalBndCoeffs
Total number of global boundary coefficients.
Definition: AssemblyMap.h:314
LibUtilities::CommSharedPtr m_comm
Communicator.
Definition: AssemblyMap.h:306
int m_successiveRHS
sucessive RHS for iterative solver
Definition: AssemblyMap.h:377
size_t m_hash
Hash for map.
Definition: AssemblyMap.h:309
int m_bndSystemBandWidth
The bandwith of the global bnd system.
Definition: AssemblyMap.h:365
GlobalSysSolnType m_solnType
The solution type of the global system.
Definition: AssemblyMap.h:363
int m_numGlobalDirBndCoeffs
Number of Global Dirichlet Boundary Coefficients.
Definition: AssemblyMap.h:318
int m_numLocalDirBndCoeffs
Number of Local Dirichlet Boundary Coefficients.
Definition: AssemblyMap.h:316
int m_numLocalBndCoeffs
Number of local boundary coefficients.
Definition: AssemblyMap.h:312
LibUtilities::SessionReaderSharedPtr m_session
Session object.
Definition: AssemblyMap.h:303
No Solution type specified.
Nektar::MultiRegions::AssemblyMap::AssemblyMap ( const LibUtilities::SessionReaderSharedPtr pSession,
const std::string  variable = "DefaultVar" 
)

Constructor with a communicator.

Definition at line 95 of file AssemblyMap.cpp.

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

97  :
98  m_session(pSession),
99  m_comm(pSession->GetComm()),
100  m_hash(0),
106  m_successiveRHS(0),
107  m_gsh(0),
108  m_bndGsh(0)
109  {
110  // Default value from Solver Info
111  m_solnType = pSession->GetSolverInfoAsEnum<GlobalSysSolnType>(
112  "GlobalSysSoln");
113  m_preconType = pSession->GetSolverInfoAsEnum<PreconditionerType>(
114  "Preconditioner");
115 
116  // Override values with data from GlobalSysSolnInfo section
117  if(pSession->DefinesGlobalSysSolnInfo(variable, "GlobalSysSoln"))
118  {
119  std::string sysSoln = pSession->GetGlobalSysSolnInfo(variable,
120  "GlobalSysSoln");
121  m_solnType = pSession->GetValueAsEnum<GlobalSysSolnType>(
122  "GlobalSysSoln", sysSoln);
123  }
124 
125  if(pSession->DefinesGlobalSysSolnInfo(variable, "Preconditioner"))
126  {
127  std::string precon = pSession->GetGlobalSysSolnInfo(variable,
128  "Preconditioner");
129  m_preconType = pSession->GetValueAsEnum<PreconditionerType>(
130  "Preconditioner", precon);
131  }
132 
133  if(pSession->DefinesGlobalSysSolnInfo(variable,
134  "IterativeSolverTolerance"))
135  {
136  m_iterativeTolerance = boost::lexical_cast<NekDouble>(
137  pSession->GetGlobalSysSolnInfo(variable,
138  "IterativeSolverTolerance").c_str());
139  }
140  else
141  {
142  pSession->LoadParameter("IterativeSolverTolerance",
145  }
146 
147 
148  if(pSession->DefinesGlobalSysSolnInfo(variable,
149  "MaxIterations"))
150  {
151  m_maxIterations = boost::lexical_cast<int>(
152  pSession->GetGlobalSysSolnInfo(variable,
153  "MaxIterations").c_str());
154  }
155  else
156  {
157  pSession->LoadParameter("MaxIterations",
159  5000);
160  }
161 
162 
163  if(pSession->DefinesGlobalSysSolnInfo(variable,"SuccessiveRHS"))
164  {
165  m_successiveRHS = boost::lexical_cast<int>(
166  pSession->GetGlobalSysSolnInfo(variable,
167  "SuccessiveRHS").c_str());
168  }
169  else
170  {
171  pSession->LoadParameter("SuccessiveRHS",
172  m_successiveRHS,0);
173  }
174 
175  }
PreconditionerType m_preconType
Type type of preconditioner to use in iterative solver.
Definition: AssemblyMap.h:368
int m_numGlobalBndCoeffs
Total number of global boundary coefficients.
Definition: AssemblyMap.h:314
LibUtilities::CommSharedPtr m_comm
Communicator.
Definition: AssemblyMap.h:306
static const NekDouble kNekIterativeTol
int m_successiveRHS
sucessive RHS for iterative solver
Definition: AssemblyMap.h:377
size_t m_hash
Hash for map.
Definition: AssemblyMap.h:309
int m_bndSystemBandWidth
The bandwith of the global bnd system.
Definition: AssemblyMap.h:365
NekDouble m_iterativeTolerance
Tolerance for iterative solver.
Definition: AssemblyMap.h:374
GlobalSysSolnType m_solnType
The solution type of the global system.
Definition: AssemblyMap.h:363
int m_numGlobalDirBndCoeffs
Number of Global Dirichlet Boundary Coefficients.
Definition: AssemblyMap.h:318
double NekDouble
int m_numLocalDirBndCoeffs
Number of Local Dirichlet Boundary Coefficients.
Definition: AssemblyMap.h:316
int m_numLocalBndCoeffs
Number of local boundary coefficients.
Definition: AssemblyMap.h:312
int m_maxIterations
Maximum iterations for iterative solver.
Definition: AssemblyMap.h:371
LibUtilities::SessionReaderSharedPtr m_session
Session object.
Definition: AssemblyMap.h:303
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 181 of file AssemblyMap.cpp.

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_nextLevelLocalToGlobalMap, m_numGlobalBndCoeffs, m_numGlobalCoeffs, m_numGlobalDirBndCoeffs, m_numLocalBndCoeffs, m_numLocalBndCoeffsPerPatch, m_numLocalDirBndCoeffs, m_numLocalIntCoeffsPerPatch, m_numPatches, m_patchMapFromPrevLevel, m_signChange, m_solnType, m_staticCondLevel, Nektar::MultiRegions::RoundNekDoubleToInt(), and sign.

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

Destructor.

Definition at line 437 of file AssemblyMap.cpp.

438  {
439  }

Member Function Documentation

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

Definition at line 750 of file AssemblyMap.cpp.

References v_Assemble().

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

753  {
754  v_Assemble(loc,global);
755  }
virtual void v_Assemble(const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global) const
void Nektar::MultiRegions::AssemblyMap::Assemble ( const NekVector< NekDouble > &  loc,
NekVector< NekDouble > &  global 
) const

Definition at line 757 of file AssemblyMap.cpp.

References v_Assemble().

760  {
761  v_Assemble(loc,global);
762  }
virtual void v_Assemble(const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global) const
void Nektar::MultiRegions::AssemblyMap::AssembleBnd ( const NekVector< NekDouble > &  loc,
NekVector< NekDouble > &  global,
int  offset 
) const

Definition at line 1070 of file AssemblyMap.cpp.

References Nektar::NekVector< DataType >::GetPtr().

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

1073  {
1074  AssembleBnd(loc.GetPtr(), global.GetPtr(), offset);
1075  }
void AssembleBnd(const NekVector< NekDouble > &loc, NekVector< NekDouble > &global, int offset) const
void Nektar::MultiRegions::AssemblyMap::AssembleBnd ( const NekVector< NekDouble > &  loc,
NekVector< NekDouble > &  global 
) const

Definition at line 1078 of file AssemblyMap.cpp.

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

1081  {
1082  AssembleBnd(loc.GetPtr(), global.GetPtr());
1083  }
void AssembleBnd(const NekVector< NekDouble > &loc, NekVector< NekDouble > &global, int offset) const
void Nektar::MultiRegions::AssemblyMap::AssembleBnd ( const Array< OneD, const NekDouble > &  loc,
Array< OneD, NekDouble > &  global,
int  offset 
) const

Definition at line 1086 of file AssemblyMap.cpp.

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

1089  {
1090  ASSERTL1(loc.num_elements() >= m_numLocalBndCoeffs,"Local array is not of correct dimension");
1091  ASSERTL1(global.num_elements() >= m_numGlobalBndCoeffs-offset,"Global array is not of correct dimension");
1092  Array<OneD,NekDouble> tmp(m_numGlobalBndCoeffs,0.0);
1093 
1094  if(m_signChange)
1095  {
1097  }
1098  else
1099  {
1100  Vmath::Assmb(m_numLocalBndCoeffs,loc.get(), m_localToGlobalBndMap.get(), tmp.get());
1101  }
1102  UniversalAssembleBnd(tmp);
1103  Vmath::Vcopy(m_numGlobalBndCoeffs-offset, tmp.get() + offset, 1, global.get(), 1);
1104  }
bool m_signChange
Flag indicating if modes require sign reversal.
Definition: AssemblyMap.h:345
int m_numGlobalBndCoeffs
Total number of global boundary coefficients.
Definition: AssemblyMap.h:314
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:695
Array< OneD, int > m_localToGlobalBndMap
Integer map of local boundary coeffs to global space.
Definition: AssemblyMap.h:348
int m_numLocalBndCoeffs
Number of local boundary coefficients.
Definition: AssemblyMap.h:312
Array< OneD, NekDouble > m_localToGlobalBndSign
Integer sign of local boundary coeffs to global space.
Definition: AssemblyMap.h:350
void UniversalAssembleBnd(Array< OneD, NekDouble > &pGlobal) const
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
void Nektar::MultiRegions::AssemblyMap::AssembleBnd ( const Array< OneD, const NekDouble > &  loc,
Array< OneD, NekDouble > &  global 
) const

Definition at line 1107 of file AssemblyMap.cpp.

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

1110  {
1111  ASSERTL1(loc.num_elements() >= m_numLocalBndCoeffs,"Local vector is not of correct dimension");
1112  ASSERTL1(global.num_elements() >= m_numGlobalBndCoeffs,"Global vector is not of correct dimension");
1113 
1114  Vmath::Zero(m_numGlobalBndCoeffs, global.get(), 1);
1115 
1116  if(m_signChange)
1117  {
1119  loc.get(), m_localToGlobalBndMap.get(), global.get());
1120  }
1121  else
1122  {
1123  Vmath::Assmb(m_numLocalBndCoeffs,loc.get(), m_localToGlobalBndMap.get(), global.get());
1124  }
1125  UniversalAssembleBnd(global);
1126  }
bool m_signChange
Flag indicating if modes require sign reversal.
Definition: AssemblyMap.h:345
int m_numGlobalBndCoeffs
Total number of global boundary coefficients.
Definition: AssemblyMap.h:314
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:695
Array< OneD, int > m_localToGlobalBndMap
Integer map of local boundary coeffs to global space.
Definition: AssemblyMap.h:348
int m_numLocalBndCoeffs
Number of local boundary coefficients.
Definition: AssemblyMap.h:312
Array< OneD, NekDouble > m_localToGlobalBndSign
Integer sign of local boundary coeffs to global space.
Definition: AssemblyMap.h:350
void UniversalAssembleBnd(Array< OneD, NekDouble > &pGlobal) const
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:359
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191
bool Nektar::MultiRegions::AssemblyMap::AtLastLevel ( ) const

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

Definition at line 1193 of file AssemblyMap.cpp.

References m_nextLevelLocalToGlobalMap.

1194  {
1195  return !( m_nextLevelLocalToGlobalMap.get() );
1196  }
AssemblyMapSharedPtr m_nextLevelLocalToGlobalMap
Map from the patches of the previous level to the patches of the current level.
Definition: AssemblyMap.h:395
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 456 of file AssemblyMap.cpp.

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

457  {
458  int i,j;
459  int cnt = 0;
460  int locSize;
461  int maxId;
462  int minId;
463  int bwidth = -1;
464  for(i = 0; i < m_numPatches; ++i)
465  {
466  locSize = m_numLocalBndCoeffsPerPatch[i];
467  maxId = -1;
468  minId = m_numLocalBndCoeffs+1;
469  for(j = 0; j < locSize; j++)
470  {
472  {
473  if(m_localToGlobalBndMap[cnt+j] > maxId)
474  {
475  maxId = m_localToGlobalBndMap[cnt+j];
476  }
477 
478  if(m_localToGlobalBndMap[cnt+j] < minId)
479  {
480  minId = m_localToGlobalBndMap[cnt+j];
481  }
482  }
483  }
484  bwidth = (bwidth>(maxId-minId))?bwidth:(maxId-minId);
485 
486  cnt+=locSize;
487  }
488 
489  m_bndSystemBandWidth = bwidth;
490  }
int m_bndSystemBandWidth
The bandwith of the global bnd system.
Definition: AssemblyMap.h:365
Array< OneD, unsigned int > m_numLocalBndCoeffsPerPatch
The number of bnd dofs per patch.
Definition: AssemblyMap.h:388
int m_numGlobalDirBndCoeffs
Number of Global Dirichlet Boundary Coefficients.
Definition: AssemblyMap.h:318
Array< OneD, int > m_localToGlobalBndMap
Integer map of local boundary coeffs to global space.
Definition: AssemblyMap.h:348
int m_numLocalBndCoeffs
Number of local boundary coefficients.
Definition: AssemblyMap.h:312
int m_numPatches
The number of patches (~elements) in the current level.
Definition: AssemblyMap.h:386
int Nektar::MultiRegions::AssemblyMap::GetBndCondCoeffsToGlobalCoeffsMap ( const int  i)

Retrieves the global index corresponding to a boundary expansion mode.

Definition at line 879 of file AssemblyMap.cpp.

References m_bndCondCoeffsToGlobalCoeffsMap.

881  {
883  }
Array< OneD, int > m_bndCondCoeffsToGlobalCoeffsMap
Integer map of bnd cond coeffs to global coefficients.
Definition: AssemblyMap.h:352
const Array< OneD, const int > & Nektar::MultiRegions::AssemblyMap::GetBndCondCoeffsToGlobalCoeffsMap ( )

Retrieves the global indices corresponding to the boundary expansion modes.

Definition at line 913 of file AssemblyMap.cpp.

References m_bndCondCoeffsToGlobalCoeffsMap.

Referenced by Nektar::CoupledAssemblyMap::CoupledAssemblyMap().

914  {
916  }
Array< OneD, int > m_bndCondCoeffsToGlobalCoeffsMap
Integer map of bnd cond coeffs to global coefficients.
Definition: AssemblyMap.h:352
NekDouble Nektar::MultiRegions::AssemblyMap::GetBndCondCoeffsToGlobalCoeffsSign ( const int  i)

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

Definition at line 900 of file AssemblyMap.cpp.

References m_bndCondCoeffsToGlobalCoeffsSign, and m_signChange.

901  {
902  if(m_signChange)
903  {
905  }
906  else
907  {
908  return 1.0;
909  }
910  }
bool m_signChange
Flag indicating if modes require sign reversal.
Definition: AssemblyMap.h:345
Array< OneD, NekDouble > m_bndCondCoeffsToGlobalCoeffsSign
Integer map of bnd cond coeffs to global coefficients.
Definition: AssemblyMap.h:354
int Nektar::MultiRegions::AssemblyMap::GetBndCondTraceToGlobalTraceMap ( const int  i)

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

Definition at line 886 of file AssemblyMap.cpp.

References ASSERTL1, and m_bndCondTraceToGlobalTraceMap.

888  {
889  ASSERTL1(i < m_bndCondTraceToGlobalTraceMap.num_elements(),
890  "Index out of range.");
892  }
Array< OneD, int > m_bndCondTraceToGlobalTraceMap
Integer map of bnd cond trace number to global trace number.
Definition: AssemblyMap.h:356
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191
const Array< OneD, const int > & Nektar::MultiRegions::AssemblyMap::GetBndCondTraceToGlobalTraceMap ( )

Definition at line 895 of file AssemblyMap.cpp.

896  {
898  }
Array< OneD, int > m_bndCondTraceToGlobalTraceMap
Integer map of bnd cond trace number to global trace number.
Definition: AssemblyMap.h:356
int Nektar::MultiRegions::AssemblyMap::GetBndSystemBandWidth ( ) const

Returns the bandwidth of the boundary system.

Definition at line 1152 of file AssemblyMap.cpp.

References m_bndSystemBandWidth.

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

1153  {
1154  return m_bndSystemBandWidth;
1155  }
int m_bndSystemBandWidth
The bandwith of the global bnd system.
Definition: AssemblyMap.h:365
LibUtilities::CommSharedPtr Nektar::MultiRegions::AssemblyMap::GetComm ( )

Retrieves the communicator.

Definition at line 672 of file AssemblyMap.cpp.

References m_comm.

673  {
674  return m_comm;
675  }
LibUtilities::CommSharedPtr m_comm
Communicator.
Definition: AssemblyMap.h:306
const Array< OneD, const int > & Nektar::MultiRegions::AssemblyMap::GetExtraDirEdges ( )

Definition at line 823 of file AssemblyMap.cpp.

References v_GetExtraDirEdges().

824  {
825  return v_GetExtraDirEdges();
826  }
virtual const Array< OneD, const int > & v_GetExtraDirEdges()
int Nektar::MultiRegions::AssemblyMap::GetFullSystemBandWidth ( ) const

Definition at line 783 of file AssemblyMap.cpp.

References v_GetFullSystemBandWidth().

784  {
785  return v_GetFullSystemBandWidth();
786  }
virtual int v_GetFullSystemBandWidth() const
GlobalSysSolnType Nektar::MultiRegions::AssemblyMap::GetGlobalSysSolnType ( ) const

Returns the method of solving global systems.

Definition at line 1199 of file AssemblyMap.cpp.

References m_solnType.

Referenced by AssemblyMap().

1200  {
1201  return m_solnType;
1202  }
GlobalSysSolnType m_solnType
The solution type of the global system.
Definition: AssemblyMap.h:363
const Array< OneD, const int > & Nektar::MultiRegions::AssemblyMap::GetGlobalToUniversalBndMap ( )

Definition at line 856 of file AssemblyMap.cpp.

References m_globalToUniversalBndMap.

Referenced by AssemblyMap().

857  {
859  }
Array< OneD, int > m_globalToUniversalBndMap
Integer map of process coeffs to universal space.
Definition: AssemblyMap.h:358
const Array< OneD, const int > & Nektar::MultiRegions::AssemblyMap::GetGlobalToUniversalBndMapUnique ( )

Definition at line 861 of file AssemblyMap.cpp.

References m_globalToUniversalBndMapUnique.

Referenced by AssemblyMap().

862  {
864  }
Array< OneD, int > m_globalToUniversalBndMapUnique
Integer map of unique process coeffs to universal space (signed)
Definition: AssemblyMap.h:360
int Nektar::MultiRegions::AssemblyMap::GetGlobalToUniversalMap ( const int  i) const

Definition at line 687 of file AssemblyMap.cpp.

References v_GetGlobalToUniversalMap().

688  {
689  return v_GetGlobalToUniversalMap(i);
690  }
virtual const Array< OneD, const int > & v_GetGlobalToUniversalMap()
const Array< OneD, const int > & Nektar::MultiRegions::AssemblyMap::GetGlobalToUniversalMap ( )

Definition at line 702 of file AssemblyMap.cpp.

References v_GetGlobalToUniversalMap().

703  {
704  return v_GetGlobalToUniversalMap();
705  }
virtual const Array< OneD, const int > & v_GetGlobalToUniversalMap()
int Nektar::MultiRegions::AssemblyMap::GetGlobalToUniversalMapUnique ( const int  i) const

Definition at line 692 of file AssemblyMap.cpp.

References v_GetGlobalToUniversalMapUnique().

693  {
695  }
virtual const Array< OneD, const int > & v_GetGlobalToUniversalMapUnique()
const Array< OneD, const int > & Nektar::MultiRegions::AssemblyMap::GetGlobalToUniversalMapUnique ( )

Definition at line 707 of file AssemblyMap.cpp.

References v_GetGlobalToUniversalMapUnique().

708  {
710  }
virtual const Array< OneD, const int > & v_GetGlobalToUniversalMapUnique()
size_t Nektar::MultiRegions::AssemblyMap::GetHash ( ) const

Retrieves the hash of this map.

Definition at line 677 of file AssemblyMap.cpp.

References m_hash.

678  {
679  return m_hash;
680  }
size_t m_hash
Hash for map.
Definition: AssemblyMap.h:309
NekDouble Nektar::MultiRegions::AssemblyMap::GetIterativeTolerance ( ) const

Definition at line 1209 of file AssemblyMap.cpp.

References m_iterativeTolerance.

1210  {
1211  return m_iterativeTolerance;
1212  }
NekDouble m_iterativeTolerance
Tolerance for iterative solver.
Definition: AssemblyMap.h:374
int Nektar::MultiRegions::AssemblyMap::GetLocalToGlobalBndMap ( const int  i) const

Retrieve the global index of a given local boundary mode.

Definition at line 833 of file AssemblyMap.cpp.

References m_localToGlobalBndMap.

Referenced by AssemblyMap().

834  {
835  return m_localToGlobalBndMap[i];
836  }
Array< OneD, int > m_localToGlobalBndMap
Integer map of local boundary coeffs to global space.
Definition: AssemblyMap.h:348
const Array< OneD, const int > & Nektar::MultiRegions::AssemblyMap::GetLocalToGlobalBndMap ( void  )

Retrieve the global indices of the local boundary modes.

Definition at line 839 of file AssemblyMap.cpp.

References m_localToGlobalBndMap.

840  {
841  return m_localToGlobalBndMap;
842  }
Array< OneD, int > m_localToGlobalBndMap
Integer map of local boundary coeffs to global space.
Definition: AssemblyMap.h:348
NekDouble Nektar::MultiRegions::AssemblyMap::GetLocalToGlobalBndSign ( const int  i) const

Retrieve the sign change of a given local boundary mode.

Definition at line 866 of file AssemblyMap.cpp.

References m_localToGlobalBndSign, and m_signChange.

Referenced by AssemblyMap().

867  {
868  if(m_signChange)
869  {
870  return m_localToGlobalBndSign[i];
871  }
872  else
873  {
874  return 1.0;
875  }
876  }
bool m_signChange
Flag indicating if modes require sign reversal.
Definition: AssemblyMap.h:345
Array< OneD, NekDouble > m_localToGlobalBndSign
Integer sign of local boundary coeffs to global space.
Definition: AssemblyMap.h:350
Array< OneD, const NekDouble > Nektar::MultiRegions::AssemblyMap::GetLocalToGlobalBndSign ( void  ) const

Retrieve the sign change for all local boundary modes.

Definition at line 851 of file AssemblyMap.cpp.

References m_localToGlobalBndSign.

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

852  {
853  return m_localToGlobalBndSign;
854  }
Array< OneD, NekDouble > m_localToGlobalBndSign
Integer sign of local boundary coeffs to global space.
Definition: AssemblyMap.h:350
int Nektar::MultiRegions::AssemblyMap::GetLocalToGlobalMap ( const int  i) const

Definition at line 682 of file AssemblyMap.cpp.

References v_GetLocalToGlobalMap().

683  {
684  return v_GetLocalToGlobalMap(i);
685  }
virtual const Array< OneD, const int > & v_GetLocalToGlobalMap()
const Array< OneD, const int > & Nektar::MultiRegions::AssemblyMap::GetLocalToGlobalMap ( )

Definition at line 697 of file AssemblyMap.cpp.

References v_GetLocalToGlobalMap().

698  {
699  return v_GetLocalToGlobalMap();
700  }
virtual const Array< OneD, const int > & v_GetLocalToGlobalMap()
NekDouble Nektar::MultiRegions::AssemblyMap::GetLocalToGlobalSign ( const int  i) const

Definition at line 712 of file AssemblyMap.cpp.

References v_GetLocalToGlobalSign().

713  {
714  return v_GetLocalToGlobalSign(i);
715  }
virtual const Array< OneD, NekDouble > & v_GetLocalToGlobalSign() const
const Array< OneD, NekDouble > & Nektar::MultiRegions::AssemblyMap::GetLocalToGlobalSign ( ) const

Definition at line 717 of file AssemblyMap.cpp.

References v_GetLocalToGlobalSign().

718  {
719  return v_GetLocalToGlobalSign();
720  }
virtual const Array< OneD, NekDouble > & v_GetLocalToGlobalSign() const
int Nektar::MultiRegions::AssemblyMap::GetLowestStaticCondLevel ( ) const
inline

Definition at line 296 of file AssemblyMap.h.

References m_lowestStaticCondLevel.

297  {
299  }
int m_lowestStaticCondLevel
Lowest static condensation level.
Definition: AssemblyMap.h:397
int Nektar::MultiRegions::AssemblyMap::GetMaxIterations ( ) const

Definition at line 1214 of file AssemblyMap.cpp.

References m_maxIterations.

1215  {
1216  return m_maxIterations;
1217  }
int m_maxIterations
Maximum iterations for iterative solver.
Definition: AssemblyMap.h:371
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 1181 of file AssemblyMap.cpp.

References m_nextLevelLocalToGlobalMap.

1182  {
1184  }
AssemblyMapSharedPtr m_nextLevelLocalToGlobalMap
Map from the patches of the previous level to the patches of the current level.
Definition: AssemblyMap.h:395
int Nektar::MultiRegions::AssemblyMap::GetNumDirEdges ( ) const

Definition at line 803 of file AssemblyMap.cpp.

References v_GetNumDirEdges().

804  {
805  return v_GetNumDirEdges();
806  }
virtual int v_GetNumDirEdges() const
int Nektar::MultiRegions::AssemblyMap::GetNumDirFaces ( ) const

Definition at line 808 of file AssemblyMap.cpp.

References v_GetNumDirFaces().

809  {
810  return v_GetNumDirFaces();
811  }
virtual int v_GetNumDirFaces() const
int Nektar::MultiRegions::AssemblyMap::GetNumGlobalBndCoeffs ( ) const

Returns the total number of global boundary coefficients.

Definition at line 935 of file AssemblyMap.cpp.

References m_numGlobalBndCoeffs.

Referenced by AssemblyMap().

936  {
937  return m_numGlobalBndCoeffs;
938  }
int m_numGlobalBndCoeffs
Total number of global boundary coefficients.
Definition: AssemblyMap.h:314
int Nektar::MultiRegions::AssemblyMap::GetNumGlobalCoeffs ( ) const

Returns the total number of global coefficients.

Definition at line 945 of file AssemblyMap.cpp.

References m_numGlobalCoeffs.

946  {
947  return m_numGlobalCoeffs;
948  }
int m_numGlobalCoeffs
Total number of global coefficients.
Definition: AssemblyMap.h:342
int Nektar::MultiRegions::AssemblyMap::GetNumGlobalDirBndCoeffs ( ) const

Returns the number of global Dirichlet boundary coefficients.

Definition at line 919 of file AssemblyMap.cpp.

References m_numGlobalDirBndCoeffs.

Referenced by AssemblyMap().

920  {
922  }
int m_numGlobalDirBndCoeffs
Number of Global Dirichlet Boundary Coefficients.
Definition: AssemblyMap.h:318
int Nektar::MultiRegions::AssemblyMap::GetNumLocalBndCoeffs ( ) const

Returns the total number of local boundary coefficients.

Definition at line 930 of file AssemblyMap.cpp.

References m_numLocalBndCoeffs.

Referenced by AssemblyMap().

931  {
932  return m_numLocalBndCoeffs;
933  }
int m_numLocalBndCoeffs
Number of local boundary coefficients.
Definition: AssemblyMap.h:312
const Array< OneD, const unsigned int > & Nektar::MultiRegions::AssemblyMap::GetNumLocalBndCoeffsPerPatch ( )

Returns the number of local boundary coefficients in each patch.

Definition at line 1168 of file AssemblyMap.cpp.

References m_numLocalBndCoeffsPerPatch.

Referenced by AssemblyMap().

1169  {
1171  }
Array< OneD, unsigned int > m_numLocalBndCoeffsPerPatch
The number of bnd dofs per patch.
Definition: AssemblyMap.h:388
int Nektar::MultiRegions::AssemblyMap::GetNumLocalCoeffs ( ) const

Returns the total number of local coefficients.

Definition at line 940 of file AssemblyMap.cpp.

References m_numLocalCoeffs.

941  {
942  return m_numLocalCoeffs;
943  }
int m_numLocalCoeffs
Total number of local coefficients.
Definition: AssemblyMap.h:331
int Nektar::MultiRegions::AssemblyMap::GetNumLocalDirBndCoeffs ( ) const

Returns the number of local Dirichlet boundary coefficients.

Definition at line 925 of file AssemblyMap.cpp.

References m_numLocalDirBndCoeffs.

Referenced by AssemblyMap().

926  {
927  return m_numLocalDirBndCoeffs;
928  }
int m_numLocalDirBndCoeffs
Number of Local Dirichlet Boundary Coefficients.
Definition: AssemblyMap.h:316
const Array< OneD, const unsigned int > & Nektar::MultiRegions::AssemblyMap::GetNumLocalIntCoeffsPerPatch ( )

Returns the number of local interior coefficients in each patch.

Definition at line 1175 of file AssemblyMap.cpp.

References m_numLocalIntCoeffsPerPatch.

1176  {
1178  }
Array< OneD, unsigned int > m_numLocalIntCoeffsPerPatch
The number of int dofs per patch.
Definition: AssemblyMap.h:390
int Nektar::MultiRegions::AssemblyMap::GetNumNonDirEdgeModes ( ) const

Definition at line 793 of file AssemblyMap.cpp.

References v_GetNumNonDirEdgeModes().

794  {
795  return v_GetNumNonDirEdgeModes();
796  }
virtual int v_GetNumNonDirEdgeModes() const
int Nektar::MultiRegions::AssemblyMap::GetNumNonDirEdges ( ) const

Definition at line 813 of file AssemblyMap.cpp.

References v_GetNumNonDirEdges().

814  {
815  return v_GetNumNonDirEdges();
816  }
virtual int v_GetNumNonDirEdges() const
int Nektar::MultiRegions::AssemblyMap::GetNumNonDirFaceModes ( ) const

Definition at line 798 of file AssemblyMap.cpp.

References v_GetNumNonDirFaceModes().

799  {
800  return v_GetNumNonDirFaceModes();
801  }
virtual int v_GetNumNonDirFaceModes() const
int Nektar::MultiRegions::AssemblyMap::GetNumNonDirFaces ( ) const

Definition at line 818 of file AssemblyMap.cpp.

References v_GetNumNonDirFaces().

819  {
820  return v_GetNumNonDirFaces();
821  }
virtual int v_GetNumNonDirFaces() const
int Nektar::MultiRegions::AssemblyMap::GetNumNonDirVertexModes ( ) const

Definition at line 788 of file AssemblyMap.cpp.

References v_GetNumNonDirVertexModes().

789  {
790  return v_GetNumNonDirVertexModes();
791  }
virtual int v_GetNumNonDirVertexModes() const
int Nektar::MultiRegions::AssemblyMap::GetNumPatches ( ) const

Returns the number of patches in this static condensation level.

Definition at line 1162 of file AssemblyMap.cpp.

References m_numPatches.

Referenced by AssemblyMap().

1163  {
1164  return m_numPatches;
1165  }
int m_numPatches
The number of patches (~elements) in the current level.
Definition: AssemblyMap.h:386
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 1187 of file AssemblyMap.cpp.

References m_patchMapFromPrevLevel.

1189  {
1190  return m_patchMapFromPrevLevel;
1191  }
PatchMapSharedPtr m_patchMapFromPrevLevel
Mapping information for previous level in MultiLevel Solver.
Definition: AssemblyMap.h:408
PreconditionerType Nektar::MultiRegions::AssemblyMap::GetPreconType ( ) const

Definition at line 1204 of file AssemblyMap.cpp.

References m_preconType.

1205  {
1206  return m_preconType;
1207  }
PreconditionerType m_preconType
Type type of preconditioner to use in iterative solver.
Definition: AssemblyMap.h:368
bool Nektar::MultiRegions::AssemblyMap::GetSignChange ( )

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

Definition at line 844 of file AssemblyMap.cpp.

References m_signChange.

Referenced by AssemblyMap().

845  {
846  return m_signChange;
847  }
bool m_signChange
Flag indicating if modes require sign reversal.
Definition: AssemblyMap.h:345
bool Nektar::MultiRegions::AssemblyMap::GetSingularSystem ( ) const

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

Definition at line 950 of file AssemblyMap.cpp.

References m_systemSingular.

951  {
952  return m_systemSingular;
953  }
bool m_systemSingular
Flag indicating if the system is singular or not.
Definition: AssemblyMap.h:320
int Nektar::MultiRegions::AssemblyMap::GetStaticCondLevel ( ) const

Returns the level of static condensation for this map.

Definition at line 1157 of file AssemblyMap.cpp.

References m_staticCondLevel.

Referenced by AssemblyMap().

1158  {
1159  return m_staticCondLevel;
1160  }
int m_staticCondLevel
The level of recursion in the case of multi-level static condensation.
Definition: AssemblyMap.h:384
int Nektar::MultiRegions::AssemblyMap::GetSuccessiveRHS ( ) const

Definition at line 1219 of file AssemblyMap.cpp.

References m_successiveRHS.

1220  {
1221  return m_successiveRHS;
1222  }
int m_successiveRHS
sucessive RHS for iterative solver
Definition: AssemblyMap.h:377
void Nektar::MultiRegions::AssemblyMap::GlobalToLocal ( const Array< OneD, const NekDouble > &  global,
Array< OneD, NekDouble > &  loc 
) const

Definition at line 736 of file AssemblyMap.cpp.

References v_GlobalToLocal().

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

739  {
740  v_GlobalToLocal(global,loc);
741  }
virtual void v_GlobalToLocal(const Array< OneD, const NekDouble > &global, Array< OneD, NekDouble > &loc) const
void Nektar::MultiRegions::AssemblyMap::GlobalToLocal ( const NekVector< NekDouble > &  global,
NekVector< NekDouble > &  loc 
) const

Definition at line 743 of file AssemblyMap.cpp.

References v_GlobalToLocal().

746  {
747  v_GlobalToLocal(global,loc);
748  }
virtual void v_GlobalToLocal(const Array< OneD, const NekDouble > &global, Array< OneD, NekDouble > &loc) const
void Nektar::MultiRegions::AssemblyMap::GlobalToLocalBnd ( const NekVector< NekDouble > &  global,
NekVector< NekDouble > &  loc,
int  offset 
) const

Definition at line 955 of file AssemblyMap.cpp.

References Nektar::NekVector< DataType >::GetPtr().

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

959  {
960  GlobalToLocalBnd(global.GetPtr(), loc.GetPtr(), offset);
961  }
void GlobalToLocalBnd(const NekVector< NekDouble > &global, NekVector< NekDouble > &loc, int offset) const
void Nektar::MultiRegions::AssemblyMap::GlobalToLocalBnd ( const NekVector< NekDouble > &  global,
NekVector< NekDouble > &  loc 
) const

Definition at line 964 of file AssemblyMap.cpp.

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

967  {
968  GlobalToLocalBnd(global.GetPtr(), loc.GetPtr());
969  }
void GlobalToLocalBnd(const NekVector< NekDouble > &global, NekVector< NekDouble > &loc, int offset) const
void Nektar::MultiRegions::AssemblyMap::GlobalToLocalBnd ( const Array< OneD, const NekDouble > &  global,
Array< OneD, NekDouble > &  loc,
int  offset 
) const

Definition at line 972 of file AssemblyMap.cpp.

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

975  {
976  ASSERTL1(loc.num_elements() >= m_numLocalBndCoeffs,"Local vector is not of correct dimension");
977  ASSERTL1(global.num_elements() >= m_numGlobalBndCoeffs-offset,"Global vector is not of correct dimension");
978 
979  // offset input data by length "offset" for Dirichlet boundary conditions.
980  Array<OneD,NekDouble> tmp(m_numGlobalBndCoeffs,0.0);
981  Vmath::Vcopy(m_numGlobalBndCoeffs-offset, global.get(), 1, tmp.get() + offset, 1);
982 
983  if(m_signChange)
984  {
986  }
987  else
988  {
989  Vmath::Gathr(m_numLocalBndCoeffs, tmp.get(), m_localToGlobalBndMap.get(), loc.get());
990  }
991  }
bool m_signChange
Flag indicating if modes require sign reversal.
Definition: AssemblyMap.h:345
int m_numGlobalBndCoeffs
Total number of global boundary coefficients.
Definition: AssemblyMap.h:314
void Gathr(int n, const T *x, const int *y, T *z)
Gather vector z[i] = x[y[i]].
Definition: Vmath.cpp:630
Array< OneD, int > m_localToGlobalBndMap
Integer map of local boundary coeffs to global space.
Definition: AssemblyMap.h:348
int m_numLocalBndCoeffs
Number of local boundary coefficients.
Definition: AssemblyMap.h:312
Array< OneD, NekDouble > m_localToGlobalBndSign
Integer sign of local boundary coeffs to global space.
Definition: AssemblyMap.h:350
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
void Nektar::MultiRegions::AssemblyMap::GlobalToLocalBnd ( const Array< OneD, const NekDouble > &  global,
Array< OneD, NekDouble > &  loc 
) const

Definition at line 994 of file AssemblyMap.cpp.

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

997  {
998  ASSERTL1(loc.num_elements() >= m_numLocalBndCoeffs,"Local vector is not of correct dimension");
999  ASSERTL1(global.num_elements() >= m_numGlobalBndCoeffs,"Global vector is not of correct dimension");
1000 
1001  if(m_signChange)
1002  {
1003  Vmath::Gathr(m_numLocalBndCoeffs, m_localToGlobalBndSign.get(), global.get(), m_localToGlobalBndMap.get(), loc.get());
1004  }
1005  else
1006  {
1007  Vmath::Gathr(m_numLocalBndCoeffs, global.get(), m_localToGlobalBndMap.get(), loc.get());
1008  }
1009  }
bool m_signChange
Flag indicating if modes require sign reversal.
Definition: AssemblyMap.h:345
int m_numGlobalBndCoeffs
Total number of global boundary coefficients.
Definition: AssemblyMap.h:314
void Gathr(int n, const T *x, const int *y, T *z)
Gather vector z[i] = x[y[i]].
Definition: Vmath.cpp:630
Array< OneD, int > m_localToGlobalBndMap
Integer map of local boundary coeffs to global space.
Definition: AssemblyMap.h:348
int m_numLocalBndCoeffs
Number of local boundary coefficients.
Definition: AssemblyMap.h:312
Array< OneD, NekDouble > m_localToGlobalBndSign
Integer sign of local boundary coeffs to global space.
Definition: AssemblyMap.h:350
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191
void Nektar::MultiRegions::AssemblyMap::GlobalToLocalBndWithoutSign ( const Array< OneD, const NekDouble > &  global,
Array< OneD, NekDouble > &  loc 
)
protected

Definition at line 1224 of file AssemblyMap.cpp.

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

Referenced by AssemblyMap().

1227  {
1228  ASSERTL1(loc.num_elements() >= m_numLocalBndCoeffs,"Local vector is not of correct dimension");
1229  ASSERTL1(global.num_elements() >= m_numGlobalBndCoeffs,"Global vector is not of correct dimension");
1230 
1231  Vmath::Gathr(m_numLocalBndCoeffs, global.get(), m_localToGlobalBndMap.get(), loc.get());
1232  }
int m_numGlobalBndCoeffs
Total number of global boundary coefficients.
Definition: AssemblyMap.h:314
void Gathr(int n, const T *x, const int *y, T *z)
Gather vector z[i] = x[y[i]].
Definition: Vmath.cpp:630
Array< OneD, int > m_localToGlobalBndMap
Integer map of local boundary coeffs to global space.
Definition: AssemblyMap.h:348
int m_numLocalBndCoeffs
Number of local boundary coefficients.
Definition: AssemblyMap.h:312
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191
boost::shared_ptr< AssemblyMap > Nektar::MultiRegions::AssemblyMap::LinearSpaceMap ( const ExpList locexp,
GlobalSysSolnType  solnType 
)

Definition at line 828 of file AssemblyMap.cpp.

References v_LinearSpaceMap().

829  {
830  return v_LinearSpaceMap(locexp, solnType);
831  }
virtual boost::shared_ptr< AssemblyMap > v_LinearSpaceMap(const ExpList &locexp, GlobalSysSolnType solnType)
Generate a linear space mapping from existing mapping.
void Nektar::MultiRegions::AssemblyMap::LocalBndToGlobal ( const NekVector< NekDouble > &  loc,
NekVector< NekDouble > &  global,
int  offset 
) const

Definition at line 1011 of file AssemblyMap.cpp.

References Nektar::NekVector< DataType >::GetPtr().

Referenced by LocalBndToGlobal().

1015  {
1016  LocalBndToGlobal(loc.GetPtr(), global.GetPtr(), offset);
1017  }
void LocalBndToGlobal(const NekVector< NekDouble > &loc, NekVector< NekDouble > &global, int offset) const
void Nektar::MultiRegions::AssemblyMap::LocalBndToGlobal ( const NekVector< NekDouble > &  loc,
NekVector< NekDouble > &  global 
) const

Definition at line 1044 of file AssemblyMap.cpp.

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

1047  {
1048  LocalBndToGlobal(loc.GetPtr(), global.GetPtr());
1049  }
void LocalBndToGlobal(const NekVector< NekDouble > &loc, NekVector< NekDouble > &global, int offset) const
void Nektar::MultiRegions::AssemblyMap::LocalBndToGlobal ( const Array< OneD, const NekDouble > &  loc,
Array< OneD, NekDouble > &  global,
int  offset 
) const

Definition at line 1020 of file AssemblyMap.cpp.

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

1024  {
1025  ASSERTL1(loc.num_elements() >= m_numLocalBndCoeffs,"Local vector is not of correct dimension");
1026  ASSERTL1(global.num_elements() >= m_numGlobalBndCoeffs-offset,"Global vector is not of correct dimension");
1027 
1028  // offset input data by length "offset" for Dirichlet boundary conditions.
1029  Array<OneD,NekDouble> tmp(m_numGlobalBndCoeffs,0.0);
1030 
1031  if(m_signChange)
1032  {
1034  }
1035  else
1036  {
1037  Vmath::Scatr(m_numLocalBndCoeffs, loc.get(), m_localToGlobalBndMap.get(), tmp.get());
1038  }
1039 
1040  UniversalAssembleBnd(tmp);
1041  Vmath::Vcopy(m_numGlobalBndCoeffs-offset, tmp.get()+offset, 1, global.get(), 1);
1042  }
bool m_signChange
Flag indicating if modes require sign reversal.
Definition: AssemblyMap.h:345
int m_numGlobalBndCoeffs
Total number of global boundary coefficients.
Definition: AssemblyMap.h:314
void Scatr(int n, const T *x, const int *y, T *z)
Scatter vector z[y[i]] = x[i].
Definition: Vmath.cpp:659
Array< OneD, int > m_localToGlobalBndMap
Integer map of local boundary coeffs to global space.
Definition: AssemblyMap.h:348
int m_numLocalBndCoeffs
Number of local boundary coefficients.
Definition: AssemblyMap.h:312
Array< OneD, NekDouble > m_localToGlobalBndSign
Integer sign of local boundary coeffs to global space.
Definition: AssemblyMap.h:350
void UniversalAssembleBnd(Array< OneD, NekDouble > &pGlobal) const
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
void Nektar::MultiRegions::AssemblyMap::LocalBndToGlobal ( const Array< OneD, const NekDouble > &  loc,
Array< OneD, NekDouble > &  global 
) const

Definition at line 1052 of file AssemblyMap.cpp.

References ASSERTL1, m_localToGlobalBndMap, m_localToGlobalBndSign, m_numGlobalBndCoeffs, m_numLocalBndCoeffs, m_signChange, and Vmath::Scatr().

1055  {
1056  ASSERTL1(loc.num_elements() >= m_numLocalBndCoeffs,"Local vector is not of correct dimension");
1057  ASSERTL1(global.num_elements() >= m_numGlobalBndCoeffs,"Global vector is not of correct dimension");
1058 
1059  if(m_signChange)
1060  {
1061  Vmath::Scatr(m_numLocalBndCoeffs, m_localToGlobalBndSign.get(), loc.get(), m_localToGlobalBndMap.get(), global.get());
1062  }
1063  else
1064  {
1065  Vmath::Scatr(m_numLocalBndCoeffs, loc.get(), m_localToGlobalBndMap.get(), global.get());
1066  }
1067  }
bool m_signChange
Flag indicating if modes require sign reversal.
Definition: AssemblyMap.h:345
int m_numGlobalBndCoeffs
Total number of global boundary coefficients.
Definition: AssemblyMap.h:314
void Scatr(int n, const T *x, const int *y, T *z)
Scatter vector z[y[i]] = x[i].
Definition: Vmath.cpp:659
Array< OneD, int > m_localToGlobalBndMap
Integer map of local boundary coeffs to global space.
Definition: AssemblyMap.h:348
int m_numLocalBndCoeffs
Number of local boundary coefficients.
Definition: AssemblyMap.h:312
Array< OneD, NekDouble > m_localToGlobalBndSign
Integer sign of local boundary coeffs to global space.
Definition: AssemblyMap.h:350
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191
void Nektar::MultiRegions::AssemblyMap::LocalToGlobal ( const Array< OneD, const NekDouble > &  loc,
Array< OneD, NekDouble > &  global 
) const

Definition at line 722 of file AssemblyMap.cpp.

References v_LocalToGlobal().

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

725  {
726  v_LocalToGlobal(loc,global);
727  }
virtual void v_LocalToGlobal(const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global) const
void Nektar::MultiRegions::AssemblyMap::LocalToGlobal ( const NekVector< NekDouble > &  loc,
NekVector< NekDouble > &  global 
) const

Definition at line 729 of file AssemblyMap.cpp.

References v_LocalToGlobal().

732  {
733  v_LocalToGlobal(loc,global);
734  }
virtual void v_LocalToGlobal(const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global) const
void Nektar::MultiRegions::AssemblyMap::PrintStats ( std::ostream &  out,
std::string  variable 
) const

Definition at line 1234 of file AssemblyMap.cpp.

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

1236  {
1238  = m_session->GetComm()->GetRowComm();
1239  bool isRoot = vRowComm->GetRank() == 0;
1240  int n = vRowComm->GetSize();
1241  int i;
1242 
1243  // Determine number of global degrees of freedom.
1244  int globBndCnt = 0, globDirCnt = 0;
1245 
1246  for (i = 0; i < m_numGlobalBndCoeffs; ++i)
1247  {
1249  {
1250  globBndCnt++;
1251 
1252  if (i < m_numGlobalDirBndCoeffs)
1253  {
1254  globDirCnt++;
1255  }
1256  }
1257  }
1258 
1259  int globCnt = m_numGlobalCoeffs - m_numGlobalBndCoeffs + globBndCnt;
1260 
1261  // Calculate maximum valency
1262  Array<OneD, NekDouble> tmpLoc (m_numLocalBndCoeffs, 1.0);
1263  Array<OneD, NekDouble> tmpGlob(m_numGlobalBndCoeffs, 0.0);
1264 
1265  Vmath::Assmb(m_numLocalBndCoeffs, tmpLoc.get(), m_localToGlobalBndMap.get(), tmpGlob.get());
1266  UniversalAssembleBnd(tmpGlob);
1267 
1268  int totGlobDof = globCnt;
1269  int totGlobBndDof = globBndCnt;
1270  int totGlobDirDof = globDirCnt;
1271  int totLocalDof = m_numLocalCoeffs;
1272  int totLocalBndDof = m_numLocalBndCoeffs;
1273  int totLocalDirDof = m_numLocalDirBndCoeffs;
1274 
1275  int meanValence = 0;
1276  int maxValence = 0;
1277  int minValence = 10000000;
1278  for (int i = 0; i < m_numGlobalBndCoeffs; ++i)
1279  {
1281  {
1282  continue;
1283  }
1284 
1285  if (tmpGlob[i] > maxValence)
1286  {
1287  maxValence = tmpGlob[i];
1288  }
1289  if (tmpGlob[i] < minValence)
1290  {
1291  minValence = tmpGlob[i];
1292  }
1293  meanValence += tmpGlob[i];
1294  }
1295 
1296  vRowComm->AllReduce(maxValence, LibUtilities::ReduceMax);
1297  vRowComm->AllReduce(minValence, LibUtilities::ReduceMin);
1298  vRowComm->AllReduce(meanValence, LibUtilities::ReduceSum);
1299  vRowComm->AllReduce(totGlobDof, LibUtilities::ReduceSum);
1300  vRowComm->AllReduce(totGlobBndDof, LibUtilities::ReduceSum);
1301  vRowComm->AllReduce(totGlobDirDof, LibUtilities::ReduceSum);
1302  vRowComm->AllReduce(totLocalDof, LibUtilities::ReduceSum);
1303  vRowComm->AllReduce(totLocalBndDof, LibUtilities::ReduceSum);
1304  vRowComm->AllReduce(totLocalDirDof, LibUtilities::ReduceSum);
1305 
1306  meanValence /= totGlobBndDof;
1307 
1308  if (isRoot)
1309  {
1310  out << "Assembly map statistics for field " << variable << ":"
1311  << endl;
1312  out << " - Number of local/global dof : "
1313  << totLocalDof << " " << totGlobDof << endl;
1314  out << " - Number of local/global boundary dof : "
1315  << totLocalBndDof << " " << totGlobBndDof << endl;
1316  out << " - Number of local/global Dirichlet dof : "
1317  << totLocalDirDof << " " << totGlobDirDof << endl;
1318  out << " - dof valency (min/max/mean) : "
1319  << minValence << " " << maxValence << " " << meanValence
1320  << endl;
1321 
1322  if (n > 1)
1323  {
1324  NekDouble mean = m_numLocalCoeffs, mean2 = mean * mean;
1325  NekDouble minval = mean, maxval = mean;
1326  Array<OneD, NekDouble> tmp(1);
1327 
1328  for (i = 1; i < n; ++i)
1329  {
1330  vRowComm->Recv(i, tmp);
1331  mean += tmp[0];
1332  mean2 += tmp[0]*tmp[0];
1333 
1334  if (tmp[0] > maxval)
1335  {
1336  maxval = tmp[0];
1337  }
1338  if (tmp[0] < minval)
1339  {
1340  minval = tmp[0];
1341  }
1342  }
1343 
1344  out << " - Local dof dist. (min/max/mean/dev) : "
1345  << minval << " " << maxval << " " << (mean / n) << " "
1346  << sqrt(mean2/n - mean*mean/n/n) << endl;
1347 
1348  vRowComm->Block();
1349 
1350  mean = minval = maxval = m_numLocalBndCoeffs;
1351  mean2 = mean * mean;
1352 
1353  for (i = 1; i < n; ++i)
1354  {
1355  vRowComm->Recv(i, tmp);
1356  mean += tmp[0];
1357  mean2 += tmp[0]*tmp[0];
1358 
1359  if (tmp[0] > maxval)
1360  {
1361  maxval = tmp[0];
1362  }
1363  if (tmp[0] < minval)
1364  {
1365  minval = tmp[0];
1366  }
1367  }
1368 
1369  out << " - Local bnd dof dist. (min/max/mean/dev) : "
1370  << minval << " " << maxval << " " << (mean / n) << " "
1371  << sqrt(mean2/n - mean*mean/n/n) << endl;
1372  }
1373  }
1374  else
1375  {
1376  Array<OneD, NekDouble> tmp(1);
1377  tmp[0] = m_numLocalCoeffs;
1378  vRowComm->Send(0, tmp);
1379  vRowComm->Block();
1380  tmp[0] = m_numLocalBndCoeffs;
1381  vRowComm->Send(0, tmp);
1382  }
1383  }
int m_numGlobalBndCoeffs
Total number of global boundary coefficients.
Definition: AssemblyMap.h:314
int m_numLocalCoeffs
Total number of local coefficients.
Definition: AssemblyMap.h:331
boost::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:53
int m_numGlobalDirBndCoeffs
Number of Global Dirichlet Boundary Coefficients.
Definition: AssemblyMap.h:318
double NekDouble
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:695
Array< OneD, int > m_localToGlobalBndMap
Integer map of local boundary coeffs to global space.
Definition: AssemblyMap.h:348
int m_numLocalDirBndCoeffs
Number of Local Dirichlet Boundary Coefficients.
Definition: AssemblyMap.h:316
int m_numLocalBndCoeffs
Number of local boundary coefficients.
Definition: AssemblyMap.h:312
LibUtilities::SessionReaderSharedPtr m_session
Session object.
Definition: AssemblyMap.h:303
Array< OneD, int > m_globalToUniversalBndMapUnique
Integer map of unique process coeffs to universal space (signed)
Definition: AssemblyMap.h:360
void UniversalAssembleBnd(Array< OneD, NekDouble > &pGlobal) const
int m_numGlobalCoeffs
Total number of global coefficients.
Definition: AssemblyMap.h:342
void Nektar::MultiRegions::AssemblyMap::SetNextLevelLocalToGlobalMap ( AssemblyMapSharedPtr  pNextLevelLocalToGlobalMap)
void Nektar::MultiRegions::AssemblyMap::UniversalAssemble ( Array< OneD, NekDouble > &  pGlobal) const
void Nektar::MultiRegions::AssemblyMap::UniversalAssemble ( NekVector< NekDouble > &  pGlobal) const

Definition at line 770 of file AssemblyMap.cpp.

References v_UniversalAssemble().

772  {
773  v_UniversalAssemble(pGlobal);
774  }
virtual void v_UniversalAssemble(Array< OneD, NekDouble > &pGlobal) const
void Nektar::MultiRegions::AssemblyMap::UniversalAssemble ( Array< OneD, NekDouble > &  pGlobal,
int  offset 
) const

Definition at line 776 of file AssemblyMap.cpp.

References v_UniversalAssemble().

779  {
780  v_UniversalAssemble(pGlobal, offset);
781  }
virtual void v_UniversalAssemble(Array< OneD, NekDouble > &pGlobal) const
void Nektar::MultiRegions::AssemblyMap::UniversalAssembleBnd ( Array< OneD, NekDouble > &  pGlobal) const

Definition at line 1128 of file AssemblyMap.cpp.

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

Referenced by AssembleBnd(), Nektar::MultiRegions::AssemblyMapCG::AssemblyMapCG(), LocalBndToGlobal(), PrintStats(), and UniversalAssembleBnd().

1130  {
1131  ASSERTL1(pGlobal.num_elements() == m_numGlobalBndCoeffs,
1132  "Wrong size.");
1133  Gs::Gather(pGlobal, Gs::gs_add, m_bndGsh);
1134  }
int m_numGlobalBndCoeffs
Total number of global boundary coefficients.
Definition: AssemblyMap.h:314
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:224
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191
void Nektar::MultiRegions::AssemblyMap::UniversalAssembleBnd ( NekVector< NekDouble > &  pGlobal) const

Definition at line 1136 of file AssemblyMap.cpp.

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

1138  {
1139  UniversalAssembleBnd(pGlobal.GetPtr());
1140  }
void UniversalAssembleBnd(Array< OneD, NekDouble > &pGlobal) const
void Nektar::MultiRegions::AssemblyMap::UniversalAssembleBnd ( Array< OneD, NekDouble > &  pGlobal,
int  offset 
) const

Definition at line 1142 of file AssemblyMap.cpp.

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

1145  {
1146  Array<OneD, NekDouble> tmp(offset);
1147  if (offset > 0) Vmath::Vcopy(offset, pGlobal, 1, tmp, 1);
1148  UniversalAssembleBnd(pGlobal);
1149  if (offset > 0) Vmath::Vcopy(offset, tmp, 1, pGlobal, 1);
1150  }
void UniversalAssembleBnd(Array< OneD, NekDouble > &pGlobal) const
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
void Nektar::MultiRegions::AssemblyMap::v_Assemble ( const Array< OneD, const NekDouble > &  loc,
Array< OneD, NekDouble > &  global 
) const
privatevirtual

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

Definition at line 573 of file AssemblyMap.cpp.

References ASSERTL0.

Referenced by Assemble().

576  {
577  ASSERTL0(false, "Not defined for this type of mapping.");
578  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
void Nektar::MultiRegions::AssemblyMap::v_Assemble ( const NekVector< NekDouble > &  loc,
NekVector< NekDouble > &  global 
) const
privatevirtual

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

Definition at line 580 of file AssemblyMap.cpp.

References ASSERTL0.

583  {
584  ASSERTL0(false, "Not defined for this type of mapping.");
585  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
const Array< OneD, const int > & Nektar::MultiRegions::AssemblyMap::v_GetExtraDirEdges ( )
privatevirtual

Reimplemented in Nektar::MultiRegions::AssemblyMapCG.

Definition at line 657 of file AssemblyMap.cpp.

References ASSERTL0.

Referenced by GetExtraDirEdges().

658  {
659  ASSERTL0(false, "Not defined for this type of mapping.");
660  static Array<OneD, const int> result;
661  return result;
662  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
int Nektar::MultiRegions::AssemblyMap::v_GetFullSystemBandWidth ( ) const
privatevirtual

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

Definition at line 609 of file AssemblyMap.cpp.

References ASSERTL0.

Referenced by GetFullSystemBandWidth().

610  {
611  ASSERTL0(false, "Not defined for this type of mapping.");
612  return 0;
613  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
int Nektar::MultiRegions::AssemblyMap::v_GetGlobalToUniversalMap ( const int  i) const
privatevirtual

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

Definition at line 499 of file AssemblyMap.cpp.

References ASSERTL0.

500  {
501  ASSERTL0(false, "Not defined for this type of mapping.");
502  return 0;
503  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
const Array< OneD, const int > & Nektar::MultiRegions::AssemblyMap::v_GetGlobalToUniversalMap ( void  )
privatevirtual

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

Definition at line 518 of file AssemblyMap.cpp.

References ASSERTL0.

Referenced by GetGlobalToUniversalMap().

519  {
520  ASSERTL0(false, "Not defined for this type of mapping.");
521  static Array<OneD, const int> result;
522  return result;
523  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
int Nektar::MultiRegions::AssemblyMap::v_GetGlobalToUniversalMapUnique ( const int  i) const
privatevirtual

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

Definition at line 505 of file AssemblyMap.cpp.

References ASSERTL0.

506  {
507  ASSERTL0(false, "Not defined for this type of mapping.");
508  return 0;
509  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
const Array< OneD, const int > & Nektar::MultiRegions::AssemblyMap::v_GetGlobalToUniversalMapUnique ( void  )
privatevirtual

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

Definition at line 525 of file AssemblyMap.cpp.

References ASSERTL0.

Referenced by GetGlobalToUniversalMapUnique().

526  {
527  ASSERTL0(false, "Not defined for this type of mapping.");
528  static Array<OneD, const int> result;
529  return result;
530  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
int Nektar::MultiRegions::AssemblyMap::v_GetLocalToGlobalMap ( const int  i) const
privatevirtual

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

Definition at line 493 of file AssemblyMap.cpp.

References ASSERTL0.

494  {
495  ASSERTL0(false, "Not defined for this type of mapping.");
496  return 0;
497  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
const Array< OneD, const int > & Nektar::MultiRegions::AssemblyMap::v_GetLocalToGlobalMap ( void  )
privatevirtual

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

Definition at line 511 of file AssemblyMap.cpp.

References ASSERTL0.

Referenced by GetLocalToGlobalMap().

512  {
513  ASSERTL0(false, "Not defined for this type of mapping.");
514  static Array<OneD,const int> result;
515  return result;
516  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
NekDouble Nektar::MultiRegions::AssemblyMap::v_GetLocalToGlobalSign ( const int  i) const
privatevirtual

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

Definition at line 532 of file AssemblyMap.cpp.

References ASSERTL0.

533  {
534  ASSERTL0(false, "Not defined for this type of mapping.");
535  return 0.0;
536  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
const Array< OneD, NekDouble > & Nektar::MultiRegions::AssemblyMap::v_GetLocalToGlobalSign ( ) const
privatevirtual

Reimplemented in Nektar::MultiRegions::AssemblyMapCG.

Definition at line 538 of file AssemblyMap.cpp.

References ASSERTL0.

Referenced by GetLocalToGlobalSign().

539  {
540  ASSERTL0(false, "Not defined for this type of mapping.");
541  static Array<OneD, NekDouble> result;
542  return result;
543  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
int Nektar::MultiRegions::AssemblyMap::v_GetNumDirEdges ( ) const
privatevirtual

Reimplemented in Nektar::MultiRegions::AssemblyMapCG.

Definition at line 633 of file AssemblyMap.cpp.

References ASSERTL0.

Referenced by GetNumDirEdges().

634  {
635  ASSERTL0(false, "Not defined for this type of mapping.");
636  return 0;
637  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
int Nektar::MultiRegions::AssemblyMap::v_GetNumDirFaces ( ) const
privatevirtual

Reimplemented in Nektar::MultiRegions::AssemblyMapCG.

Definition at line 639 of file AssemblyMap.cpp.

References ASSERTL0.

Referenced by GetNumDirFaces().

640  {
641  ASSERTL0(false, "Not defined for this type of mapping.");
642  return 0;
643  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
int Nektar::MultiRegions::AssemblyMap::v_GetNumNonDirEdgeModes ( ) const
privatevirtual

Reimplemented in Nektar::MultiRegions::AssemblyMapCG.

Definition at line 621 of file AssemblyMap.cpp.

References ASSERTL0.

Referenced by GetNumNonDirEdgeModes().

622  {
623  ASSERTL0(false, "Not defined for this type of mapping.");
624  return 0;
625  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
int Nektar::MultiRegions::AssemblyMap::v_GetNumNonDirEdges ( ) const
privatevirtual

Reimplemented in Nektar::MultiRegions::AssemblyMapCG.

Definition at line 645 of file AssemblyMap.cpp.

References ASSERTL0.

Referenced by GetNumNonDirEdges().

646  {
647  ASSERTL0(false, "Not defined for this type of mapping.");
648  return 0;
649  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
int Nektar::MultiRegions::AssemblyMap::v_GetNumNonDirFaceModes ( ) const
privatevirtual

Reimplemented in Nektar::MultiRegions::AssemblyMapCG.

Definition at line 627 of file AssemblyMap.cpp.

References ASSERTL0.

Referenced by GetNumNonDirFaceModes().

628  {
629  ASSERTL0(false, "Not defined for this type of mapping.");
630  return 0;
631  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
int Nektar::MultiRegions::AssemblyMap::v_GetNumNonDirFaces ( ) const
privatevirtual

Reimplemented in Nektar::MultiRegions::AssemblyMapCG.

Definition at line 651 of file AssemblyMap.cpp.

References ASSERTL0.

Referenced by GetNumNonDirFaces().

652  {
653  ASSERTL0(false, "Not defined for this type of mapping.");
654  return 0;
655  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
int Nektar::MultiRegions::AssemblyMap::v_GetNumNonDirVertexModes ( ) const
privatevirtual

Reimplemented in Nektar::MultiRegions::AssemblyMapCG.

Definition at line 615 of file AssemblyMap.cpp.

References ASSERTL0.

Referenced by GetNumNonDirVertexModes().

616  {
617  ASSERTL0(false, "Not defined for this type of mapping.");
618  return 0;
619  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
void Nektar::MultiRegions::AssemblyMap::v_GlobalToLocal ( const Array< OneD, const NekDouble > &  global,
Array< OneD, NekDouble > &  loc 
) const
privatevirtual

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

Definition at line 559 of file AssemblyMap.cpp.

References ASSERTL0.

Referenced by GlobalToLocal().

562  {
563  ASSERTL0(false, "Not defined for this type of mapping.");
564  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
void Nektar::MultiRegions::AssemblyMap::v_GlobalToLocal ( const NekVector< NekDouble > &  global,
NekVector< NekDouble > &  loc 
) const
privatevirtual

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

Definition at line 566 of file AssemblyMap.cpp.

References ASSERTL0.

569  {
570  ASSERTL0(false, "Not defined for this type of mapping.");
571  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
boost::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 664 of file AssemblyMap.cpp.

References ASSERTL0.

Referenced by LinearSpaceMap().

666  {
667  ASSERTL0(false, "Not defined for this sub class");
668  static boost::shared_ptr<AssemblyMap> result;
669  return result;
670  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
void Nektar::MultiRegions::AssemblyMap::v_LocalToGlobal ( const Array< OneD, const NekDouble > &  loc,
Array< OneD, NekDouble > &  global 
) const
privatevirtual

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

Definition at line 545 of file AssemblyMap.cpp.

References ASSERTL0.

Referenced by LocalToGlobal().

548  {
549  ASSERTL0(false, "Not defined for this type of mapping.");
550  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
void Nektar::MultiRegions::AssemblyMap::v_LocalToGlobal ( const NekVector< NekDouble > &  loc,
NekVector< NekDouble > &  global 
) const
privatevirtual

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

Definition at line 552 of file AssemblyMap.cpp.

References ASSERTL0.

555  {
556  ASSERTL0(false, "Not defined for this type of mapping.");
557  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
void Nektar::MultiRegions::AssemblyMap::v_UniversalAssemble ( Array< OneD, NekDouble > &  pGlobal) const
privatevirtual

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

Definition at line 587 of file AssemblyMap.cpp.

Referenced by UniversalAssemble().

589  {
590  // Do nothing here since multi-level static condensation uses a
591  // AssemblyMap and thus will call this routine in serial.
592  }
void Nektar::MultiRegions::AssemblyMap::v_UniversalAssemble ( NekVector< NekDouble > &  pGlobal) const
privatevirtual

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

Definition at line 594 of file AssemblyMap.cpp.

596  {
597  // Do nothing here since multi-level static condensation uses a
598  // AssemblyMap and thus will call this routine in serial.
599  }
void Nektar::MultiRegions::AssemblyMap::v_UniversalAssemble ( Array< OneD, NekDouble > &  pGlobal,
int  offset 
) const
privatevirtual

Reimplemented in Nektar::MultiRegions::AssemblyMapCG.

Definition at line 601 of file AssemblyMap.cpp.

604  {
605  // Do nothing here since multi-level static condensation uses a
606  // AssemblyMap and thus will call this routine in serial.
607  }

Member Data Documentation

Array<OneD,int> Nektar::MultiRegions::AssemblyMap::m_bndCondCoeffsToGlobalCoeffsMap
protected
Array<OneD,NekDouble> Nektar::MultiRegions::AssemblyMap::m_bndCondCoeffsToGlobalCoeffsSign
protected

Integer map of bnd cond coeffs to global coefficients.

Definition at line 354 of file AssemblyMap.h.

Referenced by Nektar::MultiRegions::AssemblyMapCG::AssemblyMapCG(), Nektar::CoupledAssemblyMap::CoupledAssemblyMap(), and GetBndCondCoeffsToGlobalCoeffsSign().

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

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

Definition at line 356 of file AssemblyMap.h.

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

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

The bandwith of the global bnd system.

Definition at line 365 of file AssemblyMap.h.

Referenced by CalculateBndSystemBandWidth(), and GetBndSystemBandWidth().

LibUtilities::CommSharedPtr Nektar::MultiRegions::AssemblyMap::m_comm
protected
Array<OneD,int> Nektar::MultiRegions::AssemblyMap::m_globalToUniversalBndMap
protected
Array<OneD,int> Nektar::MultiRegions::AssemblyMap::m_globalToUniversalBndMapUnique
protected
Gs::gs_data* Nektar::MultiRegions::AssemblyMap::m_gsh
protected
size_t Nektar::MultiRegions::AssemblyMap::m_hash
protected
NekDouble Nektar::MultiRegions::AssemblyMap::m_iterativeTolerance
protected

Tolerance for iterative solver.

Definition at line 374 of file AssemblyMap.h.

Referenced by AssemblyMap(), and GetIterativeTolerance().

Array<OneD,int> Nektar::MultiRegions::AssemblyMap::m_localToGlobalBndMap
protected
Array<OneD,NekDouble> Nektar::MultiRegions::AssemblyMap::m_localToGlobalBndSign
protected
int Nektar::MultiRegions::AssemblyMap::m_lowestStaticCondLevel
protected
int Nektar::MultiRegions::AssemblyMap::m_maxIterations
protected

Maximum iterations for iterative solver.

Definition at line 371 of file AssemblyMap.h.

Referenced by AssemblyMap(), and GetMaxIterations().

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

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

int Nektar::MultiRegions::AssemblyMap::m_numGlobalBndCoeffs
protected
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 342 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(), and Nektar::MultiRegions::AssemblyMapCG::v_LinearSpaceMap().

int Nektar::MultiRegions::AssemblyMap::m_numGlobalDirBndCoeffs
protected
int Nektar::MultiRegions::AssemblyMap::m_numLocalBndCoeffs
protected
Array<OneD, unsigned int> Nektar::MultiRegions::AssemblyMap::m_numLocalBndCoeffsPerPatch
protected
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 331 of file AssemblyMap.h.

Referenced by Nektar::MultiRegions::AssemblyMapCG::AssemblyMapCG(), Nektar::MultiRegions::AssemblyMapDG::AssemblyMapDG(), Nektar::MultiRegions::AssemblyMapCG::CalculateFullSystemBandWidth(), Nektar::CoupledAssemblyMap::CoupledAssemblyMap(), Nektar::CoupledLocalToGlobalC0ContMap::CoupledLocalToGlobalC0ContMap(), GetNumLocalCoeffs(), PrintStats(), Nektar::MultiRegions::AssemblyMapCG::v_Assemble(), Nektar::MultiRegions::AssemblyMapCG::v_GlobalToLocal(), and Nektar::MultiRegions::AssemblyMapCG::v_LocalToGlobal().

int Nektar::MultiRegions::AssemblyMap::m_numLocalDirBndCoeffs
protected
Array<OneD, unsigned int> Nektar::MultiRegions::AssemblyMap::m_numLocalIntCoeffsPerPatch
protected
int Nektar::MultiRegions::AssemblyMap::m_numPatches
protected
PatchMapSharedPtr Nektar::MultiRegions::AssemblyMap::m_patchMapFromPrevLevel
private

Mapping information for previous level in MultiLevel Solver.

Definition at line 408 of file AssemblyMap.h.

Referenced by AssemblyMap(), and GetPatchMapFromPrevLevel().

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

Type type of preconditioner to use in iterative solver.

Definition at line 368 of file AssemblyMap.h.

Referenced by AssemblyMap(), and GetPreconType().

LibUtilities::SessionReaderSharedPtr Nektar::MultiRegions::AssemblyMap::m_session
protected
bool Nektar::MultiRegions::AssemblyMap::m_signChange
protected
GlobalSysSolnType Nektar::MultiRegions::AssemblyMap::m_solnType
protected
int Nektar::MultiRegions::AssemblyMap::m_staticCondLevel
protected
int Nektar::MultiRegions::AssemblyMap::m_successiveRHS
protected

sucessive RHS for iterative solver

Definition at line 377 of file AssemblyMap.h.

Referenced by AssemblyMap(), and GetSuccessiveRHS().

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