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

#include <PreconditionerDiagonal.h>

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

Public Member Functions

 PreconditionerDiagonal (const std::shared_ptr< GlobalLinSys > &plinsys, const AssemblyMapSharedPtr &pLocToGloMap)
 
 ~PreconditionerDiagonal () 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

Array< OneD, NekDoublem_diagonals
 
- 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
 

Private Member Functions

void DiagonalPreconditionerSum (void)
 
void StaticCondDiagonalPreconditionerSum (void)
 

Static Private Attributes

static std::string lookupIds []
 
static std::string def
 

Detailed Description

Definition at line 46 of file PreconditionerDiagonal.h.

Constructor & Destructor Documentation

◆ PreconditionerDiagonal()

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

Definition at line 58 of file PreconditionerDiagonal.cpp.

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

◆ ~PreconditionerDiagonal()

Nektar::MultiRegions::PreconditionerDiagonal::~PreconditionerDiagonal ( )
inlineoverride

Definition at line 70 of file PreconditionerDiagonal.h.

71 {
72 }

Member Function Documentation

◆ create()

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

Creates an instance of this class.

Definition at line 50 of file PreconditionerDiagonal.h.

53 {
56 plinsys, pLocToGloMap);
57 p->InitObject();
58 return p;
59 }
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.

◆ DiagonalPreconditionerSum()

void Nektar::MultiRegions::PreconditionerDiagonal::DiagonalPreconditionerSum ( void  )
private

Diagonal preconditioner computed by summing the relevant elements of the local matrix system.

Definition at line 93 of file PreconditionerDiagonal.cpp.

94{
95 std::shared_ptr<MultiRegions::ExpList> expList =
96 ((m_linsys.lock())->GetLocMat()).lock();
97
98 auto asmMap = m_locToGloMap.lock();
99 int nGlobal = asmMap->GetNumGlobalCoeffs();
100 int nDir = asmMap->GetNumGlobalDirBndCoeffs();
101 int nInt = nGlobal - nDir;
102 int nElmt = expList->GetNumElmts();
103
104 Array<OneD, NekDouble> vOutput(nGlobal, 0.0);
105 for (int n = 0, cnt = 0; n < nElmt; ++n)
106 {
107 auto loc_mat = (m_linsys.lock())->GetBlock(n);
108 int loc_row = loc_mat->GetRows();
109 for (int i = 0; i < loc_row; ++i)
110 {
111 int gid1 = asmMap->GetLocalToGlobalMap(cnt + i);
112 vOutput[gid1] += (*loc_mat)(i, i);
113 }
114 cnt += loc_row;
115 }
116
117 // Assemble diagonal contributions across processes
118 asmMap->UniversalAssemble(vOutput);
119
120 m_diagonals = Array<OneD, NekDouble>(nInt);
121 Vmath::Sdiv(nInt, 1.0, &vOutput[nDir], 1, &m_diagonals[0], 1);
122}
const std::weak_ptr< GlobalLinSys > m_linsys
std::weak_ptr< AssemblyMap > m_locToGloMap
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 m_diagonals, Nektar::MultiRegions::Preconditioner::m_linsys, Nektar::MultiRegions::Preconditioner::m_locToGloMap, and Vmath::Sdiv().

Referenced by v_BuildPreconditioner().

◆ StaticCondDiagonalPreconditionerSum()

void Nektar::MultiRegions::PreconditionerDiagonal::StaticCondDiagonalPreconditionerSum ( void  )
private

Diagonal preconditioner defined as the inverse of the main diagonal of the Schur complement

Definition at line 129 of file PreconditionerDiagonal.cpp.

130{
131 auto asmMap = m_locToGloMap.lock();
132
133 int nGlobalBnd = asmMap->GetNumGlobalBndCoeffs();
134 int nDirBnd = asmMap->GetNumGlobalDirBndCoeffs();
135 int rows = nGlobalBnd - nDirBnd;
136
137 Array<OneD, NekDouble> vOutput(nGlobalBnd, 0.0);
138
139 // Extract diagonal contributions
140 Array<OneD, NekDouble> diagonals = AssembleStaticCondGlobalDiagonals();
141 for (unsigned int i = 0; i < rows; ++i)
142 {
143 vOutput[nDirBnd + i] = diagonals[i];
144 }
145
146 // Assemble diagonal contributions across processes
147 asmMap->UniversalAssembleBnd(vOutput);
148
149 m_diagonals = Array<OneD, NekDouble>(rows);
150 Vmath::Sdiv(rows, 1.0, &vOutput[nDirBnd], 1, &m_diagonals[0], 1);
151}
Array< OneD, NekDouble > AssembleStaticCondGlobalDiagonals()
Performs global assembly of diagonal entries to global Schur complement matrix.

References Nektar::MultiRegions::Preconditioner::AssembleStaticCondGlobalDiagonals(), m_diagonals, Nektar::MultiRegions::Preconditioner::m_locToGloMap, and Vmath::Sdiv().

Referenced by v_BuildPreconditioner().

◆ v_BuildPreconditioner()

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

Reimplemented from Nektar::MultiRegions::Preconditioner.

Reimplemented in Nektar::MultiRegions::PreconditionerJacobi.

Definition at line 69 of file PreconditionerDiagonal.cpp.

70{
71 GlobalSysSolnType solvertype = m_locToGloMap.lock()->GetGlobalSysSolnType();
72 if (solvertype == eIterativeFull)
73 {
75 }
76 else if (solvertype == eIterativeStaticCond ||
77 solvertype == eIterativeMultiLevelStaticCond ||
78 solvertype == ePETScStaticCond ||
79 solvertype == ePETScMultiLevelStaticCond)
80 {
82 }
83 else
84 {
85 ASSERTL0(0, "Unsupported solver type");
86 }
87}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208

References ASSERTL0, DiagonalPreconditionerSum(), Nektar::MultiRegions::eIterativeFull, Nektar::MultiRegions::eIterativeMultiLevelStaticCond, Nektar::MultiRegions::eIterativeStaticCond, Nektar::MultiRegions::ePETScMultiLevelStaticCond, Nektar::MultiRegions::ePETScStaticCond, Nektar::MultiRegions::Preconditioner::m_locToGloMap, and StaticCondDiagonalPreconditionerSum().

Referenced by Nektar::MultiRegions::PreconditionerJacobi::v_BuildPreconditioner().

◆ v_DoPreconditioner()

void Nektar::MultiRegions::PreconditionerDiagonal::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.

Reimplemented in Nektar::MultiRegions::PreconditionerJacobi.

Definition at line 156 of file PreconditionerDiagonal.cpp.

159{
160 auto asmMap = m_locToGloMap.lock();
161
162 GlobalSysSolnType solvertype = asmMap->GetGlobalSysSolnType();
163
164 bool isFull = solvertype == eIterativeFull ? true : false;
165
166 int nGlobal = (isFull) ? asmMap->GetNumGlobalCoeffs()
167 : asmMap->GetNumGlobalBndCoeffs();
168 int nDir = asmMap->GetNumGlobalDirBndCoeffs();
169 int nNonDir = nGlobal - nDir;
170
171 if (IsLocal)
172 {
173 Array<OneD, NekDouble> wk(nGlobal);
174 (isFull) ? asmMap->Assemble(pInput, wk)
175 : asmMap->AssembleBnd(pInput, wk);
176 Vmath::Vmul(nNonDir, wk.data() + nDir, 1, m_diagonals.data(), 1,
177 wk.data() + nDir, 1);
178 Vmath::Zero(nDir, wk, 1);
179 (isFull) ? asmMap->GlobalToLocal(wk, pOutput)
180 : asmMap->GlobalToLocalBnd(wk, pOutput);
181 }
182 else
183 {
184 Vmath::Vmul(nNonDir, &pInput[0], 1, &m_diagonals[0], 1, &pOutput[0], 1);
185 }
186}
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 Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273

References Nektar::MultiRegions::eIterativeFull, m_diagonals, Nektar::MultiRegions::Preconditioner::m_locToGloMap, Vmath::Vmul(), and Vmath::Zero().

◆ v_InitObject()

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

Reimplemented from Nektar::MultiRegions::Preconditioner.

Reimplemented in Nektar::MultiRegions::PreconditionerJacobi.

Definition at line 65 of file PreconditionerDiagonal.cpp.

66{
67}

Member Data Documentation

◆ className

string Nektar::MultiRegions::PreconditionerDiagonal::className
static
Initial value:
=
"Diagonal", PreconditionerDiagonal::create, "Diagonal Preconditioning")
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
static PreconditionerSharedPtr create(const std::shared_ptr< GlobalLinSys > &plinsys, const std::shared_ptr< AssemblyMap > &pLocToGloMap)
Creates an instance of this class.
PreconFactory & GetPreconFactory()

Name of class.

Registers the class with the Factory.

Definition at line 62 of file PreconditionerDiagonal.h.

◆ def

std::string Nektar::MultiRegions::PreconditionerDiagonal::def
staticprivate

Definition at line 91 of file PreconditionerDiagonal.h.

◆ lookupIds

std::string Nektar::MultiRegions::PreconditionerDiagonal::lookupIds[]
staticprivate

Definition at line 90 of file PreconditionerDiagonal.h.

◆ m_diagonals

Array<OneD, NekDouble> Nektar::MultiRegions::PreconditionerDiagonal::m_diagonals
protected