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)
 
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) override
 Apply a preconditioner to the conjugate gradient method. More...
 
virtual 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_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...
 

Protected Attributes

GlobalLinSysSharedPtr m_vertLinsys
 
std::shared_ptr< AssemblyMapm_vertLocToGloMap
 
- Protected Attributes inherited from Nektar::MultiRegions::Preconditioner
const std::weak_ptr< GlobalLinSysm_linsys
 
PreconditionerType 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 80 of file PreconditionerLinear.cpp.

83  : Preconditioner(plinsys, pLocToGloMap)
84 {
85 }
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.

◆ v_BuildPreconditioner()

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

Reimplemented from Nektar::MultiRegions::Preconditioner.

Definition at line 91 of file PreconditionerLinear.cpp.

92 {
93  GlobalSysSolnType sType = m_locToGloMap.lock()->GetGlobalSysSolnType();
94  ASSERTL0(sType == eIterativeStaticCond || sType == ePETScStaticCond,
95  "This type of preconditioning is not implemented "
96  "for this solver");
97 
98  std::shared_ptr<MultiRegions::ExpList> expList =
99  ((m_linsys.lock())->GetLocMat()).lock();
100 
102  expList->GetSession()->GetSolverInfoAsEnum<LinearPreconSolver>(
103  "LinearPreconSolver");
104 
105  GlobalSysSolnType linSolveType;
106 
107  switch (solveType)
108  {
109  case eLinearPreconPETSc:
110  {
111  linSolveType = ePETScFullMatrix;
112 #ifndef NEKTAR_USING_PETSC
113  NEKERROR(ErrorUtil::efatal, "Nektar++ has not been compiled with "
114  "PETSc support.");
115 #endif
116  break;
117  }
118  case eLinearPreconXxt:
119  default:
120  {
121  linSolveType = eXxtFullMatrix;
122  break;
123  }
124  }
125 
127  m_locToGloMap.lock()->LinearSpaceMap(*expList, linSolveType);
128 
129  // Generate linear solve system.
130  StdRegions::MatrixType mType =
131  m_linsys.lock()->GetKey().GetMatrixType() == StdRegions::eMass
134 
135  GlobalLinSysKey preconKey(mType, m_vertLocToGloMap,
136  (m_linsys.lock())->GetKey().GetConstFactors());
137 
138  switch (solveType)
139  {
140  case eLinearPreconXxt:
141  {
142  m_vertLinsys =
144  preconKey, expList, m_vertLocToGloMap);
145  break;
146  }
147  case eLinearPreconPETSc:
148  {
149 #ifdef NEKTAR_USING_PETSC
150  m_vertLinsys =
152  preconKey, expList, m_vertLocToGloMap);
153 #else
154  ASSERTL0(false, "Nektar++ has not been compiled with "
155  "PETSc support.");
156 #endif
157  }
158  }
159 }
#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
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, and solveType.

◆ v_DoPreconditioner()

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

Apply a preconditioner to the conjugate gradient method.

Reimplemented from Nektar::MultiRegions::Preconditioner.

Definition at line 164 of file PreconditionerLinear.cpp.

166 {
169 }
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 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 174 of file PreconditionerLinear.cpp.

178 {
179  GlobalSysSolnType solvertype = m_locToGloMap.lock()->GetGlobalSysSolnType();
180  switch (solvertype)
181  {
184  {
185  int i, val;
186  int nloc = m_vertLocToGloMap->GetNumLocalCoeffs();
187  int nglo = m_vertLocToGloMap->GetNumGlobalCoeffs();
188  // mapping from full space to vertices
189  Array<OneD, int> LocToGloBnd =
190  m_vertLocToGloMap->GetLocalToGlobalBndMap();
191 
192  // Global to local for linear solver (different from above)
193  Array<OneD, int> LocToGlo =
194  m_vertLocToGloMap->GetLocalToGlobalMap();
195 
196  // number of Dir coeffs in from full problem
197  int nDirFull = m_locToGloMap.lock()->GetNumGlobalDirBndCoeffs();
198 
199  Array<OneD, NekDouble> In(nglo, 0.0);
200  Array<OneD, NekDouble> Out(nglo, 0.0);
201 
202  // Gather rhs
203  for (i = 0; i < nloc; ++i)
204  {
205  val = LocToGloBnd[i];
206  if (val >= nDirFull)
207  {
208  In[LocToGlo[i]] = pInput[val - nDirFull];
209  }
210  }
211 
212  // Do solve without enforcing any boundary conditions.
213  m_vertLinsys->SolveLinearSystem(
214  m_vertLocToGloMap->GetNumGlobalCoeffs(), In, Out,
216  m_vertLocToGloMap->GetNumGlobalDirBndCoeffs());
217 
218  if (pNonVertOutput != NullNekDouble1DArray)
219  {
220  ASSERTL1(pNonVertOutput.size() >= pOutput.size(),
221  "Non Vert output is not of sufficient length");
222  Vmath::Vcopy(pOutput.size(), pNonVertOutput, 1, pOutput, 1);
223  }
224  else
225  {
226  // Copy input to output as a unit preconditioner on
227  // any other value
228  Vmath::Vcopy(pInput.size(), pInput, 1, pOutput, 1);
229  }
230 
231  if (pVertForce != NullNekDouble1DArray)
232  {
233  Vmath::Zero(pVertForce.size(), pVertForce, 1);
234  // Scatter back soln from linear solve
235  for (i = 0; i < nloc; ++i)
236  {
237  val = LocToGloBnd[i];
238  if (val >= nDirFull)
239  {
240  pOutput[val - nDirFull] = Out[LocToGlo[i]];
241  // copy vertex forcing into this vector
242  pVertForce[val - nDirFull] = In[LocToGlo[i]];
243  }
244  }
245  }
246  else
247  {
248  // Scatter back soln from linear solve
249  for (i = 0; i < nloc; ++i)
250  {
251  val = LocToGloBnd[i];
252  if (val >= nDirFull)
253  {
254  pOutput[val - nDirFull] = Out[LocToGlo[i]];
255  }
256  }
257  }
258  }
259  break;
260  default:
261  ASSERTL0(0, "Unsupported solver type");
262  break;
263  }
264 }
#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:492
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1255

References ASSERTL0, ASSERTL1, Nektar::MultiRegions::eIterativeStaticCond, Nektar::MultiRegions::ePETScStaticCond, 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 87 of file PreconditionerLinear.cpp.

88 {
89 }

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_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 100 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 101 of file PreconditionerLinear.h.