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

#include <Preconditioner.h>

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

Public Member Functions

 Preconditioner (const std::shared_ptr< GlobalLinSys > &plinsys, const AssemblyMapSharedPtr &pLocToGloMap)
 
virtual ~Preconditioner ()
 
void DoPreconditioner (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const bool &IsLocal=false)
 
void DoAssembleLoc (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const bool &ZeroDir)
 Apply an assembly and scatter back to lcoal array. More...
 
void DoPreconditionerWithNonVertOutput (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const Array< OneD, NekDouble > &pNonVertOutput, Array< OneD, NekDouble > &pVertForce=NullNekDouble1DArray)
 
void DoTransformBasisToLowEnergy (Array< OneD, NekDouble > &pInOut)
 
void DoTransformCoeffsFromLowEnergy (Array< OneD, NekDouble > &pInOut)
 
void DoTransformCoeffsToLowEnergy (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
 
void DoTransformBasisFromLowEnergy (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
 
void BuildPreconditioner ()
 
void InitObject ()
 
Array< OneD, NekDoubleAssembleStaticCondGlobalDiagonals ()
 Performs global assembly of diagonal entries to global Schur complement matrix. More...
 
const DNekScalBlkMatSharedPtrGetBlockTransformedSchurCompl () const
 
const DNekScalBlkMatSharedPtrGetBlockCMatrix () const
 
const DNekScalBlkMatSharedPtrGetBlockInvDMatrix () const
 
const DNekScalBlkMatSharedPtrGetBlockSchurCompl () const
 
const DNekScalBlkMatSharedPtrGetBlockTransformationMatrix () const
 
const DNekScalBlkMatSharedPtrGetBlockTransposedTransformationMatrix () const
 
DNekScalMatSharedPtr TransformedSchurCompl (int offset, int bndoffset, const std::shared_ptr< DNekScalMat > &loc_mat)
 

Protected Member Functions

virtual DNekScalMatSharedPtr v_TransformedSchurCompl (int offset, int bndoffset, const std::shared_ptr< DNekScalMat > &loc_mat)
 Get block elemental transposed transformation matrix \(\mathbf{R}^{T}\). More...
 
virtual void v_InitObject ()
 
virtual void v_DoPreconditioner (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const bool &isLocal=false)
 Apply a preconditioner to the conjugate gradient method. More...
 
virtual void v_DoPreconditionerWithNonVertOutput (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const Array< OneD, NekDouble > &pNonVertOutput, Array< OneD, NekDouble > &pVertForce)
 Apply a preconditioner to the conjugate gradient method with an output for non-vertex degrees of freedom. More...
 
virtual void v_DoTransformBasisToLowEnergy (Array< OneD, NekDouble > &pInOut)
 Transform from original basis to low energy basis. More...
 
virtual void v_DoTransformCoeffsFromLowEnergy (Array< OneD, NekDouble > &pInOut)
 Transform from low energy coeffs to orignal basis. More...
 
virtual void v_DoTransformCoeffsToLowEnergy (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
 Multiply by the block inverse transformation matrix. More...
 
virtual void v_DoTransformBasisFromLowEnergy (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
 Multiply by the block transposed inverse transformation matrix. More...
 
virtual void v_BuildPreconditioner ()
 

Protected Attributes

const std::weak_ptr< GlobalLinSysm_linsys
 
std::string m_preconType
 
DNekMatSharedPtr m_preconditioner
 
std::weak_ptr< AssemblyMapm_locToGloMap
 
LibUtilities::CommSharedPtr m_comm
 

Private Member Functions

void NullPreconditioner (void)
 

Static Private Attributes

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

Detailed Description

This class implements preconditioning for the conjugate gradient matrix solver.

This class implements diagonal preconditioning for the conjugate gradient matrix solver.

Definition at line 65 of file Preconditioner.h.

Constructor & Destructor Documentation

◆ Preconditioner()

Nektar::MultiRegions::Preconditioner::Preconditioner ( const std::shared_ptr< GlobalLinSys > &  plinsys,
const AssemblyMapSharedPtr pLocToGloMap 
)

Definition at line 58 of file Preconditioner.cpp.

60 : m_linsys(plinsys), m_preconType(pLocToGloMap->GetPreconType()),
61 m_locToGloMap(pLocToGloMap)
62{
63}
const std::weak_ptr< GlobalLinSys > m_linsys
std::weak_ptr< AssemblyMap > m_locToGloMap

◆ ~Preconditioner()

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

Definition at line 73 of file Preconditioner.h.

74 {
75 }

Member Function Documentation

◆ AssembleStaticCondGlobalDiagonals()

Array< OneD, NekDouble > Nektar::MultiRegions::Preconditioner::AssembleStaticCondGlobalDiagonals ( )

Performs global assembly of diagonal entries to global Schur complement matrix.

Definition at line 191 of file Preconditioner.cpp.

192{
193 auto asmMap = m_locToGloMap.lock();
194
195 int nGlobalBnd = asmMap->GetNumGlobalBndCoeffs();
196 int nDirBnd = asmMap->GetNumGlobalDirBndCoeffs();
197 int rows = nGlobalBnd - nDirBnd;
198
200 DNekScalMatSharedPtr bnd_mat;
201 int sign1, sign2, gid1, gid2, i, j, n, cnt;
202 Array<OneD, NekDouble> diagonals(rows, 0.0);
203
204 // Extract diagonal contributions of globally assembled
205 // schur complement matrix
206 for (cnt = n = 0; n < m_linsys.lock()->GetNumBlocks(); ++n)
207 {
208 // Get statically condensed local matrix.
209 loc_mat = (m_linsys.lock())->GetStaticCondBlock(n);
210
211 // Extract boundary block.
212 bnd_mat = loc_mat->GetBlock(0, 0);
213
214 // Offset by number of rows.
215 int bnd_row = bnd_mat->GetRows();
216
217 for (i = 0; i < bnd_row; ++i)
218 {
219 gid1 = asmMap->GetLocalToGlobalBndMap(cnt + i) - nDirBnd;
220 sign1 = asmMap->GetLocalToGlobalBndSign(cnt + i);
221
222 if (gid1 < 0)
223 {
224 continue;
225 }
226
227 for (j = 0; j < bnd_row; ++j)
228 {
229 gid2 = asmMap->GetLocalToGlobalBndMap(cnt + j) - nDirBnd;
230 sign2 = asmMap->GetLocalToGlobalBndSign(cnt + j);
231
232 if (gid2 == gid1)
233 {
234 diagonals[gid1] += sign1 * sign2 * (*bnd_mat)(i, j);
235 }
236 }
237 }
238 cnt += bnd_row;
239 }
240
241 return diagonals;
242}
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
std::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
Definition: NekTypeDefs.hpp:79

References m_linsys, and m_locToGloMap.

Referenced by Nektar::MultiRegions::PreconditionerDiagonal::StaticCondDiagonalPreconditionerSum().

◆ BuildPreconditioner()

void Nektar::MultiRegions::Preconditioner::BuildPreconditioner ( )
inline

Definition at line 241 of file Preconditioner.h.

242{
244}

References v_BuildPreconditioner().

◆ DoAssembleLoc()

void Nektar::MultiRegions::Preconditioner::DoAssembleLoc ( const Array< OneD, NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput,
const bool &  ZeroDir 
)

Apply an assembly and scatter back to lcoal array.

Definition at line 93 of file Preconditioner.cpp.

96{
97 auto assMap = m_locToGloMap.lock();
98 GlobalSysSolnType solvertype = assMap->GetGlobalSysSolnType();
99
100 if (solvertype == eIterativeFull)
101 {
102 assMap->Assemble(pInput, pOutput);
103 if (ZeroDir)
104 {
105 int nDir = assMap->GetNumGlobalDirBndCoeffs();
106 Vmath::Zero(nDir, pOutput, 1);
107 }
108 assMap->GlobalToLocal(pOutput, pOutput);
109 }
110 else // bnd version.
111 {
112 assMap->AssembleBnd(pInput, pOutput);
113 if (ZeroDir)
114 {
115 int nDir = assMap->GetNumGlobalDirBndCoeffs();
116 Vmath::Zero(nDir, pOutput, 1);
117 }
118 assMap->GlobalToLocalBnd(pOutput, pOutput);
119 }
120}
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:487

References Nektar::MultiRegions::eIterativeFull, m_locToGloMap, and Vmath::Zero().

◆ DoPreconditioner()

void Nektar::MultiRegions::Preconditioner::DoPreconditioner ( const Array< OneD, NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput,
const bool &  IsLocal = false 
)
inline

Definition at line 183 of file Preconditioner.h.

186{
187 v_DoPreconditioner(pInput, pOutput, ZeroDir);
188}
virtual void v_DoPreconditioner(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const bool &isLocal=false)
Apply a preconditioner to the conjugate gradient method.

References v_DoPreconditioner().

◆ DoPreconditionerWithNonVertOutput()

void Nektar::MultiRegions::Preconditioner::DoPreconditionerWithNonVertOutput ( const Array< OneD, NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput,
const Array< OneD, NekDouble > &  pNonVertOutput,
Array< OneD, NekDouble > &  pVertForce = NullNekDouble1DArray 
)
inline

Definition at line 193 of file Preconditioner.h.

197{
198 v_DoPreconditionerWithNonVertOutput(pInput, pOutput, pNonVertOutput,
199 pVertForce);
200}
virtual void v_DoPreconditionerWithNonVertOutput(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const Array< OneD, NekDouble > &pNonVertOutput, Array< OneD, NekDouble > &pVertForce)
Apply a preconditioner to the conjugate gradient method with an output for non-vertex degrees of free...

References v_DoPreconditionerWithNonVertOutput().

◆ DoTransformBasisFromLowEnergy()

void Nektar::MultiRegions::Preconditioner::DoTransformBasisFromLowEnergy ( const Array< OneD, NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput 
)
inline

Definition at line 232 of file Preconditioner.h.

234{
235 v_DoTransformBasisFromLowEnergy(pInput, pOutput);
236}
virtual void v_DoTransformBasisFromLowEnergy(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
Multiply by the block transposed inverse transformation matrix.

References v_DoTransformBasisFromLowEnergy().

◆ DoTransformBasisToLowEnergy()

void Nektar::MultiRegions::Preconditioner::DoTransformBasisToLowEnergy ( Array< OneD, NekDouble > &  pInOut)
inline

Definition at line 205 of file Preconditioner.h.

207{
209}
virtual void v_DoTransformBasisToLowEnergy(Array< OneD, NekDouble > &pInOut)
Transform from original basis to low energy basis.

References v_DoTransformBasisToLowEnergy().

◆ DoTransformCoeffsFromLowEnergy()

void Nektar::MultiRegions::Preconditioner::DoTransformCoeffsFromLowEnergy ( Array< OneD, NekDouble > &  pInOut)
inline

Definition at line 214 of file Preconditioner.h.

216{
218}
virtual void v_DoTransformCoeffsFromLowEnergy(Array< OneD, NekDouble > &pInOut)
Transform from low energy coeffs to orignal basis.

References v_DoTransformCoeffsFromLowEnergy().

◆ DoTransformCoeffsToLowEnergy()

void Nektar::MultiRegions::Preconditioner::DoTransformCoeffsToLowEnergy ( const Array< OneD, NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput 
)
inline

Definition at line 223 of file Preconditioner.h.

225{
226 v_DoTransformCoeffsToLowEnergy(pInput, pOutput);
227}
virtual void v_DoTransformCoeffsToLowEnergy(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
Multiply by the block inverse transformation matrix.

References v_DoTransformCoeffsToLowEnergy().

◆ GetBlockCMatrix()

const DNekScalBlkMatSharedPtr & Nektar::MultiRegions::Preconditioner::GetBlockCMatrix ( ) const
inline

◆ GetBlockInvDMatrix()

const DNekScalBlkMatSharedPtr & Nektar::MultiRegions::Preconditioner::GetBlockInvDMatrix ( ) const
inline

◆ GetBlockSchurCompl()

const DNekScalBlkMatSharedPtr & Nektar::MultiRegions::Preconditioner::GetBlockSchurCompl ( ) const
inline

◆ GetBlockTransformationMatrix()

const DNekScalBlkMatSharedPtr & Nektar::MultiRegions::Preconditioner::GetBlockTransformationMatrix ( ) const
inline

◆ GetBlockTransformedSchurCompl()

const DNekScalBlkMatSharedPtr & Nektar::MultiRegions::Preconditioner::GetBlockTransformedSchurCompl ( ) const
inline

◆ GetBlockTransposedTransformationMatrix()

const DNekScalBlkMatSharedPtr & Nektar::MultiRegions::Preconditioner::GetBlockTransposedTransformationMatrix ( ) const
inline

◆ InitObject()

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

Definition at line 166 of file Preconditioner.h.

167{
168 v_InitObject();
169}

References v_InitObject().

◆ NullPreconditioner()

void Nektar::MultiRegions::Preconditioner::NullPreconditioner ( void  )
private

◆ TransformedSchurCompl()

DNekScalMatSharedPtr Nektar::MultiRegions::Preconditioner::TransformedSchurCompl ( int  offset,
int  bndoffset,
const std::shared_ptr< DNekScalMat > &  loc_mat 
)
inline

Definition at line 174 of file Preconditioner.h.

176{
177 return v_TransformedSchurCompl(offset, bndoffset, loc_mat);
178}
virtual DNekScalMatSharedPtr v_TransformedSchurCompl(int offset, int bndoffset, const std::shared_ptr< DNekScalMat > &loc_mat)
Get block elemental transposed transformation matrix .

References v_TransformedSchurCompl().

◆ v_BuildPreconditioner()

void Nektar::MultiRegions::Preconditioner::v_BuildPreconditioner ( )
protectedvirtual

◆ v_DoPreconditioner()

void Nektar::MultiRegions::Preconditioner::v_DoPreconditioner ( const Array< OneD, NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput,
const bool &  isLocal = false 
)
protectedvirtual

Apply a preconditioner to the conjugate gradient method.

Reimplemented in Nektar::MultiRegions::PreconditionerBlock, Nektar::MultiRegions::PreconditionerDiagonal, Nektar::MultiRegions::PreconditionerNull, Nektar::MultiRegions::PreconditionerJacobi, Nektar::MultiRegions::PreconditionerLinear, Nektar::MultiRegions::PreconditionerLinearWithBlock, Nektar::MultiRegions::PreconditionerLinearWithDiag, Nektar::MultiRegions::PreconditionerLinearWithLowEnergy, and Nektar::MultiRegions::PreconditionerLowEnergy.

Definition at line 82 of file Preconditioner.cpp.

85{
86 boost::ignore_unused(pInput, pOutput, IsLocal);
87 NEKERROR(ErrorUtil::efatal, "Method does not exist");
88}
#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 DoPreconditioner().

◆ v_DoPreconditionerWithNonVertOutput()

void Nektar::MultiRegions::Preconditioner::v_DoPreconditionerWithNonVertOutput ( const Array< OneD, NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput,
const Array< OneD, NekDouble > &  pNonVertOutput,
Array< OneD, NekDouble > &  pVertForce 
)
protectedvirtual

Apply a preconditioner to the conjugate gradient method with an output for non-vertex degrees of freedom.

Reimplemented in Nektar::MultiRegions::PreconditionerLinear.

Definition at line 126 of file Preconditioner.cpp.

130{
131 boost::ignore_unused(pInput, pOutput, pNonVertOutput, pVertForce);
132 NEKERROR(ErrorUtil::efatal, "Method does not exist");
133}

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

Referenced by DoPreconditionerWithNonVertOutput().

◆ v_DoTransformBasisFromLowEnergy()

void Nektar::MultiRegions::Preconditioner::v_DoTransformBasisFromLowEnergy ( const Array< OneD, NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput 
)
protectedvirtual

Multiply by the block transposed inverse transformation matrix.

Reimplemented in Nektar::MultiRegions::PreconditionerLinearWithLowEnergy, and Nektar::MultiRegions::PreconditionerLowEnergy.

Definition at line 165 of file Preconditioner.cpp.

167{
168 boost::ignore_unused(pInput, pOutput);
169 NEKERROR(ErrorUtil::efatal, "Method does not exist");
170}

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

Referenced by DoTransformBasisFromLowEnergy().

◆ v_DoTransformBasisToLowEnergy()

void Nektar::MultiRegions::Preconditioner::v_DoTransformBasisToLowEnergy ( Array< OneD, NekDouble > &  pInOut)
protectedvirtual

Transform from original basis to low energy basis.

Reimplemented in Nektar::MultiRegions::PreconditionerLinearWithLowEnergy, and Nektar::MultiRegions::PreconditionerLowEnergy.

Definition at line 138 of file Preconditioner.cpp.

140{
141 boost::ignore_unused(pInOut);
142}

Referenced by DoTransformBasisToLowEnergy().

◆ v_DoTransformCoeffsFromLowEnergy()

void Nektar::MultiRegions::Preconditioner::v_DoTransformCoeffsFromLowEnergy ( Array< OneD, NekDouble > &  pInOut)
protectedvirtual

Transform from low energy coeffs to orignal basis.

Reimplemented in Nektar::MultiRegions::PreconditionerLinearWithLowEnergy, and Nektar::MultiRegions::PreconditionerLowEnergy.

Definition at line 147 of file Preconditioner.cpp.

149{
150 boost::ignore_unused(pInOut);
151}

Referenced by DoTransformCoeffsFromLowEnergy().

◆ v_DoTransformCoeffsToLowEnergy()

void Nektar::MultiRegions::Preconditioner::v_DoTransformCoeffsToLowEnergy ( const Array< OneD, NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput 
)
protectedvirtual

Multiply by the block inverse transformation matrix.

Reimplemented in Nektar::MultiRegions::PreconditionerLinearWithLowEnergy, and Nektar::MultiRegions::PreconditionerLowEnergy.

Definition at line 156 of file Preconditioner.cpp.

158{
159 boost::ignore_unused(pInput, pOutput);
160}

Referenced by DoTransformCoeffsToLowEnergy().

◆ v_InitObject()

void Nektar::MultiRegions::Preconditioner::v_InitObject ( )
protectedvirtual

◆ v_TransformedSchurCompl()

DNekScalMatSharedPtr Nektar::MultiRegions::Preconditioner::v_TransformedSchurCompl ( int  offset,
int  bndoffset,
const std::shared_ptr< DNekScalMat > &  loc_mat 
)
protectedvirtual

Get block elemental transposed transformation matrix \(\mathbf{R}^{T}\).

Reimplemented in Nektar::MultiRegions::PreconditionerLinearWithLowEnergy, and Nektar::MultiRegions::PreconditionerLowEnergy.

Definition at line 180 of file Preconditioner.cpp.

182{
183 boost::ignore_unused(offset, bnd_offset);
184 return loc_mat;
185}

Referenced by TransformedSchurCompl().

Member Data Documentation

◆ def

std::string Nektar::MultiRegions::Preconditioner::def
staticprivate
Initial value:
=
"Diagonal")
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 Preconditioner.h.

◆ lookupIds

std::string Nektar::MultiRegions::Preconditioner::lookupIds[]
staticprivate

Definition at line 158 of file Preconditioner.h.

◆ m_comm

LibUtilities::CommSharedPtr Nektar::MultiRegions::Preconditioner::m_comm
protected

◆ m_linsys

const std::weak_ptr<GlobalLinSys> Nektar::MultiRegions::Preconditioner::m_linsys
protected

◆ m_locToGloMap

std::weak_ptr<AssemblyMap> Nektar::MultiRegions::Preconditioner::m_locToGloMap
protected

Definition at line 125 of file Preconditioner.h.

Referenced by AssembleStaticCondGlobalDiagonals(), Nektar::MultiRegions::PreconditionerBlock::BlockPreconditionerCG(), Nektar::MultiRegions::PreconditionerBlock::BlockPreconditionerHDG(), Nektar::MultiRegions::PreconditionerLowEnergy::CreateVariablePMask(), Nektar::MultiRegions::PreconditionerDiagonal::DiagonalPreconditionerSum(), DoAssembleLoc(), Nektar::MultiRegions::PreconditionerLowEnergy::SetupBlockTransformationMatrix(), Nektar::MultiRegions::PreconditionerDiagonal::StaticCondDiagonalPreconditionerSum(), Nektar::MultiRegions::PreconditionerDiagonal::v_BuildPreconditioner(), Nektar::MultiRegions::PreconditionerJacobi::v_BuildPreconditioner(), Nektar::MultiRegions::PreconditionerLinear::v_BuildPreconditioner(), Nektar::MultiRegions::PreconditionerLowEnergy::v_BuildPreconditioner(), Nektar::MultiRegions::PreconditionerBlock::v_DoPreconditioner(), Nektar::MultiRegions::PreconditionerDiagonal::v_DoPreconditioner(), Nektar::MultiRegions::PreconditionerJacobi::v_DoPreconditioner(), Nektar::MultiRegions::PreconditionerLinearWithLowEnergy::v_DoPreconditioner(), Nektar::MultiRegions::PreconditionerLowEnergy::v_DoPreconditioner(), Nektar::MultiRegions::PreconditionerLinear::v_DoPreconditionerWithNonVertOutput(), Nektar::MultiRegions::PreconditionerLowEnergy::v_DoTransformBasisFromLowEnergy(), Nektar::MultiRegions::PreconditionerLowEnergy::v_DoTransformBasisToLowEnergy(), Nektar::MultiRegions::PreconditionerLowEnergy::v_DoTransformCoeffsFromLowEnergy(), Nektar::MultiRegions::PreconditionerLowEnergy::v_DoTransformCoeffsToLowEnergy(), Nektar::MultiRegions::PreconditionerBlock::v_InitObject(), Nektar::MultiRegions::PreconditionerLinearWithBlock::v_InitObject(), Nektar::MultiRegions::PreconditionerLinearWithDiag::v_InitObject(), Nektar::MultiRegions::PreconditionerLinearWithLowEnergy::v_InitObject(), and Nektar::MultiRegions::PreconditionerLowEnergy::v_InitObject().

◆ m_preconditioner

DNekMatSharedPtr Nektar::MultiRegions::Preconditioner::m_preconditioner
protected

Definition at line 124 of file Preconditioner.h.

◆ m_preconType

std::string Nektar::MultiRegions::Preconditioner::m_preconType
protected

Definition at line 123 of file Preconditioner.h.