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

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

: m_linsys(plinsys),
m_preconType(pLocToGloMap->GetPreconType()),
m_locToGloMap(pLocToGloMap)
{
}
virtual Nektar::MultiRegions::Preconditioner::~Preconditioner ( )
inlinevirtual

Definition at line 74 of file Preconditioner.h.

{}

Member Function Documentation

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

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

Definition at line 190 of file Preconditioner.cpp.

References m_linsys, and m_locToGloMap.

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

{
int nGlobalBnd = m_locToGloMap->GetNumGlobalBndCoeffs();
int nDirBnd = m_locToGloMap->GetNumGlobalDirBndCoeffs();
int rows = nGlobalBnd - nDirBnd;
int sign1, sign2, gid1, gid2, i, j, n, cnt;
Array<OneD, NekDouble> diagonals(rows,0.0);
// Extract diagonal contributions of globally assembled
// schur complement matrix
for (cnt = n = 0; n < m_linsys.lock()->GetNumBlocks(); ++n)
{
// Get statically condensed local matrix.
loc_mat = (m_linsys.lock())->GetStaticCondBlock(n);
// Extract boundary block.
bnd_mat = loc_mat->GetBlock(0, 0);
// Offset by number of rows.
int bnd_row = bnd_mat->GetRows();
for (i = 0; i < bnd_row; ++i)
{
gid1 = m_locToGloMap->GetLocalToGlobalBndMap (cnt + i)
- nDirBnd;
sign1 = m_locToGloMap->GetLocalToGlobalBndSign(cnt + i);
if (gid1 < 0)
{
continue;
}
for (j = 0; j < bnd_row; ++j)
{
gid2 = m_locToGloMap->GetLocalToGlobalBndMap (cnt + j)
- nDirBnd;
sign2 = m_locToGloMap->GetLocalToGlobalBndSign(cnt + j);
if (gid2 == gid1)
{
diagonals[gid1] += sign1 * sign2 * (*bnd_mat)(i, j);
}
}
}
cnt += bnd_row;
}
return diagonals;
}
void Nektar::MultiRegions::Preconditioner::BuildPreconditioner ( )
inline

Definition at line 280 of file Preconditioner.h.

References v_BuildPreconditioner().

void Nektar::MultiRegions::Preconditioner::DoMultiplybyInverseTransformationMatrix ( const Array< OneD, NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput 
)
inline
void Nektar::MultiRegions::Preconditioner::DoMultiplybyInverseTransposedTransformationMatrix ( const Array< OneD, NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput 
)
inline
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().

{
v_DoPreconditioner(pInput,pOutput);
}
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().

{
v_DoPreconditionerWithNonVertOutput(pInput,pOutput,pNonVertOutput,
pVertForce);
}
void Nektar::MultiRegions::Preconditioner::DoTransformFromLowEnergy ( Array< OneD, NekDouble > &  pInOut)
inline

Definition at line 250 of file Preconditioner.h.

References v_DoTransformFromLowEnergy().

void Nektar::MultiRegions::Preconditioner::DoTransformToLowEnergy ( Array< OneD, NekDouble > &  pInOut,
int  offset 
)
inline

Definition at line 231 of file Preconditioner.h.

References v_DoTransformToLowEnergy().

{
v_DoTransformToLowEnergy(pInOut,offset);
}
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().

{
v_DoTransformToLowEnergy(pInput,pOutput);
}
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().

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

{
return v_TransformedSchurCompl(offset,loc_mat);
}
void Nektar::MultiRegions::Preconditioner::v_BuildPreconditioner ( )
privatevirtual

Definition at line 171 of file Preconditioner.cpp.

Referenced by BuildPreconditioner().

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

Multiply by the block inverse transformation matrix.

Definition at line 154 of file Preconditioner.cpp.

References ErrorUtil::efatal, and NEKERROR.

Referenced by DoMultiplybyInverseTransformationMatrix().

{
NEKERROR(ErrorUtil::efatal,"Method does not exist" );
}
void Nektar::MultiRegions::Preconditioner::v_DoMultiplybyInverseTransposedTransformationMatrix ( const Array< OneD, NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput 
)
privatevirtual

Multiply by the block transposed inverse transformation matrix.

Definition at line 164 of file Preconditioner.cpp.

References ErrorUtil::efatal, and NEKERROR.

Referenced by DoMultiplybyInverseTransposedTransformationMatrix().

{
NEKERROR(ErrorUtil::efatal,"Method does not exist" );
}
void Nektar::MultiRegions::Preconditioner::v_DoPreconditioner ( const Array< OneD, NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput 
)
privatevirtual

Apply a preconditioner to the conjugate gradient method.

Definition at line 104 of file Preconditioner.cpp.

References ErrorUtil::efatal, and NEKERROR.

Referenced by DoPreconditioner().

{
NEKERROR(ErrorUtil::efatal,"Method does not exist" );
}
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.

Definition at line 115 of file Preconditioner.cpp.

References ErrorUtil::efatal, and NEKERROR.

Referenced by DoPreconditionerWithNonVertOutput().

{
NEKERROR(ErrorUtil::efatal, "Method does not exist");
}
void Nektar::MultiRegions::Preconditioner::v_DoTransformFromLowEnergy ( Array< OneD, NekDouble > &  pInput)
privatevirtual

Transform from low energy basis to orignal basis.

Definition at line 145 of file Preconditioner.cpp.

References Vmath::Smul().

Referenced by DoTransformFromLowEnergy().

{
Vmath::Smul(pInput.num_elements(), 1.0, pInput, 1, pInput, 1);
}
void Nektar::MultiRegions::Preconditioner::v_DoTransformToLowEnergy ( Array< OneD, NekDouble > &  pInOut,
int  offset 
)
privatevirtual

Transform from original basis to low energy basis.

Definition at line 127 of file Preconditioner.cpp.

Referenced by DoTransformToLowEnergy().

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

Transform from original basis to low energy basis.

Definition at line 136 of file Preconditioner.cpp.

{
}
void Nektar::MultiRegions::Preconditioner::v_InitObject ( )
privatevirtual

Definition at line 96 of file Preconditioner.cpp.

References ErrorUtil::efatal, and NEKERROR.

Referenced by InitObject().

{
NEKERROR(ErrorUtil::efatal,"Method does not exist" );
}
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}$.

Definition at line 179 of file Preconditioner.cpp.

Referenced by TransformedSchurCompl().

{
return loc_mat;
}

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.