Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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
< DNekSmvBsrDiagBlkMat
DNekSmvBsrDiagBlkMatSharedPtr
 

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...
 
PreconditionerSharedPtr CreatePrecon (AssemblyMapSharedPtr asmMap)
 Create a preconditioner object from the parameters defined in the supplied assembly map. 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...
 
int m_maxiter
 maximum iterations More...
 
NekDouble m_tolerance
 Tolerance of iterative solver. More...
 
NekDouble m_rhs_magnitude
 dot product of rhs to normalise stopping criterion More...
 
NekDouble m_rhs_mag_sm
 cnt to how many times rhs_magnitude is called 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 std::map< int,
RobinBCInfoSharedPtr
m_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 111 of file GlobalLinSysIterativeStaticCond.cpp.

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

115  : GlobalLinSys (pKey, pExpList, pLocToGloMap),
116  GlobalLinSysIterative (pKey, pExpList, pLocToGloMap),
117  GlobalLinSysStaticCond(pKey, pExpList, pLocToGloMap)
118  {
119  ASSERTL1((pKey.GetGlobalSysSolnType()==eIterativeStaticCond)||
120  (pKey.GetGlobalSysSolnType()==eIterativeMultiLevelStaticCond),
121  "This constructor is only valid when using static "
122  "condensation");
123  ASSERTL1(pKey.GetGlobalSysSolnType()
124  == pLocToGloMap->GetGlobalSysSolnType(),
125  "The local to global map is not set up for the requested "
126  "solution type");
127  }
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 133 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.

142  : GlobalLinSys (pKey, pExpList, pLocToGloMap),
143  GlobalLinSysIterative (pKey, pExpList, pLocToGloMap),
144  GlobalLinSysStaticCond(pKey, pExpList, pLocToGloMap)
145  {
146  m_schurCompl = pSchurCompl;
147  m_S1Blk = pSchurCompl;
148  m_BinvD = pBinvD;
149  m_C = pC;
150  m_invD = pInvD;
151  m_precon = pPrecon;
152  }
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 202 of file GlobalLinSysIterativeStaticCond.cpp.

203  {
204 
205  }

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  {
87  GlobalLinSysSharedPtr p = MemoryManager<
89  AllocateSharedPtr(pLinSysKey, pExpList, pLocToGloMap);
90  p->InitObject();
91  return p;
92  }
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:52
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 323 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().

324  {
325  LocalMatrixStorageStrategy storageStrategy =
326  m_expList.lock()->GetSession()->
327  GetSolverInfoAsEnum<LocalMatrixStorageStrategy>(
328  "LocalMatrixStorageStrategy");
329 
330  switch(storageStrategy)
331  {
334  {
335  size_t storageSize = 0;
336  int nBlk = m_schurCompl->GetNumberOfBlockRows();
337 
338  m_scale = Array<OneD, NekDouble> (nBlk, 1.0);
339  m_rows = Array<OneD, unsigned int> (nBlk, 0U);
340 
341  // Determine storage requirements for dense blocks.
342  for (int i = 0; i < nBlk; ++i)
343  {
344  m_rows[i] = m_schurCompl->GetBlock(i,i)->GetRows();
345  m_scale[i] = m_schurCompl->GetBlock(i,i)->Scale();
346  storageSize += m_rows[i] * m_rows[i];
347  }
348 
349  // Assemble dense storage blocks.
350  DNekScalMatSharedPtr loc_mat;
351  m_denseBlocks.resize(nBlk);
352  double *ptr = 0;
353 
354  if (MultiRegions::eContiguous == storageStrategy)
355  {
356  m_storage.resize (storageSize);
357  ptr = &m_storage[0];
358  }
359 
360  for (unsigned int n = 0; n < nBlk; ++n)
361  {
362  loc_mat = m_schurCompl->GetBlock(n,n);
363 
364  if (MultiRegions::eContiguous == storageStrategy)
365  {
366  int loc_lda = loc_mat->GetRows();
367  int blockSize = loc_lda * loc_lda;
368  m_denseBlocks[n] = ptr;
369  for(int i = 0; i < loc_lda; ++i)
370  {
371  for(int j = 0; j < loc_lda; ++j)
372  {
373  ptr[j*loc_lda+i] = (*loc_mat)(i,j);
374  }
375  }
376  ptr += blockSize;
378  m_expList.lock()->GetOffset_Elmt_Id(n));
379  }
380  else
381  {
382  m_denseBlocks[n] = loc_mat->GetRawPtr();
383  }
384  }
385  break;
386  }
388  {
389  DNekScalMatSharedPtr loc_mat;
390  int loc_lda;
391  int blockSize = 0;
392 
393  // First run through to split the set of local matrices into
394  // partitions of fixed block size, and count number of local
395  // matrices that belong to each partition.
396  std::vector<std::pair<int,int> > partitions;
397  for(int n = 0; n < m_schurCompl->GetNumberOfBlockRows(); ++n)
398  {
399  loc_mat = m_schurCompl->GetBlock(n,n);
400  loc_lda = loc_mat->GetRows();
401 
402  ASSERTL1(loc_lda >= 0,
403  boost::lexical_cast<std::string>(n) + "-th "
404  "matrix block in Schur complement has "
405  "rank 0!");
406 
407  if (blockSize == loc_lda)
408  {
409  partitions[partitions.size()-1].first++;
410  }
411  else
412  {
413  blockSize = loc_lda;
414  partitions.push_back(make_pair(1,loc_lda));
415  }
416  }
417 
418  MatrixStorage matStorage = eFULL;
419 
420  // Create a vector of sparse storage holders
422  sparseStorage (partitions.size());
423 
424  for (int part = 0, n = 0; part < partitions.size(); ++part)
425  {
426  BCOMatType partMat;
427 
428  for(int k = 0; k < partitions[part].first; ++k, ++n)
429  {
430  loc_mat = m_schurCompl->GetBlock(n,n);
431  loc_lda = loc_mat->GetRows();
432 
433  ASSERTL1(loc_lda == partitions[part].second,
434  boost::lexical_cast<std::string>(n) + "-th"
435  " matrix block in Schur complement has "
436  "unexpected rank");
437 
438  NekDouble scale = loc_mat->Scale();
439  if(fabs(scale-1.0) > NekConstants::kNekZeroTol)
440  {
441  Array<OneD, NekDouble> matarray(loc_lda*loc_lda);
442  Vmath::Smul(loc_lda*loc_lda,scale,
443  loc_mat->GetRawPtr(),1,&matarray[0],1);
444  partMat[make_pair(k,k)] = BCOEntryType(matarray);
445  }
446  else // scale factor is 1.0
447  {
448  partMat[make_pair(k,k)] = BCOEntryType(
449  loc_lda*loc_lda, loc_mat->GetRawPtr());
450  }
451 
453  m_expList.lock()->GetOffset_Elmt_Id(n));
454  }
455 
456  sparseStorage[part] =
459  partitions[part].first, partitions[part].first,
460  partitions[part].second, partMat, matStorage );
461  }
462 
463  // Create block diagonal matrix
465  AllocateSharedPtr(sparseStorage);
466 
467  break;
468  }
469  default:
470  ErrorUtil::NekError("Solver info property \
471  LocalMatrixStorageStrategy takes values \
472  Contiguous, Non-contiguous and Sparse");
473  }
474  }
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:129
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 230 of file GlobalLinSysIterativeStaticCond.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::convertCooToBco(), Nektar::MultiRegions::GlobalLinSys::CreatePrecon(), Nektar::eFULL, Nektar::MultiRegions::eIterativeMultiLevelStaticCond, Nektar::MultiRegions::GlobalLinSysKey::GetGlobalSysSolnType(), Nektar::MultiRegions::GlobalMatrixKey::GetMatrixType(), 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, PrepareLocalSchurComplement(), and v_UniqueMap().

232  {
233  int i,j,n,cnt,gid1,gid2;
234  NekDouble sign1,sign2;
235 
236  bool doGlobalOp = m_expList.lock()->GetGlobalOptParam()->
237  DoGlobalMatOp(m_linSysKey.GetMatrixType());
238 
239  // Set up unique map
240  v_UniqueMap();
241 
242  // Build precon again if we in multi-level static condensation (a
243  // bit of a hack)
246  {
248  m_precon->BuildPreconditioner();
249  }
250 
251  if (!doGlobalOp)
252  {
254  return;
255  }
256 
257  int nBndDofs = pLocToGloMap->GetNumGlobalBndCoeffs();
258  int NumDirBCs = pLocToGloMap->GetNumGlobalDirBndCoeffs();
259  unsigned int rows = nBndDofs - NumDirBCs;
260  unsigned int cols = nBndDofs - NumDirBCs;
261 
262  // COO sparse storage to assist in assembly
263  COOMatType gmat_coo;
264 
265  // Get the matrix storage structure
266  // (whether to store only one triangular part, if symmetric)
267  MatrixStorage matStorage = eFULL;
268 
269  // assemble globally
270  DNekScalMatSharedPtr loc_mat;
271  int loc_lda;
272  for(n = cnt = 0; n < m_schurCompl->GetNumberOfBlockRows(); ++n)
273  {
274  loc_mat = m_schurCompl->GetBlock(n,n);
275  loc_lda = loc_mat->GetRows();
276 
277  // Set up Matrix;
278  for(i = 0; i < loc_lda; ++i)
279  {
280  gid1 = pLocToGloMap->GetLocalToGlobalBndMap (cnt + i)
281  - NumDirBCs;
282  sign1 = pLocToGloMap->GetLocalToGlobalBndSign(cnt + i);
283 
284  if(gid1 >= 0)
285  {
286  for(j = 0; j < loc_lda; ++j)
287  {
288  gid2 = pLocToGloMap->GetLocalToGlobalBndMap(cnt+j)
289  - NumDirBCs;
290  sign2 = pLocToGloMap->GetLocalToGlobalBndSign(cnt+j);
291 
292  if (gid2 >= 0)
293  {
294  gmat_coo[std::make_pair(gid1,gid2)] +=
295  sign1*sign2*(*loc_mat)(i,j);
296  }
297  }
298  }
299  }
300  cnt += loc_lda;
301  }
302 
304  sparseStorage (1);
305 
306  BCOMatType partMat;
307  convertCooToBco(rows, cols, 1, gmat_coo, partMat);
308 
309  sparseStorage[0] =
311  AllocateSharedPtr(rows, cols, 1, partMat, matStorage );
312 
313  // Create block diagonal matrix
315  AllocateSharedPtr(sparseStorage);
316  }
GlobalSysSolnType GetGlobalSysSolnType() const
Return the associated solution type.
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.
boost::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
std::map< CoordType, BCOEntryType > BCOMatType
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:127
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.
PreconditionerSharedPtr CreatePrecon(AssemblyMapSharedPtr asmMap)
Create a preconditioner object from the parameters defined in the supplied assembly map...
const boost::weak_ptr< ExpList > m_expList
Local Matrix System.
Definition: GlobalLinSys.h:129
void Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::v_BasisInvTransform ( Array< OneD, NekDouble > &  pInOut)
protectedvirtual

Reimplemented from Nektar::MultiRegions::GlobalLinSysStaticCond.

Definition at line 564 of file GlobalLinSysIterativeStaticCond.cpp.

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

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

Reimplemented from Nektar::MultiRegions::GlobalLinSysStaticCond.

Definition at line 557 of file GlobalLinSysIterativeStaticCond.cpp.

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

560  {
561  m_precon->DoTransformToLowEnergy(pInOut, offset);
562  }
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 479 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.

482  {
483  int nLocal = m_locToGloMap->GetNumLocalBndCoeffs();
484  int nDir = m_locToGloMap->GetNumGlobalDirBndCoeffs();
485  bool doGlobalOp = m_expList.lock()->GetGlobalOptParam()->
486  DoGlobalMatOp(m_linSysKey.GetMatrixType());
487 
488  if(doGlobalOp)
489  {
490  // Do matrix multiply globally
491  Array<OneD, NekDouble> in = pInput + nDir;
492  Array<OneD, NekDouble> out = pOutput + nDir;
493 
494  m_sparseSchurCompl->Multiply(in,out);
495  m_locToGloMap->UniversalAssembleBnd(pOutput, nDir);
496  }
497  else if (m_sparseSchurCompl)
498  {
499  // Do matrix multiply locally using block-diagonal sparse matrix
500  Array<OneD, NekDouble> tmp = m_wsp + nLocal;
501 
502  m_locToGloMap->GlobalToLocalBnd(pInput, m_wsp);
503  m_sparseSchurCompl->Multiply(m_wsp,tmp);
504  m_locToGloMap->AssembleBnd(tmp, pOutput);
505  }
506  else
507  {
508  // Do matrix multiply locally, using direct BLAS calls
509  m_locToGloMap->GlobalToLocalBnd(pInput, m_wsp);
510  int i, cnt;
511  Array<OneD, NekDouble> tmpout = m_wsp + nLocal;
512  for (i = cnt = 0; i < m_denseBlocks.size(); cnt += m_rows[i], ++i)
513  {
514  const int rows = m_rows[i];
515  Blas::Dgemv('N', rows, rows,
516  m_scale[i], m_denseBlocks[i], rows,
517  m_wsp.get()+cnt, 1,
518  0.0, tmpout.get()+cnt, 1);
519  }
520  m_locToGloMap->AssembleBnd(tmpout, pOutput);
521  }
522  }
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:127
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:129
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 208 of file GlobalLinSysIterativeStaticCond.cpp.

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

209  {
210  DNekScalBlkMatSharedPtr schurComplBlock;
211  int scLevel = m_locToGloMap->GetStaticCondLevel();
212  DNekScalBlkMatSharedPtr sc = scLevel == 0 ? m_S1Blk : m_schurCompl;
213  DNekScalMatSharedPtr localMat = sc->GetBlock(n,n);
214  unsigned int nbdry = localMat->GetRows();
215  unsigned int nblks = 1;
216  unsigned int esize[1] = {nbdry};
217 
218  schurComplBlock = MemoryManager<DNekScalBlkMat>
219  ::AllocateSharedPtr(nblks, nblks, esize, esize);
220  schurComplBlock->SetBlock(0, 0, localMat);
221 
222  return schurComplBlock;
223  }
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 155 of file GlobalLinSysIterativeStaticCond.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::MultiRegions::GlobalLinSys::CreatePrecon(), Nektar::eDIAGONAL, Nektar::StdRegions::eHybridDGHelmBndLam, Nektar::MultiRegions::GlobalMatrixKey::GetMatrixType(), 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, and Nektar::MultiRegions::GlobalLinSysStaticCond::SetupTopLevel().

156  {
158 
159  // Allocate memory for top-level structure
161 
162  // Setup Block Matrix systems
163  int n, n_exp = m_expList.lock()->GetNumElmts();
164 
165  MatrixStorage blkmatStorage = eDIAGONAL;
166  const Array<OneD,const unsigned int>& nbdry_size
167  = m_locToGloMap->GetNumLocalBndCoeffsPerPatch();
168 
170  ::AllocateSharedPtr(nbdry_size, nbdry_size, blkmatStorage);
171 
172  // Preserve original matrix in m_S1Blk
173  for (n = 0; n < n_exp; ++n)
174  {
175  DNekScalMatSharedPtr mat = m_schurCompl->GetBlock(n, n);
176  m_S1Blk->SetBlock(n, n, mat);
177  }
178 
179  // Build preconditioner
180  m_precon->BuildPreconditioner();
181 
182  // Do transform of Schur complement matrix
183  for (n = 0; n < n_exp; ++n)
184  {
185  if (m_linSysKey.GetMatrixType() !=
187  {
188  DNekScalMatSharedPtr mat = m_S1Blk->GetBlock(n, n);
189  DNekScalMatSharedPtr t = m_precon->TransformedSchurCompl(
190  m_expList.lock()->GetOffset_Elmt_Id(n), mat);
191  m_schurCompl->SetBlock(n, n, t);
192  }
193  }
194 
195  // Construct this level
197  }
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
DNekScalBlkMatSharedPtr m_schurCompl
Block Schur complement matrix.
void Initialise(const boost::shared_ptr< AssemblyMap > &pLocToGloMap)
Definition: GlobalLinSys.h:214
boost::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
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:127
boost::shared_ptr< AssemblyMap > m_locToGloMap
Local to global map.
PreconditionerSharedPtr CreatePrecon(AssemblyMapSharedPtr asmMap)
Create a preconditioner object from the parameters defined in the supplied assembly map...
const boost::weak_ptr< ExpList > m_expList
Local Matrix System.
Definition: GlobalLinSys.h:129
DNekScalBlkMatSharedPtr Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::v_PreSolve ( int  scLevel,
NekVector< NekDouble > &  F_GlobBnd 
)
protectedvirtual

Reimplemented from Nektar::MultiRegions::GlobalLinSysStaticCond.

Definition at line 529 of file GlobalLinSysIterativeStaticCond.cpp.

References Nektar::MultiRegions::GlobalLinSys::CreatePrecon(), Nektar::NekConstants::kNekUnsetDouble, Nektar::MultiRegions::GlobalLinSysStaticCond::m_locToGloMap, Nektar::MultiRegions::GlobalLinSysIterative::m_precon, Nektar::MultiRegions::GlobalLinSysIterative::m_rhs_magnitude, m_S1Blk, Nektar::MultiRegions::GlobalLinSysStaticCond::m_schurCompl, and Nektar::MultiRegions::GlobalLinSysIterative::Set_Rhs_Magnitude().

532  {
533  if (scLevel == 0)
534  {
535  // When matrices are supplied to the constructor at the top
536  // level, the preconditioner is never set up.
537  if (!m_precon)
538  {
540  m_precon->BuildPreconditioner();
541  }
542 
543  Set_Rhs_Magnitude(F_GlobBnd);
544 
545  return m_S1Blk;
546  }
547  else
548  {
549  // for multilevel iterative solver always use rhs
550  // vector value with no weighting
552 
553  return m_schurCompl;
554  }
555  }
DNekScalBlkMatSharedPtr m_schurCompl
Block Schur complement matrix.
void Set_Rhs_Magnitude(const NekVector< NekDouble > &pIn)
static const NekDouble kNekUnsetDouble
NekDouble m_rhs_magnitude
dot product of rhs to normalise stopping criterion
boost::shared_ptr< AssemblyMap > m_locToGloMap
Local to global map.
PreconditionerSharedPtr CreatePrecon(AssemblyMapSharedPtr asmMap)
Create a preconditioner object from the parameters defined in the supplied assembly map...
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 570 of file GlobalLinSysIterativeStaticCond.cpp.

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

578  {
579  GlobalLinSysIterativeStaticCondSharedPtr sys = MemoryManager<
580  GlobalLinSysIterativeStaticCond>::AllocateSharedPtr(
581  mkey, pExpList, pSchurCompl, pBinvD, pC, pInvD, l2gMap,
582  m_precon);
583  sys->Initialise(l2gMap);
584  return sys;
585  }
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 524 of file GlobalLinSysIterativeStaticCond.cpp.

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

Referenced by v_AssembleSchurComplement().

525  {
526  m_map = m_locToGloMap->GetGlobalToUniversalBndMapUnique();
527  }
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.