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

#include <PreconditionerLinear.h>

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

Public Member Functions

 PreconditionerLinear (const std::shared_ptr< GlobalLinSys > &plinsys, const AssemblyMapSharedPtr &pLocToGloMap)
 
 ~PreconditionerLinear () 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 className1
 Name of class. More...
 

Protected Member Functions

void v_InitObject () override
 
void v_DoPreconditionerWithNonVertOutput (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const Array< OneD, NekDouble > &pNonVertOutput, Array< OneD, NekDouble > &pVertForce) override
 Apply a preconditioner to the conjugate gradient method with an output for non-vertex degrees of freedom. 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
 
void SetupInvMult (const std::shared_ptr< AssemblyMap > &pLocToGloMap)
 
- 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

GlobalLinSysSharedPtr m_vertLinsys
 
std::shared_ptr< AssemblyMapm_vertLocToGloMap
 
Array< OneD, NekDoublem_invMult
 
- 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
 

Static Private Attributes

static std::string solveType
 
static std::string solveTypeIds []
 

Detailed Description

This class implements preconditioning for the conjugate gradient matrix solver.

Definition at line 54 of file PreconditionerLinear.h.

Constructor & Destructor Documentation

◆ PreconditionerLinear()

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

Definition at line 77 of file PreconditionerLinear.cpp.

80 : Preconditioner(plinsys, pLocToGloMap)
81{
82}
Preconditioner(const std::shared_ptr< GlobalLinSys > &plinsys, const AssemblyMapSharedPtr &pLocToGloMap)

◆ ~PreconditionerLinear()

Nektar::MultiRegions::PreconditionerLinear::~PreconditionerLinear ( )
inlineoverride

Definition at line 77 of file PreconditionerLinear.h.

78 {
79 }

Member Function Documentation

◆ create()

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

Creates an instance of this class.

Definition at line 58 of file PreconditionerLinear.h.

61 {
64 plinsys, pLocToGloMap);
65 p->InitObject();
66 return p;
67 }
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.

◆ SetupInvMult()

void Nektar::MultiRegions::PreconditionerLinear::SetupInvMult ( const std::shared_ptr< AssemblyMap > &  pLocToGloMap)
protected

Create the inverse multiplicity map.

Parameters
locToGloMapLocal to global mapping information.

Definition at line 317 of file PreconditionerLinear.cpp.

319{
320 const Array<OneD, const int> &vMap = pLocToGloMap->GetLocalToGlobalMap();
321 unsigned int nGlo = pLocToGloMap->GetNumGlobalCoeffs();
322 unsigned int nEntries = pLocToGloMap->GetNumLocalCoeffs();
323 unsigned int i;
324
325 // Count the multiplicity of each global DOF on this process
326 Array<OneD, NekDouble> vCounts(nGlo, 0.0);
327 for (i = 0; i < nEntries; ++i)
328 {
329 vCounts[vMap[i]] += 1.0;
330 }
331
332 // Get universal multiplicity by globally assembling counts
333 pLocToGloMap->UniversalAssemble(vCounts);
334
335 // Construct a map of 1/multiplicity for use in XXT solve
336 m_invMult = Array<OneD, NekDouble>(nEntries);
337 for (i = 0; i < nEntries; ++i)
338 {
339 m_invMult[i] = 1.0 / vCounts[vMap[i]];
340 }
341}

References m_invMult.

Referenced by v_BuildPreconditioner().

◆ v_BuildPreconditioner()

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

Reimplemented from Nektar::MultiRegions::Preconditioner.

Definition at line 88 of file PreconditionerLinear.cpp.

89{
90 GlobalSysSolnType sType = m_locToGloMap.lock()->GetGlobalSysSolnType();
92 "This type of preconditioning is not implemented "
93 "for this solver");
94
95 std::shared_ptr<MultiRegions::ExpList> expList =
96 ((m_linsys.lock())->GetLocMat()).lock();
97
99 expList->GetSession()->GetSolverInfoAsEnum<LinearPreconSolver>(
100 "LinearPreconSolver");
101
102 GlobalSysSolnType linSolveType;
103
104 switch (solveType)
105 {
107 {
108 linSolveType = ePETScFullMatrix;
109#ifndef NEKTAR_USING_PETSC
110 NEKERROR(ErrorUtil::efatal, "Nektar++ has not been compiled with "
111 "PETSc support.");
112#endif
113 break;
114 }
115 case eLinearPreconXxt:
116 default:
117 {
118 linSolveType = eXxtFullMatrix;
119 break;
120 }
121 }
122
124 m_locToGloMap.lock()->LinearSpaceMap(*expList, linSolveType);
125
127
128 // Generate linear solve system.
130 m_linsys.lock()->GetKey().GetMatrixType() == StdRegions::eMass
133
134 GlobalLinSysKey preconKey(mType, m_vertLocToGloMap,
135 (m_linsys.lock())->GetKey().GetConstFactors());
136
137 switch (solveType)
138 {
139 case eLinearPreconXxt:
140 {
143 preconKey, expList, m_vertLocToGloMap);
144 break;
145 }
147 {
148#ifdef NEKTAR_USING_PETSC
151 preconKey, expList, m_vertLocToGloMap);
152#else
153 ASSERTL0(false, "Nektar++ has not been compiled with "
154 "PETSc support.");
155#endif
156 }
157 }
158}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:202
const std::weak_ptr< GlobalLinSys > m_linsys
std::weak_ptr< AssemblyMap > m_locToGloMap
void SetupInvMult(const std::shared_ptr< AssemblyMap > &pLocToGloMap)
std::shared_ptr< AssemblyMap > m_vertLocToGloMap

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, Nektar::ErrorUtil::efatal, Nektar::MultiRegions::eIterativeStaticCond, Nektar::MultiRegions::eLinearPreconPETSc, Nektar::MultiRegions::eLinearPreconXxt, Nektar::StdRegions::eMass, Nektar::MultiRegions::ePETScFullMatrix, Nektar::MultiRegions::ePETScStaticCond, Nektar::StdRegions::ePreconLinearSpace, Nektar::StdRegions::ePreconLinearSpaceMass, Nektar::MultiRegions::eXxtFullMatrix, Nektar::MultiRegions::Preconditioner::m_linsys, Nektar::MultiRegions::Preconditioner::m_locToGloMap, m_vertLinsys, m_vertLocToGloMap, NEKERROR, SetupInvMult(), and solveType.

◆ v_DoPreconditioner()

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

166{
167 ASSERTL0(isLocal == false, "PreconditionerLinear is only set up for Global"
168 " iteratives sovles");
171}
void v_DoPreconditionerWithNonVertOutput(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const Array< OneD, NekDouble > &pNonVertOutput, Array< OneD, NekDouble > &pVertForce) override
Apply a preconditioner to the conjugate gradient method with an output for non-vertex degrees of free...
static Array< OneD, NekDouble > NullNekDouble1DArray

References ASSERTL0, Nektar::NullNekDouble1DArray, and v_DoPreconditionerWithNonVertOutput().

◆ v_DoPreconditionerWithNonVertOutput()

void Nektar::MultiRegions::PreconditionerLinear::v_DoPreconditionerWithNonVertOutput ( const Array< OneD, NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput,
const Array< OneD, NekDouble > &  pNonVertOutput,
Array< OneD, NekDouble > &  pVertForce 
)
overrideprotectedvirtual

Apply a preconditioner to the conjugate gradient method with an output for non-vertex degrees of freedom.

Reimplemented from Nektar::MultiRegions::Preconditioner.

Definition at line 176 of file PreconditionerLinear.cpp.

180{
181 GlobalSysSolnType solvertype = m_locToGloMap.lock()->GetGlobalSysSolnType();
182 switch (solvertype)
183 {
186 {
187 int i, val;
188 int nloc = m_vertLocToGloMap->GetNumLocalCoeffs();
189 // mapping from full space to vertices
190 Array<OneD, int> LocToGloBnd =
191 m_vertLocToGloMap->GetLocalToGlobalBndMap();
192
193 // Global to local for linear solver (different from above)
194 Array<OneD, int> LocToGlo =
195 m_vertLocToGloMap->GetLocalToGlobalMap();
196
197 // number of Dir coeffs in from full problem
198 int nDirFull = m_locToGloMap.lock()->GetNumGlobalDirBndCoeffs();
199
200#if 0
201 int nglo = m_vertLocToGloMap->GetNumGlobalCoeffs();
202 Array<OneD, NekDouble> In(nglo, 0.0);
203 Array<OneD, NekDouble> Out(nglo, 0.0);
204
205 // Gather rhs to local linear space.
206 for (i = 0; i < nloc; ++i)
207 {
208 val = LocToGloBnd[i];
209 if (val >= nDirFull)
210 {
211 In[LocToGlo[i]] = pInput[val - nDirFull];
212 }
213 }
214#else
215 Array<OneD, NekDouble> In(nloc, 0.0);
216 Array<OneD, NekDouble> Out(nloc, 0.0);
217
218 // Gather rhs to local linear space.
219 for (i = 0; i < nloc; ++i)
220 {
221 val = LocToGloBnd[i];
222 if (val >= nDirFull)
223 {
224 In[i] = pInput[val - nDirFull] * m_invMult[i];
225 }
226 }
227#endif
228
229 // Do solve without enforcing any boundary conditions.
230 m_vertLinsys->SolveLinearSystem(
231 m_vertLocToGloMap->GetNumGlobalCoeffs(), In, Out,
233 m_vertLocToGloMap->GetNumGlobalDirBndCoeffs());
234
235 if (pNonVertOutput != NullNekDouble1DArray)
236 {
237 ASSERTL1(pNonVertOutput.size() >= pOutput.size(),
238 "Non Vert output is not of sufficient length");
239 Vmath::Vcopy(pOutput.size(), pNonVertOutput, 1, pOutput, 1);
240 }
241 else
242 {
243 // Copy input to output as a unit preconditioner on
244 // any other value
245 Vmath::Vcopy(pInput.size(), pInput, 1, pOutput, 1);
246 }
247
248#if 0
249 if (pVertForce != NullNekDouble1DArray)
250 {
251 Vmath::Zero(pVertForce.size(), pVertForce, 1);
252 // Scatter back soln from linear solve
253 for (i = 0; i < nloc; ++i)
254 {
255 val = LocToGloBnd[i];
256 if (val >= nDirFull)
257 {
258 pOutput[val - nDirFull] = Out[LocToGlo[i]];
259 // copy vertex forcing into this vector
260 pVertForce[val - nDirFull] = In[LocToGlo[i]];
261 }
262 }
263 }
264 else
265 {
266 // Scatter back soln from linear solve
267 for (i = 0; i < nloc; ++i)
268 {
269 val = LocToGloBnd[i];
270 if (val >= nDirFull)
271 {
272 pOutput[val - nDirFull] = Out[LocToGlo[i]];
273 }
274 }
275 }
276#else
277 if (pVertForce != NullNekDouble1DArray)
278 {
279 Vmath::Zero(pVertForce.size(), pVertForce, 1);
280 // Scatter back soln from linear solve
281 for (i = 0; i < nloc; ++i)
282 {
283 val = LocToGloBnd[i];
284 if (val >= nDirFull)
285 {
286 pOutput[val - nDirFull] = Out[i];
287 // copy vertex forcing into this vector
288 pVertForce[val - nDirFull] = In[i];
289 }
290 }
291 }
292 else
293 {
294 // Scatter back soln from linear solve
295 for (i = 0; i < nloc; ++i)
296 {
297 val = LocToGloBnd[i];
298 if (val >= nDirFull)
299 {
300 pOutput[val - nDirFull] = Out[i];
301 }
302 }
303 }
304#endif
305 }
306 break;
307 default:
308 ASSERTL0(0, "Unsupported solver type");
309 break;
310 }
311}
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825

References ASSERTL0, ASSERTL1, Nektar::MultiRegions::eIterativeStaticCond, Nektar::MultiRegions::ePETScStaticCond, m_invMult, Nektar::MultiRegions::Preconditioner::m_locToGloMap, m_vertLinsys, m_vertLocToGloMap, Nektar::NullNekDouble1DArray, Vmath::Vcopy(), and Vmath::Zero().

Referenced by v_DoPreconditioner().

◆ v_InitObject()

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

Reimplemented from Nektar::MultiRegions::Preconditioner.

Definition at line 84 of file PreconditionerLinear.cpp.

85{
86}

Member Data Documentation

◆ className1

string Nektar::MultiRegions::PreconditionerLinear::className1
static
Initial value:
=
"FullLinearSpace", PreconditionerLinear::create,
"Full Linear space inverse 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.
PreconFactory & GetPreconFactory()

Name of class.

Registers the class with the Factory.

Definition at line 70 of file PreconditionerLinear.h.

◆ m_invMult

Array<OneD, NekDouble> Nektar::MultiRegions::PreconditionerLinear::m_invMult
protected

Definition at line 85 of file PreconditionerLinear.h.

Referenced by SetupInvMult(), and v_DoPreconditionerWithNonVertOutput().

◆ m_vertLinsys

GlobalLinSysSharedPtr Nektar::MultiRegions::PreconditionerLinear::m_vertLinsys
protected

◆ m_vertLocToGloMap

std::shared_ptr<AssemblyMap> Nektar::MultiRegions::PreconditionerLinear::m_vertLocToGloMap
protected

◆ solveType

std::string Nektar::MultiRegions::PreconditionerLinear::solveType
staticprivate
Initial value:
=
"Xxt")
static std::string RegisterDefaultSolverInfo(const std::string &pName, const std::string &pValue)
Registers the default string value of a solver info property.

Definition at line 103 of file PreconditionerLinear.h.

Referenced by v_BuildPreconditioner().

◆ solveTypeIds

std::string Nektar::MultiRegions::PreconditionerLinear::solveTypeIds
staticprivate
Initial value:
= {
"LinearPreconSolver", "PETSc", MultiRegions::eLinearPreconPETSc),
"LinearPreconSolver", "Xxt", MultiRegions::eLinearPreconXxt)}
static std::string RegisterEnumValue(std::string pEnum, std::string pString, int pEnumValue)
Registers an enumeration value.

Definition at line 104 of file PreconditionerLinear.h.