35 #ifndef NEKTAR_LIB_MULTIREGIONS_GLOBALLINSYSPETSC_H
36 #define NEKTAR_LIB_MULTIREGIONS_GLOBALLINSYSPETSC_H
46 namespace MultiRegions
65 const boost::weak_ptr<ExpList> &pExp,
66 const boost::shared_ptr<AssemblyMap> &pLocToGloMap);
140 Vec &in, Vec &out,
ShellCtx *ctx,
bool precon);
VecScatter m_ctx
PETSc scatter context that takes us between Nektar++ global ordering and PETSc vector ordering...
GlobalLinSysPETSc * linSys
Pointer to the original calling object.
PreconditionerSharedPtr m_precon
boost::shared_ptr< AssemblyMap > AssemblyMapSharedPtr
void CalculateReordering(const Array< OneD, const int > &glo2uniMap, const Array< OneD, const int > &glo2unique, const AssemblyMapSharedPtr &pLocToGloMap)
Calculate a reordering of universal IDs for PETSc.
PC m_pc
PCShell for preconditioner.
#define MULTI_REGIONS_EXPORT
static std::string matMult
boost::shared_ptr< Preconditioner > PreconditionerSharedPtr
void SetUpSolver(NekDouble tolerance)
Set up KSP solver object.
Internal struct for MatShell and PCShell calls to store current context for callback.
vector< int > m_reorderedMap
Reordering that takes universal IDs to a unique row in the PETSc matrix.
Vec m_x
PETSc vector objects used for local storage.
static PetscErrorCode DoPreconditioner(PC pc, Vec in, Vec out)
Apply preconditioning using Nektar++ routines.
int nGlobal
Number of global degrees of freedom.
GlobalLinSysPETSc(const GlobalLinSysKey &pKey, const boost::weak_ptr< ExpList > &pExp, const boost::shared_ptr< AssemblyMap > &pLocToGloMap)
Constructor for full direct matrix solve.
Mat m_matrix
PETSc matrix object.
int m_nLocal
Number of unique degrees of freedom on this process.
KSP m_ksp
KSP object that represents solver system.
virtual void v_DoMatrixMultiply(const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)=0
virtual ~GlobalLinSysPETSc()
Clean up PETSc objects.
static PetscErrorCode DoMatrixMultiply(Mat M, Vec in, Vec out)
Perform matrix multiplication using Nektar++ routines.
Describe a linear system.
A PETSc global linear system.
void SetUpScatter()
Set up PETSc local (equivalent to Nektar++ global) and global (equivalent to universal) scatter maps...
int nDir
Number of Dirichlet degrees of freedom.
virtual void v_SolveLinearSystem(const int pNumRows, const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const AssemblyMapSharedPtr &locToGloMap, const int pNumDir)
Solve linear system using PETSc.
static std::string matMultIds[]
static PetscErrorCode DoDestroyMatCtx(Mat M)
Destroy matrix shell context object.
void SetUpMatVec(int nGlobal, int nDir)
Construct PETSc matrix and vector handles.
static void DoNekppOperation(Vec &in, Vec &out, ShellCtx *ctx, bool precon)
Perform either matrix multiplication or preconditioning using Nektar++ routines.
static PetscErrorCode DoDestroyPCCtx(PC pc)
Destroy preconditioner context object.
PETScMatMult m_matMult
Enumerator to select matrix multiplication type.