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)
 
void DoPreconditionerWithNonVertOutput (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const Array< OneD, NekDouble > &pNonVertOutput, Array< OneD, NekDouble > &pVertForce=NullNekDouble1DArray)
 
void DoTransformToLowEnergy (Array< OneD, NekDouble > &pInOut, int offset)
 
void DoTransformToLowEnergy (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
 
void DoTransformFromLowEnergy (Array< OneD, NekDouble > &pInOut)
 
void DoMultiplybyInverseTransformationMatrix (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
 
void DoMultiplybyInverseTransposedTransformationMatrix (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...
 

Protected Attributes

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

Private Member Functions

void NullPreconditioner (void)
 
virtual void v_InitObject ()
 
virtual void v_DoPreconditioner (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
 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_DoTransformToLowEnergy (Array< OneD, NekDouble > &pInOut, int offset)
 Transform from original basis to low energy basis. More...
 
virtual void v_DoTransformToLowEnergy (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
 Transform from original basis to low energy basis. More...
 
virtual void v_DoTransformFromLowEnergy (Array< OneD, NekDouble > &pInput)
 Transform from low energy basis to orignal basis. More...
 
virtual void v_DoMultiplybyInverseTransformationMatrix (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
 Multiply by the block inverse transformation matrix. More...
 
virtual void v_DoMultiplybyInverseTransposedTransformationMatrix (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
 Multiply by the block transposed inverse transformation matrix. More...
 
virtual void v_BuildPreconditioner ()
 

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 64 of file Preconditioner.h.

Constructor & Destructor Documentation

◆ Preconditioner()

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

Definition at line 77 of file Preconditioner.cpp.

80  : m_linsys(plinsys),
81  m_preconType(pLocToGloMap->GetPreconType()),
82  m_locToGloMap(pLocToGloMap)
83  {
84  }
const std::weak_ptr< GlobalLinSys > m_linsys
std::weak_ptr< AssemblyMap > m_locToGloMap

◆ ~Preconditioner()

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

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 197 of file Preconditioner.cpp.

References m_linsys, and m_locToGloMap.

Referenced by Nektar::MultiRegions::PreconditionerDiagonal::StaticCondDiagonalPreconditionerSum(), and ~Preconditioner().

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

◆ BuildPreconditioner()

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

Definition at line 281 of file Preconditioner.h.

References v_BuildPreconditioner().

Referenced by ~Preconditioner().

282  {
284  }

◆ DoMultiplybyInverseTransformationMatrix()

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

Definition at line 260 of file Preconditioner.h.

References DoMultiplybyInverseTransposedTransformationMatrix(), and v_DoMultiplybyInverseTransformationMatrix().

Referenced by ~Preconditioner().

263  {
265  }
virtual void v_DoMultiplybyInverseTransformationMatrix(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
Multiply by the block inverse transformation matrix.

◆ DoMultiplybyInverseTransposedTransformationMatrix()

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

Definition at line 271 of file Preconditioner.h.

References v_DoMultiplybyInverseTransposedTransformationMatrix().

Referenced by DoMultiplybyInverseTransformationMatrix(), and ~Preconditioner().

274  {
276  }
virtual void v_DoMultiplybyInverseTransposedTransformationMatrix(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
Multiply by the block transposed inverse transformation matrix.

◆ DoPreconditioner()

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

Definition at line 208 of file Preconditioner.h.

References v_DoPreconditioner().

Referenced by ~Preconditioner().

211  {
212  v_DoPreconditioner(pInput,pOutput);
213  }
virtual void v_DoPreconditioner(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
Apply a preconditioner to the conjugate gradient method.

◆ 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 219 of file Preconditioner.h.

References v_DoPreconditionerWithNonVertOutput().

Referenced by ~Preconditioner().

224  {
225  v_DoPreconditionerWithNonVertOutput(pInput,pOutput,pNonVertOutput,
226  pVertForce);
227  }
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...

◆ DoTransformFromLowEnergy()

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

Definition at line 251 of file Preconditioner.h.

References v_DoTransformFromLowEnergy().

Referenced by ~Preconditioner().

253  {
255  }
virtual void v_DoTransformFromLowEnergy(Array< OneD, NekDouble > &pInput)
Transform from low energy basis to orignal basis.

◆ DoTransformToLowEnergy() [1/2]

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

Definition at line 232 of file Preconditioner.h.

References v_DoTransformToLowEnergy().

Referenced by ~Preconditioner().

234  {
235  v_DoTransformToLowEnergy(pInOut,offset);
236  }
virtual void v_DoTransformToLowEnergy(Array< OneD, NekDouble > &pInOut, int offset)
Transform from original basis to low energy basis.

◆ DoTransformToLowEnergy() [2/2]

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

Definition at line 241 of file Preconditioner.h.

References v_DoTransformToLowEnergy().

244  {
245  v_DoTransformToLowEnergy(pInput,pOutput);
246  }
virtual void v_DoTransformToLowEnergy(Array< OneD, NekDouble > &pInOut, int offset)
Transform from original basis to low energy basis.

◆ GetBlockCMatrix()

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

Referenced by ~Preconditioner().

◆ GetBlockInvDMatrix()

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

Referenced by ~Preconditioner().

◆ GetBlockSchurCompl()

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

Referenced by ~Preconditioner().

◆ GetBlockTransformationMatrix()

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

Referenced by ~Preconditioner().

◆ GetBlockTransformedSchurCompl()

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

Referenced by ~Preconditioner().

◆ GetBlockTransposedTransformationMatrix()

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

Referenced by ~Preconditioner().

◆ InitObject()

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

Definition at line 190 of file Preconditioner.h.

References v_InitObject().

Referenced by ~Preconditioner().

191  {
192  v_InitObject();
193  }

◆ 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 198 of file Preconditioner.h.

References v_TransformedSchurCompl().

Referenced by ~Preconditioner().

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

◆ v_BuildPreconditioner()

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

◆ v_DoMultiplybyInverseTransformationMatrix()

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

Multiply by the block inverse transformation matrix.

Reimplemented in Nektar::MultiRegions::PreconditionerLowEnergy.

Definition at line 157 of file Preconditioner.cpp.

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

Referenced by DoMultiplybyInverseTransformationMatrix().

160  {
161  boost::ignore_unused(pInput, pOutput);
162  NEKERROR(ErrorUtil::efatal,"Method does not exist" );
163  }
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mod...
Definition: ErrorUtil.hpp:209

◆ v_DoMultiplybyInverseTransposedTransformationMatrix()

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

Multiply by the block transposed inverse transformation matrix.

Reimplemented in Nektar::MultiRegions::PreconditionerLowEnergy.

Definition at line 168 of file Preconditioner.cpp.

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

Referenced by DoMultiplybyInverseTransposedTransformationMatrix().

171  {
172  boost::ignore_unused(pInput, pOutput);
173  NEKERROR(ErrorUtil::efatal,"Method does not exist" );
174  }
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mod...
Definition: ErrorUtil.hpp:209

◆ v_DoPreconditioner()

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

Apply a preconditioner to the conjugate gradient method.

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

Definition at line 103 of file Preconditioner.cpp.

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

Referenced by DoPreconditioner().

106  {
107  boost::ignore_unused(pInput, pOutput);
108  NEKERROR(ErrorUtil::efatal,"Method does not exist" );
109  }
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mod...
Definition: ErrorUtil.hpp:209

◆ 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 
)
privatevirtual

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 115 of file Preconditioner.cpp.

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

Referenced by DoPreconditionerWithNonVertOutput().

120  {
121  boost::ignore_unused(pInput, pOutput, pNonVertOutput, pVertForce);
122  NEKERROR(ErrorUtil::efatal, "Method does not exist");
123  }
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mod...
Definition: ErrorUtil.hpp:209

◆ v_DoTransformFromLowEnergy()

void Nektar::MultiRegions::Preconditioner::v_DoTransformFromLowEnergy ( Array< OneD, NekDouble > &  pInput)
privatevirtual

Transform from low energy basis to orignal basis.

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

Definition at line 148 of file Preconditioner.cpp.

References Vmath::Smul().

Referenced by DoTransformFromLowEnergy().

150  {
151  Vmath::Smul(pInput.num_elements(), 1.0, pInput, 1, pInput, 1);
152  }
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:216

◆ v_DoTransformToLowEnergy() [1/2]

void Nektar::MultiRegions::Preconditioner::v_DoTransformToLowEnergy ( Array< OneD, NekDouble > &  pInOut,
int  offset 
)
privatevirtual

Transform from original basis to low energy basis.

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

Definition at line 128 of file Preconditioner.cpp.

Referenced by DoTransformToLowEnergy().

131  {
132  boost::ignore_unused(pInOut, offset);
133  }

◆ v_DoTransformToLowEnergy() [2/2]

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

Transform from original basis to low energy basis.

Reimplemented in Nektar::MultiRegions::PreconditionerLowEnergy.

Definition at line 138 of file Preconditioner.cpp.

141  {
142  boost::ignore_unused(pInOut, pOutput);
143  }

◆ v_InitObject()

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

◆ 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::PreconditionerLowEnergy, and Nektar::MultiRegions::PreconditionerLinearWithLowEnergy.

Definition at line 184 of file Preconditioner.cpp.

Referenced by TransformedSchurCompl().

187  {
188  boost::ignore_unused(offset, bnd_offset);
189  return loc_mat;
190  }

Member Data Documentation

◆ def

std::string Nektar::MultiRegions::Preconditioner::def
staticprivate
Initial value:

Definition at line 183 of file Preconditioner.h.

◆ lookupIds

std::string Nektar::MultiRegions::Preconditioner::lookupIds
staticprivate
Initial value:
= {
"Preconditioner", "Null", eNull),
"Preconditioner", "Diagonal", eDiagonal),
"Preconditioner", "FullLinearSpaceWithDiagonal",
"Preconditioner", "FullLinearSpace",eLinear),
"Preconditioner", "LowEnergyBlock",eLowEnergy),
"Preconditioner", "FullLinearSpaceWithLowEnergyBlock",
"Preconditioner", "Block",eBlock),
"Preconditioner", "FullLinearSpaceWithBlock",eLinearWithBlock),
}

Definition at line 182 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 136 of file Preconditioner.h.

Referenced by AssembleStaticCondGlobalDiagonals(), Nektar::MultiRegions::PreconditionerBlock::BlockPreconditionerCG(), Nektar::MultiRegions::PreconditionerBlock::BlockPreconditionerHDG(), Nektar::MultiRegions::PreconditionerLowEnergy::CreateMultiplicityMap(), Nektar::MultiRegions::PreconditionerDiagonal::DiagonalPreconditionerSum(), Nektar::MultiRegions::PreconditionerLowEnergy::SetupBlockTransformationMatrix(), Nektar::MultiRegions::PreconditionerDiagonal::StaticCondDiagonalPreconditionerSum(), Nektar::MultiRegions::PreconditionerDiagonal::v_BuildPreconditioner(), Nektar::MultiRegions::PreconditionerLinear::v_BuildPreconditioner(), Nektar::MultiRegions::PreconditionerLowEnergy::v_BuildPreconditioner(), Nektar::MultiRegions::PreconditionerLowEnergy::v_DoMultiplybyInverseTransformationMatrix(), Nektar::MultiRegions::PreconditionerLowEnergy::v_DoMultiplybyInverseTransposedTransformationMatrix(), Nektar::MultiRegions::PreconditionerBlock::v_DoPreconditioner(), Nektar::MultiRegions::PreconditionerDiagonal::v_DoPreconditioner(), Nektar::MultiRegions::PreconditionerLowEnergy::v_DoPreconditioner(), Nektar::MultiRegions::PreconditionerLinear::v_DoPreconditionerWithNonVertOutput(), Nektar::MultiRegions::PreconditionerLowEnergy::v_DoTransformFromLowEnergy(), Nektar::MultiRegions::PreconditionerLowEnergy::v_DoTransformToLowEnergy(), Nektar::MultiRegions::PreconditionerLinearWithLowEnergy::v_InitObject(), Nektar::MultiRegions::PreconditionerLinearWithBlock::v_InitObject(), Nektar::MultiRegions::PreconditionerLinearWithDiag::v_InitObject(), Nektar::MultiRegions::PreconditionerBlock::v_InitObject(), and Nektar::MultiRegions::PreconditionerLowEnergy::v_InitObject().

◆ m_preconditioner

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

Definition at line 135 of file Preconditioner.h.

◆ m_preconType

PreconditionerType Nektar::MultiRegions::Preconditioner::m_preconType
protected

Definition at line 134 of file Preconditioner.h.