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)
 
virtual ~PreconditionerDiagonal ()
 
- 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 className
 Name of class. More...
 

Protected Member Functions

virtual void v_InitObject () override
 
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_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...
 

Protected Attributes

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

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 48 of file PreconditionerDiagonal.h.

Constructor & Destructor Documentation

◆ PreconditionerDiagonal()

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

Definition at line 60 of file PreconditionerDiagonal.cpp.

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

◆ ~PreconditionerDiagonal()

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

Definition at line 72 of file PreconditionerDiagonal.h.

73  {
74  }

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 52 of file PreconditionerDiagonal.h.

55  {
58  plinsys, pLocToGloMap);
59  p->InitObject();
60  return p;
61  }
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.

◆ DiagonalPreconditionerSum()

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

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

Definition at line 95 of file PreconditionerDiagonal.cpp.

96 {
97  std::shared_ptr<MultiRegions::ExpList> expList =
98  ((m_linsys.lock())->GetLocMat()).lock();
99 
101 
102  auto asmMap = m_locToGloMap.lock();
103 
104  int i, j, n, cnt, gid1, gid2;
105  NekDouble sign1, sign2, value;
106  int nGlobal = asmMap->GetNumGlobalCoeffs();
107  int nDir = asmMap->GetNumGlobalDirBndCoeffs();
108  int nInt = nGlobal - nDir;
109 
110  // fill global matrix
111  DNekScalMatSharedPtr loc_mat;
112  Array<OneD, NekDouble> vOutput(nGlobal, 0.0);
113 
114  int loc_row;
115  int nElmt = expList->GetNumElmts();
116  for (n = cnt = 0; n < nElmt; ++n)
117  {
118  loc_mat = (m_linsys.lock())->GetBlock(n);
119  loc_row = loc_mat->GetRows();
120 
121  for (i = 0; i < loc_row; ++i)
122  {
123  gid1 = asmMap->GetLocalToGlobalMap(cnt + i) - nDir;
124  sign1 = asmMap->GetLocalToGlobalSign(cnt + i);
125 
126  if (gid1 >= 0)
127  {
128  for (j = 0; j < loc_row; ++j)
129  {
130  gid2 = asmMap->GetLocalToGlobalMap(cnt + j) - nDir;
131  sign2 = asmMap->GetLocalToGlobalSign(cnt + j);
132  if (gid2 == gid1)
133  {
134  // When global matrix is symmetric,
135  // only add the value for the upper
136  // triangular part in order to avoid
137  // entries to be entered twice
138  value = vOutput[gid1 + nDir] +
139  sign1 * sign2 * (*loc_mat)(i, j);
140  vOutput[gid1 + nDir] = value;
141  }
142  }
143  }
144  }
145  cnt += loc_row;
146  }
147 
148  // Assemble diagonal contributions across processes
149  asmMap->UniversalAssemble(vOutput);
150 
151  m_diagonals = Array<OneD, NekDouble>(nInt);
152  Vmath::Sdiv(nInt, 1.0, &vOutput[nDir], 1, &m_diagonals[0], 1);
153 }
const std::weak_ptr< GlobalLinSys > m_linsys
std::weak_ptr< AssemblyMap > m_locToGloMap
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:68
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
double NekDouble
void Sdiv(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha/y.
Definition: Vmath.cpp:324

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 160 of file PreconditionerDiagonal.cpp.

161 {
162  auto asmMap = m_locToGloMap.lock();
163 
164  int nGlobalBnd = asmMap->GetNumGlobalBndCoeffs();
165  int nDirBnd = asmMap->GetNumGlobalDirBndCoeffs();
166  int rows = nGlobalBnd - nDirBnd;
167 
168  Array<OneD, NekDouble> vOutput(nGlobalBnd, 0.0);
169 
170  // Extract diagonal contributions
171  Array<OneD, NekDouble> diagonals = AssembleStaticCondGlobalDiagonals();
172  for (unsigned int i = 0; i < rows; ++i)
173  {
174  vOutput[nDirBnd + i] = diagonals[i];
175  }
176 
177  // Assemble diagonal contributions across processes
178  asmMap->UniversalAssembleBnd(vOutput);
179 
180  m_diagonals = Array<OneD, NekDouble>(rows);
181  Vmath::Sdiv(rows, 1.0, &vOutput[nDirBnd], 1, &m_diagonals[0], 1);
182 }
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.

Definition at line 71 of file PreconditionerDiagonal.cpp.

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

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().

◆ v_DoPreconditioner()

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

189 {
190  auto asmMap = m_locToGloMap.lock();
191 
192  GlobalSysSolnType solvertype = asmMap->GetGlobalSysSolnType();
193 
194  int nGlobal = solvertype == eIterativeFull
195  ? asmMap->GetNumGlobalCoeffs()
196  : asmMap->GetNumGlobalBndCoeffs();
197  int nDir = asmMap->GetNumGlobalDirBndCoeffs();
198  int nNonDir = nGlobal - nDir;
199  Vmath::Vmul(nNonDir, &pInput[0], 1, &m_diagonals[0], 1, &pOutput[0], 1);
200 }
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.cpp:209

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

◆ v_InitObject()

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

Reimplemented from Nektar::MultiRegions::Preconditioner.

Definition at line 67 of file PreconditionerDiagonal.cpp.

68 {
69 }

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.
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 64 of file PreconditionerDiagonal.h.

◆ def

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

Definition at line 92 of file PreconditionerDiagonal.h.

◆ lookupIds

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

Definition at line 91 of file PreconditionerDiagonal.h.

◆ m_diagonals

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