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 Member Functions | |
void | v_InitObject () override |
void | v_DoPreconditioner (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const bool &isLocal=false) override |
Apply preconditioner to pInput and store the result in pOutput . More... | |
void | v_BuildPreconditioner () override |
Protected Member Functions inherited from Nektar::MultiRegions::Preconditioner | |
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 | |
bool | m_isFull |
DNekBlkMatSharedPtr | m_blkMat |
Gs::gs_data * | m_gs_mapping |
Protected Attributes inherited from Nektar::MultiRegions::Preconditioner | |
const std::weak_ptr< GlobalLinSys > | m_linsys |
std::string | m_preconType |
DNekMatSharedPtr | m_preconditioner |
std::weak_ptr< AssemblyMap > | m_locToGloMap |
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... | |
Definition at line 48 of file PreconditionerBlock.h.
Nektar::MultiRegions::PreconditionerBlock::PreconditionerBlock | ( | const std::shared_ptr< GlobalLinSys > & | plinsys, |
const AssemblyMapSharedPtr & | pLocToGloMap | ||
) |
Definition at line 61 of file PreconditionerBlock.cpp.
|
inlineoverride |
Definition at line 71 of file PreconditionerBlock.h.
References Gs::Free(), and m_gs_mapping.
Construct a block preconditioner from \(\mathbf{S}_{1}\) for the continuous Galerkin system.
The preconditioner for the statically-condensed system 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 (vv), edge (eb) and face blocks (fb) respectively.
In case of the full system \(\mathbf{A}\), the preconditioner includes the interior blocks as well. The full block preconditioner is defined as:
\[\mathbf{M}^{-1}=\left[\begin{array}{cccc} \mathrm{Diag}[(\mathbf{A})_{vv}] & & & \\ & (\mathbf{A})_{eb} & & \\ & & (\mathbf{A})_{fb} & \\ & & & (\mathbf{A})_{int} \end{array}\right] \]
where the interior portion is added at the end following the ordering inside data structures.
Process boundary matrices
Definition at line 123 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, m_gs_mapping, m_isFull, 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 657 of file PreconditionerBlock.cpp.
References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::eDIAGONAL, Gs::Gather(), Gs::gs_add, Gs::Init(), m_blkMat, 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 52 of file PreconditionerBlock.h.
References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), and CellMLToNektar.cellml_metadata::p.
|
overrideprotectedvirtual |
Reimplemented from Nektar::MultiRegions::Preconditioner.
Definition at line 79 of file PreconditionerBlock.cpp.
References BlockPreconditionerCG(), BlockPreconditionerHDG(), Nektar::StdRegions::eHybridDGHelmBndLam, Nektar::MultiRegions::GlobalMatrixKey::GetMatrixType(), and Nektar::MultiRegions::Preconditioner::m_linsys.
|
overrideprotectedvirtual |
Apply preconditioner to pInput
and store the result in pOutput
.
Reimplemented from Nektar::MultiRegions::Preconditioner.
Definition at line 814 of file PreconditionerBlock.cpp.
References Nektar::eWrapper, m_isFull, Nektar::MultiRegions::Preconditioner::m_locToGloMap, Nektar::UnitTests::z(), and Vmath::Zero().
|
overrideprotectedvirtual |
Reimplemented from Nektar::MultiRegions::Preconditioner.
Definition at line 68 of file PreconditionerBlock.cpp.
References ASSERTL0, Nektar::MultiRegions::eIterativeFull, Nektar::MultiRegions::eIterativeStaticCond, Nektar::MultiRegions::ePETScFullMatrix, Nektar::MultiRegions::ePETScStaticCond, m_isFull, and Nektar::MultiRegions::Preconditioner::m_locToGloMap.
|
static |
Name of class.
Registers the class with the Factory.
Definition at line 64 of file PreconditionerBlock.h.
|
protected |
Definition at line 78 of file PreconditionerBlock.h.
Referenced by BlockPreconditionerCG(), and BlockPreconditionerHDG().
|
protected |
Definition at line 81 of file PreconditionerBlock.h.
Referenced by BlockPreconditionerCG(), and ~PreconditionerBlock().
|
protected |
Definition at line 77 of file PreconditionerBlock.h.
Referenced by BlockPreconditionerCG(), v_DoPreconditioner(), and v_InitObject().