Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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:
Inheritance graph
[legend]
Collaboration diagram for Nektar::MultiRegions::PreconditionerDiagonal:
Collaboration graph
[legend]

Public Member Functions

 PreconditionerDiagonal (const boost::shared_ptr< GlobalLinSys > &plinsys, const AssemblyMapSharedPtr &pLocToGloMap)
 
virtual ~PreconditionerDiagonal ()
 
- Public Member Functions inherited from Nektar::MultiRegions::Preconditioner
 Preconditioner (const boost::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 DoTransformToLowEnergy (Array< OneD, NekDouble > &pInOut, int offset)
 
void DoTransformToLowEnergy (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
 
void DoTransformFromLowEnergy (Array< OneD, NekDouble > &pInOut)
 
void DoMultiplybyInverseTransformationMatrix (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
 
void DoMultiplybyInverseTransposedTransformationMatrix (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, const boost::shared_ptr< DNekScalMat > &loc_mat)
 

Static Public Member Functions

static PreconditionerSharedPtr create (const boost::shared_ptr< GlobalLinSys > &plinsys, const boost::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
 
PreconditionerType m_preconType
 
- Protected Attributes inherited from Nektar::MultiRegions::Preconditioner
const boost::weak_ptr
< GlobalLinSys
m_linsys
 
PreconditionerType m_preconType
 
DNekMatSharedPtr m_preconditioner
 
boost::shared_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, const boost::shared_ptr< DNekScalMat > &loc_mat)
 Get block elemental transposed transformation matrix $\mathbf{R}^{T}$. More...
 

Detailed Description

Definition at line 50 of file PreconditionerDiagonal.h.

Constructor & Destructor Documentation

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

Definition at line 61 of file PreconditionerDiagonal.cpp.

64  : Preconditioner(plinsys, pLocToGloMap),
65  m_preconType(pLocToGloMap->GetPreconType())
66  {
67  }
Preconditioner(const boost::shared_ptr< GlobalLinSys > &plinsys, const AssemblyMapSharedPtr &pLocToGloMap)
virtual Nektar::MultiRegions::PreconditionerDiagonal::~PreconditionerDiagonal ( )
inlinevirtual

Definition at line 73 of file PreconditionerDiagonal.h.

73 {}

Member Function Documentation

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

Creates an instance of this class.

Definition at line 54 of file PreconditionerDiagonal.h.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr().

58  {
60  p->InitObject();
61  return p;
62  }
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
boost::shared_ptr< Preconditioner > PreconditionerSharedPtr
Definition: GlobalLinSys.h:62
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.

References m_diagonals, Nektar::MultiRegions::Preconditioner::m_linsys, Nektar::MultiRegions::Preconditioner::m_locToGloMap, and Vmath::Sdiv().

Referenced by v_BuildPreconditioner().

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

Diagonal preconditioner defined as the inverse of the main diagonal of the Schur complement

Definition at line 161 of file PreconditionerDiagonal.cpp.

References Nektar::MultiRegions::Preconditioner::AssembleStaticCondGlobalDiagonals(), m_diagonals, Nektar::MultiRegions::Preconditioner::m_locToGloMap, and Vmath::Sdiv().

Referenced by v_BuildPreconditioner().

162  {
163  int nGlobalBnd = m_locToGloMap->GetNumGlobalBndCoeffs();
164  int nDirBnd = m_locToGloMap->GetNumGlobalDirBndCoeffs();
165  int rows = nGlobalBnd - nDirBnd;
166 
167  Array<OneD, NekDouble> vOutput(nGlobalBnd,0.0);
168 
169  // Extract diagonal contributions
171  for (unsigned int i = 0; i < rows; ++i)
172  {
173  vOutput[nDirBnd + i] = diagonals[i];
174  }
175 
176  // Assemble diagonal contributions across processes
177  m_locToGloMap->UniversalAssembleBnd(vOutput);
178 
180  Vmath::Sdiv(rows, 1.0, &vOutput[nDirBnd], 1, &m_diagonals[0], 1);
181  }
boost::shared_ptr< AssemblyMap > m_locToGloMap
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:257
Array< OneD, NekDouble > AssembleStaticCondGlobalDiagonals()
Performs global assembly of diagonal entries to global Schur complement matrix.
void Nektar::MultiRegions::PreconditionerDiagonal::v_BuildPreconditioner ( )
privatevirtual

Reimplemented from Nektar::MultiRegions::Preconditioner.

Definition at line 73 of file PreconditionerDiagonal.cpp.

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

74  {
75  GlobalSysSolnType solvertype =
76  m_locToGloMap->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  }
boost::shared_ptr< AssemblyMap > m_locToGloMap
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
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 186 of file PreconditionerDiagonal.cpp.

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

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

string Nektar::MultiRegions::PreconditionerDiagonal::className
static
Initial value:
"Diagonal",
"Diagonal Preconditioning")

Name of class.

Registers the class with the Factory.

Definition at line 65 of file PreconditionerDiagonal.h.

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

Definition at line 96 of file PreconditionerDiagonal.h.

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

Definition at line 95 of file PreconditionerDiagonal.h.

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

Definition at line 79 of file PreconditionerDiagonal.h.