Nektar++
Public Member Functions | Static Public Member Functions | Static Public Attributes | 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 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)
 
virtual void v_InitObject ()
 
virtual void v_DoPreconditioner (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
 Apply a preconditioner to the conjugate gradient method. More...
 
virtual void v_BuildPreconditioner ()
 

Static Private Attributes

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

Additional Inherited Members

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

Detailed Description

Definition at line 49 of file PreconditionerDiagonal.h.

Constructor & Destructor Documentation

◆ PreconditionerDiagonal()

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

This class implements diagonal preconditioning for the conjugate

gradient matrix solver.

Definition at line 62 of file PreconditionerDiagonal.cpp.

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

◆ ~PreconditionerDiagonal()

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

Definition at line 71 of file PreconditionerDiagonal.h.

71 {}

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

56  {
58  p->InitObject();
59  return p;
60  }
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 98 of file PreconditionerDiagonal.cpp.

99  {
100  std::shared_ptr<MultiRegions::ExpList> expList =
101  ((m_linsys.lock())->GetLocMat()).lock();
102 
104 
105  auto asmMap = m_locToGloMap.lock();
106 
107  int i,j,n,cnt,gid1,gid2;
108  NekDouble sign1,sign2,value;
109  int nGlobal = asmMap->GetNumGlobalCoeffs();
110  int nDir = asmMap->GetNumGlobalDirBndCoeffs();
111  int nInt = nGlobal - nDir;
112 
113  // fill global matrix
114  DNekScalMatSharedPtr loc_mat;
115  Array<OneD, NekDouble> vOutput(nGlobal,0.0);
116 
117  int loc_row;
118  int nElmt = expList->GetNumElmts();
119  for(n = cnt = 0; n < nElmt; ++n)
120  {
121  loc_mat = (m_linsys.lock())->GetBlock(n);
122  loc_row = loc_mat->GetRows();
123 
124  for(i = 0; i < loc_row; ++i)
125  {
126  gid1 = asmMap->GetLocalToGlobalMap(cnt+i)-nDir;
127  sign1 = asmMap->GetLocalToGlobalSign(cnt+i);
128 
129  if(gid1 >= 0)
130  {
131  for(j = 0; j < loc_row; ++j)
132  {
133  gid2 = asmMap->GetLocalToGlobalMap(cnt+j)
134  - nDir;
135  sign2 = asmMap->GetLocalToGlobalSign(cnt+j);
136  if(gid2 == gid1)
137  {
138  // When global matrix is symmetric,
139  // only add the value for the upper
140  // triangular part in order to avoid
141  // entries to be entered twice
142  value = vOutput[gid1 + nDir]
143  + sign1*sign2*(*loc_mat)(i,j);
144  vOutput[gid1 + nDir] = value;
145  }
146  }
147  }
148  }
149  cnt += loc_row;
150  }
151 
152  // Assemble diagonal contributions across processes
153  asmMap->UniversalAssemble(vOutput);
154 
155  m_diagonals = Array<OneD, NekDouble> (nInt);
156  Vmath::Sdiv(nInt, 1.0, &vOutput[nDir], 1, &m_diagonals[0], 1);
157  }
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:291

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

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

Reimplemented from Nektar::MultiRegions::Preconditioner.

Definition at line 73 of file PreconditionerDiagonal.cpp.

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

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 
)
privatevirtual

Apply a preconditioner to the conjugate gradient method.

Reimplemented from Nektar::MultiRegions::Preconditioner.

Definition at line 191 of file PreconditionerDiagonal.cpp.

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

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

◆ v_InitObject()

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

Reimplemented from Nektar::MultiRegions::Preconditioner.

Definition at line 69 of file PreconditionerDiagonal.cpp.

70  {
71  }

Member Data Documentation

◆ className

string Nektar::MultiRegions::PreconditionerDiagonal::className
static
Initial value:
"Diagonal",
"Diagonal Preconditioning")
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:200
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 63 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