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 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)
 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
 
PreconditionerType 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 76 of file Preconditioner.cpp.

78  : m_linsys(plinsys), m_preconType(pLocToGloMap->GetPreconType()),
79  m_locToGloMap(pLocToGloMap)
80 {
81 }
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 176 of file Preconditioner.cpp.

177 {
178  auto asmMap = m_locToGloMap.lock();
179 
180  int nGlobalBnd = asmMap->GetNumGlobalBndCoeffs();
181  int nDirBnd = asmMap->GetNumGlobalDirBndCoeffs();
182  int rows = nGlobalBnd - nDirBnd;
183 
184  DNekScalBlkMatSharedPtr loc_mat;
185  DNekScalMatSharedPtr bnd_mat;
186  int sign1, sign2, gid1, gid2, i, j, n, cnt;
187  Array<OneD, NekDouble> diagonals(rows, 0.0);
188 
189  // Extract diagonal contributions of globally assembled
190  // schur complement matrix
191  for (cnt = n = 0; n < m_linsys.lock()->GetNumBlocks(); ++n)
192  {
193  // Get statically condensed local matrix.
194  loc_mat = (m_linsys.lock())->GetStaticCondBlock(n);
195 
196  // Extract boundary block.
197  bnd_mat = loc_mat->GetBlock(0, 0);
198 
199  // Offset by number of rows.
200  int bnd_row = bnd_mat->GetRows();
201 
202  for (i = 0; i < bnd_row; ++i)
203  {
204  gid1 = asmMap->GetLocalToGlobalBndMap(cnt + i) - nDirBnd;
205  sign1 = asmMap->GetLocalToGlobalBndSign(cnt + i);
206 
207  if (gid1 < 0)
208  {
209  continue;
210  }
211 
212  for (j = 0; j < bnd_row; ++j)
213  {
214  gid2 = asmMap->GetLocalToGlobalBndMap(cnt + j) - nDirBnd;
215  sign2 = asmMap->GetLocalToGlobalBndSign(cnt + j);
216 
217  if (gid2 == gid1)
218  {
219  diagonals[gid1] += sign1 * sign2 * (*bnd_mat)(i, j);
220  }
221  }
222  }
223  cnt += bnd_row;
224  }
225 
226  return diagonals;
227 }
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 235 of file Preconditioner.h.

236 {
238 }

References v_BuildPreconditioner().

◆ DoPreconditioner()

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

Definition at line 178 of file Preconditioner.h.

180 {
181  v_DoPreconditioner(pInput, pOutput);
182 }
virtual void v_DoPreconditioner(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
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 187 of file Preconditioner.h.

191 {
192  v_DoPreconditionerWithNonVertOutput(pInput, pOutput, pNonVertOutput,
193  pVertForce);
194 }
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 226 of file Preconditioner.h.

228 {
229  v_DoTransformBasisFromLowEnergy(pInput, pOutput);
230 }
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 199 of file Preconditioner.h.

201 {
203 }
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 208 of file Preconditioner.h.

210 {
212 }
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 217 of file Preconditioner.h.

219 {
220  v_DoTransformCoeffsToLowEnergy(pInput, pOutput);
221 }
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 161 of file Preconditioner.h.

162 {
163  v_InitObject();
164 }

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

171 {
172  return v_TransformedSchurCompl(offset, bndoffset, loc_mat);
173 }
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 
)
protectedvirtual

Apply a preconditioner to the conjugate gradient method.

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

Definition at line 100 of file Preconditioner.cpp.

102 {
103  boost::ignore_unused(pInput, pOutput);
104  NEKERROR(ErrorUtil::efatal, "Method does not exist");
105 }
#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 111 of file Preconditioner.cpp.

115 {
116  boost::ignore_unused(pInput, pOutput, pNonVertOutput, pVertForce);
117  NEKERROR(ErrorUtil::efatal, "Method does not exist");
118 }

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

Definition at line 150 of file Preconditioner.cpp.

152 {
153  boost::ignore_unused(pInput, pOutput);
154  NEKERROR(ErrorUtil::efatal, "Method does not exist");
155 }

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

Definition at line 123 of file Preconditioner.cpp.

125 {
126  boost::ignore_unused(pInOut);
127 }

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

Definition at line 132 of file Preconditioner.cpp.

134 {
135  boost::ignore_unused(pInOut);
136 }

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

Definition at line 141 of file Preconditioner.cpp.

143 {
144  boost::ignore_unused(pInput, pOutput);
145 }

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

Definition at line 165 of file Preconditioner.cpp.

167 {
168  boost::ignore_unused(offset, bnd_offset);
169  return loc_mat;
170 }

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

◆ lookupIds

std::string Nektar::MultiRegions::Preconditioner::lookupIds
staticprivate
Initial value:
= {
LibUtilities::SessionReader::RegisterEnumValue("Preconditioner", "Diagonal",
"Preconditioner", "FullLinearSpaceWithDiagonal", eLinearWithDiagonal),
"FullLinearSpace", eLinear),
"Preconditioner", "LowEnergyBlock", eLowEnergy),
"Preconditioner", "FullLinearSpaceWithLowEnergyBlock",
"Preconditioner", "FullLinearSpaceWithBlock", eLinearWithBlock),
}
static std::string RegisterEnumValue(std::string pEnum, std::string pString, int pEnumValue)
Registers an enumeration value.
@ eNull
No Solution type specified.

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

Referenced by AssembleStaticCondGlobalDiagonals(), Nektar::MultiRegions::PreconditionerBlock::BlockPreconditionerCG(), Nektar::MultiRegions::PreconditionerBlock::BlockPreconditionerHDG(), Nektar::MultiRegions::PreconditionerLowEnergy::CreateVariablePMask(), 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::PreconditionerBlock::v_DoPreconditioner(), Nektar::MultiRegions::PreconditionerDiagonal::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 120 of file Preconditioner.h.

◆ m_preconType

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

Definition at line 119 of file Preconditioner.h.