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)
 
virtual ~PreconditionerLinear ()
 
- 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

virtual void v_InitObject () override
 
virtual 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...
 
virtual 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...
 
virtual 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 56 of file PreconditionerLinear.h.

Constructor & Destructor Documentation

◆ PreconditionerLinear()

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

Definition at line 79 of file PreconditionerLinear.cpp.

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

◆ ~PreconditionerLinear()

virtual Nektar::MultiRegions::PreconditionerLinear::~PreconditionerLinear ( )
inlinevirtual

Definition at line 79 of file PreconditionerLinear.h.

80 {
81 }

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 60 of file PreconditionerLinear.h.

63 {
66 plinsys, pLocToGloMap);
67 p->InitObject();
68 return p;
69 }
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:60

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 318 of file PreconditionerLinear.cpp.

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

References m_invMult.

Referenced by v_BuildPreconditioner().

◆ v_BuildPreconditioner()

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

Reimplemented from Nektar::MultiRegions::Preconditioner.

Definition at line 90 of file PreconditionerLinear.cpp.

91{
92 GlobalSysSolnType sType = m_locToGloMap.lock()->GetGlobalSysSolnType();
94 "This type of preconditioning is not implemented "
95 "for this solver");
96
97 std::shared_ptr<MultiRegions::ExpList> expList =
98 ((m_linsys.lock())->GetLocMat()).lock();
99
101 expList->GetSession()->GetSolverInfoAsEnum<LinearPreconSolver>(
102 "LinearPreconSolver");
103
104 GlobalSysSolnType linSolveType;
105
106 switch (solveType)
107 {
109 {
110 linSolveType = ePETScFullMatrix;
111#ifndef NEKTAR_USING_PETSC
112 NEKERROR(ErrorUtil::efatal, "Nektar++ has not been compiled with "
113 "PETSc support.");
114#endif
115 break;
116 }
117 case eLinearPreconXxt:
118 default:
119 {
120 linSolveType = eXxtFullMatrix;
121 break;
122 }
123 }
124
126 m_locToGloMap.lock()->LinearSpaceMap(*expList, linSolveType);
127
129
130 // Generate linear solve system.
132 m_linsys.lock()->GetKey().GetMatrixType() == StdRegions::eMass
135
136 GlobalLinSysKey preconKey(mType, m_vertLocToGloMap,
137 (m_linsys.lock())->GetKey().GetConstFactors());
138
139 switch (solveType)
140 {
141 case eLinearPreconXxt:
142 {
145 preconKey, expList, m_vertLocToGloMap);
146 break;
147 }
149 {
150#ifdef NEKTAR_USING_PETSC
153 preconKey, expList, m_vertLocToGloMap);
154#else
155 ASSERTL0(false, "Nektar++ has not been compiled with "
156 "PETSc support.");
157#endif
158 }
159 }
160}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:209
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 165 of file PreconditionerLinear.cpp.

168{
169 ASSERTL0(isLocal == false, "PreconditionerLinear is only set up for Global"
170 " iteratives sovles");
173}
virtual 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 178 of file PreconditionerLinear.cpp.

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

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 86 of file PreconditionerLinear.cpp.

87{
88}

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:198
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 72 of file PreconditionerLinear.h.

◆ m_invMult

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

Definition at line 87 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 105 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 106 of file PreconditionerLinear.h.