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::PreconditionerLinear Class Reference

#include <PreconditionerLinear.h>

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

Public Member Functions

 PreconditionerLinear (const boost::shared_ptr< GlobalLinSys > &plinsys, const AssemblyMapSharedPtr &pLocToGloMap)
virtual ~PreconditionerLinear ()
- 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.
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.

Static Public Attributes

static std::string className1
 Name of class.

Protected Attributes

GlobalLinSysSharedPtr m_vertLinsys
boost::shared_ptr< AssemblyMapm_vertLocToGloMap
- 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

virtual void v_InitObject ()
virtual void v_DoPreconditionerWithNonVertOutput (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const Array< OneD, NekDouble > &pNonVertOutput, Array< OneD, NekDouble > &pVertForce)
virtual void v_DoPreconditioner (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
virtual void v_BuildPreconditioner ()

Static Private Attributes

static std::string solveType
static std::string solveTypeIds []

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

Detailed Description

This class implements preconditioning for the conjugate 

gradient matrix solver.

Definition at line 58 of file PreconditionerLinear.h.

Constructor & Destructor Documentation

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

Definition at line 86 of file PreconditionerLinear.cpp.

: Preconditioner(plinsys, pLocToGloMap)
{
}
virtual Nektar::MultiRegions::PreconditionerLinear::~PreconditionerLinear ( )
inlinevirtual

Definition at line 80 of file PreconditionerLinear.h.

{}

Member Function Documentation

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

Creates an instance of this class.

Definition at line 62 of file PreconditionerLinear.h.

{
p->InitObject();
return p;
}
void Nektar::MultiRegions::PreconditionerLinear::v_BuildPreconditioner ( )
privatevirtual

Definition at line 97 of file PreconditionerLinear.cpp.

References ASSERTL0, Nektar::MultiRegions::eIterativeStaticCond, Nektar::MultiRegions::eLinearPreconPETSc, Nektar::MultiRegions::eLinearPreconXxt, Nektar::StdRegions::eMass, Nektar::MultiRegions::ePETScFullMatrix, Nektar::StdRegions::ePreconLinearSpace, Nektar::StdRegions::ePreconLinearSpaceMass, Nektar::MultiRegions::eXxtFullMatrix, Nektar::MultiRegions::Preconditioner::m_linsys, Nektar::MultiRegions::Preconditioner::m_locToGloMap, m_vertLinsys, m_vertLocToGloMap, and solveType.

{
GlobalSysSolnType sType = m_locToGloMap->GetGlobalSysSolnType();
"This type of preconditioning is not implemented "
"for this solver");
boost::shared_ptr<MultiRegions::ExpList>
expList=((m_linsys.lock())->GetLocMat()).lock();
expList->GetSession()->GetSolverInfoAsEnum<LinearPreconSolver>(
"LinearPreconSolver");
GlobalSysSolnType linSolveType;
switch(solveType)
{
{
linSolveType = eXxtFullMatrix;
break;
}
{
#ifdef NEKTAR_USING_PETSC
linSolveType = ePETScFullMatrix;
#else
ASSERTL0(false, "Nektar++ has not been compiled with "
"PETSc support.");
#endif
}
}
m_vertLocToGloMap = m_locToGloMap->LinearSpaceMap(*expList, linSolveType);
// Generate linear solve system.
m_linsys.lock()->GetKey().GetMatrixType() == StdRegions::eMass ?
GlobalLinSysKey preconKey(mType, m_vertLocToGloMap, (m_linsys.lock())->GetKey().GetConstFactors());
switch(solveType)
{
{
break;
}
{
#ifdef NEKTAR_USING_PETSC
#else
ASSERTL0(false, "Nektar++ has not been compiled with "
"PETSc support.");
#endif
}
}
}
void Nektar::MultiRegions::PreconditionerLinear::v_DoPreconditioner ( const Array< OneD, NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput 
)
privatevirtual
void Nektar::MultiRegions::PreconditionerLinear::v_DoPreconditionerWithNonVertOutput ( const Array< OneD, NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput,
const Array< OneD, NekDouble > &  pNonVertOutput,
Array< OneD, NekDouble > &  pVertForce 
)
privatevirtual

Definition at line 176 of file PreconditionerLinear.cpp.

References ASSERTL0, ASSERTL1, Nektar::MultiRegions::eIterativeStaticCond, Nektar::MultiRegions::Preconditioner::m_locToGloMap, m_vertLinsys, m_vertLocToGloMap, Nektar::NullNekDouble1DArray, Vmath::Vcopy(), and Vmath::Zero().

Referenced by v_DoPreconditioner().

{
GlobalSysSolnType solvertype=m_locToGloMap->GetGlobalSysSolnType();
switch(solvertype)
{
{
int i,val;
int nloc = m_vertLocToGloMap->GetNumLocalCoeffs();
int nglo = m_vertLocToGloMap->GetNumGlobalCoeffs();
// mapping from full space to vertices
Array<OneD, int> LocToGloBnd = m_vertLocToGloMap->GetLocalToGlobalBndMap();
// Global to local for linear solver (different from above)
Array<OneD, int> LocToGlo = m_vertLocToGloMap->GetLocalToGlobalMap();
// number of Dir coeffs in from full problem
int nDirFull = m_locToGloMap->GetNumGlobalDirBndCoeffs();
Array<OneD,NekDouble> In(nglo,0.0);
Array<OneD,NekDouble> Out(nglo,0.0);
// Gather rhs
for(i = 0; i < nloc; ++i)
{
val = LocToGloBnd[i];
if(val >= nDirFull)
{
In[LocToGlo[i]] = pInput[val-nDirFull];
}
}
// Do solve without enforcing any boundary conditions.
m_vertLinsys->SolveLinearSystem(m_vertLocToGloMap->GetNumLocalCoeffs(),
if(pNonVertOutput != NullNekDouble1DArray)
{
ASSERTL1(pNonVertOutput.num_elements() >= pOutput.num_elements(),"Non Vert output is not of sufficient length");
Vmath::Vcopy(pOutput.num_elements(),pNonVertOutput,1,pOutput,1);
}
else
{
//Copy input to output as a unit preconditioner on
//any other value
Vmath::Vcopy(pInput.num_elements(),pInput,1,pOutput,1);
}
if(pVertForce != NullNekDouble1DArray)
{
Vmath::Zero(pVertForce.num_elements(),pVertForce,1);
// Scatter back soln from linear solve
for(i = 0; i < nloc; ++i)
{
val = LocToGloBnd[i];
if(val >= nDirFull)
{
pOutput[val-nDirFull] = Out[LocToGlo[i]];
// copy vertex forcing into this vector
pVertForce[val-nDirFull] = In[LocToGlo[i]];
}
}
}
else
{
// Scatter back soln from linear solve
for(i = 0; i < nloc; ++i)
{
val = LocToGloBnd[i];
if(val >= nDirFull)
{
pOutput[val-nDirFull] = Out[LocToGlo[i]];
}
}
}
}
break;
default:
ASSERTL0(0,"Unsupported solver type");
break;
}
}
void Nektar::MultiRegions::PreconditionerLinear::v_InitObject ( )
privatevirtual

Definition at line 93 of file PreconditionerLinear.cpp.

{
}

Member Data Documentation

string Nektar::MultiRegions::PreconditionerLinear::className1
static
Initial value:
"FullLinearSpace",
"Full Linear space inverse Preconditioning")

Name of class.

Registers the class with the Factory.

Definition at line 73 of file PreconditionerLinear.h.

GlobalLinSysSharedPtr Nektar::MultiRegions::PreconditionerLinear::m_vertLinsys
protected
boost::shared_ptr<AssemblyMap> Nektar::MultiRegions::PreconditionerLinear::m_vertLocToGloMap
protected
std::string Nektar::MultiRegions::PreconditionerLinear::solveType
staticprivate
Initial value:

Definition at line 87 of file PreconditionerLinear.h.

Referenced by v_BuildPreconditioner().

std::string Nektar::MultiRegions::PreconditionerLinear::solveTypeIds
staticprivate
Initial value:

Definition at line 88 of file PreconditionerLinear.h.