Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
AssemblyMap.h
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File AssemblyMap.h
4 //
5 // For more information, please see: http://www.nektar.info
6 //
7 // The MIT License
8 //
9 // Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10 // Department of Aeronautics, Imperial College London (UK), and Scientific
11 // Computing and Imaging Institute, University of Utah (USA).
12 //
13 // License for the specific language governing rights and limitations under
14 // Permission is hereby granted, free of charge, to any person obtaining a
15 // copy of this software and associated documentation files (the "Software"),
16 // to deal in the Software without restriction, including without limitation
17 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
18 // and/or sell copies of the Software, and to permit persons to whom the
19 // Software is furnished to do so, subject to the following conditions:
20 //
21 // The above copyright notice and this permission notice shall be included
22 // in all copies or substantial portions of the Software.
23 //
24 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30 // DEALINGS IN THE SOFTWARE.
31 //
32 // Description: Assembly (e.g. local to global) base mapping routines
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
36 #ifndef MULTIREGIONS_ASSEMBLY_MAP_H
37 #define MULTIREGIONS_ASSEMBLY_MAP_H
38 
42 #include <vector>
45 
46 
47 namespace Nektar
48 {
49  namespace MultiRegions
50  {
51  // Forward declarations
52  class AssemblyMap;
53  class ExpList;
54  typedef boost::shared_ptr<AssemblyMap> AssemblyMapSharedPtr;
55  static AssemblyMapSharedPtr NullAssemblyMapSharedPtr;
56 
57  /// Base class for constructing local to global mapping of degrees of
58  /// freedom.
60  {
61  public:
62  /// Default constructor.
64  /// Constructor with a communicator
67  const std::string variable = "DefaultVar");
68 
69  /// Constructor for next level in multi-level static condensation.
70  MULTI_REGIONS_EXPORT AssemblyMap(AssemblyMap* oldLevelMap,
71  const BottomUpSubStructuredGraphSharedPtr& multiLevelGraph);
72  /// Destructor.
74 
75  /// Retrieves the communicator
77 
78  /// Retrieves the hash of this map
79  MULTI_REGIONS_EXPORT size_t GetHash() const;
80 
81  MULTI_REGIONS_EXPORT int GetLocalToGlobalMap(const int i) const;
82 
83  MULTI_REGIONS_EXPORT int GetGlobalToUniversalMap(const int i) const;
84 
86 
88 
90 
92 
94 
96 
99  Array<OneD, NekDouble>& global,
100  bool useComm = true) const;
101 
103  const NekVector<NekDouble>& loc,
104  NekVector< NekDouble>& global,
105  bool useComm = true) const;
106 
108  const Array<OneD, const NekDouble>& global,
109  Array<OneD, NekDouble>& loc) const;
110 
112  const NekVector<NekDouble>& global,
113  NekVector< NekDouble>& loc) const;
114 
116  const Array<OneD, const NekDouble> &loc,
117  Array<OneD, NekDouble> &global) const;
118 
120  const NekVector<NekDouble>& loc,
121  NekVector< NekDouble>& global) const;
122 
124  Array<OneD, NekDouble>& pGlobal) const;
125 
127  NekVector< NekDouble>& pGlobal) const;
128 
130  Array<OneD, NekDouble>& pGlobal,
131  int offset) const;
132 
133  /// Retrieve the global index of a given local boundary mode.
134  MULTI_REGIONS_EXPORT int GetLocalToGlobalBndMap(const int i) const;
135  /// Retrieve the global indices of the local boundary modes.
137 
139 
141 
142  /// Returns true if using a modal expansion requiring a change of
143  /// sign of some modes.
145 
146  /// Retrieve the sign change of a given local boundary mode.
148  /// Retrieve the sign change for all local boundary modes.
150  /// Retrieves the global index corresponding to a boundary expansion
151  /// mode.
153  /// Retrieves the global indices corresponding to the boundary
154  /// expansion modes.
157  /// Returns the modal sign associated with a given boundary
158  /// expansion mode.
160 
161  /// Returns the global index of the boundary trace giving the
162  /// index on the boundary expansion
166 
167  /// Returns the number of global Dirichlet boundary coefficients.
169  /// Returns the number of local Dirichlet boundary coefficients.
171  /// Returns the total number of global boundary coefficients.
173  /// Returns the total number of local boundary coefficients.
175  /// Returns the total number of local coefficients.
177  /// Returns the total number of global coefficients.
179  /// Retrieves if the system is singular (true) or not (false)
181 
182  ///
184  const NekVector<NekDouble>& global,
186  int offset) const;
187 
189  const NekVector<NekDouble>& global,
190  NekVector<NekDouble>& loc) const;
191 
193  const Array<OneD, const NekDouble>& global,
195  int offset) const;
196 
198  const Array<OneD, const NekDouble>& global,
199  Array<OneD,NekDouble>& loc) const;
200 
202  const NekVector<NekDouble>& loc,
203  NekVector<NekDouble>& global,
204  int offset) const;
205 
207  const NekVector<NekDouble>& loc,
208  NekVector<NekDouble>& global) const;
209 
211  const Array<OneD, const NekDouble>& loc,
212  Array<OneD,NekDouble>& global,
213  int offset) const;
214 
216  const Array<OneD, const NekDouble>& loc,
217  Array<OneD,NekDouble>& global) const;
218 
220  NekVector<NekDouble>& global, int offset) const;
221 
223  NekVector<NekDouble>& global) const;
224 
226  Array<OneD, NekDouble>& global, int offset) const;
227 
229  Array<OneD, NekDouble>& global) const;
230 
232  Array<OneD, NekDouble>& pGlobal) const;
233 
235  NekVector< NekDouble>& pGlobal) const;
236 
238  Array<OneD, NekDouble>& pGlobal,
239  int offset) const;
240 
242 
244 
246 
248 
250 
252 
254 
256 
257  MULTI_REGIONS_EXPORT void PrintStats(std::ostream &out, std::string variable, bool printHeader = true) const;
258 
261 
262  MULTI_REGIONS_EXPORT boost::shared_ptr<AssemblyMap> LinearSpaceMap(const ExpList &locexp, GlobalSysSolnType solnType);
263 
264  /// Returns the bandwidth of the boundary system.
266  /// Returns the level of static condensation for this map.
268  /// Returns the number of patches in this static condensation level.
270  /// Returns the number of local boundary coefficients in each patch.
273  /// Returns the number of local interior coefficients in each patch.
276  /// Returns the local to global mapping for the next level in the
277  /// multi-level static condensation.
278  MULTI_REGIONS_EXPORT const AssemblyMapSharedPtr
280 
281  MULTI_REGIONS_EXPORT void SetNextLevelLocalToGlobalMap( AssemblyMapSharedPtr pNextLevelLocalToGlobalMap);
282 
283  /// Returns the patch map from the previous level
284  /// of the multi-level static condensation.
286  GetPatchMapFromPrevLevel(void) const;
287 
288  /// Returns true if this is the last level in the multi-level
289  /// static condensation.
290  MULTI_REGIONS_EXPORT bool AtLastLevel() const;
291  /// Returns the method of solving global systems.
297 
299  {
301  }
302 
303  protected:
304  /// Session object
306 
307  /// Communicator
309 
310  /// Hash for map
311  size_t m_hash;
312 
313  /// Number of local boundary coefficients
315  /// Total number of global boundary coefficients
317  /// Number of Local Dirichlet Boundary Coefficients
319  /// Number of Global Dirichlet Boundary Coefficients
321  /// Flag indicating if the system is singular or not
323 
324  /// Total number of local coefficients
325  /** This corresponds to the number of total number of coefficients
326  * - For CG this corresponds to the total of bnd + int DOFs
327  * - For DG this corresponds to the number of bnd DOFs.
328  * This means that #m_numLocalCoeffs = #m_numLocalBndCoeffs
329  * This way, we can consider the trace-system solve as a
330  * statically condensed solve without interior DOFs. This allows
331  * us to use the same global system solver for both cases.
332  */
334 
335  /// Total number of global coefficients
336  /** This corresponds to the number of total number of coefficients
337  * - For CG this corresponds to the total of bnd + int DOFs.
338  * - For DG this corresponds to the number of bnd DOFs.
339  * This means that #m_numGlobalCoeffs = #m_numGlobalBndCoeffs
340  * This way, we can consider the trace-system solve as a
341  * statically condensed solve without interior DOFs. This allows
342  * us to use the same global system solver for both cases.
343  */
345 
346  /// Flag indicating if modes require sign reversal.
348 
349  /// Integer map of local boundary coeffs to global space
351  /// Integer sign of local boundary coeffs to global space
353  /// Integer map of bnd cond coeffs to global coefficients
355  /// Integer map of bnd cond coeffs to global coefficients
357  /// Integer map of bnd cond trace number to global trace number
359  /// Integer map of process coeffs to universal space
361  /// Integer map of unique process coeffs to universal space (signed)
363 
364  /// The solution type of the global system
366  /// The bandwith of the global bnd system
368 
369  /// Type type of preconditioner to use in iterative solver.
371 
372  /// Maximum iterations for iterative solver
374 
375  /// Tolerance for iterative solver
377 
378  /// sucessive RHS for iterative solver
380 
383 
384  /// The level of recursion in the case of multi-level static
385  /// condensation.
387  /// The number of patches (~elements) in the current level
389  /// The number of bnd dofs per patch
391  /// The number of int dofs per patch
393  /// Map from the patches of the previous level to the patches of
394  /// the current level
395 
396  /// The local to global mapping of the next level of recursion
397  AssemblyMapSharedPtr m_nextLevelLocalToGlobalMap;
398  /// Lowest static condensation level.
400 
401  /// Calculates the bandwidth of the boundary system.
403 
405  const Array<OneD, const NekDouble>& global,
406  Array<OneD,NekDouble>& loc);
407 
408  private:
409  /// Mapping information for previous level in MultiLevel Solver
411 
412  virtual int v_GetLocalToGlobalMap(const int i) const;
413 
414  virtual int v_GetGlobalToUniversalMap(const int i) const;
415 
416  virtual int v_GetGlobalToUniversalMapUnique(const int i) const;
417 
419 
421 
423 
424  virtual NekDouble v_GetLocalToGlobalSign(const int i) const;
425 
426  virtual const Array<OneD, NekDouble>& v_GetLocalToGlobalSign() const;
427 
428  virtual void v_LocalToGlobal(
429  const Array<OneD, const NekDouble>& loc,
430  Array<OneD, NekDouble>& global,
431  bool useComm) const;
432 
433  virtual void v_LocalToGlobal(
434  const NekVector<NekDouble>& loc,
435  NekVector< NekDouble>& global,
436  bool useComm) const;
437 
438  virtual void v_GlobalToLocal(
439  const Array<OneD, const NekDouble>& global,
440  Array<OneD, NekDouble>& loc) const;
441 
442  virtual void v_GlobalToLocal(
443  const NekVector<NekDouble>& global,
444  NekVector< NekDouble>& loc) const;
445 
446  virtual void v_Assemble(
447  const Array<OneD, const NekDouble> &loc,
448  Array<OneD, NekDouble> &global) const;
449 
450  virtual void v_Assemble(
451  const NekVector<NekDouble>& loc,
452  NekVector< NekDouble>& global) const;
453 
454  virtual void v_UniversalAssemble(
455  Array<OneD, NekDouble>& pGlobal) const;
456 
457  virtual void v_UniversalAssemble(
458  NekVector< NekDouble>& pGlobal) const;
459 
460  virtual void v_UniversalAssemble(
461  Array<OneD, NekDouble>& pGlobal,
462  int offset) const;
463 
464  virtual int v_GetFullSystemBandWidth() const;
465 
466  virtual int v_GetNumNonDirVertexModes() const;
467 
468  virtual int v_GetNumNonDirEdgeModes() const;
469 
470  virtual int v_GetNumNonDirFaceModes() const;
471 
472  virtual int v_GetNumDirEdges() const;
473 
474  virtual int v_GetNumDirFaces() const;
475 
476  virtual int v_GetNumNonDirEdges() const;
477 
478  virtual int v_GetNumNonDirFaces() const;
479 
480  virtual const Array<OneD, const int>&
482 
483  /// Generate a linear space mapping from existing mapping
484  virtual boost::shared_ptr<AssemblyMap> v_LinearSpaceMap(
485  const ExpList &locexp, GlobalSysSolnType solnType);
486  };
487 
488 
489  } // end of namespace
490 } // end of namespace
491 
492 #endif //MULTIREGIONS_ASSEMBLY_MAP_H
493 
494 
const Array< OneD, const int > & GetBndCondTraceToGlobalTraceMap()
void PrintStats(std::ostream &out, std::string variable, bool printHeader=true) const
PreconditionerType GetPreconType() const
bool m_systemSingular
Flag indicating if the system is singular or not.
Definition: AssemblyMap.h:322
const Array< OneD, const int > & GetGlobalToUniversalMap()
bool m_signChange
Flag indicating if modes require sign reversal.
Definition: AssemblyMap.h:347
virtual int v_GetNumNonDirFaceModes() const
const Array< OneD, NekDouble > & GetLocalToGlobalSign() const
PreconditionerType m_preconType
Type type of preconditioner to use in iterative solver.
Definition: AssemblyMap.h:370
int m_numGlobalBndCoeffs
Total number of global boundary coefficients.
Definition: AssemblyMap.h:316
LibUtilities::CommSharedPtr GetComm()
Retrieves the communicator.
LibUtilities::CommSharedPtr m_comm
Communicator.
Definition: AssemblyMap.h:308
virtual ~AssemblyMap()
Destructor.
virtual void v_LocalToGlobal(const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global, bool useComm) const
boost::shared_ptr< AssemblyMap > AssemblyMapSharedPtr
Definition: AssemblyMap.h:53
int GetNumLocalCoeffs() const
Returns the total number of local coefficients.
virtual int v_GetNumNonDirVertexModes() const
Array< OneD, int > m_bndCondTraceToGlobalTraceMap
Integer map of bnd cond trace number to global trace number.
Definition: AssemblyMap.h:358
int GetBndSystemBandWidth() const
Returns the bandwidth of the boundary system.
virtual int v_GetNumNonDirFaces() const
size_t GetHash() const
Retrieves the hash of this map.
boost::shared_ptr< BottomUpSubStructuredGraph > BottomUpSubStructuredGraphSharedPtr
#define MULTI_REGIONS_EXPORT
bool AtLastLevel() const
Returns true if this is the last level in the multi-level static condensation.
int GetNumPatches() const
Returns the number of patches in this static condensation level.
Array< OneD, const NekDouble > GetLocalToGlobalBndSign() const
Retrieve the sign change for all local boundary modes.
NekDouble GetIterativeTolerance() const
int m_numLocalCoeffs
Total number of local coefficients.
Definition: AssemblyMap.h:333
const Array< OneD, const int > & GetExtraDirEdges()
boost::shared_ptr< SessionReader > SessionReaderSharedPtr
Definition: MeshPartition.h:51
PatchMapSharedPtr m_patchMapFromPrevLevel
Mapping information for previous level in MultiLevel Solver.
Definition: AssemblyMap.h:410
virtual void v_GlobalToLocal(const Array< OneD, const NekDouble > &global, Array< OneD, NekDouble > &loc) const
boost::shared_ptr< PatchMap > PatchMapSharedPtr
const Array< OneD, const int > & GetGlobalToUniversalBndMapUnique()
virtual void v_UniversalAssemble(Array< OneD, NekDouble > &pGlobal) const
const Array< OneD, const int > & GetLocalToGlobalBndMap()
Retrieve the global indices of the local boundary modes.
AssemblyMapSharedPtr m_nextLevelLocalToGlobalMap
Map from the patches of the previous level to the patches of the current level.
Definition: AssemblyMap.h:397
Base class for constructing local to global mapping of degrees of freedom.
Definition: AssemblyMap.h:59
int m_successiveRHS
sucessive RHS for iterative solver
Definition: AssemblyMap.h:379
size_t m_hash
Hash for map.
Definition: AssemblyMap.h:311
void UniversalAssemble(Array< OneD, NekDouble > &pGlobal) const
Base class for all multi-elemental spectral/hp expansions.
Definition: ExpList.h:101
boost::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:55
void LocalToGlobal(const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global, bool useComm=true) const
int GetNumGlobalDirBndCoeffs() const
Returns the number of global Dirichlet boundary coefficients.
virtual boost::shared_ptr< AssemblyMap > v_LinearSpaceMap(const ExpList &locexp, GlobalSysSolnType solnType)
Generate a linear space mapping from existing mapping.
GlobalSysSolnType GetGlobalSysSolnType() const
Returns the method of solving global systems.
int m_bndSystemBandWidth
The bandwith of the global bnd system.
Definition: AssemblyMap.h:367
NekDouble m_iterativeTolerance
Tolerance for iterative solver.
Definition: AssemblyMap.h:376
Array< OneD, unsigned int > m_numLocalBndCoeffsPerPatch
The number of bnd dofs per patch.
Definition: AssemblyMap.h:390
boost::shared_ptr< AssemblyMap > LinearSpaceMap(const ExpList &locexp, GlobalSysSolnType solnType)
GlobalSysSolnType m_solnType
The solution type of the global system.
Definition: AssemblyMap.h:365
Array< OneD, NekDouble > m_bndCondCoeffsToGlobalCoeffsSign
Integer map of bnd cond coeffs to global coefficients.
Definition: AssemblyMap.h:356
int m_numGlobalDirBndCoeffs
Number of Global Dirichlet Boundary Coefficients.
Definition: AssemblyMap.h:320
NekDouble GetBndCondCoeffsToGlobalCoeffsSign(const int i)
Returns the modal sign associated with a given boundary expansion mode.
int GetStaticCondLevel() const
Returns the level of static condensation for this map.
void GlobalToLocal(const Array< OneD, const NekDouble > &global, Array< OneD, NekDouble > &loc) const
virtual int v_GetNumDirFaces() const
virtual int v_GetFullSystemBandWidth() const
void SetNextLevelLocalToGlobalMap(AssemblyMapSharedPtr pNextLevelLocalToGlobalMap)
const PatchMapSharedPtr & GetPatchMapFromPrevLevel(void) const
Returns the patch map from the previous level of the multi-level static condensation.
double NekDouble
Array< OneD, unsigned int > m_numLocalIntCoeffsPerPatch
The number of int dofs per patch.
Definition: AssemblyMap.h:392
int m_lowestStaticCondLevel
Lowest static condensation level.
Definition: AssemblyMap.h:399
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:350
const Array< OneD, const int > & GetBndCondCoeffsToGlobalCoeffsMap()
Retrieves the global indices corresponding to the boundary expansion modes.
int m_numLocalDirBndCoeffs
Number of Local Dirichlet Boundary Coefficients.
Definition: AssemblyMap.h:318
bool GetSignChange()
Returns true if using a modal expansion requiring a change of sign of some modes. ...
const Array< OneD, const int > & GetLocalToGlobalMap()
virtual void v_Assemble(const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global) const
const Array< OneD, const int > & GetGlobalToUniversalBndMap()
virtual int v_GetNumNonDirEdgeModes() const
Array< OneD, int > m_bndCondCoeffsToGlobalCoeffsMap
Integer map of bnd cond coeffs to global coefficients.
Definition: AssemblyMap.h:354
int m_numLocalBndCoeffs
Number of local boundary coefficients.
Definition: AssemblyMap.h:314
static AssemblyMapSharedPtr NullAssemblyMapSharedPtr
Definition: AssemblyMap.h:55
int m_staticCondLevel
The level of recursion in the case of multi-level static condensation.
Definition: AssemblyMap.h:386
Array< OneD, NekDouble > m_localToGlobalBndSign
Integer sign of local boundary coeffs to global space.
Definition: AssemblyMap.h:352
int GetNumGlobalBndCoeffs() const
Returns the total number of global boundary coefficients.
AssemblyMap()
Default constructor.
Definition: AssemblyMap.cpp:79
int m_maxIterations
Maximum iterations for iterative solver.
Definition: AssemblyMap.h:373
const Array< OneD, const unsigned int > & GetNumLocalBndCoeffsPerPatch()
Returns the number of local boundary coefficients in each patch.
void GlobalToLocalBndWithoutSign(const Array< OneD, const NekDouble > &global, Array< OneD, NekDouble > &loc)
const Array< OneD, const int > & GetGlobalToUniversalMapUnique()
int GetNumLocalDirBndCoeffs() const
Returns the number of local Dirichlet boundary coefficients.
LibUtilities::SessionReaderSharedPtr m_session
Session object.
Definition: AssemblyMap.h:305
void LocalBndToGlobal(const NekVector< NekDouble > &loc, NekVector< NekDouble > &global, int offset) const
void AssembleBnd(const NekVector< NekDouble > &loc, NekVector< NekDouble > &global, int offset) const
virtual const Array< OneD, NekDouble > & v_GetLocalToGlobalSign() const
virtual const Array< OneD, const int > & v_GetExtraDirEdges()
int GetNumGlobalCoeffs() const
Returns the total number of global coefficients.
int GetNumLocalBndCoeffs() const
Returns the total number of local boundary coefficients.
Array< OneD, int > m_globalToUniversalBndMap
Integer map of process coeffs to universal space.
Definition: AssemblyMap.h:360
Array< OneD, int > m_globalToUniversalBndMapUnique
Integer map of unique process coeffs to universal space (signed)
Definition: AssemblyMap.h:362
virtual const Array< OneD, const int > & v_GetGlobalToUniversalMapUnique()
void UniversalAssembleBnd(Array< OneD, NekDouble > &pGlobal) const
int m_numGlobalCoeffs
Total number of global coefficients.
Definition: AssemblyMap.h:344
void GlobalToLocalBnd(const NekVector< NekDouble > &global, NekVector< NekDouble > &loc, int offset) const
const Array< OneD, const unsigned int > & GetNumLocalIntCoeffsPerPatch()
Returns the number of local interior coefficients in each patch.
virtual int v_GetNumNonDirEdges() const
virtual const Array< OneD, const int > & v_GetGlobalToUniversalMap()
const AssemblyMapSharedPtr GetNextLevelLocalToGlobalMap() const
Returns the local to global mapping for the next level in the multi-level static condensation.
void Assemble(const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global) const
virtual const Array< OneD, const int > & v_GetLocalToGlobalMap()
virtual int v_GetNumDirEdges() const
int m_numPatches
The number of patches (~elements) in the current level.
Definition: AssemblyMap.h:388
bool GetSingularSystem() const
Retrieves if the system is singular (true) or not (false)