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

#include <PreconditionerLinearWithBlock.h>

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

Public Member Functions

 PreconditionerLinearWithBlock (const std::shared_ptr< GlobalLinSys > &plinsys, const AssemblyMapSharedPtr &pLocToGloMap)
 
 ~PreconditionerLinearWithBlock () override
 
- Public Member Functions inherited from Nektar::MultiRegions::Preconditioner
 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)
 

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 a preconditioner to the conjugate gradient method. 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

PreconditionerSharedPtr m_linSpacePrecon
 
PreconditionerSharedPtr m_blockPrecon
 
- Protected Attributes inherited from Nektar::MultiRegions::Preconditioner
const std::weak_ptr< GlobalLinSysm_linsys
 
std::string m_preconType
 
DNekMatSharedPtr m_preconditioner
 
std::weak_ptr< AssemblyMapm_locToGloMap
 
LibUtilities::CommSharedPtr m_comm
 

Detailed Description

This class implements preconditioning for the conjugate gradient matrix solver.

Definition at line 47 of file PreconditionerLinearWithBlock.h.

Constructor & Destructor Documentation

◆ PreconditionerLinearWithBlock()

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

Definition at line 60 of file PreconditionerLinearWithBlock.cpp.

63 : Preconditioner(plinsys, pLocToGloMap)
64{
65}
Preconditioner(const std::shared_ptr< GlobalLinSys > &plinsys, const AssemblyMapSharedPtr &pLocToGloMap)

◆ ~PreconditionerLinearWithBlock()

Nektar::MultiRegions::PreconditionerLinearWithBlock::~PreconditionerLinearWithBlock ( )
inlineoverride

Definition at line 70 of file PreconditionerLinearWithBlock.h.

71 {
72 }

Member Function Documentation

◆ create()

static PreconditionerSharedPtr Nektar::MultiRegions::PreconditionerLinearWithBlock::create ( const std::shared_ptr< GlobalLinSys > &  plinsys,
const std::shared_ptr< AssemblyMap > &  pLocToGloMap 
)
inlinestatic

Creates an instance of this class.

Definition at line 51 of file PreconditionerLinearWithBlock.h.

54 {
57 plinsys, pLocToGloMap);
58 p->InitObject();
59 return p;
60 }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::shared_ptr< Preconditioner > PreconditionerSharedPtr
Definition: GlobalLinSys.h:58

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), and CellMLToNektar.cellml_metadata::p.

◆ v_BuildPreconditioner()

void Nektar::MultiRegions::PreconditionerLinearWithBlock::v_BuildPreconditioner ( )
overrideprotectedvirtual

Reimplemented from Nektar::MultiRegions::Preconditioner.

Definition at line 81 of file PreconditionerLinearWithBlock.cpp.

82{
83 m_linSpacePrecon->BuildPreconditioner();
84 m_blockPrecon->BuildPreconditioner();
85}

References m_blockPrecon, and m_linSpacePrecon.

◆ v_DoPreconditioner()

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

Apply a preconditioner to the conjugate gradient method.

Reimplemented from Nektar::MultiRegions::Preconditioner.

Definition at line 90 of file PreconditionerLinearWithBlock.cpp.

93{
94 ASSERTL0(isLocal == false, "PreconditionerLinearWithBlock is only set up "
95 "for Global iteratives sovles");
96 int nGlobal = pInput.size();
97
98 Array<OneD, NekDouble> OutputBlock(nGlobal, 0.0);
99 Array<OneD, NekDouble> OutputLinear(nGlobal, 0.0);
100 Array<OneD, NekDouble> InputLinear(nGlobal, 0.0);
101
102 // Apply Low Energy preconditioner
103 m_blockPrecon->DoPreconditioner(pInput, OutputBlock);
104
105 // Apply linear space preconditioner
106 m_linSpacePrecon->DoPreconditionerWithNonVertOutput(pInput, pOutput,
107 OutputBlock);
108}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208

References ASSERTL0, m_blockPrecon, and m_linSpacePrecon.

◆ v_InitObject()

void Nektar::MultiRegions::PreconditionerLinearWithBlock::v_InitObject ( )
overrideprotectedvirtual

Member Data Documentation

◆ className

string Nektar::MultiRegions::PreconditionerLinearWithBlock::className
static
Initial value:
=
"FullLinearSpaceWithBlock", PreconditionerLinearWithBlock::create,
"Full Linear space and block preconditioning")
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:197
static PreconditionerSharedPtr create(const std::shared_ptr< GlobalLinSys > &plinsys, const std::shared_ptr< AssemblyMap > &pLocToGloMap)
Creates an instance of this class.

Name of class.

Registers the class with the Factory.

Definition at line 63 of file PreconditionerLinearWithBlock.h.

◆ m_blockPrecon

PreconditionerSharedPtr Nektar::MultiRegions::PreconditionerLinearWithBlock::m_blockPrecon
protected

◆ m_linSpacePrecon

PreconditionerSharedPtr Nektar::MultiRegions::PreconditionerLinearWithBlock::m_linSpacePrecon
protected