Nektar++
|
#include <PreconditionerBlock.h>
Static Public Member Functions | |
static PreconditionerSharedPtr | create (const std::shared_ptr< GlobalLinSys > &plinsys, const std::shared_ptr< AssemblyMap > &pLocToGloMap) |
Creates an instance of this class. More... | |
Static Public Attributes | |
static std::string | className |
Name of class. More... | |
Protected Attributes | |
DNekBlkMatSharedPtr | m_blkMat |
![]() | |
const std::weak_ptr< GlobalLinSys > | m_linsys |
PreconditionerType | m_preconType |
DNekMatSharedPtr | m_preconditioner |
std::weak_ptr< AssemblyMap > | m_locToGloMap |
LibUtilities::CommSharedPtr | m_comm |
Private Member Functions | |
void | BlockPreconditionerCG (void) |
Construct a block preconditioner from \(\mathbf{S}_{1}\) for the continuous Galerkin system. More... | |
void | BlockPreconditionerHDG (void) |
Construct a block preconditioner for the hybridized discontinuous Galerkin system. More... | |
virtual void | v_InitObject () |
virtual void | v_DoPreconditioner (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput) |
Apply preconditioner to pInput and store the result in pOutput . More... | |
virtual void | v_BuildPreconditioner () |
Additional Inherited Members | |
![]() | |
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... | |
Definition at line 51 of file PreconditionerBlock.h.
Nektar::MultiRegions::PreconditionerBlock::PreconditionerBlock | ( | const std::shared_ptr< GlobalLinSys > & | plinsys, |
const AssemblyMapSharedPtr & | pLocToGloMap | ||
) |
Definition at line 66 of file PreconditionerBlock.cpp.
|
inlinevirtual |
Definition at line 73 of file PreconditionerBlock.h.
Construct a block preconditioner from \(\mathbf{S}_{1}\) for the continuous Galerkin system.
The preconditioner is defined as:
\[\mathbf{M}^{-1}=\left[\begin{array}{ccc} \mathrm{Diag}[(\mathbf{S_{1}})_{vv}] & & \\ & (\mathbf{S}_{1})_{eb} & \\ & & (\mathbf{S}_{1})_{fb} \end{array}\right] \]
where \(\mathbf{S}_{1}\) is the local Schur complement matrix for each element and the subscript denotes the portion of the Schur complement associated with the vertex, edge and face blocks respectively.
Definition at line 113 of file PreconditionerBlock.cpp.
References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL1, Nektar::MultiRegions::DeterminePeriodicEdgeOrientId(), Nektar::MultiRegions::DeterminePeriodicFaceOrient(), Nektar::eDIAGONAL, Gs::Gather(), Gs::gs_add, Gs::Init(), m_blkMat, Nektar::MultiRegions::Preconditioner::m_comm, Nektar::MultiRegions::Preconditioner::m_linsys, Nektar::MultiRegions::Preconditioner::m_locToGloMap, Nektar::LibUtilities::ReduceMax, sign, and Vmath::Vadd().
Referenced by v_BuildPreconditioner().
Construct a block preconditioner for the hybridized discontinuous Galerkin system.
This system is built in a similar fashion to its continuous variant found in PreconditionerBlock::BlockPreconditionerCG. In this setting however, the matrix is constructed as:
\[ M^{-1} = \mathrm{Diag}[ (\mathbf{S_{1}})_{f}^{-1} ] \]
where each matrix is the Schur complement system restricted to a single face of the trace system.
Definition at line 539 of file PreconditionerBlock.cpp.
References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::eDIAGONAL, Gs::Gather(), Nektar::MultiRegions::AssemblyMap::GetNumGlobalDirBndCoeffs(), Gs::gs_add, Gs::Init(), m_blkMat, Nektar::MultiRegions::Preconditioner::m_comm, Nektar::MultiRegions::Preconditioner::m_linsys, Nektar::MultiRegions::Preconditioner::m_locToGloMap, and Nektar::LibUtilities::ReduceMax.
Referenced by v_BuildPreconditioner().
|
inlinestatic |
Creates an instance of this class.
Definition at line 55 of file PreconditionerBlock.h.
References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), and CellMLToNektar.cellml_metadata::p.
|
privatevirtual |
Reimplemented from Nektar::MultiRegions::Preconditioner.
Definition at line 82 of file PreconditionerBlock.cpp.
References BlockPreconditionerCG(), BlockPreconditionerHDG(), Nektar::StdRegions::eHybridDGHelmBndLam, Nektar::MultiRegions::GlobalMatrixKey::GetMatrixType(), and Nektar::MultiRegions::Preconditioner::m_linsys.
|
privatevirtual |
Apply preconditioner to pInput
and store the result in pOutput
.
Reimplemented from Nektar::MultiRegions::Preconditioner.
Definition at line 698 of file PreconditionerBlock.cpp.
References Nektar::eWrapper, and Nektar::MultiRegions::Preconditioner::m_locToGloMap.
|
privatevirtual |
Reimplemented from Nektar::MultiRegions::Preconditioner.
Definition at line 73 of file PreconditionerBlock.cpp.
References ASSERTL0, Nektar::MultiRegions::eIterativeStaticCond, Nektar::MultiRegions::ePETScStaticCond, and Nektar::MultiRegions::Preconditioner::m_locToGloMap.
|
static |
Name of class.
Registers the class with the Factory.
Definition at line 66 of file PreconditionerBlock.h.
|
protected |
Definition at line 76 of file PreconditionerBlock.h.
Referenced by BlockPreconditionerCG(), and BlockPreconditionerHDG().