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

A global linear system. More...

#include <GlobalLinSysIterativeStaticCond.h>

Inheritance diagram for Nektar::MultiRegions::GlobalLinSysIterativeStaticCond:
Inheritance graph
[legend]
Collaboration diagram for Nektar::MultiRegions::GlobalLinSysIterativeStaticCond:
Collaboration graph
[legend]

Public Types

typedef NekSparseDiagBlkMatrix< StorageSmvBsr< NekDouble > > DNekSmvBsrDiagBlkMat
 
typedef boost::shared_ptr< DNekSmvBsrDiagBlkMatDNekSmvBsrDiagBlkMatSharedPtr
 

Public Member Functions

 GlobalLinSysIterativeStaticCond (const GlobalLinSysKey &mkey, const boost::weak_ptr< ExpList > &pExpList, const boost::shared_ptr< AssemblyMap > &locToGloMap)
 Constructor for full direct matrix solve. More...
 
 GlobalLinSysIterativeStaticCond (const GlobalLinSysKey &mkey, const boost::weak_ptr< ExpList > &pExpList, const DNekScalBlkMatSharedPtr pSchurCompl, const DNekScalBlkMatSharedPtr pBinvD, const DNekScalBlkMatSharedPtr pC, const DNekScalBlkMatSharedPtr pInvD, const boost::shared_ptr< AssemblyMap > &locToGloMap, const PreconditionerSharedPtr pPrecon)
 Constructor for full direct matrix solve. More...
 
virtual ~GlobalLinSysIterativeStaticCond ()
 
- Public Member Functions inherited from Nektar::MultiRegions::GlobalLinSysIterative
 GlobalLinSysIterative (const GlobalLinSysKey &pKey, const boost::weak_ptr< ExpList > &pExpList, const boost::shared_ptr< AssemblyMap > &pLocToGloMap)
 Constructor for full direct matrix solve. More...
 
virtual ~GlobalLinSysIterative ()
 
- Public Member Functions inherited from Nektar::MultiRegions::GlobalLinSys
 GlobalLinSys (const GlobalLinSysKey &pKey, const boost::weak_ptr< ExpList > &pExpList, const boost::shared_ptr< AssemblyMap > &pLocToGloMap)
 Constructor for full direct matrix solve. More...
 
virtual ~GlobalLinSys ()
 
const GlobalLinSysKeyGetKey (void) const
 Returns the key associated with the system. More...
 
const boost::weak_ptr< ExpList > & GetLocMat (void) const
 
void InitObject ()
 
void Initialise (const boost::shared_ptr< AssemblyMap > &pLocToGloMap)
 
void Solve (const Array< OneD, const NekDouble > &in, Array< OneD, NekDouble > &out, const AssemblyMapSharedPtr &locToGloMap, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
 Solve the linear system for given input and output vectors using a specified local to global map. More...
 
boost::shared_ptr< GlobalLinSysGetSharedThisPtr ()
 Returns a shared pointer to the current object. More...
 
int GetNumBlocks ()
 
DNekScalMatSharedPtr GetBlock (unsigned int n)
 
DNekScalBlkMatSharedPtr GetStaticCondBlock (unsigned int n)
 
void DropStaticCondBlock (unsigned int n)
 
void SolveLinearSystem (const int pNumRows, const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const AssemblyMapSharedPtr &locToGloMap, const int pNumDir=0)
 Solve the linear system for given input and output vectors. More...
 
- Public Member Functions inherited from Nektar::MultiRegions::GlobalLinSysStaticCond
 GlobalLinSysStaticCond (const GlobalLinSysKey &mkey, const boost::weak_ptr< ExpList > &pExpList, const boost::shared_ptr< AssemblyMap > &locToGloMap)
 Constructor for full direct matrix solve. More...
 
virtual ~GlobalLinSysStaticCond ()
 

Static Public Member Functions

static GlobalLinSysSharedPtr create (const GlobalLinSysKey &pLinSysKey, const boost::weak_ptr< ExpList > &pExpList, const boost::shared_ptr< AssemblyMap > &pLocToGloMap)
 Creates an instance of this class. More...
 

Static Public Attributes

static std::string className
 Name of class. More...
 
static std::string className2
 

Protected Member Functions

virtual DNekScalBlkMatSharedPtr v_GetStaticCondBlock (unsigned int n)
 Retrieves a the static condensation block matrices from n-th expansion using the matrix key provided by the m_linSysKey. More...
 
virtual GlobalLinSysStaticCondSharedPtr v_Recurse (const GlobalLinSysKey &mkey, const boost::weak_ptr< ExpList > &pExpList, const DNekScalBlkMatSharedPtr pSchurCompl, const DNekScalBlkMatSharedPtr pBinvD, const DNekScalBlkMatSharedPtr pC, const DNekScalBlkMatSharedPtr pInvD, const boost::shared_ptr< AssemblyMap > &locToGloMap)
 
virtual DNekScalBlkMatSharedPtr v_PreSolve (int scLevel, NekVector< NekDouble > &F_GlobBnd)
 
virtual void v_BasisTransform (Array< OneD, NekDouble > &pInOut, int offset)
 
virtual void v_BasisInvTransform (Array< OneD, NekDouble > &pInOut)
 
- Protected Member Functions inherited from Nektar::MultiRegions::GlobalLinSysIterative
void DoAconjugateProjection (const int pNumRows, const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const AssemblyMapSharedPtr &locToGloMap, const int pNumDir)
 A-conjugate projection technique. More...
 
void DoConjugateGradient (const int pNumRows, const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const AssemblyMapSharedPtr &locToGloMap, const int pNumDir)
 Actual iterative solve. More...
 
void Set_Rhs_Magnitude (const NekVector< NekDouble > &pIn)
 
- Protected Member Functions inherited from Nektar::MultiRegions::GlobalLinSys
virtual DNekScalMatSharedPtr v_GetBlock (unsigned int n)
 Retrieves the block matrix from n-th expansion using the matrix key provided by the m_linSysKey. More...
 
virtual void v_DropStaticCondBlock (unsigned int n)
 Releases the static condensation block matrices from NekManager of n-th expansion using the matrix key provided by the m_linSysKey. More...
 
- Protected Member Functions inherited from Nektar::MultiRegions::GlobalLinSysStaticCond
virtual int v_GetNumBlocks ()
 Get the number of blocks in this system. More...
 
virtual void v_Solve (const Array< OneD, const NekDouble > &in, Array< OneD, NekDouble > &out, const AssemblyMapSharedPtr &locToGloMap, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
 Solve the linear system for given input and output vectors using a specified local to global map. More...
 
virtual void v_Initialise (const boost::shared_ptr< AssemblyMap > &locToGloMap)
 Initialise this object. More...
 
void SetupTopLevel (const boost::shared_ptr< AssemblyMap > &locToGloMap)
 Set up the storage for the Schur complement or the top level of the multi-level Schur complement. More...
 
void ConstructNextLevelCondensedSystem (const boost::shared_ptr< AssemblyMap > &locToGloMap)
 

Private Member Functions

virtual void v_InitObject ()
 
void v_AssembleSchurComplement (const boost::shared_ptr< AssemblyMap > locToGloMap)
 Assemble the Schur complement matrix. More...
 
void PrepareLocalSchurComplement ()
 Prepares local representation of Schur complement stored as a sparse block-diagonal matrix. More...
 
virtual void v_DoMatrixMultiply (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
 Perform a Shur-complement matrix multiply operation. More...
 
virtual void v_UniqueMap ()
 

Private Attributes

DNekScalBlkMatSharedPtr m_S1Blk
 
std::vector< double > m_storage
 Dense storage for block Schur complement matrix. More...
 
std::vector< const double * > m_denseBlocks
 Vector of pointers to local matrix data. More...
 
Array< OneD, unsigned int > m_rows
 Ranks of local matrices. More...
 
Array< OneD, NekDoublem_scale
 Scaling factors for local matrices. More...
 
DNekSmvBsrDiagBlkMatSharedPtr m_sparseSchurCompl
 Sparse representation of Schur complement matrix at this level. More...
 

Static Private Attributes

static std::string storagedef
 Utility strings. More...
 
static std::string storagelookupIds []
 

Additional Inherited Members

- Protected Attributes inherited from Nektar::MultiRegions::GlobalLinSysIterative
Array< OneD, int > m_map
 Global to universal unique map. More...
 
NekDouble m_tolerance
 Tolerance of iterative solver. More...
 
NekDouble m_rhs_magnitude
 dot product of rhs to normalise stopping criterion More...
 
PreconditionerSharedPtr m_precon
 
MultiRegions::PreconditionerType m_precontype
 
int m_totalIterations
 
bool m_useProjection
 Whether to apply projection technique. More...
 
bool m_root
 Provide verbose output and root if parallel. More...
 
bool m_verbose
 
boost::circular_buffer< Array< OneD, NekDouble > > m_prevLinSol
 Storage for solutions to previous linear problems. More...
 
int m_numPrevSols
 Total counter of previous solutions. More...
 
- Protected Attributes inherited from Nektar::MultiRegions::GlobalLinSys
const GlobalLinSysKey m_linSysKey
 Key associated with this linear system. More...
 
const boost::weak_ptr< ExpListm_expList
 Local Matrix System. More...
 
const map< int, RobinBCInfoSharedPtrm_robinBCInfo
 Robin boundary info. More...
 
- Protected Attributes inherited from Nektar::MultiRegions::GlobalLinSysStaticCond
GlobalLinSysStaticCondSharedPtr m_recursiveSchurCompl
 Schur complement for Direct Static Condensation. More...
 
DNekScalBlkMatSharedPtr m_schurCompl
 Block Schur complement matrix. More...
 
DNekScalBlkMatSharedPtr m_BinvD
 Block $ BD^{-1} $ matrix. More...
 
DNekScalBlkMatSharedPtr m_C
 Block $ C $ matrix. More...
 
DNekScalBlkMatSharedPtr m_invD
 Block $ D^{-1} $ matrix. More...
 
boost::shared_ptr< AssemblyMapm_locToGloMap
 Local to global map. More...
 
Array< OneD, NekDoublem_wsp
 Workspace array for matrix multiplication. More...
 

Detailed Description

A global linear system.

Solves a linear system iteratively using single- or multi-level static condensation.

Definition at line 72 of file GlobalLinSysIterativeStaticCond.h.

Member Typedef Documentation

Definition at line 77 of file GlobalLinSysIterativeStaticCond.h.

Definition at line 79 of file GlobalLinSysIterativeStaticCond.h.

Constructor & Destructor Documentation

Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::GlobalLinSysIterativeStaticCond ( const GlobalLinSysKey pKey,
const boost::weak_ptr< ExpList > &  pExpList,
const boost::shared_ptr< AssemblyMap > &  pLocToGloMap 
)

Constructor for full direct matrix solve.

For a matrix system of the form

\[ \left[ \begin{array}{cc} \boldsymbol{A} & \boldsymbol{B}\\ \boldsymbol{C} & \boldsymbol{D} \end{array} \right] \left[ \begin{array}{c} \boldsymbol{x_1}\\ \boldsymbol{x_2} \end{array}\right] = \left[ \begin{array}{c} \boldsymbol{y_1}\\ \boldsymbol{y_2} \end{array}\right], \]

where $\boldsymbol{D}$ and $(\boldsymbol{A-BD^{-1}C})$ are invertible, store and assemble a static condensation system, according to a given local to global mapping. #m_linSys is constructed by AssembleSchurComplement().

Parameters
mKeyAssociated matrix key.
pLocMatSysLocalMatrixSystem
locToGloMapLocal to global mapping.

Definition at line 109 of file GlobalLinSysIterativeStaticCond.cpp.

References ASSERTL1, Nektar::MultiRegions::eIterativeMultiLevelStaticCond, Nektar::MultiRegions::eIterativeStaticCond, and Nektar::MultiRegions::GlobalLinSysKey::GetGlobalSysSolnType().

113  : GlobalLinSys (pKey, pExpList, pLocToGloMap),
114  GlobalLinSysIterative (pKey, pExpList, pLocToGloMap),
115  GlobalLinSysStaticCond(pKey, pExpList, pLocToGloMap)
116  {
117  ASSERTL1((pKey.GetGlobalSysSolnType()==eIterativeStaticCond)||
118  (pKey.GetGlobalSysSolnType()==eIterativeMultiLevelStaticCond),
119  "This constructor is only valid when using static "
120  "condensation");
121  ASSERTL1(pKey.GetGlobalSysSolnType()
122  == pLocToGloMap->GetGlobalSysSolnType(),
123  "The local to global map is not set up for the requested "
124  "solution type");
125  }
GlobalLinSysIterative(const GlobalLinSysKey &pKey, const boost::weak_ptr< ExpList > &pExpList, const boost::shared_ptr< AssemblyMap > &pLocToGloMap)
Constructor for full direct matrix solve.
GlobalLinSys(const GlobalLinSysKey &pKey, const boost::weak_ptr< ExpList > &pExpList, const boost::shared_ptr< AssemblyMap > &pLocToGloMap)
Constructor for full direct matrix solve.
GlobalLinSysStaticCond(const GlobalLinSysKey &mkey, const boost::weak_ptr< ExpList > &pExpList, const boost::shared_ptr< AssemblyMap > &locToGloMap)
Constructor for full direct matrix solve.
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191
Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::GlobalLinSysIterativeStaticCond ( const GlobalLinSysKey mkey,
const boost::weak_ptr< ExpList > &  pExpList,
const DNekScalBlkMatSharedPtr  pSchurCompl,
const DNekScalBlkMatSharedPtr  pBinvD,
const DNekScalBlkMatSharedPtr  pC,
const DNekScalBlkMatSharedPtr  pInvD,
const boost::shared_ptr< AssemblyMap > &  locToGloMap,
const PreconditionerSharedPtr  pPrecon 
)

Constructor for full direct matrix solve.

Definition at line 131 of file GlobalLinSysIterativeStaticCond.cpp.

References Nektar::MultiRegions::GlobalLinSysStaticCond::m_BinvD, Nektar::MultiRegions::GlobalLinSysStaticCond::m_C, Nektar::MultiRegions::GlobalLinSysStaticCond::m_invD, Nektar::MultiRegions::GlobalLinSysIterative::m_precon, m_S1Blk, and Nektar::MultiRegions::GlobalLinSysStaticCond::m_schurCompl.

140  : GlobalLinSys (pKey, pExpList, pLocToGloMap),
141  GlobalLinSysIterative (pKey, pExpList, pLocToGloMap),
142  GlobalLinSysStaticCond(pKey, pExpList, pLocToGloMap)
143  {
144  m_schurCompl = pSchurCompl;
145  m_S1Blk = pSchurCompl;
146  m_BinvD = pBinvD;
147  m_C = pC;
148  m_invD = pInvD;
149  m_precon = pPrecon;
150  }
DNekScalBlkMatSharedPtr m_schurCompl
Block Schur complement matrix.
DNekScalBlkMatSharedPtr m_C
Block matrix.
DNekScalBlkMatSharedPtr m_invD
Block matrix.
GlobalLinSysIterative(const GlobalLinSysKey &pKey, const boost::weak_ptr< ExpList > &pExpList, const boost::shared_ptr< AssemblyMap > &pLocToGloMap)
Constructor for full direct matrix solve.
GlobalLinSys(const GlobalLinSysKey &pKey, const boost::weak_ptr< ExpList > &pExpList, const boost::shared_ptr< AssemblyMap > &pLocToGloMap)
Constructor for full direct matrix solve.
GlobalLinSysStaticCond(const GlobalLinSysKey &mkey, const boost::weak_ptr< ExpList > &pExpList, const boost::shared_ptr< AssemblyMap > &locToGloMap)
Constructor for full direct matrix solve.
DNekScalBlkMatSharedPtr m_BinvD
Block matrix.
Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::~GlobalLinSysIterativeStaticCond ( )
virtual

Definition at line 205 of file GlobalLinSysIterativeStaticCond.cpp.

206  {
207 
208  }

Member Function Documentation

static GlobalLinSysSharedPtr Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::create ( const GlobalLinSysKey pLinSysKey,
const boost::weak_ptr< ExpList > &  pExpList,
const boost::shared_ptr< AssemblyMap > &  pLocToGloMap 
)
inlinestatic

Creates an instance of this class.

Definition at line 82 of file GlobalLinSysIterativeStaticCond.h.

86  {
89  AllocateSharedPtr(pLinSysKey, pExpList, pLocToGloMap);
90  p->InitObject();
91  return p;
92  }
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
GlobalLinSysIterativeStaticCond(const GlobalLinSysKey &mkey, const boost::weak_ptr< ExpList > &pExpList, const boost::shared_ptr< AssemblyMap > &locToGloMap)
Constructor for full direct matrix solve.
boost::shared_ptr< GlobalLinSys > GlobalLinSysSharedPtr
Pointer to a GlobalLinSys object.
Definition: GlobalLinSys.h:51
void Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::PrepareLocalSchurComplement ( )
private

Prepares local representation of Schur complement stored as a sparse block-diagonal matrix.

Populates sparse block-diagonal schur complement matrix from the block matrices stored in #m_blkMatrices.

Definition at line 330 of file GlobalLinSysIterativeStaticCond.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL1, Nektar::MultiRegions::eContiguous, Nektar::eFULL, Nektar::MultiRegions::eNonContiguous, Nektar::MultiRegions::eSparse, Nektar::NekConstants::kNekZeroTol, m_denseBlocks, Nektar::MultiRegions::GlobalLinSys::m_expList, m_rows, m_scale, Nektar::MultiRegions::GlobalLinSysStaticCond::m_schurCompl, m_sparseSchurCompl, m_storage, Vmath::Smul(), and Nektar::MultiRegions::GlobalLinSys::v_DropStaticCondBlock().

Referenced by v_AssembleSchurComplement().

331  {
332  LocalMatrixStorageStrategy storageStrategy =
333  m_expList.lock()->GetSession()->
334  GetSolverInfoAsEnum<LocalMatrixStorageStrategy>(
335  "LocalMatrixStorageStrategy");
336 
337  switch(storageStrategy)
338  {
341  {
342  size_t storageSize = 0;
343  int nBlk = m_schurCompl->GetNumberOfBlockRows();
344 
345  m_scale = Array<OneD, NekDouble> (nBlk, 1.0);
346  m_rows = Array<OneD, unsigned int> (nBlk, 0U);
347 
348  // Determine storage requirements for dense blocks.
349  for (int i = 0; i < nBlk; ++i)
350  {
351  m_rows[i] = m_schurCompl->GetBlock(i,i)->GetRows();
352  m_scale[i] = m_schurCompl->GetBlock(i,i)->Scale();
353  storageSize += m_rows[i] * m_rows[i];
354  }
355 
356  // Assemble dense storage blocks.
357  DNekScalMatSharedPtr loc_mat;
358  m_denseBlocks.resize(nBlk);
359  double *ptr = 0;
360 
361  if (MultiRegions::eContiguous == storageStrategy)
362  {
363  m_storage.resize (storageSize);
364  ptr = &m_storage[0];
365  }
366 
367  for (unsigned int n = 0; n < nBlk; ++n)
368  {
369  loc_mat = m_schurCompl->GetBlock(n,n);
370 
371  if (MultiRegions::eContiguous == storageStrategy)
372  {
373  int loc_lda = loc_mat->GetRows();
374  int blockSize = loc_lda * loc_lda;
375  m_denseBlocks[n] = ptr;
376  for(int i = 0; i < loc_lda; ++i)
377  {
378  for(int j = 0; j < loc_lda; ++j)
379  {
380  ptr[j*loc_lda+i] = (*loc_mat)(i,j);
381  }
382  }
383  ptr += blockSize;
385  m_expList.lock()->GetOffset_Elmt_Id(n));
386  }
387  else
388  {
389  m_denseBlocks[n] = loc_mat->GetRawPtr();
390  }
391  }
392  break;
393  }
395  {
396  DNekScalMatSharedPtr loc_mat;
397  int loc_lda;
398  int blockSize = 0;
399 
400  // First run through to split the set of local matrices into
401  // partitions of fixed block size, and count number of local
402  // matrices that belong to each partition.
403  std::vector<std::pair<int,int> > partitions;
404  for(int n = 0; n < m_schurCompl->GetNumberOfBlockRows(); ++n)
405  {
406  loc_mat = m_schurCompl->GetBlock(n,n);
407  loc_lda = loc_mat->GetRows();
408 
409  ASSERTL1(loc_lda >= 0,
410  boost::lexical_cast<std::string>(n) + "-th "
411  "matrix block in Schur complement has "
412  "rank 0!");
413 
414  if (blockSize == loc_lda)
415  {
416  partitions[partitions.size()-1].first++;
417  }
418  else
419  {
420  blockSize = loc_lda;
421  partitions.push_back(make_pair(1,loc_lda));
422  }
423  }
424 
425  MatrixStorage matStorage = eFULL;
426 
427  // Create a vector of sparse storage holders
429  sparseStorage (partitions.size());
430 
431  for (int part = 0, n = 0; part < partitions.size(); ++part)
432  {
433  BCOMatType partMat;
434 
435  for(int k = 0; k < partitions[part].first; ++k, ++n)
436  {
437  loc_mat = m_schurCompl->GetBlock(n,n);
438  loc_lda = loc_mat->GetRows();
439 
440  ASSERTL1(loc_lda == partitions[part].second,
441  boost::lexical_cast<std::string>(n) + "-th"
442  " matrix block in Schur complement has "
443  "unexpected rank");
444 
445  NekDouble scale = loc_mat->Scale();
446  if(fabs(scale-1.0) > NekConstants::kNekZeroTol)
447  {
448  Array<OneD, NekDouble> matarray(loc_lda*loc_lda);
449  Vmath::Smul(loc_lda*loc_lda,scale,
450  loc_mat->GetRawPtr(),1,&matarray[0],1);
451  partMat[make_pair(k,k)] = BCOEntryType(matarray);
452  }
453  else // scale factor is 1.0
454  {
455  partMat[make_pair(k,k)] = BCOEntryType(
456  loc_lda*loc_lda, loc_mat->GetRawPtr());
457  }
458 
460  m_expList.lock()->GetOffset_Elmt_Id(n));
461  }
462 
463  sparseStorage[part] =
466  partitions[part].first, partitions[part].first,
467  partitions[part].second, partMat, matStorage );
468  }
469 
470  // Create block diagonal matrix
472  AllocateSharedPtr(sparseStorage);
473 
474  break;
475  }
476  default:
477  ErrorUtil::NekError("Solver info property \
478  LocalMatrixStorageStrategy takes values \
479  Contiguous, Non-contiguous and Sparse");
480  }
481  }
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
DNekScalBlkMatSharedPtr m_schurCompl
Block Schur complement matrix.
std::vector< double > m_storage
Dense storage for block Schur complement matrix.
boost::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
Array< OneD, NekDouble > BCOEntryType
static const NekDouble kNekZeroTol
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:199
Array< OneD, NekDouble > m_scale
Scaling factors for local matrices.
std::map< CoordType, BCOEntryType > BCOMatType
Array< OneD, unsigned int > m_rows
Ranks of local matrices.
double NekDouble
std::vector< const double * > m_denseBlocks
Vector of pointers to local matrix data.
DNekSmvBsrDiagBlkMatSharedPtr m_sparseSchurCompl
Sparse representation of Schur complement matrix at this level.
virtual void v_DropStaticCondBlock(unsigned int n)
Releases the static condensation block matrices from NekManager of n-th expansion using the matrix ke...
Array< OneD, SparseStorageSharedPtr > SparseStorageSharedPtrVector
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191
const boost::weak_ptr< ExpList > m_expList
Local Matrix System.
Definition: GlobalLinSys.h:125
void Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::v_AssembleSchurComplement ( const boost::shared_ptr< AssemblyMap locToGloMap)
privatevirtual

Assemble the Schur complement matrix.

Assemble the schur complement matrix from the block matrices stored in #m_blkMatrices and the given local to global mapping information.

Parameters
locToGloMapLocal to global mapping information.

Reimplemented from Nektar::MultiRegions::GlobalLinSysStaticCond.

Definition at line 233 of file GlobalLinSysIterativeStaticCond.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::convertCooToBco(), Nektar::LibUtilities::NekFactory< tKey, tBase, >::CreateInstance(), Nektar::eFULL, Nektar::MultiRegions::eIterativeMultiLevelStaticCond, Nektar::MultiRegions::GlobalLinSysKey::GetGlobalSysSolnType(), Nektar::MultiRegions::GlobalMatrixKey::GetMatrixType(), Nektar::MultiRegions::GetPreconFactory(), Nektar::MultiRegions::GlobalLinSys::GetSharedThisPtr(), Nektar::MultiRegions::GlobalLinSys::m_expList, Nektar::MultiRegions::GlobalLinSys::m_linSysKey, Nektar::MultiRegions::GlobalLinSysStaticCond::m_locToGloMap, Nektar::MultiRegions::GlobalLinSysIterative::m_precon, Nektar::MultiRegions::GlobalLinSysStaticCond::m_schurCompl, m_sparseSchurCompl, Nektar::MultiRegions::PreconditionerTypeMap, PrepareLocalSchurComplement(), and v_UniqueMap().

235  {
236  int i,j,n,cnt,gid1,gid2;
237  NekDouble sign1,sign2;
238 
239  bool doGlobalOp = m_expList.lock()->GetGlobalOptParam()->
240  DoGlobalMatOp(m_linSysKey.GetMatrixType());
241 
242  // Set up unique map
243  v_UniqueMap();
244 
245  // Build precon again if we in multi-level static condensation (a
246  // bit of a hack)
248  {
250  = m_locToGloMap->GetPreconType();
251  std::string PreconType
254  PreconType,GetSharedThisPtr(),m_locToGloMap);
255  m_precon->BuildPreconditioner();
256  }
257 
258  if (!doGlobalOp)
259  {
261  return;
262  }
263 
264  int nBndDofs = pLocToGloMap->GetNumGlobalBndCoeffs();
265  int NumDirBCs = pLocToGloMap->GetNumGlobalDirBndCoeffs();
266  unsigned int rows = nBndDofs - NumDirBCs;
267  unsigned int cols = nBndDofs - NumDirBCs;
268 
269  // COO sparse storage to assist in assembly
270  COOMatType gmat_coo;
271 
272  // Get the matrix storage structure
273  // (whether to store only one triangular part, if symmetric)
274  MatrixStorage matStorage = eFULL;
275 
276  // assemble globally
277  DNekScalMatSharedPtr loc_mat;
278  int loc_lda;
279  for(n = cnt = 0; n < m_schurCompl->GetNumberOfBlockRows(); ++n)
280  {
281  loc_mat = m_schurCompl->GetBlock(n,n);
282  loc_lda = loc_mat->GetRows();
283 
284  // Set up Matrix;
285  for(i = 0; i < loc_lda; ++i)
286  {
287  gid1 = pLocToGloMap->GetLocalToGlobalBndMap (cnt + i)
288  - NumDirBCs;
289  sign1 = pLocToGloMap->GetLocalToGlobalBndSign(cnt + i);
290 
291  if(gid1 >= 0)
292  {
293  for(j = 0; j < loc_lda; ++j)
294  {
295  gid2 = pLocToGloMap->GetLocalToGlobalBndMap(cnt+j)
296  - NumDirBCs;
297  sign2 = pLocToGloMap->GetLocalToGlobalBndSign(cnt+j);
298 
299  if (gid2 >= 0)
300  {
301  gmat_coo[std::make_pair(gid1,gid2)] +=
302  sign1*sign2*(*loc_mat)(i,j);
303  }
304  }
305  }
306  }
307  cnt += loc_lda;
308  }
309 
311  sparseStorage (1);
312 
313  BCOMatType partMat;
314  convertCooToBco(rows, cols, 1, gmat_coo, partMat);
315 
316  sparseStorage[0] =
318  AllocateSharedPtr(rows, cols, 1, partMat, matStorage );
319 
320  // Create block diagonal matrix
322  AllocateSharedPtr(sparseStorage);
323  }
GlobalSysSolnType GetGlobalSysSolnType() const
Return the associated solution type.
tBaseSharedPtr CreateInstance(tKey idKey BOOST_PP_COMMA_IF(MAX_PARAM) BOOST_PP_ENUM_BINARY_PARAMS(MAX_PARAM, tParam, x))
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:162
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
std::map< CoordType, NekDouble > COOMatType
DNekScalBlkMatSharedPtr m_schurCompl
Block Schur complement matrix.
PreconFactory & GetPreconFactory()
boost::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
std::map< CoordType, BCOEntryType > BCOMatType
boost::shared_ptr< GlobalLinSys > GetSharedThisPtr()
Returns a shared pointer to the current object.
Definition: GlobalLinSys.h:103
double NekDouble
void PrepareLocalSchurComplement()
Prepares local representation of Schur complement stored as a sparse block-diagonal matrix...
StdRegions::MatrixType GetMatrixType() const
Return the matrix type.
const GlobalLinSysKey m_linSysKey
Key associated with this linear system.
Definition: GlobalLinSys.h:123
DNekSmvBsrDiagBlkMatSharedPtr m_sparseSchurCompl
Sparse representation of Schur complement matrix at this level.
Array< OneD, SparseStorageSharedPtr > SparseStorageSharedPtrVector
void convertCooToBco(const unsigned int blkRows, const unsigned int blkColumns, const unsigned int blkDim, const COOMatType &cooMat, BCOMatType &bcoMat)
Definition: SparseUtils.cpp:48
boost::shared_ptr< AssemblyMap > m_locToGloMap
Local to global map.
const char *const PreconditionerTypeMap[]
const boost::weak_ptr< ExpList > m_expList
Local Matrix System.
Definition: GlobalLinSys.h:125
void Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::v_BasisInvTransform ( Array< OneD, NekDouble > &  pInOut)
protectedvirtual

Reimplemented from Nektar::MultiRegions::GlobalLinSysStaticCond.

Definition at line 571 of file GlobalLinSysIterativeStaticCond.cpp.

References Nektar::MultiRegions::GlobalLinSysIterative::m_precon.

573  {
574  m_precon->DoTransformFromLowEnergy(pInOut);
575  }
void Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::v_BasisTransform ( Array< OneD, NekDouble > &  pInOut,
int  offset 
)
protectedvirtual

Reimplemented from Nektar::MultiRegions::GlobalLinSysStaticCond.

Definition at line 564 of file GlobalLinSysIterativeStaticCond.cpp.

References Nektar::MultiRegions::GlobalLinSysIterative::m_precon.

567  {
568  m_precon->DoTransformToLowEnergy(pInOut, offset);
569  }
void Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::v_DoMatrixMultiply ( const Array< OneD, NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput 
)
privatevirtual

Perform a Shur-complement matrix multiply operation.

Implements Nektar::MultiRegions::GlobalLinSysIterative.

Definition at line 486 of file GlobalLinSysIterativeStaticCond.cpp.

References Nektar::MultiRegions::GlobalMatrixKey::GetMatrixType(), m_denseBlocks, Nektar::MultiRegions::GlobalLinSys::m_expList, Nektar::MultiRegions::GlobalLinSys::m_linSysKey, Nektar::MultiRegions::GlobalLinSysStaticCond::m_locToGloMap, m_rows, m_scale, m_sparseSchurCompl, and Nektar::MultiRegions::GlobalLinSysStaticCond::m_wsp.

489  {
490  int nLocal = m_locToGloMap->GetNumLocalBndCoeffs();
491  int nDir = m_locToGloMap->GetNumGlobalDirBndCoeffs();
492  bool doGlobalOp = m_expList.lock()->GetGlobalOptParam()->
493  DoGlobalMatOp(m_linSysKey.GetMatrixType());
494 
495  if(doGlobalOp)
496  {
497  // Do matrix multiply globally
498  Array<OneD, NekDouble> in = pInput + nDir;
499  Array<OneD, NekDouble> out = pOutput + nDir;
500 
501  m_sparseSchurCompl->Multiply(in,out);
502  m_locToGloMap->UniversalAssembleBnd(pOutput, nDir);
503  }
504  else if (m_sparseSchurCompl)
505  {
506  // Do matrix multiply locally using block-diagonal sparse matrix
507  Array<OneD, NekDouble> tmp = m_wsp + nLocal;
508 
509  m_locToGloMap->GlobalToLocalBnd(pInput, m_wsp);
510  m_sparseSchurCompl->Multiply(m_wsp,tmp);
511  m_locToGloMap->AssembleBnd(tmp, pOutput);
512  }
513  else
514  {
515  // Do matrix multiply locally, using direct BLAS calls
516  m_locToGloMap->GlobalToLocalBnd(pInput, m_wsp);
517  int i, cnt;
518  Array<OneD, NekDouble> tmpout = m_wsp + nLocal;
519  for (i = cnt = 0; i < m_denseBlocks.size(); cnt += m_rows[i], ++i)
520  {
521  const int rows = m_rows[i];
522  Blas::Dgemv('N', rows, rows,
523  m_scale[i], m_denseBlocks[i], rows,
524  m_wsp.get()+cnt, 1,
525  0.0, tmpout.get()+cnt, 1);
526  }
527  m_locToGloMap->AssembleBnd(tmpout, pOutput);
528  }
529  }
Array< OneD, NekDouble > m_wsp
Workspace array for matrix multiplication.
Array< OneD, NekDouble > m_scale
Scaling factors for local matrices.
Array< OneD, unsigned int > m_rows
Ranks of local matrices.
std::vector< const double * > m_denseBlocks
Vector of pointers to local matrix data.
StdRegions::MatrixType GetMatrixType() const
Return the matrix type.
const GlobalLinSysKey m_linSysKey
Key associated with this linear system.
Definition: GlobalLinSys.h:123
DNekSmvBsrDiagBlkMatSharedPtr m_sparseSchurCompl
Sparse representation of Schur complement matrix at this level.
boost::shared_ptr< AssemblyMap > m_locToGloMap
Local to global map.
const boost::weak_ptr< ExpList > m_expList
Local Matrix System.
Definition: GlobalLinSys.h:125
DNekScalBlkMatSharedPtr Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::v_GetStaticCondBlock ( unsigned int  n)
protectedvirtual

Retrieves a the static condensation block matrices from n-th expansion using the matrix key provided by the m_linSysKey.

Parameters
nNumber of the expansion
Returns
2x2 Block matrix holding the static condensation matrices for the n-th expansion.

Reimplemented from Nektar::MultiRegions::GlobalLinSys.

Definition at line 211 of file GlobalLinSysIterativeStaticCond.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::MultiRegions::GlobalLinSysStaticCond::m_locToGloMap, m_S1Blk, and Nektar::MultiRegions::GlobalLinSysStaticCond::m_schurCompl.

212  {
213  DNekScalBlkMatSharedPtr schurComplBlock;
214  int scLevel = m_locToGloMap->GetStaticCondLevel();
215  DNekScalBlkMatSharedPtr sc = scLevel == 0 ? m_S1Blk : m_schurCompl;
216  DNekScalMatSharedPtr localMat = sc->GetBlock(n,n);
217  unsigned int nbdry = localMat->GetRows();
218  unsigned int nblks = 1;
219  unsigned int esize[1] = {nbdry};
220 
221  schurComplBlock = MemoryManager<DNekScalBlkMat>
222  ::AllocateSharedPtr(nblks, nblks, esize, esize);
223  schurComplBlock->SetBlock(0, 0, localMat);
224 
225  return schurComplBlock;
226  }
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
DNekScalBlkMatSharedPtr m_schurCompl
Block Schur complement matrix.
boost::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
boost::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
Definition: NekTypeDefs.hpp:74
boost::shared_ptr< AssemblyMap > m_locToGloMap
Local to global map.
void Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::v_InitObject ( )
privatevirtual

Reimplemented from Nektar::MultiRegions::GlobalLinSysStaticCond.

Definition at line 153 of file GlobalLinSysIterativeStaticCond.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::LibUtilities::NekFactory< tKey, tBase, >::CreateInstance(), Nektar::eDIAGONAL, Nektar::StdRegions::eHybridDGHelmBndLam, Nektar::MultiRegions::GlobalMatrixKey::GetMatrixType(), Nektar::MultiRegions::GetPreconFactory(), Nektar::MultiRegions::GlobalLinSys::GetSharedThisPtr(), Nektar::MultiRegions::GlobalLinSys::Initialise(), Nektar::MultiRegions::GlobalLinSys::m_expList, Nektar::MultiRegions::GlobalLinSys::m_linSysKey, Nektar::MultiRegions::GlobalLinSysStaticCond::m_locToGloMap, Nektar::MultiRegions::GlobalLinSysIterative::m_precon, m_S1Blk, Nektar::MultiRegions::GlobalLinSysStaticCond::m_schurCompl, Nektar::MultiRegions::PreconditionerTypeMap, and Nektar::MultiRegions::GlobalLinSysStaticCond::SetupTopLevel().

154  {
156  = m_locToGloMap->GetPreconType();
157  std::string PreconType
160  PreconType,GetSharedThisPtr(),m_locToGloMap);
161 
162  // Allocate memory for top-level structure
164 
165  // Setup Block Matrix systems
166  int n, n_exp = m_expList.lock()->GetNumElmts();
167 
168  MatrixStorage blkmatStorage = eDIAGONAL;
169  const Array<OneD,const unsigned int>& nbdry_size
170  = m_locToGloMap->GetNumLocalBndCoeffsPerPatch();
171 
173  ::AllocateSharedPtr(nbdry_size, nbdry_size , blkmatStorage);
174 
175  // Preserve original matrix in m_S1Blk
176  for (n = 0; n < n_exp; ++n)
177  {
178  DNekScalMatSharedPtr mat = m_schurCompl->GetBlock(n, n);
179  m_S1Blk->SetBlock(n, n, mat);
180  }
181 
182  // Build preconditioner
183  m_precon->BuildPreconditioner();
184 
185  // Do transform of Schur complement matrix
186  for (n = 0; n < n_exp; ++n)
187  {
188  if (m_linSysKey.GetMatrixType() !=
190  {
191  DNekScalMatSharedPtr mat = m_S1Blk->GetBlock(n, n);
192  DNekScalMatSharedPtr t = m_precon->TransformedSchurCompl(
193  m_expList.lock()->GetOffset_Elmt_Id(n), mat);
194  m_schurCompl->SetBlock(n, n, t);
195  }
196  }
197 
198  // Construct this level
200  }
tBaseSharedPtr CreateInstance(tKey idKey BOOST_PP_COMMA_IF(MAX_PARAM) BOOST_PP_ENUM_BINARY_PARAMS(MAX_PARAM, tParam, x))
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:162
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
DNekScalBlkMatSharedPtr m_schurCompl
Block Schur complement matrix.
PreconFactory & GetPreconFactory()
void Initialise(const boost::shared_ptr< AssemblyMap > &pLocToGloMap)
Definition: GlobalLinSys.h:208
boost::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
boost::shared_ptr< GlobalLinSys > GetSharedThisPtr()
Returns a shared pointer to the current object.
Definition: GlobalLinSys.h:103
void SetupTopLevel(const boost::shared_ptr< AssemblyMap > &locToGloMap)
Set up the storage for the Schur complement or the top level of the multi-level Schur complement...
StdRegions::MatrixType GetMatrixType() const
Return the matrix type.
const GlobalLinSysKey m_linSysKey
Key associated with this linear system.
Definition: GlobalLinSys.h:123
boost::shared_ptr< AssemblyMap > m_locToGloMap
Local to global map.
const char *const PreconditionerTypeMap[]
const boost::weak_ptr< ExpList > m_expList
Local Matrix System.
Definition: GlobalLinSys.h:125
DNekScalBlkMatSharedPtr Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::v_PreSolve ( int  scLevel,
NekVector< NekDouble > &  F_GlobBnd 
)
protectedvirtual

Reimplemented from Nektar::MultiRegions::GlobalLinSysStaticCond.

Definition at line 536 of file GlobalLinSysIterativeStaticCond.cpp.

References Nektar::LibUtilities::NekFactory< tKey, tBase, >::CreateInstance(), Nektar::MultiRegions::GetPreconFactory(), Nektar::MultiRegions::GlobalLinSys::GetSharedThisPtr(), Nektar::MultiRegions::GlobalLinSysStaticCond::m_locToGloMap, Nektar::MultiRegions::GlobalLinSysIterative::m_precon, m_S1Blk, Nektar::MultiRegions::GlobalLinSysStaticCond::m_schurCompl, Nektar::MultiRegions::PreconditionerTypeMap, and Nektar::MultiRegions::GlobalLinSysIterative::Set_Rhs_Magnitude().

539  {
540  if (scLevel == 0)
541  {
542  // When matrices are supplied to the constructor at the top
543  // level, the preconditioner is never set up.
544  if (!m_precon)
545  {
547  = m_locToGloMap->GetPreconType();
548  std::string PreconType
551  PreconType, GetSharedThisPtr(), m_locToGloMap);
552  m_precon->BuildPreconditioner();
553  }
554 
555  Set_Rhs_Magnitude(F_GlobBnd);
556  return m_S1Blk;
557  }
558  else
559  {
560  return m_schurCompl;
561  }
562  }
tBaseSharedPtr CreateInstance(tKey idKey BOOST_PP_COMMA_IF(MAX_PARAM) BOOST_PP_ENUM_BINARY_PARAMS(MAX_PARAM, tParam, x))
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:162
DNekScalBlkMatSharedPtr m_schurCompl
Block Schur complement matrix.
void Set_Rhs_Magnitude(const NekVector< NekDouble > &pIn)
PreconFactory & GetPreconFactory()
boost::shared_ptr< GlobalLinSys > GetSharedThisPtr()
Returns a shared pointer to the current object.
Definition: GlobalLinSys.h:103
boost::shared_ptr< AssemblyMap > m_locToGloMap
Local to global map.
const char *const PreconditionerTypeMap[]
GlobalLinSysStaticCondSharedPtr Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::v_Recurse ( const GlobalLinSysKey mkey,
const boost::weak_ptr< ExpList > &  pExpList,
const DNekScalBlkMatSharedPtr  pSchurCompl,
const DNekScalBlkMatSharedPtr  pBinvD,
const DNekScalBlkMatSharedPtr  pC,
const DNekScalBlkMatSharedPtr  pInvD,
const boost::shared_ptr< AssemblyMap > &  locToGloMap 
)
protectedvirtual

Implements Nektar::MultiRegions::GlobalLinSysStaticCond.

Definition at line 577 of file GlobalLinSysIterativeStaticCond.cpp.

References Nektar::MultiRegions::GlobalLinSysIterative::m_precon.

585  {
587  GlobalLinSysIterativeStaticCond>::AllocateSharedPtr(
588  mkey, pExpList, pSchurCompl, pBinvD, pC, pInvD, l2gMap,
589  m_precon);
590  sys->Initialise(l2gMap);
591  return sys;
592  }
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
GlobalLinSysIterativeStaticCond(const GlobalLinSysKey &mkey, const boost::weak_ptr< ExpList > &pExpList, const boost::shared_ptr< AssemblyMap > &locToGloMap)
Constructor for full direct matrix solve.
boost::shared_ptr< GlobalLinSysIterativeStaticCond > GlobalLinSysIterativeStaticCondSharedPtr
void Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::v_UniqueMap ( )
privatevirtual

Implements Nektar::MultiRegions::GlobalLinSysIterative.

Definition at line 531 of file GlobalLinSysIterativeStaticCond.cpp.

References Nektar::MultiRegions::GlobalLinSysStaticCond::m_locToGloMap, and Nektar::MultiRegions::GlobalLinSysIterative::m_map.

Referenced by v_AssembleSchurComplement().

532  {
533  m_map = m_locToGloMap->GetGlobalToUniversalBndMapUnique();
534  }
Array< OneD, int > m_map
Global to universal unique map.
boost::shared_ptr< AssemblyMap > m_locToGloMap
Local to global map.

Member Data Documentation

string Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::className
static
Initial value:
"IterativeStaticCond",
"Iterative static condensation.")

Name of class.

Registers the class with the Factory.

Definition at line 95 of file GlobalLinSysIterativeStaticCond.h.

string Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::className2
static
Initial value:
"IterativeMultiLevelStaticCond",
"Iterative multi-level static condensation.")

Definition at line 96 of file GlobalLinSysIterativeStaticCond.h.

std::vector<const double*> Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::m_denseBlocks
private

Vector of pointers to local matrix data.

Definition at line 142 of file GlobalLinSysIterativeStaticCond.h.

Referenced by PrepareLocalSchurComplement(), and v_DoMatrixMultiply().

Array<OneD, unsigned int> Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::m_rows
private

Ranks of local matrices.

Definition at line 144 of file GlobalLinSysIterativeStaticCond.h.

Referenced by PrepareLocalSchurComplement(), and v_DoMatrixMultiply().

DNekScalBlkMatSharedPtr Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::m_S1Blk
private
Array<OneD, NekDouble> Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::m_scale
private

Scaling factors for local matrices.

Definition at line 146 of file GlobalLinSysIterativeStaticCond.h.

Referenced by PrepareLocalSchurComplement(), and v_DoMatrixMultiply().

DNekSmvBsrDiagBlkMatSharedPtr Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::m_sparseSchurCompl
private

Sparse representation of Schur complement matrix at this level.

Definition at line 148 of file GlobalLinSysIterativeStaticCond.h.

Referenced by PrepareLocalSchurComplement(), v_AssembleSchurComplement(), and v_DoMatrixMultiply().

std::vector<double> Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::m_storage
private

Dense storage for block Schur complement matrix.

Definition at line 140 of file GlobalLinSysIterativeStaticCond.h.

Referenced by PrepareLocalSchurComplement().

std::string Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::storagedef
staticprivate
Initial value:
=
"LocalMatrixStorageStrategy",
"Sparse")

Utility strings.

Definition at line 150 of file GlobalLinSysIterativeStaticCond.h.

std::string Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::storagelookupIds
staticprivate
Initial value:
= {
"LocalMatrixStorageStrategy",
"Contiguous",
"LocalMatrixStorageStrategy",
"Non-contiguous",
"LocalMatrixStorageStrategy",
"Sparse",
}

Definition at line 151 of file GlobalLinSysIterativeStaticCond.h.