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)
 
void DropBlock (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 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)
 
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 void v_DropBlock (unsigned int n)
 Releases the local block matrix from NekManager of 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)
 

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 71 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 188 of file GlobalLinSys.cpp.

191  : m_linSysKey(pKey), m_expList(pExpList),
192  m_robinBCInfo(m_expList.lock()->GetRobinBCInfo()),
193  m_verbose(
194  m_expList.lock()->GetSession()->DefinesCmdLineArgument("verbose"))
195 {
196  boost::ignore_unused(pLocToGloMap);
197 }
const std::weak_ptr< ExpList > m_expList
Local Matrix System.
Definition: GlobalLinSys.h:124
const std::map< int, RobinBCInfoSharedPtr > m_robinBCInfo
Robin boundary info.
Definition: GlobalLinSys.h:126
const GlobalLinSysKey m_linSysKey
Key associated with this linear system.
Definition: GlobalLinSys.h:122

◆ ~GlobalLinSys()

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

Definition at line 80 of file GlobalLinSys.h.

81  {
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 214 of file GlobalLinSys.cpp.

215 {
216  PreconditionerType pType = asmMap->GetPreconType();
217  std::string PreconType = MultiRegions::PreconditionerTypeMap[pType];
218  return GetPreconFactory().CreateInstance(PreconType, GetSharedThisPtr(),
219  asmMap);
220 }
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:144
std::shared_ptr< GlobalLinSys > GetSharedThisPtr()
Returns a shared pointer to the current object.
Definition: GlobalLinSys.h:102
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().

◆ DropBlock()

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

Definition at line 216 of file GlobalLinSys.h.

217 {
218  return v_DropBlock(n);
219 }
virtual void v_DropBlock(unsigned int n)
Releases the local block matrix from NekManager of n-th expansion using the matrix key provided by th...

References v_DropBlock().

◆ DropStaticCondBlock()

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

Definition at line 226 of file GlobalLinSys.h.

227 {
228  return v_DropStaticCondBlock(n);
229 }
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 211 of file GlobalLinSys.h.

212 {
213  return v_GetBlock(n);
214 }
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 237 of file GlobalLinSys.cpp.

238 {
239 
240  std::shared_ptr<MultiRegions::ExpList> expList = m_expList.lock();
241 
242  LocalRegions::ExpansionSharedPtr vExp = expList->GetExp(n);
243 
244  // need to be initialised with zero size for non variable
245  // coefficient case
246  StdRegions::VarCoeffMap vVarCoeffMap;
247 
249 
250  // setup variable factors
251  if (m_linSysKey.GetNVarFactors() > 0)
252  {
253  if (m_linSysKey.GetVarFactors().count(
255  {
256  vConstFactorMap[StdRegions::eFactorSVVDiffCoeff] =
258 
261  "VarCoeffSVVCuroffRatio is set but "
262  " not FactorSVVCutoffRatio");
263 
264  vConstFactorMap[StdRegions::eFactorSVVCutoffRatio] =
266  }
267 
268  if (m_linSysKey.GetVarFactors().count(
270  {
271  vConstFactorMap[StdRegions::eFactorSVVPowerKerDiffCoeff] =
274  }
275 
276  if (m_linSysKey.GetVarFactors().count(
278  {
279  vConstFactorMap[StdRegions::eFactorSVVDGKerDiffCoeff] =
282  }
283  }
284 
285  // retrieve variable coefficients
286  if (m_linSysKey.GetNVarCoeffs() > 0)
287  {
289  expList->GetPhys_Offset(n),
290  vExp->GetTotPoints());
291  }
292 
293  LocalRegions::MatrixKey matkey(m_linSysKey.GetMatrixType(),
294  vExp->DetShapeType(), *vExp, vConstFactorMap,
295  vVarCoeffMap);
296  return matkey;
297 }
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
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
VarCoeffMap RestrictCoeffMap(const VarCoeffMap &m, size_t offset, size_t cnt)
Definition: StdRegions.hpp:346
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:399
std::map< StdRegions::VarCoeffType, VarCoeffEntry > VarCoeffMap
Definition: StdRegions.hpp:343

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, m_linSysKey, and Nektar::StdRegions::RestrictCoeffMap().

Referenced by v_DropBlock(), 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 165 of file GlobalLinSys.h.

166 {
167  return m_linSysKey;
168 }

References m_linSysKey.

◆ GetLocMat()

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

Definition at line 173 of file GlobalLinSys.h.

174 {
175  return m_expList;
176 }

References m_expList.

◆ GetNumBlocks()

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

Definition at line 231 of file GlobalLinSys.h.

232 {
233  return v_GetNumBlocks();
234 }
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 102 of file GlobalLinSys.h.

103  {
104  return shared_from_this();
105  }

Referenced by CreatePrecon().

◆ GetStaticCondBlock()

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

Definition at line 221 of file GlobalLinSys.h.

222 {
223  return v_GetStaticCondBlock(n);
224 }
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 200 of file GlobalLinSys.h.

201 {
202  v_InitObject();
203 }

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 181 of file GlobalLinSys.h.

185 {
186  v_Solve(in, out, locToGloMap, dirForcing);
187 }
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 192 of file GlobalLinSys.h.

196 {
197  v_SolveLinearSystem(pNumRows, pInput, pOutput, locToGloMap, pNumDir);
198 }
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_DropBlock()

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

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

Parameters
nNumber of the expansion

Definition at line 414 of file GlobalLinSys.cpp.

415 {
416  LocalRegions::ExpansionSharedPtr vExp = m_expList.lock()->GetExp(n);
417  vExp->DropLocMatrix(GetBlockMatrixKey(n));
418 }
LocalRegions::MatrixKey GetBlockMatrixKey(unsigned int n)

References GetBlockMatrixKey(), and m_expList.

Referenced by DropBlock().

◆ 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 402 of file GlobalLinSys.cpp.

403 {
404  LocalRegions::ExpansionSharedPtr vExp = m_expList.lock()->GetExp(n);
405  vExp->DropLocStaticCondMatrix(GetBlockMatrixKey(n));
406 }

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 306 of file GlobalLinSys.cpp.

307 {
308  LocalRegions::ExpansionSharedPtr vExp = m_expList.lock()->GetExp(n);
309  DNekScalMatSharedPtr loc_mat;
310  loc_mat = vExp->GetLocMatrix(GetBlockMatrixKey(n));
311 
312  // apply robin boundary conditions to the matrix.
313  if (m_robinBCInfo.count(n) != 0) // add robin mass matrix
314  {
316 
317  // declare local matrix from scaled matrix.
318  int rows = loc_mat->GetRows();
319  int cols = loc_mat->GetColumns();
320  const NekDouble *dat = loc_mat->GetRawPtr();
321  DNekMatSharedPtr new_mat =
323  Blas::Dscal(rows * cols, loc_mat->Scale(), new_mat->GetRawPtr(), 1);
324 
325  // add local matrix contribution
326  for (rBC = m_robinBCInfo.find(n)->second; rBC; rBC = rBC->next)
327  {
328  vExp->AddRobinMassMatrix(rBC->m_robinID,
329  rBC->m_robinPrimitiveCoeffs, new_mat);
330  }
331 
332  // redeclare loc_mat to point to new_mat plus the scalar.
333  loc_mat = MemoryManager<DNekScalMat>::AllocateSharedPtr(1.0, new_mat);
334  }
335 
336  // finally return the matrix.
337  return loc_mat;
338 }
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:168
std::shared_ptr< RobinBCInfo > RobinBCInfoSharedPtr
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75
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 228 of file GlobalLinSys.cpp.

229 {
230  return m_expList.lock()->GetExpSize();
231 }

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 348 of file GlobalLinSys.cpp.

349 {
350 
351  LocalRegions::ExpansionSharedPtr vExp = m_expList.lock()->GetExp(n);
352  DNekScalBlkMatSharedPtr loc_mat;
353  loc_mat = vExp->GetLocStaticCondMatrix(GetBlockMatrixKey(n));
354 
355  if (m_robinBCInfo.count(n) != 0) // add robin mass matrix
356  {
357  DNekScalMatSharedPtr tmp_mat;
359 
360  tmp_mat = loc_mat->GetBlock(0, 0);
361 
362  // declare local matrix from scaled matrix.
363  int rows = tmp_mat->GetRows();
364  int cols = tmp_mat->GetColumns();
365  const NekDouble *dat = tmp_mat->GetRawPtr();
366  DNekMatSharedPtr new_mat =
368  Blas::Dscal(rows * cols, tmp_mat->Scale(), new_mat->GetRawPtr(), 1);
369 
370  // add local matrix contribution
371  for (rBC = m_robinBCInfo.find(n)->second; rBC; rBC = rBC->next)
372  {
373  vExp->AddRobinMassMatrix(rBC->m_robinID,
374  rBC->m_robinPrimitiveCoeffs, new_mat);
375  }
376 
377  // redeclare loc_mat to point to new_mat plus the scalar.
378  tmp_mat = MemoryManager<DNekScalMat>::AllocateSharedPtr(1.0, new_mat);
379  DNekScalBlkMatSharedPtr new_loc_mat;
380  unsigned int exp_size[] = {tmp_mat->GetRows(),
381  loc_mat->GetBlock(1, 1)->GetRows()};
382  unsigned int nblks = 2;
384  nblks, nblks, exp_size, exp_size);
385 
386  new_loc_mat->SetBlock(0, 0, tmp_mat);
387  new_loc_mat->SetBlock(0, 1, loc_mat->GetBlock(0, 1));
388  new_loc_mat->SetBlock(1, 0, loc_mat->GetBlock(1, 0));
389  new_loc_mat->SetBlock(1, 1, loc_mat->GetBlock(1, 1));
390  loc_mat = new_loc_mat;
391  }
392 
393  return loc_mat;
394 }
std::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
Definition: NekTypeDefs.hpp:79

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)
protectedvirtual

Reimplemented in Nektar::MultiRegions::GlobalLinSysStaticCond.

Definition at line 425 of file GlobalLinSys.cpp.

427 {
428  boost::ignore_unused(pLocToGloMap);
429  NEKERROR(ErrorUtil::efatal, "Method does not exist");
430 }
#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 ( )
protectedvirtual

◆ 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 
)
protectedpure 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 
)
protectedpure 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 159 of file GlobalLinSys.h.

◆ lookupIds

std::string Nektar::MultiRegions::GlobalLinSys::lookupIds
staticprivate
Initial value:
= {
"GlobalSysSoln", "DirectFull", MultiRegions::eDirectFullMatrix),
"GlobalSysSoln", "DirectStaticCond", MultiRegions::eDirectStaticCond),
"GlobalSysSoln", "DirectMultiLevelStaticCond",
"GlobalSysSoln", "IterativeFull", MultiRegions::eIterativeFull),
"GlobalSysSoln", "IterativeStaticCond",
"GlobalSysSoln", "IterativeMultiLevelStaticCond",
"GlobalSysSoln", "XxtFull", MultiRegions::eXxtFullMatrix),
"GlobalSysSoln", "XxtStaticCond", MultiRegions::eXxtStaticCond),
"GlobalSysSoln", "XxtMultiLevelStaticCond",
"GlobalSysSoln", "PETScFull", MultiRegions::ePETScFullMatrix),
"GlobalSysSoln", "PETScStaticCond", MultiRegions::ePETScStaticCond),
"GlobalSysSoln", "PETScMultiLevelStaticCond",
static std::string RegisterEnumValue(std::string pEnum, std::string pString, int pEnumValue)
Registers an enumeration value.

Definition at line 158 of file GlobalLinSys.h.

◆ m_expList

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

Local Matrix System.

Definition at line 124 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_DropBlock(), 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