16#ifndef dealii_petsc_solver_h
17# define dealii_petsc_solver_h
22# ifdef DEAL_II_WITH_PETSC
31# ifdef DEAL_II_WITH_SLEPC
39# ifdef DEAL_II_WITH_SLEPC
43 class TransformationBase;
54 class PreconditionBase;
210 static PetscErrorCode
212 const PetscInt iteration,
213 const PetscReal residual_norm,
214 KSPConvergedReason *reason,
253# ifdef DEAL_II_WITH_SLEPC
968 static PetscErrorCode
970 const PetscInt iteration,
971 const PetscReal residual_norm,
972 KSPConvergedReason *reason,
virtual void set_solver_type(KSP &ksp) const =0
SolverControl & control() const
SolverControl & solver_control
void initialize(const PreconditionBase &preconditioner)
static PetscErrorCode convergence_test(KSP ksp, const PetscInt iteration, const PetscReal residual_norm, KSPConvergedReason *reason, void *solver_control)
virtual ~SolverBase()=default
const MPI_Comm mpi_communicator
void set_prefix(const std::string &prefix)
void solve(const MatrixBase &A, VectorBase &x, const VectorBase &b, const PreconditionBase &preconditioner)
std::unique_ptr< SolverData > solver_data
SolverBase(SolverControl &cn, const MPI_Comm &mpi_communicator)
const AdditionalData additional_data
SolverBiCG(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
virtual void set_solver_type(KSP &ksp) const override
SolverBicgstab(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
virtual void set_solver_type(KSP &ksp) const override
const AdditionalData additional_data
SolverCGS(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
const AdditionalData additional_data
virtual void set_solver_type(KSP &ksp) const override
SolverCG(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
const AdditionalData additional_data
virtual void set_solver_type(KSP &ksp) const override
virtual void set_solver_type(KSP &ksp) const override
SolverCR(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
const AdditionalData additional_data
const AdditionalData additional_data
virtual void set_solver_type(KSP &ksp) const override
SolverChebychev(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
const AdditionalData additional_data
SolverGMRES(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
virtual void set_solver_type(KSP &ksp) const override
virtual void set_solver_type(KSP &ksp) const override
SolverLSQR(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
const AdditionalData additional_data
SolverPreOnly(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
const AdditionalData additional_data
virtual void set_solver_type(KSP &ksp) const override
SolverRichardson(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
virtual void set_solver_type(KSP &ksp) const override
const AdditionalData additional_data
const AdditionalData additional_data
virtual void set_solver_type(KSP &ksp) const override
SolverTCQMR(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
SolverTFQMR(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
const AdditionalData additional_data
virtual void set_solver_type(KSP &ksp) const override
void solve(const MatrixBase &A, VectorBase &x, const VectorBase &b)
virtual void set_solver_type(KSP &ksp) const override
const AdditionalData additional_data
void set_symmetric_mode(const bool flag)
static PetscErrorCode convergence_test(KSP ksp, const PetscInt iteration, const PetscReal residual_norm, KSPConvergedReason *reason, void *solver_control)
SparseDirectMUMPS(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
std::unique_ptr< SolverDataMUMPS > solver_data
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
AdditionalData(const unsigned int restart_parameter=30, const bool right_preconditioning=false)
bool right_preconditioning
unsigned int restart_parameter
AdditionalData(const double omega=1)