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

#include <AssemblyMapDG.h>

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

Public Member Functions

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

Static Public Member Functions

static void RealignTraceElement (Array< OneD, int > &toAlign, StdRegions::Orientation orient, int nquad1, int nquad2=0)
 

Protected Member Functions

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

Protected Attributes

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

Detailed Description

Definition at line 49 of file AssemblyMapDG.h.

Constructor & Destructor Documentation

◆ AssemblyMapDG() [1/2]

Nektar::MultiRegions::AssemblyMapDG::AssemblyMapDG ( )

Default constructor.

Definition at line 51 of file AssemblyMapDG.cpp.

52{
53}
int m_numDirichletBndPhys
Number of physical dirichlet boundary values in trace.
Definition: AssemblyMapDG.h:98

◆ AssemblyMapDG() [2/2]

Nektar::MultiRegions::AssemblyMapDG::AssemblyMapDG ( const LibUtilities::SessionReaderSharedPtr pSession,
const SpatialDomains::MeshGraphSharedPtr graph1D,
const ExpListSharedPtr trace,
const ExpList locExp,
const Array< OneD, const MultiRegions::ExpListSharedPtr > &  bndConstraint,
const Array< OneD, const SpatialDomains::BoundaryConditionShPtr > &  bndCond,
const PeriodicMap periodicTrace,
const std::string  variable = "DefaultVar" 
)

Constructor for trace map for one-dimensional expansion.

Definition at line 59 of file AssemblyMapDG.cpp.

66 : AssemblyMap(pSession, locExp.GetComm(), variable)
67{
68 int i, j, k, cnt, id, id1, gid;
69 int order_e = 0;
70 int nTraceExp = trace->GetExpSize();
71 int nbnd = bndCondExp.size();
72
76
77 const LocalRegions::ExpansionVector expList = *(locExp.GetExp());
78 int nel = expList.size();
79
80 map<int, int> meshTraceId;
81
82 m_signChange = true;
83
84 // determine mapping from geometry edges to trace
85 for (i = 0; i < nTraceExp; ++i)
86 {
87 meshTraceId[trace->GetExp(i)->GetGeom()->GetGlobalID()] = i;
88 }
89
90 // Count total number of trace elements
91 cnt = 0;
92 for (i = 0; i < nel; ++i)
93 {
94 cnt += expList[i]->GetNtraces();
95 }
96
97 Array<OneD, LocalRegions::ExpansionSharedPtr> tracemap(cnt);
99 Array<OneD, Array<OneD, LocalRegions::ExpansionSharedPtr>>(nel);
100
101 // set up trace expansions links;
102 for (cnt = i = 0; i < nel; ++i)
103 {
104 m_elmtToTrace[i] = tracemap + cnt;
105 exp = expList[i];
106
107 for (j = 0; j < exp->GetNtraces(); ++j)
108 {
109 id = exp->GetGeom()->GetTid(j);
110
111 if (meshTraceId.count(id) > 0)
112 {
113 m_elmtToTrace[i][j] =
114 trace->GetExp(meshTraceId.find(id)->second);
115 }
116 else
117 {
118 ASSERTL0(false, "Failed to find trace map");
119 }
120 }
121
122 cnt += exp->GetNtraces();
123 }
124
125 // Set up boundary mapping
126 cnt = 0;
127 for (i = 0; i < nbnd; ++i)
128 {
129 if (bndCond[i]->GetBoundaryConditionType() == SpatialDomains::ePeriodic)
130 {
131 continue;
132 }
133 cnt += bndCondExp[i]->GetExpSize();
134 }
135
136 set<int> dirTrace;
137
140
141 for (i = 0; i < bndCondExp.size(); ++i)
142 {
143 for (j = 0; j < bndCondExp[i]->GetExpSize(); ++j)
144 {
145 bndExp = bndCondExp[i]->GetExp(j);
146 traceGeom = bndExp->GetGeom();
147 id = traceGeom->GetGlobalID();
148
149 if (bndCond[i]->GetBoundaryConditionType() ==
151 {
152 m_numLocalDirBndCoeffs += bndExp->GetNcoeffs();
153 m_numDirichletBndPhys += bndExp->GetTotPoints();
154 dirTrace.insert(id);
155 }
156 }
157 }
158
159 // Set up integer mapping array and sign change for each degree of
160 // freedom + initialise some more data members.
163 m_numPatches = nel;
164 m_numLocalBndCoeffsPerPatch = Array<OneD, unsigned int>(nel);
165 m_numLocalIntCoeffsPerPatch = Array<OneD, unsigned int>(nel);
166
167 int nbndry = 0;
168 for (i = 0; i < nel; ++i) // count number of elements in array
169 {
170 int BndCoeffs = expList[i]->NumDGBndryCoeffs();
171 nbndry += BndCoeffs;
173 m_numLocalBndCoeffsPerPatch[i] = (unsigned int)BndCoeffs;
174 }
175
177 m_numLocalBndCoeffs = nbndry;
178 m_numLocalCoeffs = nbndry;
179 m_localToGlobalBndMap = Array<OneD, int>(nbndry);
180 m_localToGlobalBndSign = Array<OneD, NekDouble>(nbndry, 1);
181
182 // Set up array for potential mesh optimsation
183 Array<OneD, int> traceElmtGid(nTraceExp, -1);
184 int nDir = 0;
185 cnt = 0;
186
187 // We are now going to construct a graph of the mesh which can be
188 // reordered depending on the type of solver we would like to use.
189 typedef boost::adjacency_list<boost::setS, boost::vecS, boost::undirectedS>
190 BoostGraph;
191
192 BoostGraph boostGraphObj;
193 int trace_id, trace_id1;
194 int dirOffset = 0;
195
196 // make trace renumbering map where first solved trace starts
197 // at 0 so we can set up graph.
198 for (i = 0; i < nTraceExp; ++i)
199 {
200 id = trace->GetExp(i)->GetGeom()->GetGlobalID();
201
202 if (dirTrace.count(id) == 0)
203 {
204 // Initial put in element ordering (starting from zero) into
205 // traceElmtGid
206 boost::add_vertex(boostGraphObj);
207 traceElmtGid[i] = cnt++;
208 }
209 else
210 {
211 // Use existing offset for Dirichlet edges
212 traceElmtGid[i] = dirOffset;
213 dirOffset += trace->GetExp(i)->GetNcoeffs();
214 nDir++;
215 }
216 }
217
218 // Set up boost Graph
219 for (i = 0; i < nel; ++i)
220 {
221 exp = expList[i];
222
223 for (j = 0; j < exp->GetNtraces(); ++j)
224 {
225 // Add trace to boost graph for non-Dirichlet Boundary
226 traceGeom = m_elmtToTrace[i][j]->GetGeom();
227 id = traceGeom->GetGlobalID();
228 trace_id = meshTraceId.find(id)->second;
229
230 if (dirTrace.count(id) == 0)
231 {
232 for (k = j + 1; k < exp->GetNtraces(); ++k)
233 {
234 traceGeom = m_elmtToTrace[i][k]->GetGeom();
235 id1 = traceGeom->GetGlobalID();
236 trace_id1 = meshTraceId.find(id1)->second;
237
238 if (dirTrace.count(id1) == 0)
239 {
240 boost::add_edge((size_t)traceElmtGid[trace_id],
241 (size_t)traceElmtGid[trace_id1],
242 boostGraphObj);
243 }
244 }
245 }
246 }
247 }
248
249 int nGraphVerts = nTraceExp - nDir;
250 Array<OneD, int> perm(nGraphVerts);
251 Array<OneD, int> iperm(nGraphVerts);
252 Array<OneD, int> vwgts(nGraphVerts);
254
255 for (i = 0; i < nGraphVerts; ++i)
256 {
257 vwgts[i] = trace->GetExp(i + nDir)->GetNcoeffs();
258 }
259
260 if (nGraphVerts)
261 {
262 switch (m_solnType)
263 {
265 case eIterativeFull:
267 case eXxtFullMatrix:
268 case eXxtStaticCond:
269 case ePETScFullMatrix:
270 case ePETScStaticCond:
271 {
272 NoReordering(boostGraphObj, perm, iperm);
273 break;
274 }
276 {
277 CuthillMckeeReordering(boostGraphObj, perm, iperm);
278 break;
279 }
284 {
285 MultiLevelBisectionReordering(boostGraphObj, perm, iperm,
286 bottomUpGraph);
287 break;
288 }
289 default:
290 {
291 ASSERTL0(false, "Unrecognised solution type");
292 }
293 }
294 }
295
296 // Recast the permutation so that it can be used as a map from old
297 // trace ID to new trace ID
299 for (i = 0; i < nTraceExp - nDir; ++i)
300 {
301 traceElmtGid[perm[i] + nDir] = cnt;
302 cnt += trace->GetExp(perm[i] + nDir)->GetNcoeffs();
303 }
304
305 // Now have trace edges Gid position
306 cnt = 0;
307 for (i = 0; i < nel; ++i)
308 {
309 exp = expList[i];
310
311 for (j = 0; j < exp->GetNtraces(); ++j)
312 {
313 traceGeom = m_elmtToTrace[i][j]->GetGeom();
314 id = traceGeom->GetGlobalID();
315 gid = traceElmtGid[meshTraceId.find(id)->second];
316
317 const int nDim = exp->GetNumBases();
318
319 if (nDim == 1)
320 {
321 order_e = 1;
322 m_localToGlobalBndMap[cnt] = gid;
323 }
324 else if (nDim == 2)
325 {
326 order_e = exp->GetTraceNcoeffs(j);
327
328 if (exp->GetTraceOrient(j) == StdRegions::eForwards)
329 {
330 for (k = 0; k < order_e; ++k)
331 {
332 m_localToGlobalBndMap[k + cnt] = gid + k;
333 }
334 }
335 else
336 {
337 switch (m_elmtToTrace[i][j]->GetBasisType(0))
338 {
340 {
341 // reverse vertex order
342 m_localToGlobalBndMap[cnt] = gid + 1;
343 m_localToGlobalBndMap[cnt + 1] = gid;
344 for (k = 2; k < order_e; ++k)
345 {
346 m_localToGlobalBndMap[k + cnt] = gid + k;
347 }
348
349 // negate odd modes
350 for (k = 3; k < order_e; k += 2)
351 {
352 m_localToGlobalBndSign[cnt + k] = -1.0;
353 }
354 break;
355 }
357 {
358 // reverse order
359 for (k = 0; k < order_e; ++k)
360 {
361 m_localToGlobalBndMap[cnt + order_e - k - 1] =
362 gid + k;
363 }
364 break;
365 }
367 {
368 // reverse order
369 for (k = 0; k < order_e; ++k)
370 {
371 m_localToGlobalBndMap[cnt + order_e - k - 1] =
372 gid + k;
373 }
374 break;
375 }
376 default:
377 {
378 ASSERTL0(false, "Boundary type not permitted");
379 }
380 }
381 }
382 }
383 else if (nDim == 3)
384 {
385 order_e = exp->GetTraceNcoeffs(j);
386
387 std::map<int, int> orientMap;
388
389 Array<OneD, unsigned int> elmMap1(order_e);
390 Array<OneD, int> elmSign1(order_e);
391 Array<OneD, unsigned int> elmMap2(order_e);
392 Array<OneD, int> elmSign2(order_e);
393
394 StdRegions::Orientation fo = exp->GetTraceOrient(j);
395
396 // Construct mapping which will permute global IDs
397 // according to face orientations.
398 exp->GetTraceToElementMap(j, elmMap1, elmSign1, fo);
399 exp->GetTraceToElementMap(j, elmMap2, elmSign2,
401
402 for (k = 0; k < elmMap1.size(); ++k)
403 {
404 // Find the elemental co-efficient in the original
405 // mapping.
406 int idx = -1;
407 for (int l = 0; l < elmMap2.size(); ++l)
408 {
409 if (elmMap1[k] == elmMap2[l])
410 {
411 idx = l;
412 break;
413 }
414 }
415
416 ASSERTL2(idx != -1, "Problem with face to element map!");
417 orientMap[k] = idx;
418 }
419
420 for (k = 0; k < order_e; ++k)
421 {
422 m_localToGlobalBndMap[k + cnt] = gid + orientMap[k];
423 m_localToGlobalBndSign[k + cnt] = elmSign2[orientMap[k]];
424 }
425 }
426
427 cnt += order_e;
428 }
429 }
430
431 // set up identify map for lcoal to lcoal
432 m_localToLocalBndMap = Array<OneD, int>(m_numLocalBndCoeffs);
433
434 // local to bnd map is just a copy
435 for (i = 0; i < m_numLocalBndCoeffs; ++i)
436 {
438 }
439
440 m_numGlobalBndCoeffs = trace->GetNcoeffs();
442
444
445 // set up m_bndCondCoeffsToLocalTraceMap
446 // Number of boundary expansions
447 int nbndexp = 0;
448 int bndTotal = 0;
449 int bndOffset;
450 int traceOffset;
451
452 cnt = 0;
453 for (i = 0; i < nbnd; ++i)
454 {
455 if (bndCond[i]->GetBoundaryConditionType() == SpatialDomains::ePeriodic)
456 {
457 continue;
458 }
459 cnt += bndCondExp[i]->GetNcoeffs();
460 nbndexp += bndCondExp[i]->GetExpSize();
461 }
462
463 m_bndCondCoeffsToLocalTraceMap = Array<OneD, int>(cnt);
464 m_bndCondIDToGlobalTraceID = Array<OneD, int>(nbndexp);
465
466 cnt = 0;
467 for (i = 0; i < bndCondExp.size(); ++i)
468 {
469 if (bndCond[i]->GetBoundaryConditionType() == SpatialDomains::ePeriodic)
470 {
471 continue;
472 }
473
474 for (j = 0; j < bndCondExp[i]->GetExpSize(); ++j)
475 {
476 bndExp = bndCondExp[i]->GetExp(j);
477 id = bndExp->GetGeom()->GetGlobalID();
478
479 int meshId = meshTraceId.find(id)->second;
480 m_bndCondIDToGlobalTraceID[cnt++] = meshId;
481
482 // initialy set up map with global bnd location
483 // and then use the localToGlobalBndMap to invert
484 // since this is a one to one mapping on boundaries
485 traceOffset = traceElmtGid[meshId];
486 bndOffset = bndCondExp[i]->GetCoeff_Offset(j) + bndTotal;
487
488 for (k = 0; k < bndExp->GetNcoeffs(); ++k)
489 {
490 m_bndCondCoeffsToLocalTraceMap[bndOffset + k] = traceOffset + k;
491 }
492 }
493 bndTotal += bndCondExp[i]->GetNcoeffs();
494 }
495
496 // generate an inverse local to global bnd map;
497 map<int, int> invLocToGloMap;
498 for (i = 0; i < nbndry; ++i)
499 {
500 invLocToGloMap[m_localToGlobalBndMap[i]] = i;
501 }
502
503 // reset bndCondCoeffToLocalTraceMap to hold local rather
504 // than global reference
505 for (i = 0; i < m_bndCondCoeffsToLocalTraceMap.size(); ++i)
506 {
508 invLocToGloMap[m_bndCondCoeffsToLocalTraceMap[i]];
509 }
510
511 // Now set up mapping from global coefficients to universal.
512 ExpListSharedPtr tr = std::dynamic_pointer_cast<ExpList>(trace);
513 SetUpUniversalDGMap(locExp);
514
517 locExp, tr, m_elmtToTrace, bndCondExp, bndCond, periodicTrace));
518
523 nGraphVerts)
524 {
525 if (m_staticCondLevel < (bottomUpGraph->GetNlevels() - 1))
526 {
527 Array<OneD, int> vwgts_perm(nGraphVerts);
528
529 for (int i = 0; i < nGraphVerts; i++)
530 {
531 vwgts_perm[i] = vwgts[perm[i]];
532 }
533
534 bottomUpGraph->ExpandGraphWithVertexWeights(vwgts_perm);
537 bottomUpGraph);
538 }
539 }
540
541 m_hash =
543}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
Definition: ErrorUtil.hpp:265
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
Array< OneD, Array< OneD, LocalRegions::ExpansionSharedPtr > > m_elmtToTrace
list of edge expansions for a given element
AssemblyCommDGSharedPtr m_assemblyComm
void SetUpUniversalDGMap(const ExpList &locExp)
int m_lowestStaticCondLevel
Lowest static condensation level.
Definition: AssemblyMap.h:439
GlobalSysSolnType m_solnType
The solution type of the global system.
Definition: AssemblyMap.h:402
int m_numLocalCoeffs
Total number of local coefficients.
Definition: AssemblyMap.h:364
Array< OneD, int > m_bndCondCoeffsToLocalTraceMap
Integer map of bnd cond coeff to local trace coeff.
Definition: AssemblyMap.h:393
bool m_signChange
Flag indicating if modes require sign reversal.
Definition: AssemblyMap.h:378
int m_numGlobalCoeffs
Total number of global coefficients.
Definition: AssemblyMap.h:375
void CalculateBndSystemBandWidth()
Calculates the bandwidth of the boundary system.
Array< OneD, NekDouble > m_localToGlobalBndSign
Integer sign of local boundary coeffs to global space.
Definition: AssemblyMap.h:383
Array< OneD, unsigned int > m_numLocalBndCoeffsPerPatch
The number of bnd dofs per patch.
Definition: AssemblyMap.h:430
int m_numLocalBndCoeffs
Number of local boundary coefficients.
Definition: AssemblyMap.h:345
AssemblyMapSharedPtr m_nextLevelLocalToGlobalMap
Map from the patches of the previous level to the patches of the current level.
Definition: AssemblyMap.h:437
int m_staticCondLevel
The level of recursion in the case of multi-level static condensation.
Definition: AssemblyMap.h:426
int m_numLocalDirBndCoeffs
Number of Local Dirichlet Boundary Coefficients.
Definition: AssemblyMap.h:349
int m_numGlobalDirBndCoeffs
Number of Global Dirichlet Boundary Coefficients.
Definition: AssemblyMap.h:351
Array< OneD, unsigned int > m_numLocalIntCoeffsPerPatch
The number of int dofs per patch.
Definition: AssemblyMap.h:432
Array< OneD, int > m_localToGlobalBndMap
Integer map of local coeffs to global Boundary Dofs.
Definition: AssemblyMap.h:381
AssemblyMap()
Default constructor.
Definition: AssemblyMap.cpp:79
Array< OneD, int > m_localToLocalBndMap
Integer map of local boundary coeffs to local boundary system numbering.
Definition: AssemblyMap.h:385
int m_numPatches
The number of patches (~elements) in the current level.
Definition: AssemblyMap.h:428
int m_numGlobalBndCoeffs
Total number of global boundary coefficients.
Definition: AssemblyMap.h:347
Array< OneD, int > m_bndCondIDToGlobalTraceID
Integer map of bnd cond trace number to global trace number.
Definition: AssemblyMap.h:395
@ eGauss_Lagrange
Lagrange Polynomials using the Gauss points.
Definition: BasisType.h:57
@ eGLL_Lagrange
Lagrange for SEM basis .
Definition: BasisType.h:56
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:48
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:66
std::vector< ExpansionSharedPtr > ExpansionVector
Definition: Expansion.h:68
void NoReordering(const BoostGraph &graph, Array< OneD, int > &perm, Array< OneD, int > &iperm)
std::shared_ptr< AssemblyCommDG > AssemblyCommDGSharedPtr
void CuthillMckeeReordering(const BoostGraph &graph, Array< OneD, int > &perm, Array< OneD, int > &iperm)
void MultiLevelBisectionReordering(const BoostGraph &graph, Array< OneD, int > &perm, Array< OneD, int > &iperm, BottomUpSubStructuredGraphSharedPtr &substructgraph, std::set< int > partVerts, int mdswitch)
std::shared_ptr< BottomUpSubStructuredGraph > BottomUpSubStructuredGraphSharedPtr
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
std::shared_ptr< Geometry > GeometrySharedPtr
Definition: Geometry.h:51
std::size_t hash_range(Iter first, Iter last)
Definition: HashUtils.hpp:64

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, ASSERTL2, Nektar::MultiRegions::AssemblyMap::CalculateBndSystemBandWidth(), Nektar::MultiRegions::CuthillMckeeReordering(), Nektar::StdRegions::eDir1FwdDir1_Dir2FwdDir2, Nektar::MultiRegions::eDirectFullMatrix, Nektar::MultiRegions::eDirectMultiLevelStaticCond, Nektar::MultiRegions::eDirectStaticCond, Nektar::SpatialDomains::eDirichlet, Nektar::StdRegions::eForwards, Nektar::LibUtilities::eGauss_Lagrange, Nektar::LibUtilities::eGLL_Lagrange, Nektar::MultiRegions::eIterativeFull, Nektar::MultiRegions::eIterativeMultiLevelStaticCond, Nektar::MultiRegions::eIterativeStaticCond, Nektar::LibUtilities::eModified_A, Nektar::SpatialDomains::ePeriodic, Nektar::MultiRegions::ePETScFullMatrix, Nektar::MultiRegions::ePETScMultiLevelStaticCond, Nektar::MultiRegions::ePETScStaticCond, Nektar::MultiRegions::eXxtFullMatrix, Nektar::MultiRegions::eXxtMultiLevelStaticCond, Nektar::MultiRegions::eXxtStaticCond, Nektar::MultiRegions::ExpList::GetExp(), Nektar::hash_range(), m_assemblyComm, Nektar::MultiRegions::AssemblyMap::m_bndCondCoeffsToLocalTraceMap, Nektar::MultiRegions::AssemblyMap::m_bndCondIDToGlobalTraceID, m_elmtToTrace, Nektar::MultiRegions::AssemblyMap::m_hash, Nektar::MultiRegions::AssemblyMap::m_localToGlobalBndMap, Nektar::MultiRegions::AssemblyMap::m_localToGlobalBndSign, Nektar::MultiRegions::AssemblyMap::m_localToLocalBndMap, Nektar::MultiRegions::AssemblyMap::m_lowestStaticCondLevel, Nektar::MultiRegions::AssemblyMap::m_nextLevelLocalToGlobalMap, m_numDirichletBndPhys, Nektar::MultiRegions::AssemblyMap::m_numGlobalBndCoeffs, Nektar::MultiRegions::AssemblyMap::m_numGlobalCoeffs, Nektar::MultiRegions::AssemblyMap::m_numGlobalDirBndCoeffs, Nektar::MultiRegions::AssemblyMap::m_numLocalBndCoeffs, Nektar::MultiRegions::AssemblyMap::m_numLocalBndCoeffsPerPatch, Nektar::MultiRegions::AssemblyMap::m_numLocalCoeffs, Nektar::MultiRegions::AssemblyMap::m_numLocalDirBndCoeffs, Nektar::MultiRegions::AssemblyMap::m_numLocalIntCoeffsPerPatch, Nektar::MultiRegions::AssemblyMap::m_numPatches, Nektar::MultiRegions::AssemblyMap::m_signChange, Nektar::MultiRegions::AssemblyMap::m_solnType, Nektar::MultiRegions::AssemblyMap::m_staticCondLevel, Nektar::MultiRegions::MultiLevelBisectionReordering(), Nektar::MultiRegions::NoReordering(), and SetUpUniversalDGMap().

◆ ~AssemblyMapDG()

Nektar::MultiRegions::AssemblyMapDG::~AssemblyMapDG ( )
override

Destructor.

Definition at line 55 of file AssemblyMapDG.cpp.

56{
57}

Member Function Documentation

◆ GetAssemblyCommDG()

AssemblyCommDGSharedPtr Nektar::MultiRegions::AssemblyMapDG::GetAssemblyCommDG ( )

Definition at line 908 of file AssemblyMapDG.cpp.

909{
910 return m_assemblyComm;
911}

References m_assemblyComm.

◆ GetElmtToTrace() [1/2]

Array< OneD, Array< OneD, LocalRegions::ExpansionSharedPtr > > & Nektar::MultiRegions::AssemblyMapDG::GetElmtToTrace ( )

Definition at line 903 of file AssemblyMapDG.cpp.

904{
905 return m_elmtToTrace;
906}

References m_elmtToTrace.

◆ GetElmtToTrace() [2/2]

Array< OneD, LocalRegions::ExpansionSharedPtr > & Nektar::MultiRegions::AssemblyMapDG::GetElmtToTrace ( const int  i)

Definition at line 895 of file AssemblyMapDG.cpp.

897{
898 ASSERTL1(i >= 0 && i < m_elmtToTrace.size(), "i is out of range");
899 return m_elmtToTrace[i];
900}
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242

References ASSERTL1, and m_elmtToTrace.

◆ GetNumDirichletBndPhys()

int Nektar::MultiRegions::AssemblyMapDG::GetNumDirichletBndPhys ( )

Return the number of boundary segments on which Dirichlet boundary conditions are imposed.

Definition at line 890 of file AssemblyMapDG.cpp.

891{
893}

References m_numDirichletBndPhys.

◆ RealignTraceElement()

void Nektar::MultiRegions::AssemblyMapDG::RealignTraceElement ( Array< OneD, int > &  toAlign,
StdRegions::Orientation  orient,
int  nquad1,
int  nquad2 = 0 
)
static

Changes toAlign quadrature point order, where the realignment is given by orient, which defines the mapping needed to go between the original ordering and the new desired ordering.

Parameters
[in,out]toAlignData to reorder
[in]orientThe transformation to perform
[in]nquad1Quadrature points in direction 1
[in]nquad2Quadrature points in direction 2

Definition at line 732 of file AssemblyMapDG.cpp.

735{
736 if (orient == StdRegions::eBackwards)
737 {
738 ASSERTL1(nquad2 == 0, "nquad2 != 0 for reorienation");
739 for (int i = 0; i < nquad1 / 2; ++i)
740 {
741 swap(toAlign[i], toAlign[nquad1 - 1 - i]);
742 }
743 }
744 else if (orient != StdRegions::eForwards)
745 {
746 ASSERTL1(nquad2 != 0, "nquad2 == 0 for reorienation");
747
748 Array<OneD, int> tmp(nquad1 * nquad2);
749
750 // Copy transpose.
755 {
756 for (int i = 0; i < nquad2; ++i)
757 {
758 for (int j = 0; j < nquad1; ++j)
759 {
760 tmp[i * nquad1 + j] = toAlign[j * nquad2 + i];
761 }
762 }
763 }
764 else
765 {
766 for (int i = 0; i < nquad2; ++i)
767 {
768 for (int j = 0; j < nquad1; ++j)
769 {
770 tmp[i * nquad1 + j] = toAlign[i * nquad1 + j];
771 }
772 }
773 }
774
779 {
780 // Reverse x direction
781 for (int i = 0; i < nquad2; ++i)
782 {
783 for (int j = 0; j < nquad1 / 2; ++j)
784 {
785 swap(tmp[i * nquad1 + j], tmp[i * nquad1 + nquad1 - j - 1]);
786 }
787 }
788 }
789
794 {
795 // Reverse y direction
796 for (int j = 0; j < nquad1; ++j)
797 {
798 for (int i = 0; i < nquad2 / 2; ++i)
799 {
800 swap(tmp[i * nquad1 + j],
801 tmp[(nquad2 - i - 1) * nquad1 + j]);
802 }
803 }
804 }
805 Vmath::Vcopy(nquad1 * nquad2, tmp, 1, toAlign, 1);
806 }
807}
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825

References ASSERTL1, Nektar::StdRegions::eBackwards, Nektar::StdRegions::eDir1BwdDir1_Dir2BwdDir2, Nektar::StdRegions::eDir1BwdDir1_Dir2FwdDir2, Nektar::StdRegions::eDir1BwdDir2_Dir2BwdDir1, Nektar::StdRegions::eDir1BwdDir2_Dir2FwdDir1, Nektar::StdRegions::eDir1FwdDir1_Dir2BwdDir2, Nektar::StdRegions::eDir1FwdDir2_Dir2BwdDir1, Nektar::StdRegions::eDir1FwdDir2_Dir2FwdDir1, Nektar::StdRegions::eForwards, and Vmath::Vcopy().

Referenced by Nektar::MultiRegions::AssemblyCommDG::InitialiseStructure().

◆ SetUpUniversalDGMap()

void Nektar::MultiRegions::AssemblyMapDG::SetUpUniversalDGMap ( const ExpList locExp)
protected

Constructs a mapping between the process-local global numbering and a universal numbering of the trace space expansion. The universal numbering is defined by the mesh edge IDs to enforce consistency across processes.

Parameters
locExpList of local elemental expansions.

Definition at line 553 of file AssemblyMapDG.cpp.

554{
556 int cnt = 0;
557 int id = 0;
558 int order_e = 0;
559 int vGlobalId = 0;
560 int maxDof = 0;
561 int dof = 0;
562 int nDim = 0;
563 int i, j, k;
564
565 const LocalRegions::ExpansionVector &locExpVector = *(locExp.GetExp());
566
567 // Initialise the global to universal maps.
572
573 // Loop over all the elements in the domain and compute max
574 // DOF. Reduce across all processes to get universal maximum.
575 for (i = 0; i < locExpVector.size(); ++i)
576 {
577 locExpansion = locExpVector[i];
578 nDim = locExpansion->GetShapeDimension();
579
580 // Loop over all edges of element i
581 if (nDim == 1)
582 {
583 maxDof = (1 > maxDof ? 1 : maxDof);
584 }
585 else
586 {
587 for (j = 0; j < locExpansion->GetNtraces(); ++j)
588 {
589 dof = locExpansion->GetTraceNcoeffs(j);
590 maxDof = (dof > maxDof ? dof : maxDof);
591 }
592 }
593 }
594 m_comm->GetRowComm()->AllReduce(maxDof, LibUtilities::ReduceMax);
595
596 // Now have trace edges Gid position
597 cnt = 0;
598 for (i = 0; i < locExpVector.size(); ++i)
599 {
600 locExpansion = locExpVector[i];
601 nDim = locExpansion->GetShapeDimension();
602
603 // Populate mapping for each edge of the element.
604 if (nDim == 1)
605 {
606 int nverts = locExpansion->GetNverts();
607 for (j = 0; j < nverts; ++j)
608 {
610 m_elmtToTrace[i][j]->as<LocalRegions::PointExp>();
611 id = locPointExp->GetGeom()->GetGlobalID();
612 vGlobalId = m_localToGlobalBndMap[cnt + j];
613 m_globalToUniversalBndMap[vGlobalId] = id * maxDof + j + 1;
614 }
615 cnt += nverts;
616 }
617 else if (nDim == 2)
618 {
619 for (j = 0; j < locExpansion->GetNtraces(); ++j)
620 {
622 m_elmtToTrace[i][j]->as<LocalRegions::SegExp>();
623
624 id = locSegExp->GetGeom()->GetGlobalID();
625 order_e = locExpansion->GetTraceNcoeffs(j);
626
627 map<int, int> orientMap;
628 Array<OneD, unsigned int> map1(order_e), map2(order_e);
629 Array<OneD, int> sign1(order_e), sign2(order_e);
630
631 locExpansion->GetTraceToElementMap(j, map1, sign1,
633 locExpansion->GetTraceToElementMap(
634 j, map2, sign2, locExpansion->GetTraceOrient(j));
635
636 for (k = 0; k < map1.size(); ++k)
637 {
638 // Find the elemental co-efficient in the original
639 // mapping.
640 int idx = -1;
641 for (int l = 0; l < map2.size(); ++l)
642 {
643 if (map1[k] == map2[l])
644 {
645 idx = l;
646 break;
647 }
648 }
649
650 ASSERTL2(idx != -1, "Problem with face to"
651 " element map!");
652 orientMap[k] = idx;
653 }
654
655 for (k = 0; k < order_e; ++k)
656 {
657 vGlobalId = m_localToGlobalBndMap[k + cnt];
658 m_globalToUniversalBndMap[vGlobalId] =
659 id * maxDof + orientMap[k] + 1;
660 }
661 cnt += order_e;
662 }
663 }
664 else if (nDim == 3) // This could likely be combined with nDim == 2
665 {
666 for (j = 0; j < locExpansion->GetNtraces(); ++j)
667 {
669 m_elmtToTrace[i][j]->as<LocalRegions::Expansion2D>();
670
671 id = locFaceExp->GetGeom()->GetGlobalID();
672 order_e = locExpansion->GetTraceNcoeffs(j);
673
674 map<int, int> orientMap;
675 Array<OneD, unsigned int> map1(order_e), map2(order_e);
676 Array<OneD, int> sign1(order_e), sign2(order_e);
677
678 locExpansion->GetTraceToElementMap(
680 locExpansion->GetTraceToElementMap(
681 j, map2, sign2, locExpansion->GetTraceOrient(j));
682
683 for (k = 0; k < map1.size(); ++k)
684 {
685 // Find the elemental co-efficient in the original
686 // mapping.
687 int idx = -1;
688 for (int l = 0; l < map2.size(); ++l)
689 {
690 if (map1[k] == map2[l])
691 {
692 idx = l;
693 break;
694 }
695 }
696
697 ASSERTL2(idx != -1, "Problem with face to "
698 "element map!");
699 orientMap[k] = idx;
700 }
701
702 for (k = 0; k < order_e; ++k)
703 {
704 vGlobalId = m_localToGlobalBndMap[k + cnt];
705 m_globalToUniversalBndMap[vGlobalId] =
706 id * maxDof + orientMap[k] + 1;
707 }
708 cnt += order_e;
709 }
710 }
711 }
712
713 // Initialise GSlib and populate the unique map.
714 Array<OneD, long> tmp(m_globalToUniversalBndMap.size());
715 for (i = 0; i < m_globalToUniversalBndMap.size(); ++i)
716 {
717 tmp[i] = m_globalToUniversalBndMap[i];
718 }
719
720 bool verbose = m_comm->IsParallelInTime()
721 ? m_comm->GetTimeComm()->GetRank() == 0
722 : true;
723
724 m_bndGsh = m_gsh = Gs::Init(tmp, m_comm->GetRowComm(), verbose);
725 Gs::Unique(tmp, m_comm->GetRowComm());
726 for (i = 0; i < m_globalToUniversalBndMap.size(); ++i)
727 {
728 m_globalToUniversalBndMapUnique[i] = (tmp[i] >= 0 ? 1 : 0);
729 }
730}
Array< OneD, int > m_globalToUniversalBndMap
Integer map of process coeffs to universal space.
Definition: AssemblyMap.h:397
Array< OneD, int > m_globalToUniversalBndMapUnique
Integer map of unique process coeffs to universal space (signed)
Definition: AssemblyMap.h:399
LibUtilities::CommSharedPtr m_comm
Communicator.
Definition: AssemblyMap.h:336
static gs_data * Init(const Nektar::Array< OneD, long > pId, const LibUtilities::CommSharedPtr &pComm, bool verbose=true)
Initialise Gather-Scatter map.
Definition: GsLib.hpp:190
static void Unique(const Nektar::Array< OneD, long > pId, const LibUtilities::CommSharedPtr &pComm)
Updates pId to negate all-but-one references to each universal ID.
Definition: GsLib.hpp:225
std::shared_ptr< SegExp > SegExpSharedPtr
Definition: SegExp.h:248
std::shared_ptr< PointExp > PointExpSharedPtr
Definition: PointExp.h:132
std::shared_ptr< Expansion2D > Expansion2DSharedPtr
Definition: Expansion1D.h:46

References ASSERTL2, Nektar::StdRegions::eDir1FwdDir1_Dir2FwdDir2, Nektar::StdRegions::eForwards, Nektar::MultiRegions::ExpList::GetExp(), Nektar::LocalRegions::Expansion::GetGeom(), Nektar::LocalRegions::PointExp::GetGeom(), Gs::Init(), Nektar::MultiRegions::AssemblyMap::m_bndGsh, Nektar::MultiRegions::AssemblyMap::m_comm, m_elmtToTrace, Nektar::MultiRegions::AssemblyMap::m_globalToUniversalBndMap, Nektar::MultiRegions::AssemblyMap::m_globalToUniversalBndMapUnique, Nektar::MultiRegions::AssemblyMap::m_gsh, Nektar::MultiRegions::AssemblyMap::m_localToGlobalBndMap, Nektar::MultiRegions::AssemblyMap::m_numGlobalBndCoeffs, Nektar::LibUtilities::ReduceMax, and Gs::Unique().

Referenced by AssemblyMapDG().

◆ v_Assemble() [1/2]

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

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 863 of file AssemblyMapDG.cpp.

865{
866 AssembleBnd(loc, global);
867}
void AssembleBnd(const NekVector< NekDouble > &loc, NekVector< NekDouble > &global, int offset) const

References Nektar::MultiRegions::AssemblyMap::AssembleBnd(), and CG_Iterations::loc.

◆ v_Assemble() [2/2]

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

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 869 of file AssemblyMapDG.cpp.

871{
872 AssembleBnd(loc, global);
873}

References Nektar::MultiRegions::AssemblyMap::AssembleBnd(), and CG_Iterations::loc.

◆ v_GetFullSystemBandWidth()

int Nektar::MultiRegions::AssemblyMapDG::v_GetFullSystemBandWidth ( ) const
overrideprotectedvirtual

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 885 of file AssemblyMapDG.cpp.

886{
887 return GetBndSystemBandWidth();
888}
int GetBndSystemBandWidth() const
Returns the bandwidth of the boundary system.

References Nektar::MultiRegions::AssemblyMap::GetBndSystemBandWidth().

◆ v_GetGlobalToUniversalMap() [1/2]

const Array< OneD, const int > & Nektar::MultiRegions::AssemblyMapDG::v_GetGlobalToUniversalMap ( void  )
overrideprotectedvirtual

◆ v_GetGlobalToUniversalMap() [2/2]

int Nektar::MultiRegions::AssemblyMapDG::v_GetGlobalToUniversalMap ( const int  i) const
overrideprotectedvirtual

◆ v_GetGlobalToUniversalMapUnique() [1/2]

const Array< OneD, const int > & Nektar::MultiRegions::AssemblyMapDG::v_GetGlobalToUniversalMapUnique ( void  )
overrideprotectedvirtual

◆ v_GetGlobalToUniversalMapUnique() [2/2]

int Nektar::MultiRegions::AssemblyMapDG::v_GetGlobalToUniversalMapUnique ( const int  i) const
overrideprotectedvirtual

◆ v_GetLocalToGlobalMap() [1/2]

const Array< OneD, const int > & Nektar::MultiRegions::AssemblyMapDG::v_GetLocalToGlobalMap ( void  )
overrideprotectedvirtual

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 824 of file AssemblyMapDG.cpp.

825{
827}

References Nektar::MultiRegions::AssemblyMap::m_localToGlobalBndMap.

◆ v_GetLocalToGlobalMap() [2/2]

int Nektar::MultiRegions::AssemblyMapDG::v_GetLocalToGlobalMap ( const int  i) const
overrideprotectedvirtual

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 809 of file AssemblyMapDG.cpp.

810{
811 return m_localToGlobalBndMap[i];
812}

References Nektar::MultiRegions::AssemblyMap::m_localToGlobalBndMap.

◆ v_GetLocalToGlobalSign()

NekDouble Nektar::MultiRegions::AssemblyMapDG::v_GetLocalToGlobalSign ( const int  i) const
overrideprotectedvirtual

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 839 of file AssemblyMapDG.cpp.

840{
841 return GetLocalToGlobalBndSign(i);
842}
Array< OneD, const NekDouble > GetLocalToGlobalBndSign() const
Retrieve the sign change for all local boundary modes.

References Nektar::MultiRegions::AssemblyMap::GetLocalToGlobalBndSign().

◆ v_GlobalToLocal() [1/2]

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

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 851 of file AssemblyMapDG.cpp.

853{
854 GlobalToLocalBnd(global, loc);
855}
void GlobalToLocalBnd(const NekVector< NekDouble > &global, NekVector< NekDouble > &loc, int offset) const

References Nektar::MultiRegions::AssemblyMap::GlobalToLocalBnd(), and CG_Iterations::loc.

◆ v_GlobalToLocal() [2/2]

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

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 857 of file AssemblyMapDG.cpp.

859{
860 GlobalToLocalBnd(global, loc);
861}

References Nektar::MultiRegions::AssemblyMap::GlobalToLocalBnd(), and CG_Iterations::loc.

◆ v_LocalToGlobal()

void Nektar::MultiRegions::AssemblyMapDG::v_LocalToGlobal ( const Array< OneD, const NekDouble > &  loc,
Array< OneD, NekDouble > &  global,
bool  useComm = false 
) const
overrideprotectedvirtual

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 844 of file AssemblyMapDG.cpp.

847{
848 LocalBndToGlobal(loc, global, useComm);
849}
void LocalBndToGlobal(const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global, int offset, bool UseComm=true) const

References CG_Iterations::loc, and Nektar::MultiRegions::AssemblyMap::LocalBndToGlobal().

◆ v_UniversalAssemble() [1/2]

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

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 875 of file AssemblyMapDG.cpp.

876{
877 Gs::Gather(pGlobal, Gs::gs_add, m_gsh);
878}
static void Gather(Nektar::Array< OneD, NekDouble > pU, gs_op pOp, gs_data *pGsh, Nektar::Array< OneD, NekDouble > pBuffer=NullNekDouble1DArray)
Performs a gather-scatter operation of the provided values.
Definition: GsLib.hpp:265
@ gs_add
Definition: GsLib.hpp:60

References Gs::Gather(), Gs::gs_add, and Nektar::MultiRegions::AssemblyMap::m_gsh.

◆ v_UniversalAssemble() [2/2]

void Nektar::MultiRegions::AssemblyMapDG::v_UniversalAssemble ( NekVector< NekDouble > &  pGlobal) const
overrideprotectedvirtual

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 880 of file AssemblyMapDG.cpp.

881{
882 UniversalAssemble(pGlobal.GetPtr());
883}
void UniversalAssemble(Array< OneD, NekDouble > &pGlobal) const
Array< OneD, DataType > & GetPtr()
Definition: NekVector.cpp:217

References Nektar::NekVector< DataType >::GetPtr(), and Nektar::MultiRegions::AssemblyMap::UniversalAssemble().

Member Data Documentation

◆ m_assemblyComm

AssemblyCommDGSharedPtr Nektar::MultiRegions::AssemblyMapDG::m_assemblyComm
protected

Definition at line 100 of file AssemblyMapDG.h.

Referenced by AssemblyMapDG(), and GetAssemblyCommDG().

◆ m_elmtToTrace

Array<OneD, Array<OneD, LocalRegions::ExpansionSharedPtr> > Nektar::MultiRegions::AssemblyMapDG::m_elmtToTrace
protected

list of edge expansions for a given element

Definition at line 103 of file AssemblyMapDG.h.

Referenced by AssemblyMapDG(), GetElmtToTrace(), and SetUpUniversalDGMap().

◆ m_numDirichletBndPhys

int Nektar::MultiRegions::AssemblyMapDG::m_numDirichletBndPhys
protected

Number of physical dirichlet boundary values in trace.

Definition at line 98 of file AssemblyMapDG.h.

Referenced by AssemblyMapDG(), and GetNumDirichletBndPhys().