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 | List of all members
Nektar::MultiRegions::PreconditionerLinearWithLowEnergy Class Reference

#include <PreconditionerLinearWithLowEnergy.h>

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

Public Member Functions

 PreconditionerLinearWithLowEnergy (const boost::shared_ptr< GlobalLinSys > &plinsys, const AssemblyMapSharedPtr &pLocToGloMap)
virtual ~PreconditionerLinearWithLowEnergy ()
- 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 className
 Name of class.

Protected Attributes

PreconditionerSharedPtr m_linSpacePrecon
PreconditionerSharedPtr m_lowEnergyPrecon
- 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_DoTransformToLowEnergy (Array< OneD, NekDouble > &pInOut, int offset)
virtual void v_DoTransformFromLowEnergy (Array< OneD, NekDouble > &pInput)
virtual DNekScalMatSharedPtr v_TransformedSchurCompl (int offset, const boost::shared_ptr< DNekScalMat > &loc_mat)
virtual void v_DoPreconditioner (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
virtual void v_BuildPreconditioner ()

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 49 of file PreconditionerLinearWithLowEnergy.h.

Constructor & Destructor Documentation

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

Definition at line 63 of file PreconditionerLinearWithLowEnergy.cpp.

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

Definition at line 71 of file PreconditionerLinearWithLowEnergy.h.

{}

Member Function Documentation

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

Creates an instance of this class.

Definition at line 53 of file PreconditionerLinearWithLowEnergy.h.

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

Definition at line 110 of file PreconditionerLinearWithLowEnergy.cpp.

References m_linSpacePrecon, and m_lowEnergyPrecon.

{
m_linSpacePrecon->BuildPreconditioner();
m_lowEnergyPrecon->BuildPreconditioner();
}
void Nektar::MultiRegions::PreconditionerLinearWithLowEnergy::v_DoPreconditioner ( const Array< OneD, NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput 
)
privatevirtual

Definition at line 120 of file PreconditionerLinearWithLowEnergy.cpp.

References m_linSpacePrecon, m_lowEnergyPrecon, and Vmath::Vadd().

{
int nGlobal = pInput.num_elements();
Array<OneD, NekDouble> OutputLowEnergy(nGlobal, 0.0);
Array<OneD, NekDouble> tmp(nGlobal, 0.0);
Array<OneD, NekDouble> OutputLinear(nGlobal, 0.0);
Array<OneD, NekDouble> InputLinear(nGlobal, 0.0);
//Apply Low Energy preconditioner
m_lowEnergyPrecon->DoPreconditioner(pInput, OutputLowEnergy);
//Transform input from low energy to original basis
m_lowEnergyPrecon->DoMultiplybyInverseTransformationMatrix(pInput, InputLinear);
//Apply linear space preconditioner
m_linSpacePrecon->DoPreconditionerWithNonVertOutput(InputLinear, OutputLinear, tmp);
m_lowEnergyPrecon->DoMultiplybyInverseTransposedTransformationMatrix(OutputLinear,pOutput);
Vmath::Vadd(nGlobal,pOutput,1,OutputLowEnergy,1,pOutput,1);
}
void Nektar::MultiRegions::PreconditionerLinearWithLowEnergy::v_DoTransformFromLowEnergy ( Array< OneD, NekDouble > &  pInput)
privatevirtual

Definition at line 92 of file PreconditionerLinearWithLowEnergy.cpp.

References m_lowEnergyPrecon.

{
m_lowEnergyPrecon->DoTransformFromLowEnergy(pInput);
}
void Nektar::MultiRegions::PreconditionerLinearWithLowEnergy::v_DoTransformToLowEnergy ( Array< OneD, NekDouble > &  pInOut,
int  offset 
)
privatevirtual

Definition at line 82 of file PreconditionerLinearWithLowEnergy.cpp.

References m_lowEnergyPrecon.

{
m_lowEnergyPrecon->DoTransformToLowEnergy(pInOut,offset);
}
void Nektar::MultiRegions::PreconditionerLinearWithLowEnergy::v_InitObject ( )
privatevirtual
DNekScalMatSharedPtr Nektar::MultiRegions::PreconditionerLinearWithLowEnergy::v_TransformedSchurCompl ( int  offset,
const boost::shared_ptr< DNekScalMat > &  loc_mat 
)
privatevirtual

Definition at line 100 of file PreconditionerLinearWithLowEnergy.cpp.

References m_lowEnergyPrecon.

{
returnval=m_lowEnergyPrecon->TransformedSchurCompl(offset,loc_mat);
return returnval;
}

Member Data Documentation

string Nektar::MultiRegions::PreconditionerLinearWithLowEnergy::className
static
Initial value:
"FullLinearSpaceWithLowEnergyBlock",
"Full Linear space and low energy block preconditioning")

Name of class.

Registers the class with the Factory.

Definition at line 64 of file PreconditionerLinearWithLowEnergy.h.

PreconditionerSharedPtr Nektar::MultiRegions::PreconditionerLinearWithLowEnergy::m_linSpacePrecon
protected
PreconditionerSharedPtr Nektar::MultiRegions::PreconditionerLinearWithLowEnergy::m_lowEnergyPrecon
protected