Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | Static Private Attributes | List of all members
Nektar::MultiRegions::GlobalLinSys Class Reference

A global linear system. More...

#include <GlobalLinSys.h>

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

Public Member Functions

 GlobalLinSys (const GlobalLinSysKey &pKey, const boost::weak_ptr< ExpList > &pExpList, const boost::shared_ptr< AssemblyMap > &pLocToGloMap)
 Constructor for full direct matrix solve.
virtual ~GlobalLinSys ()
const GlobalLinSysKeyGetKey (void) const
 Returns the key associated with the system.
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.
boost::shared_ptr< GlobalLinSysGetSharedThisPtr ()
 Returns a shared pointer to the current object.
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.

Protected Member Functions

virtual int v_GetNumBlocks ()
 Get the number of blocks in this system.
virtual DNekScalMatSharedPtr v_GetBlock (unsigned int n)
 Retrieves the block matrix from n-th expansion using the matrix key provided by the m_linSysKey.
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.
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.

Protected Attributes

const GlobalLinSysKey m_linSysKey
 Key associated with this linear system.
const boost::weak_ptr< ExpListm_expList
 Local Matrix System.
const map< int,
RobinBCInfoSharedPtr
m_robinBCInfo
 Robin boundary info.

Private 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.
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.
virtual void v_InitObject ()
virtual void v_Initialise (const boost::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 70 of file GlobalLinSys.h.

Constructor & Destructor Documentation

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.

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.

:
m_linSysKey(pKey),
m_expList(pExpList),
m_robinBCInfo(m_expList.lock()->GetRobinBCInfo())
{
}
virtual Nektar::MultiRegions::GlobalLinSys::~GlobalLinSys ( )
inlinevirtual

Definition at line 80 of file GlobalLinSys.h.

{}

Member Function Documentation

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

Definition at line 224 of file GlobalLinSys.h.

References v_DropStaticCondBlock().

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

Definition at line 214 of file GlobalLinSys.h.

References v_GetBlock().

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

{
return v_GetBlock(n);
}
const GlobalLinSysKey & Nektar::MultiRegions::GlobalLinSys::GetKey ( void  ) const
inline

Returns the key associated with the system.

Definition at line 163 of file GlobalLinSys.h.

References m_linSysKey.

{
return m_linSysKey;
}
const boost::weak_ptr< ExpList > & Nektar::MultiRegions::GlobalLinSys::GetLocMat ( void  ) const
inline

Definition at line 171 of file GlobalLinSys.h.

References m_expList.

{
return m_expList;
}
int Nektar::MultiRegions::GlobalLinSys::GetNumBlocks ( )
inline

Definition at line 229 of file GlobalLinSys.h.

References v_GetNumBlocks().

{
return v_GetNumBlocks();
}
boost::shared_ptr<GlobalLinSys> Nektar::MultiRegions::GlobalLinSys::GetSharedThisPtr ( )
inline

Returns a shared pointer to the current object.

Definition at line 103 of file GlobalLinSys.h.

Referenced by Nektar::MultiRegions::GlobalLinSysIterative::DoConjugateGradient(), Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::v_AssembleSchurComplement(), and Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::v_InitObject().

{
return shared_from_this();
}
DNekScalBlkMatSharedPtr Nektar::MultiRegions::GlobalLinSys::GetStaticCondBlock ( unsigned int  n)
inline

Definition at line 219 of file GlobalLinSys.h.

References v_GetStaticCondBlock().

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

Definition at line 208 of file GlobalLinSys.h.

References v_Initialise().

Referenced by Nektar::MultiRegions::GlobalLinSysStaticCond::v_InitObject(), and Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::v_InitObject().

{
v_Initialise(pLocToGloMap);
}
void Nektar::MultiRegions::GlobalLinSys::InitObject ( )
inline

Definition at line 203 of file GlobalLinSys.h.

References v_InitObject().

{
}
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 180 of file GlobalLinSys.h.

References v_Solve().

Referenced by Nektar::MultiRegions::GlobalLinSysXxt::v_SolveLinearSystem().

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

References v_SolveLinearSystem().

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

{
v_SolveLinearSystem(pNumRows, pInput, pOutput, locToGloMap, pNumDir);
}
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 367 of file GlobalLinSys.cpp.

References Nektar::MultiRegions::GlobalMatrixKey::GetConstFactors(), Nektar::MultiRegions::GlobalMatrixKey::GetMatrixType(), Nektar::MultiRegions::GlobalMatrixKey::GetNVarCoeffs(), Nektar::MultiRegions::GlobalMatrixKey::GetVarCoeffs(), m_expList, and m_linSysKey.

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

{
boost::shared_ptr<MultiRegions::ExpList> expList = m_expList.lock();
StdRegions::StdExpansionSharedPtr vExp = expList->GetExp(n);
// need to be initialised with zero size for non variable
// coefficient case
// retrieve variable coefficients
{
StdRegions::VarCoeffMap::const_iterator x;
int cnt = expList->GetPhys_Offset(n);
for (x = m_linSysKey.GetVarCoeffs().begin();
x != m_linSysKey.GetVarCoeffs().end (); ++x)
{
vVarCoeffMap[x->first] = x->second + cnt;
}
}
vExp->DetShapeType(),
*vExp,
vVarCoeffMap);
vExp->DropLocStaticCondMatrix(matkey);
}
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 227 of file GlobalLinSys.cpp.

References Nektar::MultiRegions::GlobalMatrixKey::GetConstFactors(), Nektar::MultiRegions::GlobalMatrixKey::GetMatrixType(), Nektar::MultiRegions::GlobalMatrixKey::GetNVarCoeffs(), Nektar::MultiRegions::GlobalMatrixKey::GetVarCoeffs(), m_expList, m_linSysKey, and m_robinBCInfo.

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

{
boost::shared_ptr<MultiRegions::ExpList> expList = m_expList.lock();
int cnt = 0;
boost::dynamic_pointer_cast<LocalRegions::Expansion>(
expList->GetExp(n));
// need to be initialised with zero size for non variable
// coefficient case
// retrieve variable coefficients
{
StdRegions::VarCoeffMap::const_iterator x;
cnt = expList->GetPhys_Offset(n);
for (x = m_linSysKey.GetVarCoeffs().begin();
x != m_linSysKey.GetVarCoeffs().end(); ++x)
{
vVarCoeffMap[x->first] = x->second + cnt;
}
}
vExp->DetShapeType(),
vVarCoeffMap);
loc_mat = vExp->GetLocMatrix(matkey);
// apply robin boundary conditions to the matrix.
if(m_robinBCInfo.count(n) != 0) // add robin mass matrix
{
// declare local matrix from scaled matrix.
int rows = loc_mat->GetRows();
int cols = loc_mat->GetColumns();
const NekDouble *dat = loc_mat->GetRawPtr();
AllocateSharedPtr(rows,cols,dat);
Blas::Dscal(rows*cols,loc_mat->Scale(),new_mat->GetRawPtr(),1);
// add local matrix contribution
for(rBC = m_robinBCInfo.find(n)->second;rBC; rBC = rBC->next)
{
vExp->AddRobinMassMatrix(
rBC->m_robinID, rBC->m_robinPrimitiveCoeffs, new_mat);
}
// redeclare loc_mat to point to new_mat plus the scalar.
1.0, new_mat);
}
// finally return the matrix.
return loc_mat;
}
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 215 of file GlobalLinSys.cpp.

References m_expList.

Referenced by GetNumBlocks().

{
return m_expList.lock()->GetExpSize();
}
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::GlobalLinSysIterativeStaticCond.

Definition at line 297 of file GlobalLinSys.cpp.

References Nektar::MultiRegions::GlobalMatrixKey::GetConstFactors(), Nektar::MultiRegions::GlobalMatrixKey::GetMatrixType(), Nektar::MultiRegions::GlobalMatrixKey::GetNVarCoeffs(), Nektar::MultiRegions::GlobalMatrixKey::GetVarCoeffs(), m_expList, m_linSysKey, and m_robinBCInfo.

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

{
boost::shared_ptr<MultiRegions::ExpList> expList = m_expList.lock();
int cnt = 0;
StdRegions::StdExpansionSharedPtr vExp = expList->GetExp(n);
// need to be initialised with zero size for non variable
// coefficient case
// retrieve variable coefficients
{
StdRegions::VarCoeffMap::const_iterator x;
cnt = expList->GetPhys_Offset(n);
for (x = m_linSysKey.GetVarCoeffs().begin();
x != m_linSysKey.GetVarCoeffs().end (); ++x)
{
vVarCoeffMap[x->first] = x->second + cnt;
}
}
vExp->DetShapeType(),
*vExp,
vVarCoeffMap);
loc_mat = vExp->GetLocStaticCondMatrix(matkey);
if(m_robinBCInfo.count(n) != 0) // add robin mass matrix
{
tmp_mat = loc_mat->GetBlock(0,0);
// declare local matrix from scaled matrix.
int rows = tmp_mat->GetRows();
int cols = tmp_mat->GetColumns();
const NekDouble *dat = tmp_mat->GetRawPtr();
AllocateSharedPtr(rows, cols, dat);
Blas::Dscal(rows*cols,tmp_mat->Scale(),new_mat->GetRawPtr(),1);
// add local matrix contribution
for(rBC = m_robinBCInfo.find(n)->second;rBC; rBC = rBC->next)
{
vExp->AddRobinMassMatrix(
rBC->m_robinID, rBC->m_robinPrimitiveCoeffs, new_mat);
}
// redeclare loc_mat to point to new_mat plus the scalar.
1.0, new_mat);
loc_mat->SetBlock(0,0,tmp_mat);
}
return loc_mat;
}
void Nektar::MultiRegions::GlobalLinSys::v_Initialise ( const boost::shared_ptr< AssemblyMap > &  pLocToGloMap)
privatevirtual

Reimplemented in Nektar::MultiRegions::GlobalLinSysStaticCond.

Definition at line 403 of file GlobalLinSys.cpp.

References ErrorUtil::efatal, and NEKERROR.

Referenced by Initialise().

{
NEKERROR(ErrorUtil::efatal, "Method does not exist" );
}
void Nektar::MultiRegions::GlobalLinSys::v_InitObject ( )
privatevirtual

Reimplemented in Nektar::MultiRegions::GlobalLinSysStaticCond.

Definition at line 398 of file GlobalLinSys.cpp.

References ErrorUtil::efatal, and NEKERROR.

Referenced by InitObject().

{
NEKERROR(ErrorUtil::efatal, "Method does not exist" );
}
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

Solve a linear system based on mapping.

Implemented in Nektar::MultiRegions::GlobalLinSysStaticCond.

Referenced by Solve().

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

Solve a basic matrix system.

Implemented in Nektar::MultiRegions::GlobalLinSysDirect, Nektar::MultiRegions::GlobalLinSysPETSc, and Nektar::MultiRegions::GlobalLinSysXxt.

Referenced by SolveLinearSystem().

Member Data Documentation

std::string Nektar::MultiRegions::GlobalLinSys::def
staticprivate
Initial value:
"DirectMultiLevelStaticCond")

Definition at line 156 of file GlobalLinSys.h.

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

Definition at line 155 of file GlobalLinSys.h.

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

Local Matrix System.

Definition at line 125 of file GlobalLinSys.h.

Referenced by Nektar::MultiRegions::GlobalLinSysDirectFull::AssembleFullMatrix(), Nektar::MultiRegions::GlobalLinSysXxtFull::AssembleMatrixArrays(), Nektar::MultiRegions::GlobalLinSysIterative::CalculateAnorm(), Nektar::MultiRegions::GlobalLinSysPETSc::CalculateReordering(), Nektar::MultiRegions::GlobalLinSysStaticCond::ConstructNextLevelCondensedSystem(), Nektar::MultiRegions::GlobalLinSysIterative::DoAconjugateProjection(), Nektar::MultiRegions::GlobalLinSysIterative::DoConjugateGradient(), GetLocMat(), Nektar::MultiRegions::GlobalLinSysIterative::GlobalLinSysIterative(), 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::GlobalLinSysIterativeStaticCond::v_AssembleSchurComplement(), Nektar::MultiRegions::GlobalLinSysIterativeFull::v_DoMatrixMultiply(), Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::v_DoMatrixMultiply(), v_DropStaticCondBlock(), v_GetBlock(), v_GetNumBlocks(), v_GetStaticCondBlock(), Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::v_InitObject(), Nektar::MultiRegions::GlobalLinSysDirectFull::v_Solve(), Nektar::MultiRegions::GlobalLinSysXxtFull::v_Solve(), Nektar::MultiRegions::GlobalLinSysIterativeFull::v_Solve(), and Nektar::MultiRegions::GlobalLinSysPETScFull::v_Solve().

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

Key associated with this linear system.

Definition at line 123 of file GlobalLinSys.h.

Referenced by Nektar::MultiRegions::GlobalLinSysDirectFull::AssembleFullMatrix(), Nektar::MultiRegions::GlobalLinSysXxtFull::AssembleMatrixArrays(), Nektar::MultiRegions::GlobalLinSysStaticCond::ConstructNextLevelCondensedSystem(), Nektar::MultiRegions::GlobalLinSysDirectStaticCond::DetermineMatrixStorage(), GetKey(), Nektar::MultiRegions::GlobalLinSysDirectFull::GlobalLinSysDirectFull(), Nektar::MultiRegions::GlobalLinSysIterativeFull::GlobalLinSysIterativeFull(), Nektar::MultiRegions::GlobalLinSysXxtFull::GlobalLinSysXxtFull(), Nektar::MultiRegions::GlobalLinSysStaticCond::SetupTopLevel(), Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::v_AssembleSchurComplement(), Nektar::MultiRegions::GlobalLinSysIterativeFull::v_DoMatrixMultiply(), Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::v_DoMatrixMultiply(), v_DropStaticCondBlock(), v_GetBlock(), v_GetStaticCondBlock(), Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::v_InitObject(), Nektar::MultiRegions::GlobalLinSysDirectFull::v_Solve(), Nektar::MultiRegions::GlobalLinSysXxtFull::v_Solve(), Nektar::MultiRegions::GlobalLinSysIterativeFull::v_Solve(), and Nektar::MultiRegions::GlobalLinSysPETScFull::v_Solve().

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

Robin boundary info.

Definition at line 127 of file GlobalLinSys.h.

Referenced by v_GetBlock(), and v_GetStaticCondBlock().