Nektar++
Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | Static Private Attributes | List of all members
Nektar::MultiRegions::GlobalLinSys Class Referenceabstract

A global linear system. More...

#include <GlobalLinSys.h>

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

Public Member Functions

 GlobalLinSys (const GlobalLinSysKey &pKey, const std::weak_ptr< ExpList > &pExpList, const std::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 std::weak_ptr< ExpList > & GetLocMat (void) const
 
void InitObject ()
 
void Initialise (const std::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...
 
std::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...
 

Protected Member Functions

virtual int v_GetNumBlocks ()
 Get the number of blocks in this system. More...
 
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 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 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 Attributes

const GlobalLinSysKey m_linSysKey
 Key associated with this linear system. More...
 
const std::weak_ptr< ExpListm_expList
 Local Matrix System. More...
 
const std::map< int, RobinBCInfoSharedPtrm_robinBCInfo
 Robin boundary info. More...
 
bool m_verbose
 

Private Member Functions

LocalRegions::MatrixKey GetBlockMatrixKey (unsigned int n)
 
virtual void v_Solve (const Array< OneD, const NekDouble > &in, Array< OneD, NekDouble > &out, const AssemblyMapSharedPtr &locToGloMap, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)=0
 Solve a linear system based on mapping. More...
 
virtual void v_SolveLinearSystem (const int pNumRows, const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const AssemblyMapSharedPtr &locToGloMap, const int pNumDir)=0
 Solve a basic matrix system. More...
 
virtual void v_InitObject ()
 
virtual void v_Initialise (const std::shared_ptr< AssemblyMap > &pLocToGloMap)
 

Static Private Attributes

static std::string lookupIds []
 
static std::string def
 

Detailed Description

A global linear system.

Consider the linear system \(\boldsymbol{M\hat{u}}_g=\boldsymbol{\hat{f}}\). Distinguishing between the boundary and interior components of \(\boldsymbol{\hat{u}}_g\) and \(\boldsymbol{\hat{f}}\) using \(\boldsymbol{\hat{u}}_b\), \(\boldsymbol{\hat{u}}_i\) and \(\boldsymbol{\hat{f}}_b\), \(\boldsymbol{\hat{f}}_i\) respectively, this system can be split into its constituent parts as

\[\left[\begin{array}{cc} \boldsymbol{M}_b&\boldsymbol{M}_{c1}\\ \boldsymbol{M}_{c2}&\boldsymbol{M}_i\\ \end{array}\right] \left[\begin{array}{c} \boldsymbol{\hat{u}_b}\\ \boldsymbol{\hat{u}_i}\\ \end{array}\right]= \left[\begin{array}{c} \boldsymbol{\hat{f}_b}\\ \boldsymbol{\hat{f}_i}\\ \end{array}\right]\]

where \(\boldsymbol{M}_b\) represents the components of \(\boldsymbol{M}\) resulting from boundary-boundary mode interactions, \(\boldsymbol{M}_{c1}\) and \(\boldsymbol{M}_{c2}\) represent the components resulting from coupling between the boundary-interior modes, and \(\boldsymbol{M}_i\) represents the components of \(\boldsymbol{M}\) resulting from interior-interior mode interactions.

The solution of the linear system can now be determined in two steps:

\begin{eqnarray*} \mathrm{step 1:}&\quad&(\boldsymbol{M}_b-\boldsymbol{M}_{c1} \boldsymbol{M}_i^{-1}\boldsymbol{M}_{c2}) \boldsymbol{\hat{u}_b} = \boldsymbol{\hat{f}}_b - \boldsymbol{M}_{c1}\boldsymbol{M}_i^{-1} \boldsymbol{\hat{f}}_i,\nonumber \\ \mathrm{step 2:}&\quad&\boldsymbol{\hat{u}_i}=\boldsymbol{M}_i^{-1} \left( \boldsymbol{\hat{f}}_i - \boldsymbol{M}_{c2}\boldsymbol{\hat{u}_b} \right). \nonumber \\ \end{eqnarray*}

As the inverse of \(\boldsymbol{M}_i^{-1}\) is

\[ \boldsymbol{M}_i^{-1} = \left [\underline{\boldsymbol{M}^e_i} \right ]^{-1} = \underline{[\boldsymbol{M}^e_i]}^{-1} \]

and the following operations can be evaluated as,

\begin{eqnarray*} \boldsymbol{M}_{c1}\boldsymbol{M}_i^{-1}\boldsymbol{\hat{f}}_i & =& \boldsymbol{\mathcal{A}}_b^T \underline{\boldsymbol{M}^e_{c1}} \underline{[\boldsymbol{M}^e_i]}^{-1} \boldsymbol{\hat{f}}_i \\ \boldsymbol{M}_{c2} \boldsymbol{\hat{u}_b} &=& \underline{\boldsymbol{M}^e_{c2}} \boldsymbol{\mathcal{A}}_b \boldsymbol{\hat{u}_b}.\end{eqnarray*}

where \(\boldsymbol{\mathcal{A}}_b \) is the permutation matrix which scatters from global to local degrees of freedom, only the following four matrices should be constructed:

The first three matrices are just a concatenation of the corresponding local matrices and they can be created as such. They also allow for an elemental evaluation of the operations concerned.

The global Schur complement however should be assembled from the concatenation of the local elemental Schur complements, that is,

\[ \boldsymbol{M}_{\mathrm{Schur}}=\boldsymbol{M}_b - \boldsymbol{M}_{c1} \boldsymbol{M}_i^{-1} \boldsymbol{M}_{c2} = \boldsymbol{\mathcal{A}}_b^T \left [\underline{\boldsymbol{M}^e_b - \boldsymbol{M}^e_{c1} [\boldsymbol{M}^e_i]^{-1} (\boldsymbol{M}^e_{c2})} \right ] \boldsymbol{\mathcal{A}}_b \]

and it is the only matrix operation that need to be evaluated on a global level when using static condensation. However, due to the size and sparsity of the matrix \(\boldsymbol{\mathcal{A}}_b\), it is more efficient to assemble the global Schur matrix using the mapping array bmap \([e][i]\) contained in the input argument locToGloMap. The global Schur complement is then constructed as:

\[\boldsymbol{M}_{\mathrm{Schur}}\left[\mathrm{\a bmap}[e][i]\right] \left[\mathrm{\a bmap}[e][j]\right]=\mathrm{\a bsign}[e][i]\cdot \mathrm{\a bsign}[e][j] \cdot\boldsymbol{M}^e_{\mathrm{Schur}}[i][j]\]

All four matrices are stored in the GlobalLinSys returned by this function.

Definition at line 72 of file GlobalLinSys.h.

Constructor & Destructor Documentation

◆ GlobalLinSys()

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

Constructor for full direct matrix solve.

Given a block matrix, construct a global matrix system according to a local to global mapping. #m_linSys is constructed by AssembleFullMatrix().

Parameters
pkeyAssociated linear system key.
locToGloMapLocal to global mapping.

Definition at line 196 of file GlobalLinSys.cpp.

199  :
200  m_linSysKey(pKey),
201  m_expList(pExpList),
202  m_robinBCInfo(m_expList.lock()->GetRobinBCInfo()),
203  m_verbose(m_expList.lock()->GetSession()->
204  DefinesCmdLineArgument("verbose"))
205  {
206  boost::ignore_unused(pLocToGloMap);
207  }
const std::weak_ptr< ExpList > m_expList
Local Matrix System.
Definition: GlobalLinSys.h:127
const std::map< int, RobinBCInfoSharedPtr > m_robinBCInfo
Robin boundary info.
Definition: GlobalLinSys.h:129
const GlobalLinSysKey m_linSysKey
Key associated with this linear system.
Definition: GlobalLinSys.h:125

◆ ~GlobalLinSys()

virtual Nektar::MultiRegions::GlobalLinSys::~GlobalLinSys ( )
inlinevirtual

Definition at line 82 of file GlobalLinSys.h.

82 {}

Member Function Documentation

◆ CreatePrecon()

PreconditionerSharedPtr Nektar::MultiRegions::GlobalLinSys::CreatePrecon ( AssemblyMapSharedPtr  asmMap)
protected

Create a preconditioner object from the parameters defined in the supplied assembly map.

Parameters
asmMapAssembly map used to construct the global system.

Definition at line 224 of file GlobalLinSys.cpp.

226  {
227  PreconditionerType pType = asmMap->GetPreconType();
228  std::string PreconType = MultiRegions::PreconditionerTypeMap[pType];
230  PreconType, GetSharedThisPtr(), asmMap);
231  }
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:145
std::shared_ptr< GlobalLinSys > GetSharedThisPtr()
Returns a shared pointer to the current object.
Definition: GlobalLinSys.h:105
const char *const PreconditionerTypeMap[]
PreconFactory & GetPreconFactory()

References Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), Nektar::MultiRegions::GetPreconFactory(), GetSharedThisPtr(), and Nektar::MultiRegions::PreconditionerTypeMap.

Referenced by Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::v_AssembleSchurComplement(), Nektar::MultiRegions::GlobalLinSysPETScStaticCond::v_AssembleSchurComplement(), Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::v_InitObject(), Nektar::MultiRegions::GlobalLinSysPETScStaticCond::v_InitObject(), Nektar::MultiRegions::GlobalLinSysPETScStaticCond::v_PreSolve(), Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::v_PreSolve(), Nektar::MultiRegions::GlobalLinSysIterative::v_SolveLinearSystem(), and Nektar::MultiRegions::GlobalLinSysPETSc::v_SolveLinearSystem().

◆ DropStaticCondBlock()

void Nektar::MultiRegions::GlobalLinSys::DropStaticCondBlock ( unsigned int  n)
inline

Definition at line 232 of file GlobalLinSys.h.

233  {
234  return v_DropStaticCondBlock(n);
235  }
virtual void v_DropStaticCondBlock(unsigned int n)
Releases the static condensation block matrices from NekManager of n-th expansion using the matrix ke...

References v_DropStaticCondBlock().

◆ GetBlock()

DNekScalMatSharedPtr Nektar::MultiRegions::GlobalLinSys::GetBlock ( unsigned int  n)
inline

Definition at line 222 of file GlobalLinSys.h.

223  {
224  return v_GetBlock(n);
225  }
virtual DNekScalMatSharedPtr v_GetBlock(unsigned int n)
Retrieves the block matrix from n-th expansion using the matrix key provided by the m_linSysKey.

References v_GetBlock().

Referenced by Nektar::MultiRegions::GlobalLinSysDirectFull::AssembleFullMatrix(), Nektar::MultiRegions::GlobalLinSysXxtFull::AssembleMatrixArrays(), and Nektar::MultiRegions::GlobalLinSysPETScFull::GlobalLinSysPETScFull().

◆ GetBlockMatrixKey()

LocalRegions::MatrixKey Nektar::MultiRegions::GlobalLinSys::GetBlockMatrixKey ( unsigned int  n)
private

Assemble the matrix key for each block n

Definition at line 250 of file GlobalLinSys.cpp.

251  {
252 
253  std::shared_ptr<MultiRegions::ExpList> expList = m_expList.lock();
254  int cnt = 0;
255 
256  LocalRegions::ExpansionSharedPtr vExp = expList->GetExp( n );
257 
258  // need to be initialised with zero size for non variable
259  // coefficient case
260  StdRegions::VarCoeffMap vVarCoeffMap;
261 
262  StdRegions::ConstFactorMap vConstFactorMap =
264 
265  // setup variable factors
266  if(m_linSysKey.GetNVarFactors() > 0)
267  {
270  {
271  vConstFactorMap[StdRegions::eFactorSVVDiffCoeff] =
274 
277  "VarCoeffSVVCuroffRatio is set but "
278  " not FactorSVVCutoffRatio");
279 
280  vConstFactorMap[StdRegions::eFactorSVVCutoffRatio] =
283 
284  }
285 
288  {
289  vConstFactorMap[StdRegions::eFactorSVVPowerKerDiffCoeff] =
292  }
293 
296  {
297  vConstFactorMap[StdRegions::eFactorSVVDGKerDiffCoeff] =
300  }
301  }
302 
303  // retrieve variable coefficients
304  if(m_linSysKey.GetNVarCoeffs() > 0)
305  {
306  cnt = expList->GetPhys_Offset(n);
307 
308  for (auto &x : m_linSysKey.GetVarCoeffs())
309  {
310  vVarCoeffMap[x.first] = x.second + cnt;
311  }
312  }
313 
314 
315  LocalRegions::MatrixKey matkey(m_linSysKey.GetMatrixType(),
316  vExp->DetShapeType(),
317  *vExp,
318  vConstFactorMap,
319  vVarCoeffMap);
320  return matkey;
321  }
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:250
const Array< OneD, const NekDouble > & GetVarFactors(const StdRegions::ConstFactorType &coeff) const
const StdRegions::ConstFactorMap & GetConstFactors() const
Returns all the constants.
const StdRegions::VarCoeffMap & GetVarCoeffs() const
StdRegions::MatrixType GetMatrixType() const
Return the matrix type.
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:68
std::map< StdRegions::VarCoeffType, Array< OneD, NekDouble > > VarCoeffMap
Definition: StdRegions.hpp:272
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:314

References ASSERTL1, Nektar::StdRegions::eFactorSVVCutoffRatio, Nektar::StdRegions::eFactorSVVDGKerDiffCoeff, Nektar::StdRegions::eFactorSVVDiffCoeff, Nektar::StdRegions::eFactorSVVPowerKerDiffCoeff, Nektar::MultiRegions::GlobalMatrixKey::GetConstFactors(), Nektar::MultiRegions::GlobalMatrixKey::GetMatrixType(), Nektar::MultiRegions::GlobalMatrixKey::GetNVarCoeffs(), Nektar::MultiRegions::GlobalLinSysKey::GetNVarFactors(), Nektar::MultiRegions::GlobalMatrixKey::GetVarCoeffs(), Nektar::MultiRegions::GlobalLinSysKey::GetVarFactors(), m_expList, and m_linSysKey.

Referenced by v_DropStaticCondBlock(), v_GetBlock(), and v_GetStaticCondBlock().

◆ GetKey()

const GlobalLinSysKey & Nektar::MultiRegions::GlobalLinSys::GetKey ( void  ) const
inline

Returns the key associated with the system.

Definition at line 171 of file GlobalLinSys.h.

172  {
173  return m_linSysKey;
174  }

References m_linSysKey.

◆ GetLocMat()

const std::weak_ptr< ExpList > & Nektar::MultiRegions::GlobalLinSys::GetLocMat ( void  ) const
inline

Definition at line 179 of file GlobalLinSys.h.

180  {
181  return m_expList;
182  }

References m_expList.

◆ GetNumBlocks()

int Nektar::MultiRegions::GlobalLinSys::GetNumBlocks ( )
inline

Definition at line 237 of file GlobalLinSys.h.

238  {
239  return v_GetNumBlocks();
240  }
virtual int v_GetNumBlocks()
Get the number of blocks in this system.

References v_GetNumBlocks().

◆ GetSharedThisPtr()

std::shared_ptr<GlobalLinSys> Nektar::MultiRegions::GlobalLinSys::GetSharedThisPtr ( )
inline

Returns a shared pointer to the current object.

Definition at line 105 of file GlobalLinSys.h.

106  {
107  return shared_from_this();
108  }

Referenced by CreatePrecon().

◆ GetStaticCondBlock()

DNekScalBlkMatSharedPtr Nektar::MultiRegions::GlobalLinSys::GetStaticCondBlock ( unsigned int  n)
inline

Definition at line 227 of file GlobalLinSys.h.

228  {
229  return v_GetStaticCondBlock(n);
230  }
virtual DNekScalBlkMatSharedPtr v_GetStaticCondBlock(unsigned int n)
Retrieves a the static condensation block matrices from n-th expansion using the matrix key provided ...

References v_GetStaticCondBlock().

◆ Initialise()

void Nektar::MultiRegions::GlobalLinSys::Initialise ( const std::shared_ptr< AssemblyMap > &  pLocToGloMap)
inline

◆ InitObject()

void Nektar::MultiRegions::GlobalLinSys::InitObject ( )
inline

Definition at line 211 of file GlobalLinSys.h.

212  {
213  v_InitObject();
214  }

References v_InitObject().

◆ Solve()

void Nektar::MultiRegions::GlobalLinSys::Solve ( const Array< OneD, const NekDouble > &  in,
Array< OneD, NekDouble > &  out,
const AssemblyMapSharedPtr locToGloMap,
const Array< OneD, const NekDouble > &  dirForcing = NullNekDouble1DArray 
)
inline

Solve the linear system for given input and output vectors using a specified local to global map.

Definition at line 188 of file GlobalLinSys.h.

193  {
194  v_Solve(in,out,locToGloMap,dirForcing);
195  }
virtual void v_Solve(const Array< OneD, const NekDouble > &in, Array< OneD, NekDouble > &out, const AssemblyMapSharedPtr &locToGloMap, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)=0
Solve a linear system based on mapping.

References v_Solve().

◆ SolveLinearSystem()

void Nektar::MultiRegions::GlobalLinSys::SolveLinearSystem ( const int  pNumRows,
const Array< OneD, const NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput,
const AssemblyMapSharedPtr locToGloMap,
const int  pNumDir = 0 
)
inline

Solve the linear system for given input and output vectors.

Definition at line 201 of file GlobalLinSys.h.

207  {
208  v_SolveLinearSystem(pNumRows, pInput, pOutput, locToGloMap, pNumDir);
209  }
virtual void v_SolveLinearSystem(const int pNumRows, const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const AssemblyMapSharedPtr &locToGloMap, const int pNumDir)=0
Solve a basic matrix system.

References v_SolveLinearSystem().

Referenced by Nektar::MultiRegions::GlobalLinSysIterativeFull::v_Solve(), Nektar::MultiRegions::GlobalLinSysPETScFull::v_Solve(), Nektar::MultiRegions::GlobalLinSysStaticCond::v_Solve(), Nektar::MultiRegions::GlobalLinSysXxtFull::v_Solve(), and Nektar::MultiRegions::GlobalLinSysDirectFull::v_Solve().

◆ v_DropStaticCondBlock()

void Nektar::MultiRegions::GlobalLinSys::v_DropStaticCondBlock ( unsigned int  n)
protectedvirtual

Releases the static condensation block matrices from NekManager of n-th expansion using the matrix key provided by the m_linSysKey.

Parameters
nNumber of the expansion

Definition at line 428 of file GlobalLinSys.cpp.

429  {
430  LocalRegions::ExpansionSharedPtr vExp = m_expList.lock()->GetExp( n );
431  vExp->DropLocStaticCondMatrix(GetBlockMatrixKey(n));
432  }
LocalRegions::MatrixKey GetBlockMatrixKey(unsigned int n)

References GetBlockMatrixKey(), and m_expList.

Referenced by DropStaticCondBlock(), and Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::PrepareLocalSchurComplement().

◆ v_GetBlock()

DNekScalMatSharedPtr Nektar::MultiRegions::GlobalLinSys::v_GetBlock ( unsigned int  n)
protectedvirtual

Retrieves the block matrix from n-th expansion using the matrix key provided by the m_linSysKey.

Parameters
nNumber of the expansion.
Returns
Block matrix for the specified expansion.

Definition at line 330 of file GlobalLinSys.cpp.

331  {
332  LocalRegions::ExpansionSharedPtr vExp = m_expList.lock()->GetExp( n );
333  DNekScalMatSharedPtr loc_mat;
334  loc_mat = vExp->GetLocMatrix(GetBlockMatrixKey(n));
335 
336  // apply robin boundary conditions to the matrix.
337  if(m_robinBCInfo.count(n) != 0) // add robin mass matrix
338  {
340 
341  // declare local matrix from scaled matrix.
342  int rows = loc_mat->GetRows();
343  int cols = loc_mat->GetColumns();
344  const NekDouble *dat = loc_mat->GetRawPtr();
346  AllocateSharedPtr(rows,cols,dat);
347  Blas::Dscal(rows*cols,loc_mat->Scale(),new_mat->GetRawPtr(),1);
348 
349  // add local matrix contribution
350  for(rBC = m_robinBCInfo.find(n)->second;rBC; rBC = rBC->next)
351  {
352  vExp->AddRobinMassMatrix(
353  rBC->m_robinID, rBC->m_robinPrimitiveCoeffs, new_mat);
354  }
355 
356  // redeclare loc_mat to point to new_mat plus the scalar.
358  1.0, new_mat);
359  }
360 
361  // finally return the matrix.
362  return loc_mat;
363  }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
static void Dscal(const int &n, const double &alpha, double *x, const int &incx)
BLAS level 1: x = alpha x.
Definition: Blas.hpp:182
std::shared_ptr< RobinBCInfo > RobinBCInfoSharedPtr
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:69
double NekDouble

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Blas::Dscal(), GetBlockMatrixKey(), m_expList, and m_robinBCInfo.

Referenced by GetBlock(), and Nektar::MultiRegions::GlobalLinSysStaticCond::SetupTopLevel().

◆ v_GetNumBlocks()

int Nektar::MultiRegions::GlobalLinSys::v_GetNumBlocks ( )
protectedvirtual

Get the number of blocks in this system.

At the top level this corresponds to the number of elements in the expansion list.

Reimplemented in Nektar::MultiRegions::GlobalLinSysStaticCond.

Definition at line 240 of file GlobalLinSys.cpp.

241  {
242  return m_expList.lock()->GetExpSize();
243  }

References m_expList.

Referenced by GetNumBlocks().

◆ v_GetStaticCondBlock()

DNekScalBlkMatSharedPtr Nektar::MultiRegions::GlobalLinSys::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 in Nektar::MultiRegions::GlobalLinSysPETScStaticCond, and Nektar::MultiRegions::GlobalLinSysIterativeStaticCond.

Definition at line 373 of file GlobalLinSys.cpp.

375  {
376 
377  LocalRegions::ExpansionSharedPtr vExp = m_expList.lock()->GetExp( n );
378  DNekScalBlkMatSharedPtr loc_mat;
379  loc_mat = vExp->GetLocStaticCondMatrix(GetBlockMatrixKey(n));
380 
381  if(m_robinBCInfo.count(n) != 0) // add robin mass matrix
382  {
383  DNekScalMatSharedPtr tmp_mat;
385 
386  tmp_mat = loc_mat->GetBlock(0,0);
387 
388  // declare local matrix from scaled matrix.
389  int rows = tmp_mat->GetRows();
390  int cols = tmp_mat->GetColumns();
391  const NekDouble *dat = tmp_mat->GetRawPtr();
393  AllocateSharedPtr(rows, cols, dat);
394  Blas::Dscal(rows*cols,tmp_mat->Scale(),new_mat->GetRawPtr(),1);
395 
396  // add local matrix contribution
397  for(rBC = m_robinBCInfo.find(n)->second;rBC; rBC = rBC->next)
398  {
399  vExp->AddRobinMassMatrix(
400  rBC->m_robinID, rBC->m_robinPrimitiveCoeffs, new_mat);
401  }
402 
403  // redeclare loc_mat to point to new_mat plus the scalar.
405  1.0, new_mat);
406  DNekScalBlkMatSharedPtr new_loc_mat;
407  unsigned int exp_size[] = {tmp_mat->GetRows(), loc_mat->GetBlock(1,1)->GetRows()};
408  unsigned int nblks = 2;
409  new_loc_mat = MemoryManager<DNekScalBlkMat>::AllocateSharedPtr(nblks, nblks, exp_size, exp_size);
410 
411 
412  new_loc_mat->SetBlock(0,0,tmp_mat);
413  new_loc_mat->SetBlock(0,1,loc_mat->GetBlock(0,1));
414  new_loc_mat->SetBlock(1,0,loc_mat->GetBlock(1,0));
415  new_loc_mat->SetBlock(1,1,loc_mat->GetBlock(1,1));
416  loc_mat = new_loc_mat;
417  }
418 
419  return loc_mat;
420  }
std::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
Definition: NekTypeDefs.hpp:73

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Blas::Dscal(), GetBlockMatrixKey(), m_expList, and m_robinBCInfo.

Referenced by GetStaticCondBlock(), and Nektar::MultiRegions::GlobalLinSysStaticCond::SetupTopLevel().

◆ v_Initialise()

void Nektar::MultiRegions::GlobalLinSys::v_Initialise ( const std::shared_ptr< AssemblyMap > &  pLocToGloMap)
privatevirtual

Reimplemented in Nektar::MultiRegions::GlobalLinSysStaticCond.

Definition at line 439 of file GlobalLinSys.cpp.

441  {
442  boost::ignore_unused(pLocToGloMap);
443  NEKERROR(ErrorUtil::efatal, "Method does not exist" );
444  }
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:209

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by Initialise().

◆ v_InitObject()

void Nektar::MultiRegions::GlobalLinSys::v_InitObject ( )
privatevirtual

◆ v_Solve()

virtual void Nektar::MultiRegions::GlobalLinSys::v_Solve ( const Array< OneD, const NekDouble > &  in,
Array< OneD, NekDouble > &  out,
const AssemblyMapSharedPtr locToGloMap,
const Array< OneD, const NekDouble > &  dirForcing = NullNekDouble1DArray 
)
privatepure virtual

◆ v_SolveLinearSystem()

virtual void Nektar::MultiRegions::GlobalLinSys::v_SolveLinearSystem ( const int  pNumRows,
const Array< OneD, const NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput,
const AssemblyMapSharedPtr locToGloMap,
const int  pNumDir 
)
privatepure virtual

Member Data Documentation

◆ def

std::string Nektar::MultiRegions::GlobalLinSys::def
staticprivate
Initial value:
"DirectStaticCond")
static std::string RegisterDefaultSolverInfo(const std::string &pName, const std::string &pValue)
Registers the default string value of a solver info property.

Definition at line 164 of file GlobalLinSys.h.

◆ lookupIds

std::string Nektar::MultiRegions::GlobalLinSys::lookupIds
staticprivate

Definition at line 163 of file GlobalLinSys.h.

◆ m_expList

const std::weak_ptr<ExpList> Nektar::MultiRegions::GlobalLinSys::m_expList
protected

Local Matrix System.

Definition at line 127 of file GlobalLinSys.h.

Referenced by Nektar::MultiRegions::GlobalLinSysDirectFull::AssembleFullMatrix(), Nektar::MultiRegions::GlobalLinSysXxtFull::AssembleMatrixArrays(), Nektar::MultiRegions::GlobalLinSysPETSc::CalculateReordering(), Nektar::MultiRegions::GlobalLinSysStaticCond::ConstructNextLevelCondensedSystem(), Nektar::MultiRegions::GlobalLinSysIterative::DoProjection(), GetBlockMatrixKey(), GetLocMat(), Nektar::MultiRegions::GlobalLinSysIterative::GlobalLinSysIterative(), Nektar::MultiRegions::GlobalLinSysPETSc::GlobalLinSysPETSc(), Nektar::MultiRegions::GlobalLinSysPETScFull::GlobalLinSysPETScFull(), Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::PrepareLocalSchurComplement(), Nektar::MultiRegions::GlobalLinSysIterative::Set_Rhs_Magnitude(), Nektar::MultiRegions::GlobalLinSysStaticCond::SetupTopLevel(), Nektar::MultiRegions::GlobalLinSysIterative::UpdateKnownSolutions(), Nektar::MultiRegions::GlobalLinSysXxtStaticCond::v_AssembleSchurComplement(), Nektar::MultiRegions::GlobalLinSysPETScFull::v_DoMatrixMultiply(), Nektar::MultiRegions::GlobalLinSysIterativeFull::v_DoMatrixMultiply(), v_DropStaticCondBlock(), v_GetBlock(), v_GetNumBlocks(), v_GetStaticCondBlock(), Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::v_InitObject(), Nektar::MultiRegions::GlobalLinSysPETScStaticCond::v_InitObject(), Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::v_PreSolve(), Nektar::MultiRegions::GlobalLinSysIterativeFull::v_Solve(), Nektar::MultiRegions::GlobalLinSysPETScFull::v_Solve(), Nektar::MultiRegions::GlobalLinSysXxtFull::v_Solve(), Nektar::MultiRegions::GlobalLinSysDirectFull::v_Solve(), and Nektar::MultiRegions::GlobalLinSysIterative::v_SolveLinearSystem().

◆ m_linSysKey

const GlobalLinSysKey Nektar::MultiRegions::GlobalLinSys::m_linSysKey
protected

◆ m_robinBCInfo

const std::map<int, RobinBCInfoSharedPtr> Nektar::MultiRegions::GlobalLinSys::m_robinBCInfo
protected

◆ m_verbose

bool Nektar::MultiRegions::GlobalLinSys::m_verbose
protected