Chombo + EB  3.2
Public Member Functions | Static Public Member Functions | Public Attributes | Protected Member Functions | Protected Attributes | List of all members
PetscSolver< T > Class Template Referenceabstract

#include <PetscSolver.H>

Inheritance diagram for PetscSolver< T >:
Inheritance graph
[legend]

Public Member Functions

 PetscSolver ()
 
virtual ~PetscSolver ()
 
void destroy ()
 
virtual void setHomogeneous (bool a_homogeneous)
 
virtual void setFunctionAndJacobian (PetscErrorCode(*f)(SNES, Vec, Vec, void *), PetscErrorCode(*j)(SNES, Vec, Mat *, Mat *, MatStructure *, void *))
 
virtual void define (Real a_dx, bool a_homogeneous=true)
 
virtual void define (LinearOp< T > *, bool a_homogeneous=true)
 
virtual void solve (T &a_phi, const T &a_rhs)
 
int solve_private (T &a_phi, const T &a_rhs)
 
virtual void solve_mfree (T &a_phi, const T &a_rhs, LinearOp< T > *a_op)
 
int solve_mfree_private (T &a_phi, const T &a_rhs, LinearOp< T > *a_op)
 
Real computeResidual ()
 
virtual int applyOp (T &a_phi, const T &a_rhs)
 
void setInitialGuessNonzero (bool b=true)
 
void setOptionsPrefix (const char prefix[])
 
KSP getKSP ()
 
virtual int getNNZPerRow () const
 
virtual bool supportNNZExact () const
 
virtual void rhsOperation (const T &a_rhs)
 
virtual Real addBCdiagValue (const IntVect &a_iv, const IntVect &a_jv, const T &a_rhs, const DataIndex &a_datInd, const Real coeff=1)
 
int resetOperator ()
 
Real normInfinity (const T &a_phi) const
 
int setup_solver (const T &a_phi)
 
int create_mat_vec (const T &a_phi)
 
PetscErrorCode putPetscInChombo (T &a_phi, const Vec xx)
 
PetscErrorCode putChomboInPetsc (Vec out, const T &a_phi)
 
virtual int formMatrix (Mat, const T *=0, PetscInt my0=0, PetscInt nloc=0, PetscInt *d=0, PetscInt *o=0)=0
 
virtual BaseFab< Real > & getRegFab (T &a_fab, const DataIndex &a_datInd)=0
 
virtual const BaseFab< Real > & getRegFab (const T &a_fab, const DataIndex &a_datInd) const =0
 
virtual const BaseFab< Real > & getRegFab (const T &a_fab, const DataIndex &a_datInd, Box &a_box) const =0
 
void setNull (bool n=true)
 
- Public Member Functions inherited from LinearSolver< T >
virtual ~LinearSolver ()
 
virtual void setConvergenceMetrics (Real a_metric, Real a_tolerance)
 Set a convergence metric, along with solver tolerance, if desired. More...
 

Static Public Member Functions

static PetscErrorCode ksp_monitor_pout (KSP ksp, PetscInt it, PetscReal rnorm, void *ctx)
 
static PetscErrorCode apply_mfree (Mat A, Vec x, Vec f)
 

Public Attributes

LinearOp< T > * m_op_mfree
 
m_phi_mfree
 
m_Lphi_mfree
 
bool m_mfree_homogeneous
 
bool m_homogeneous
 
Real m_dx
 
Mat m_mat
 
void * m_ctx
 
Vec m_xx
 
Vec m_rr
 
Vec m_bb
 
SNES m_snes
 
KSP m_ksp
 
LevelData< BaseFab< PetscInt > > m_gids
 
PetscInt m_gid0
 

Protected Member Functions

virtual Real addBCrhsValue (const IntVect &a_iv, const T &a_phi, const DataIndex &a_datInd, const Real &coeff=1)
 

Protected Attributes

PetscInt m_defined
 
PetscErrorCode(* m_function )(SNES, Vec, Vec, void *)
 
PetscErrorCode(* m_jacobian )(SNES, Vec, Mat *, Mat *, MatStructure *, void *)
 
char m_prestring [32]
 
bool m_null
 
bool m_nz_init_guess
 
LevelData< BaseFab< bool > > m_bccode
 

Detailed Description

template<class T>
class PetscSolver< T >

Interface to PETSC linear solvers BC implementation hooks: addBCdiagValue (formMatrix) and addBCrhsValue (solveprivate) m_bccode should be encoded true for cells that are at the domain boundaries at formMatrix time.

Constructor & Destructor Documentation

◆ PetscSolver()

template<class T >
PetscSolver< T >::PetscSolver ( )

◆ ~PetscSolver()

template<class T >
PetscSolver< T >::~PetscSolver ( )
virtual

Member Function Documentation

◆ destroy()

template<class T >
void PetscSolver< T >::destroy ( )

◆ setHomogeneous()

template<class T>
virtual void PetscSolver< T >::setHomogeneous ( bool  a_homogeneous)
inlinevirtual

Set whether the solver uses only homogeneous boundary conditions

Implements LinearSolver< T >.

◆ setFunctionAndJacobian()

template<class T>
virtual void PetscSolver< T >::setFunctionAndJacobian ( PetscErrorCode(*)(SNES, Vec, Vec, void *)  f,
PetscErrorCode(*)(SNES, Vec, Mat *, Mat *, MatStructure *, void *)  j 
)
inlinevirtual

Set Function F(u) = 0 and Jacobian dF(u)/du for nonlinear solves

◆ ksp_monitor_pout()

template<class T>
static PetscErrorCode PetscSolver< T >::ksp_monitor_pout ( KSP  ksp,
PetscInt  it,
PetscReal  rnorm,
void *  ctx 
)
inlinestatic

◆ define() [1/2]

template<class T >
void PetscSolver< T >::define ( Real  a_dx,
bool  a_homogeneous = true 
)
virtual

◆ define() [2/2]

template<class T>
void PetscSolver< T >::define ( LinearOp< T > *  a_operator,
bool  a_homogeneous = true 
)
virtual

Define the operator and whether it is a homogeneous solver or not. The LinearSolver does not take over ownership of this a_operator object. It does not call delete on it when the LinearSolver is deleted. It is meant to be like a late-binding reference. If you created a_operator with new, you should call delete on it after LinearSolver is deleted if you want to avoid memory leaks.

Implements LinearSolver< T >.

Reimplemented in PetscSolverPoisson< T >.

◆ solve()

template<class T>
void PetscSolver< T >::solve ( T &  a_phi,
const T &  a_rhs 
)
virtual

◆ solve_private()

template<class T>
int PetscSolver< T >::solve_private ( T &  a_phi,
const T &  a_rhs 
)

◆ solve_mfree()

template<class T>
void PetscSolver< T >::solve_mfree ( T &  a_phi,
const T &  a_rhs,
LinearOp< T > *  a_op 
)
virtual

Solve a_op*a_phi = a_rhs using the PETSC matrix free functions The preconditioner (for which a matrix is formed) need not be the same as the actual operator.

Referenced by PetscSolver< LevelData< FArrayBox > >::ksp_monitor_pout().

◆ solve_mfree_private()

template<class T>
int PetscSolver< T >::solve_mfree_private ( T &  a_phi,
const T &  a_rhs,
LinearOp< T > *  a_op 
)

◆ apply_mfree()

template<class T >
PetscErrorCode PetscSolver< T >::apply_mfree ( Mat  A,
Vec  x,
Vec  f 
)
static

◆ computeResidual()

template<class T >
Real PetscSolver< T >::computeResidual ( )

viewResidual

◆ applyOp()

template<class T>
int PetscSolver< T >::applyOp ( T &  a_phi,
const T &  a_rhs 
)
virtual

apply

◆ setInitialGuessNonzero()

template<class T>
void PetscSolver< T >::setInitialGuessNonzero ( bool  b = true)
inline

set initial guess non-zero

◆ setOptionsPrefix()

template<class T>
void PetscSolver< T >::setOptionsPrefix ( const char  prefix[])
inline

set initial guess non-zero

◆ getKSP()

template<class T>
KSP PetscSolver< T >::getKSP ( )
inline

◆ getNNZPerRow()

template<class T>
virtual int PetscSolver< T >::getNNZPerRow ( ) const
inlinevirtual

get an estimate of the number of nnz/row for matrix allocation

Reimplemented in PetscSolverViscousTensor< T >.

Referenced by PetscSolver< LevelData< FArrayBox > >::create_mat_vec().

◆ supportNNZExact()

template<class T>
virtual bool PetscSolver< T >::supportNNZExact ( ) const
inlinevirtual

do I support having formMatrix precompute exact NNZ/row

Referenced by PetscSolver< LevelData< FArrayBox > >::create_mat_vec().

◆ rhsOperation()

template<class T>
virtual void PetscSolver< T >::rhsOperation ( const T &  a_rhs)
inlinevirtual

◆ addBCdiagValue()

template<class T>
virtual Real PetscSolver< T >::addBCdiagValue ( const IntVect a_iv,
const IntVect a_jv,
const T &  a_rhs,
const DataIndex a_datInd,
const Real  coeff = 1 
)
inlinevirtual

handling boundary conditions, turn it into a term to be added to the diag term this function coordinates with addBCrhsValue for Dirichlet BC for Neumann BC no RHS modification is required

◆ resetOperator()

template<class T>
int PetscSolver< T >::resetOperator ( )
inline

◆ normInfinity()

template<class T>
Real PetscSolver< T >::normInfinity ( const T &  a_phi) const

◆ addBCrhsValue()

template<class T>
virtual Real PetscSolver< T >::addBCrhsValue ( const IntVect a_iv,
const T &  a_phi,
const DataIndex a_datInd,
const Real coeff = 1 
)
inlineprotectedvirtual

handling boundary conditions, turn it into a term to be added to the RHS this function coordinates with addBCdiagValue for Dirichlet BC, should return a term that is to be added to RHS

Referenced by PetscSolver< LevelData< FArrayBox > >::putChomboInPetsc(), PetscSolver< LevelData< FArrayBox > >::solve_mfree_private(), and PetscSolver< LevelData< FArrayBox > >::solve_private().

◆ setup_solver()

template<class T>
int PetscSolver< T >::setup_solver ( const T &  a_phi)

◆ create_mat_vec()

template<class T>
int PetscSolver< T >::create_mat_vec ( const T &  a_phi)

◆ putPetscInChombo()

template<class T>
PetscErrorCode PetscSolver< T >::putPetscInChombo ( T &  a_phi,
const Vec  xx 
)

◆ putChomboInPetsc()

template<class T>
PetscErrorCode PetscSolver< T >::putChomboInPetsc ( Vec  out,
const T &  a_phi 
)

◆ formMatrix()

template<class T>
virtual int PetscSolver< T >::formMatrix ( Mat  ,
const T *  = 0,
PetscInt  my0 = 0,
PetscInt  nloc = 0,
PetscInt *  d = 0,
PetscInt *  o = 0 
)
pure virtual

◆ getRegFab() [1/3]

template<class T>
virtual BaseFab<Real>& PetscSolver< T >::getRegFab ( T &  a_fab,
const DataIndex a_datInd 
)
pure virtual

◆ getRegFab() [2/3]

template<class T>
virtual const BaseFab<Real>& PetscSolver< T >::getRegFab ( const T &  a_fab,
const DataIndex a_datInd 
) const
pure virtual

◆ getRegFab() [3/3]

template<class T>
virtual const BaseFab<Real>& PetscSolver< T >::getRegFab ( const T &  a_fab,
const DataIndex a_datInd,
Box a_box 
) const
pure virtual

◆ setNull()

template<class T >
void PetscSolver< T >::setNull ( bool  n = true)

Member Data Documentation

◆ m_op_mfree

template<class T>
LinearOp<T>* PetscSolver< T >::m_op_mfree

◆ m_phi_mfree

template<class T>
T PetscSolver< T >::m_phi_mfree

◆ m_Lphi_mfree

template<class T>
T PetscSolver< T >::m_Lphi_mfree

◆ m_mfree_homogeneous

template<class T>
bool PetscSolver< T >::m_mfree_homogeneous

◆ m_homogeneous

template<class T>
bool PetscSolver< T >::m_homogeneous

◆ m_dx

template<class T>
Real PetscSolver< T >::m_dx

◆ m_mat

template<class T>
Mat PetscSolver< T >::m_mat

◆ m_ctx

template<class T>
void* PetscSolver< T >::m_ctx

◆ m_xx

template<class T>
Vec PetscSolver< T >::m_xx

◆ m_rr

template<class T>
Vec PetscSolver< T >::m_rr

◆ m_bb

template<class T>
Vec PetscSolver< T >::m_bb

◆ m_snes

template<class T>
SNES PetscSolver< T >::m_snes

◆ m_ksp

template<class T>
KSP PetscSolver< T >::m_ksp

◆ m_defined

template<class T>
PetscInt PetscSolver< T >::m_defined
protected

◆ m_function

template<class T>
PetscErrorCode(* PetscSolver< T >::m_function) (SNES, Vec, Vec, void *)
protected

◆ m_jacobian

template<class T>
PetscErrorCode(* PetscSolver< T >::m_jacobian) (SNES, Vec, Mat *, Mat *, MatStructure *, void *)
protected

◆ m_prestring

template<class T>
char PetscSolver< T >::m_prestring[32]
protected

◆ m_null

template<class T>
bool PetscSolver< T >::m_null
protected

◆ m_nz_init_guess

template<class T>
bool PetscSolver< T >::m_nz_init_guess
protected

◆ m_bccode

template<class T>
LevelData<BaseFab<bool> > PetscSolver< T >::m_bccode
protected

◆ m_gids

template<class T>
LevelData<BaseFab<PetscInt> > PetscSolver< T >::m_gids

◆ m_gid0

template<class T>
PetscInt PetscSolver< T >::m_gid0

The documentation for this class was generated from the following files: