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

#include <PreconditionerLinearWithLowEnergy.h>

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

Public Member Functions

 PreconditionerLinearWithLowEnergy (const std::shared_ptr< GlobalLinSys > &plinsys, const AssemblyMapSharedPtr &pLocToGloMap)
 
 ~PreconditionerLinearWithLowEnergy () 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_DoTransformBasisToLowEnergy (Array< OneD, NekDouble > &pInOut) override
 Transform from original basis to low energy basis. More...
 
void v_DoTransformCoeffsFromLowEnergy (Array< OneD, NekDouble > &pInOut) override
 Transform from low energy coeffs to orignal basis. More...
 
void v_DoTransformCoeffsToLowEnergy (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput) override
 Multiply by the block inverse transformation matrix. More...
 
void v_DoTransformBasisFromLowEnergy (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput) override
 Multiply by the block transposed inverse transformation matrix. More...
 
DNekScalMatSharedPtr v_TransformedSchurCompl (int n, int offset, const std::shared_ptr< DNekScalMat > &loc_mat) override
 Get block elemental transposed transformation matrix \(\mathbf{R}^{T}\). More...
 
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_lowEnergyPrecon
 
Array< OneD, NekDoublem_invMultiplicity
 
- 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 PreconditionerLinearWithLowEnergy.h.

Constructor & Destructor Documentation

◆ PreconditionerLinearWithLowEnergy()

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

Definition at line 62 of file PreconditionerLinearWithLowEnergy.cpp.

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

◆ ~PreconditionerLinearWithLowEnergy()

Nektar::MultiRegions::PreconditionerLinearWithLowEnergy::~PreconditionerLinearWithLowEnergy ( )
inlineoverride

Definition at line 70 of file PreconditionerLinearWithLowEnergy.h.

71 {
72 }

Member Function Documentation

◆ create()

static PreconditionerSharedPtr Nektar::MultiRegions::PreconditionerLinearWithLowEnergy::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 PreconditionerLinearWithLowEnergy.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::PreconditionerLinearWithLowEnergy::v_BuildPreconditioner ( )
overrideprotectedvirtual

Reimplemented from Nektar::MultiRegions::Preconditioner.

Definition at line 142 of file PreconditionerLinearWithLowEnergy.cpp.

143{
144 m_linSpacePrecon->BuildPreconditioner();
145 m_lowEnergyPrecon->BuildPreconditioner();
146}

References m_linSpacePrecon, and m_lowEnergyPrecon.

◆ v_DoPreconditioner()

void Nektar::MultiRegions::PreconditionerLinearWithLowEnergy::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 151 of file PreconditionerLinearWithLowEnergy.cpp.

154{
155 ASSERTL0(isLocal == false, "PreconditionerLinearWithLowEnergy"
156 " is only set up for Global iteratives sovles");
157 int nDirBndDofs = m_locToGloMap.lock()->GetNumGlobalDirBndCoeffs();
158 int nGlobHomBndDofs =
159 m_locToGloMap.lock()->GetNumGlobalBndCoeffs() - nDirBndDofs;
160 int nLocBndDofs = m_locToGloMap.lock()->GetNumLocalBndCoeffs();
161
162 Array<OneD, NekDouble> tmp(nGlobHomBndDofs, 0.0);
163
164 Array<OneD, NekDouble> InputLinear(nGlobHomBndDofs);
165 Array<OneD, NekDouble> OutputLowEnergy(nGlobHomBndDofs);
166 Array<OneD, NekDouble> OutputLinear(nGlobHomBndDofs);
167 Array<OneD, NekDouble> local(nLocBndDofs);
168
169 // Transform input from low energy to original basis
170 Vmath::Vmul(nGlobHomBndDofs, m_invMultiplicity, 1, pInput, 1, OutputLinear,
171 1);
172 m_locToGloMap.lock()->GlobalToLocalBnd(OutputLinear, local, nDirBndDofs);
173 m_lowEnergyPrecon->DoTransformBasisFromLowEnergy(local, local);
174 m_locToGloMap.lock()->AssembleBnd(local, InputLinear, nDirBndDofs);
175
176 // Apply linear space preconditioner
177 m_linSpacePrecon->DoPreconditionerWithNonVertOutput(InputLinear,
178 OutputLinear, tmp);
179
180 // transform coefficients back to low energy space
181 m_locToGloMap.lock()->GlobalToLocalBnd(OutputLinear, local, nDirBndDofs);
182 m_lowEnergyPrecon->DoTransformCoeffsToLowEnergy(local, local);
183 m_locToGloMap.lock()->LocalBndToGlobal(local, pOutput, nDirBndDofs, false);
184
185 // Apply Low Energy preconditioner
186 m_lowEnergyPrecon->DoPreconditioner(pInput, OutputLowEnergy);
187
188 ASSERTL1(pOutput.size() >= nGlobHomBndDofs, "Output array is not correct");
189 Vmath::Vadd(nGlobHomBndDofs, pOutput, 1, OutputLowEnergy, 1, pOutput, 1);
190}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
std::weak_ptr< AssemblyMap > m_locToGloMap
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.hpp:72
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.hpp:180

References ASSERTL0, ASSERTL1, m_invMultiplicity, m_linSpacePrecon, Nektar::MultiRegions::Preconditioner::m_locToGloMap, m_lowEnergyPrecon, Vmath::Vadd(), and Vmath::Vmul().

◆ v_DoTransformBasisFromLowEnergy()

void Nektar::MultiRegions::PreconditionerLinearWithLowEnergy::v_DoTransformBasisFromLowEnergy ( const Array< OneD, NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput 
)
overrideprotectedvirtual

Multiply by the block transposed inverse transformation matrix.

Reimplemented from Nektar::MultiRegions::Preconditioner.

Definition at line 125 of file PreconditionerLinearWithLowEnergy.cpp.

127{
128 m_lowEnergyPrecon->DoTransformBasisFromLowEnergy(pInput, pOutput);
129}

References m_lowEnergyPrecon.

◆ v_DoTransformBasisToLowEnergy()

void Nektar::MultiRegions::PreconditionerLinearWithLowEnergy::v_DoTransformBasisToLowEnergy ( Array< OneD, NekDouble > &  pInOut)
overrideprotectedvirtual

Transform from original basis to low energy basis.

Reimplemented from Nektar::MultiRegions::Preconditioner.

Definition at line 98 of file PreconditionerLinearWithLowEnergy.cpp.

100{
101 m_lowEnergyPrecon->DoTransformBasisToLowEnergy(pInOut);
102}

References m_lowEnergyPrecon.

◆ v_DoTransformCoeffsFromLowEnergy()

void Nektar::MultiRegions::PreconditionerLinearWithLowEnergy::v_DoTransformCoeffsFromLowEnergy ( Array< OneD, NekDouble > &  pInOut)
overrideprotectedvirtual

Transform from low energy coeffs to orignal basis.

Reimplemented from Nektar::MultiRegions::Preconditioner.

Definition at line 107 of file PreconditionerLinearWithLowEnergy.cpp.

109{
110 m_lowEnergyPrecon->DoTransformCoeffsFromLowEnergy(pInput);
111}

References m_lowEnergyPrecon.

◆ v_DoTransformCoeffsToLowEnergy()

void Nektar::MultiRegions::PreconditionerLinearWithLowEnergy::v_DoTransformCoeffsToLowEnergy ( const Array< OneD, NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput 
)
overrideprotectedvirtual

Multiply by the block inverse transformation matrix.

Reimplemented from Nektar::MultiRegions::Preconditioner.

Definition at line 116 of file PreconditionerLinearWithLowEnergy.cpp.

118{
119 m_lowEnergyPrecon->DoTransformCoeffsToLowEnergy(pInput, pOutput);
120}

References m_lowEnergyPrecon.

◆ v_InitObject()

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

Reimplemented from Nektar::MultiRegions::Preconditioner.

Definition at line 72 of file PreconditionerLinearWithLowEnergy.cpp.

73{
75 "FullLinearSpace", m_linsys.lock(), m_locToGloMap.lock());
77 "LowEnergyBlock", m_linsys.lock(), m_locToGloMap.lock());
78
79 // Set up multiplicity array for inverse transposed transformation matrix
80 int nDirBnd = m_locToGloMap.lock()->GetNumGlobalDirBndCoeffs();
81 int nGlobHomBnd = m_locToGloMap.lock()->GetNumGlobalBndCoeffs() - nDirBnd;
82 int nLocBnd = m_locToGloMap.lock()->GetNumLocalBndCoeffs();
83
84 m_invMultiplicity = Array<OneD, NekDouble>(nGlobHomBnd, 1.0);
85 Array<OneD, NekDouble> loc(nLocBnd);
86
87 // need to scatter from global array to handle sign changes
88 m_locToGloMap.lock()->GlobalToLocalBnd(m_invMultiplicity, loc, nDirBnd);
89
90 // Now assemble values back together to get multiplicity
91 m_locToGloMap.lock()->AssembleBnd(loc, m_invMultiplicity, nDirBnd);
92 Vmath::Sdiv(nGlobHomBnd, 1.0, m_invMultiplicity, 1, m_invMultiplicity, 1);
93}
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:143
const std::weak_ptr< GlobalLinSys > m_linsys
PreconFactory & GetPreconFactory()
void Sdiv(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha/x.
Definition: Vmath.hpp:154

References Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), Nektar::MultiRegions::GetPreconFactory(), CG_Iterations::loc, m_invMultiplicity, m_linSpacePrecon, Nektar::MultiRegions::Preconditioner::m_linsys, Nektar::MultiRegions::Preconditioner::m_locToGloMap, m_lowEnergyPrecon, and Vmath::Sdiv().

◆ v_TransformedSchurCompl()

DNekScalMatSharedPtr Nektar::MultiRegions::PreconditionerLinearWithLowEnergy::v_TransformedSchurCompl ( int  offset,
int  bndoffset,
const std::shared_ptr< DNekScalMat > &  loc_mat 
)
overrideprotectedvirtual

Get block elemental transposed transformation matrix \(\mathbf{R}^{T}\).

Reimplemented from Nektar::MultiRegions::Preconditioner.

Definition at line 131 of file PreconditionerLinearWithLowEnergy.cpp.

133{
134 DNekScalMatSharedPtr returnval;
135 returnval = m_lowEnergyPrecon->TransformedSchurCompl(n, offset, loc_mat);
136 return returnval;
137}
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr

References m_lowEnergyPrecon.

Member Data Documentation

◆ className

string Nektar::MultiRegions::PreconditionerLinearWithLowEnergy::className
static
Initial value:
=
"FullLinearSpaceWithLowEnergyBlock",
"Full Linear space and low energy 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 PreconditionerLinearWithLowEnergy.h.

◆ m_invMultiplicity

Array<OneD, NekDouble> Nektar::MultiRegions::PreconditionerLinearWithLowEnergy::m_invMultiplicity
protected

Definition at line 78 of file PreconditionerLinearWithLowEnergy.h.

Referenced by v_DoPreconditioner(), and v_InitObject().

◆ m_linSpacePrecon

PreconditionerSharedPtr Nektar::MultiRegions::PreconditionerLinearWithLowEnergy::m_linSpacePrecon
protected

◆ m_lowEnergyPrecon

PreconditionerSharedPtr Nektar::MultiRegions::PreconditionerLinearWithLowEnergy::m_lowEnergyPrecon
protected