Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | Static Private Attributes | List of all members
Nektar::MultiRegions::Preconditioner Class Reference

#include <Preconditioner.h>

Inheritance diagram for Nektar::MultiRegions::Preconditioner:
Inheritance graph
[legend]
Collaboration diagram for Nektar::MultiRegions::Preconditioner:
Collaboration graph
[legend]

Public Member Functions

 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)
 

Protected Member Functions

virtual DNekScalMatSharedPtr v_TransformedSchurCompl (int offset, const boost::shared_ptr< DNekScalMat > &loc_mat)
 Get block elemental transposed transformation matrix $\mathbf{R}^{T}$. More...
 

Protected Attributes

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 NullPreconditioner (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_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_DoTransformToLowEnergy (Array< OneD, NekDouble > &pInOut, int offset)
 Transform from original basis to low energy basis. More...
 
virtual void v_DoTransformToLowEnergy (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
 Transform from original basis to low energy basis. More...
 
virtual void v_DoTransformFromLowEnergy (Array< OneD, NekDouble > &pInput)
 Transform from low energy basis to orignal basis. More...
 
virtual void v_DoMultiplybyInverseTransformationMatrix (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
 Multiply by the block inverse transformation matrix. More...
 
virtual void v_DoMultiplybyInverseTransposedTransformationMatrix (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
 Multiply by the block transposed inverse transformation matrix. More...
 
virtual void v_BuildPreconditioner ()
 

Static Private Attributes

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

Detailed Description

This class implements preconditioning for the conjugate gradient matrix solver.

This class implements diagonal preconditioning for the conjugate gradient matrix solver.

Definition at line 66 of file Preconditioner.h.

Constructor & Destructor Documentation

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

Definition at line 76 of file Preconditioner.cpp.

79  : m_linsys(plinsys),
80  m_preconType(pLocToGloMap->GetPreconType()),
81  m_locToGloMap(pLocToGloMap)
82  {
83  }
boost::shared_ptr< AssemblyMap > m_locToGloMap
const boost::weak_ptr< GlobalLinSys > m_linsys
virtual Nektar::MultiRegions::Preconditioner::~Preconditioner ( )
inlinevirtual

Definition at line 74 of file Preconditioner.h.

74 {}

Member Function Documentation

Array< OneD, NekDouble > Nektar::MultiRegions::Preconditioner::AssembleStaticCondGlobalDiagonals ( )

Performs global assembly of diagonal entries to global Schur complement matrix.

Definition at line 191 of file Preconditioner.cpp.

References m_linsys, and m_locToGloMap.

Referenced by Nektar::MultiRegions::PreconditionerDiagonal::StaticCondDiagonalPreconditionerSum().

192  {
193  int nGlobalBnd = m_locToGloMap->GetNumGlobalBndCoeffs();
194  int nDirBnd = m_locToGloMap->GetNumGlobalDirBndCoeffs();
195  int rows = nGlobalBnd - nDirBnd;
196 
197  DNekScalBlkMatSharedPtr loc_mat;
198  DNekScalMatSharedPtr bnd_mat;
199  int sign1, sign2, gid1, gid2, i, j, n, cnt;
200  Array<OneD, NekDouble> diagonals(rows,0.0);
201 
202  // Extract diagonal contributions of globally assembled
203  // schur complement matrix
204  for (cnt = n = 0; n < m_linsys.lock()->GetNumBlocks(); ++n)
205  {
206  // Get statically condensed local matrix.
207  loc_mat = (m_linsys.lock())->GetStaticCondBlock(n);
208 
209  // Extract boundary block.
210  bnd_mat = loc_mat->GetBlock(0, 0);
211 
212  // Offset by number of rows.
213  int bnd_row = bnd_mat->GetRows();
214 
215  for (i = 0; i < bnd_row; ++i)
216  {
217  gid1 = m_locToGloMap->GetLocalToGlobalBndMap (cnt + i)
218  - nDirBnd;
219  sign1 = m_locToGloMap->GetLocalToGlobalBndSign(cnt + i);
220 
221  if (gid1 < 0)
222  {
223  continue;
224  }
225 
226  for (j = 0; j < bnd_row; ++j)
227  {
228  gid2 = m_locToGloMap->GetLocalToGlobalBndMap (cnt + j)
229  - nDirBnd;
230  sign2 = m_locToGloMap->GetLocalToGlobalBndSign(cnt + j);
231 
232  if (gid2 == gid1)
233  {
234  diagonals[gid1] += sign1 * sign2 * (*bnd_mat)(i, j);
235  }
236  }
237  }
238  cnt += bnd_row;
239  }
240 
241  return diagonals;
242  }
boost::shared_ptr< AssemblyMap > m_locToGloMap
boost::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
const boost::weak_ptr< GlobalLinSys > m_linsys
boost::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
Definition: NekTypeDefs.hpp:74
void Nektar::MultiRegions::Preconditioner::BuildPreconditioner ( )
inline

Definition at line 280 of file Preconditioner.h.

References v_BuildPreconditioner().

281  {
283  }
void Nektar::MultiRegions::Preconditioner::DoMultiplybyInverseTransformationMatrix ( const Array< OneD, NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput 
)
inline

Definition at line 259 of file Preconditioner.h.

References v_DoMultiplybyInverseTransformationMatrix().

262  {
264  }
virtual void v_DoMultiplybyInverseTransformationMatrix(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
Multiply by the block inverse transformation matrix.
void Nektar::MultiRegions::Preconditioner::DoMultiplybyInverseTransposedTransformationMatrix ( const Array< OneD, NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput 
)
inline

Definition at line 270 of file Preconditioner.h.

References v_DoMultiplybyInverseTransposedTransformationMatrix().

273  {
275  }
virtual void v_DoMultiplybyInverseTransposedTransformationMatrix(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
Multiply by the block transposed inverse transformation matrix.
void Nektar::MultiRegions::Preconditioner::DoPreconditioner ( const Array< OneD, NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput 
)
inline

Definition at line 207 of file Preconditioner.h.

References v_DoPreconditioner().

210  {
211  v_DoPreconditioner(pInput,pOutput);
212  }
virtual void v_DoPreconditioner(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
Apply a preconditioner to the conjugate gradient method.
void Nektar::MultiRegions::Preconditioner::DoPreconditionerWithNonVertOutput ( const Array< OneD, NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput,
const Array< OneD, NekDouble > &  pNonVertOutput,
Array< OneD, NekDouble > &  pVertForce = NullNekDouble1DArray 
)
inline

Definition at line 218 of file Preconditioner.h.

References v_DoPreconditionerWithNonVertOutput().

223  {
224  v_DoPreconditionerWithNonVertOutput(pInput,pOutput,pNonVertOutput,
225  pVertForce);
226  }
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 free...
void Nektar::MultiRegions::Preconditioner::DoTransformFromLowEnergy ( Array< OneD, NekDouble > &  pInOut)
inline

Definition at line 250 of file Preconditioner.h.

References v_DoTransformFromLowEnergy().

252  {
254  }
virtual void v_DoTransformFromLowEnergy(Array< OneD, NekDouble > &pInput)
Transform from low energy basis to orignal basis.
void Nektar::MultiRegions::Preconditioner::DoTransformToLowEnergy ( Array< OneD, NekDouble > &  pInOut,
int  offset 
)
inline

Definition at line 231 of file Preconditioner.h.

References v_DoTransformToLowEnergy().

233  {
234  v_DoTransformToLowEnergy(pInOut,offset);
235  }
virtual void v_DoTransformToLowEnergy(Array< OneD, NekDouble > &pInOut, int offset)
Transform from original basis to low energy basis.
void Nektar::MultiRegions::Preconditioner::DoTransformToLowEnergy ( const Array< OneD, NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput 
)
inline

Definition at line 240 of file Preconditioner.h.

References v_DoTransformToLowEnergy().

243  {
244  v_DoTransformToLowEnergy(pInput,pOutput);
245  }
virtual void v_DoTransformToLowEnergy(Array< OneD, NekDouble > &pInOut, int offset)
Transform from original basis to low energy basis.
const DNekScalBlkMatSharedPtr& Nektar::MultiRegions::Preconditioner::GetBlockCMatrix ( ) const
inline
const DNekScalBlkMatSharedPtr& Nektar::MultiRegions::Preconditioner::GetBlockInvDMatrix ( ) const
inline
const DNekScalBlkMatSharedPtr& Nektar::MultiRegions::Preconditioner::GetBlockSchurCompl ( ) const
inline
const DNekScalBlkMatSharedPtr& Nektar::MultiRegions::Preconditioner::GetBlockTransformationMatrix ( ) const
inline
const DNekScalBlkMatSharedPtr& Nektar::MultiRegions::Preconditioner::GetBlockTransformedSchurCompl ( ) const
inline
const DNekScalBlkMatSharedPtr& Nektar::MultiRegions::Preconditioner::GetBlockTransposedTransformationMatrix ( ) const
inline
void Nektar::MultiRegions::Preconditioner::InitObject ( )
inline

Definition at line 190 of file Preconditioner.h.

References v_InitObject().

191  {
192  v_InitObject();
193  }
void Nektar::MultiRegions::Preconditioner::NullPreconditioner ( void  )
private
DNekScalMatSharedPtr Nektar::MultiRegions::Preconditioner::TransformedSchurCompl ( int  offset,
const boost::shared_ptr< DNekScalMat > &  loc_mat 
)
inline

Definition at line 198 of file Preconditioner.h.

References v_TransformedSchurCompl().

200  {
201  return v_TransformedSchurCompl(offset,loc_mat);
202  }
virtual DNekScalMatSharedPtr v_TransformedSchurCompl(int offset, const boost::shared_ptr< DNekScalMat > &loc_mat)
Get block elemental transposed transformation matrix .
void Nektar::MultiRegions::Preconditioner::v_BuildPreconditioner ( )
privatevirtual
void Nektar::MultiRegions::Preconditioner::v_DoMultiplybyInverseTransformationMatrix ( const Array< OneD, NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput 
)
privatevirtual

Multiply by the block inverse transformation matrix.

Reimplemented in Nektar::MultiRegions::PreconditionerLowEnergy.

Definition at line 155 of file Preconditioner.cpp.

References ErrorUtil::efatal, and NEKERROR.

Referenced by DoMultiplybyInverseTransformationMatrix().

158  {
159  NEKERROR(ErrorUtil::efatal,"Method does not exist" );
160  }
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mod...
Definition: ErrorUtil.hpp:158
void Nektar::MultiRegions::Preconditioner::v_DoMultiplybyInverseTransposedTransformationMatrix ( const Array< OneD, NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput 
)
privatevirtual

Multiply by the block transposed inverse transformation matrix.

Reimplemented in Nektar::MultiRegions::PreconditionerLowEnergy.

Definition at line 165 of file Preconditioner.cpp.

References ErrorUtil::efatal, and NEKERROR.

Referenced by DoMultiplybyInverseTransposedTransformationMatrix().

168  {
169  NEKERROR(ErrorUtil::efatal,"Method does not exist" );
170  }
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mod...
Definition: ErrorUtil.hpp:158
void Nektar::MultiRegions::Preconditioner::v_DoPreconditioner ( const Array< OneD, NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput 
)
privatevirtual

Apply a preconditioner to the conjugate gradient method.

Reimplemented in Nektar::MultiRegions::PreconditionerNull, Nektar::MultiRegions::PreconditionerLowEnergy, Nektar::MultiRegions::PreconditionerLinear, Nektar::MultiRegions::PreconditionerLinearWithLowEnergy, Nektar::MultiRegions::PreconditionerDiagonal, Nektar::MultiRegions::PreconditionerBlock, Nektar::MultiRegions::PreconditionerLinearWithDiag, and Nektar::MultiRegions::PreconditionerLinearWithBlock.

Definition at line 105 of file Preconditioner.cpp.

References ErrorUtil::efatal, and NEKERROR.

Referenced by DoPreconditioner().

108  {
109  NEKERROR(ErrorUtil::efatal,"Method does not exist" );
110  }
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mod...
Definition: ErrorUtil.hpp:158
void Nektar::MultiRegions::Preconditioner::v_DoPreconditionerWithNonVertOutput ( const Array< OneD, NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput,
const Array< OneD, NekDouble > &  pNonVertOutput,
Array< OneD, NekDouble > &  pVertForce 
)
privatevirtual

Apply a preconditioner to the conjugate gradient method with an output for non-vertex degrees of freedom.

Reimplemented in Nektar::MultiRegions::PreconditionerLinear.

Definition at line 116 of file Preconditioner.cpp.

References ErrorUtil::efatal, and NEKERROR.

Referenced by DoPreconditionerWithNonVertOutput().

121  {
122  NEKERROR(ErrorUtil::efatal, "Method does not exist");
123  }
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mod...
Definition: ErrorUtil.hpp:158
void Nektar::MultiRegions::Preconditioner::v_DoTransformFromLowEnergy ( Array< OneD, NekDouble > &  pInput)
privatevirtual

Transform from low energy basis to orignal basis.

Reimplemented in Nektar::MultiRegions::PreconditionerLowEnergy, and Nektar::MultiRegions::PreconditionerLinearWithLowEnergy.

Definition at line 146 of file Preconditioner.cpp.

References Vmath::Smul().

Referenced by DoTransformFromLowEnergy().

148  {
149  Vmath::Smul(pInput.num_elements(), 1.0, pInput, 1, pInput, 1);
150  }
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:199
void Nektar::MultiRegions::Preconditioner::v_DoTransformToLowEnergy ( Array< OneD, NekDouble > &  pInOut,
int  offset 
)
privatevirtual

Transform from original basis to low energy basis.

Reimplemented in Nektar::MultiRegions::PreconditionerLowEnergy, and Nektar::MultiRegions::PreconditionerLinearWithLowEnergy.

Definition at line 128 of file Preconditioner.cpp.

Referenced by DoTransformToLowEnergy().

131  {
132  }
void Nektar::MultiRegions::Preconditioner::v_DoTransformToLowEnergy ( const Array< OneD, NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput 
)
privatevirtual

Transform from original basis to low energy basis.

Reimplemented in Nektar::MultiRegions::PreconditionerLowEnergy.

Definition at line 137 of file Preconditioner.cpp.

140  {
141  }
void Nektar::MultiRegions::Preconditioner::v_InitObject ( )
privatevirtual
DNekScalMatSharedPtr Nektar::MultiRegions::Preconditioner::v_TransformedSchurCompl ( int  offset,
const boost::shared_ptr< DNekScalMat > &  loc_mat 
)
protectedvirtual

Get block elemental transposed transformation matrix $\mathbf{R}^{T}$.

Reimplemented in Nektar::MultiRegions::PreconditionerLowEnergy, and Nektar::MultiRegions::PreconditionerLinearWithLowEnergy.

Definition at line 180 of file Preconditioner.cpp.

Referenced by TransformedSchurCompl().

182  {
183  return loc_mat;
184  }

Member Data Documentation

std::string Nektar::MultiRegions::Preconditioner::def
staticprivate
Initial value:

Definition at line 183 of file Preconditioner.h.

std::string Nektar::MultiRegions::Preconditioner::lookupIds
staticprivate
Initial value:
= {
"Preconditioner", "Null", eNull),
"Preconditioner", "Diagonal", eDiagonal),
"Preconditioner", "FullLinearSpaceWithDiagonal",
"Preconditioner", "FullLinearSpace",eLinear),
"Preconditioner", "LowEnergyBlock",eLowEnergy),
"Preconditioner", "FullLinearSpaceWithLowEnergyBlock",
"Preconditioner", "Block",eBlock),
"Preconditioner", "FullLinearSpaceWithBlock",eLinearWithBlock),
}

Definition at line 182 of file Preconditioner.h.

LibUtilities::CommSharedPtr Nektar::MultiRegions::Preconditioner::m_comm
protected
const boost::weak_ptr<GlobalLinSys> Nektar::MultiRegions::Preconditioner::m_linsys
protected
boost::shared_ptr<AssemblyMap> Nektar::MultiRegions::Preconditioner::m_locToGloMap
protected
DNekMatSharedPtr Nektar::MultiRegions::Preconditioner::m_preconditioner
protected

Definition at line 136 of file Preconditioner.h.

PreconditionerType Nektar::MultiRegions::Preconditioner::m_preconType
protected

Definition at line 135 of file Preconditioner.h.