17 #include <deal.II/base/logstream.h> 18 #include <deal.II/lac/petsc_solver.h> 20 #ifdef DEAL_II_WITH_PETSC 22 # include <deal.II/lac/exceptions.h> 23 # include <deal.II/lac/petsc_compatibility.h> 24 # include <deal.II/lac/petsc_matrix_base.h> 25 # include <deal.II/lac/petsc_vector_base.h> 26 # include <deal.II/lac/petsc_precondition.h> 29 #include <petscversion.h> 31 DEAL_II_NAMESPACE_OPEN
89 const Mat B = preconditioner;
91 ExcMessage(
"PETSc preconditioner should have an" 92 "associated matrix set to be used in solver."));
96 #if DEAL_II_PETSC_VERSION_LT(3, 5, 0) 99 ierr = KSPSetOperators (
solver_data->ksp, A, preconditioner,
100 SAME_PRECONDITIONER);
102 ierr = KSPSetOperators (
solver_data->ksp, A, preconditioner);
170 const PetscInt iteration,
171 const PetscReal residual_norm,
172 KSPConvergedReason *reason,
173 void *solver_control_x)
182 case ::SolverControl::iterate:
183 *reason = KSP_CONVERGED_ITERATING;
186 case ::SolverControl::success:
187 *reason =
static_cast<KSPConvergedReason
>(1);
190 case ::SolverControl::failure:
192 *reason = KSP_DIVERGED_ITS;
194 *reason = KSP_DIVERGED_DTOL;
264 PetscErrorCode ierr = KSPSetType (ksp, KSPRICHARDSON);
274 ierr = KSPSetInitialGuessNonzero (ksp, PETSC_TRUE);
299 const MPI_Comm &mpi_communicator,
303 additional_data (data)
314 #if DEAL_II_PETSC_VERSION_LT(3,3,0) 315 PetscErrorCode ierr = KSPSetType (ksp, KSPCHEBYCHEV);
317 PetscErrorCode ierr = KSPSetType (ksp, KSPCHEBYSHEV);
324 ierr = KSPSetInitialGuessNonzero (ksp, PETSC_TRUE);
332 const MPI_Comm &mpi_communicator,
336 additional_data (data)
343 PetscErrorCode ierr = KSPSetType (ksp, KSPCG);
349 ierr = KSPSetInitialGuessNonzero (ksp, PETSC_TRUE);
357 const MPI_Comm &mpi_communicator,
361 additional_data (data)
368 PetscErrorCode ierr = KSPSetType (ksp, KSPBICG);
374 ierr = KSPSetInitialGuessNonzero (ksp, PETSC_TRUE);
383 const bool right_preconditioning)
385 restart_parameter (restart_parameter),
386 right_preconditioning (right_preconditioning)
403 PetscErrorCode ierr = KSPSetType (ksp, KSPGMRES);
426 int (*fun_ptr)(KSP,int);
427 ierr = PetscObjectQueryFunction((PetscObject)(ksp),
428 "KSPGMRESSetRestart_C",
429 (
void (* *)())&fun_ptr);
439 #if DEAL_II_PETSC_VERSION_LT(3,2,0) 440 ierr = KSPSetPreconditionerSide(ksp, PC_RIGHT);
442 ierr = KSPSetPCSide(ksp, PC_RIGHT);
451 ierr = KSPSetInitialGuessNonzero (ksp, PETSC_TRUE);
459 const MPI_Comm &mpi_communicator,
463 additional_data (data)
470 PetscErrorCode ierr = KSPSetType (ksp, KSPBCGS);
476 ierr = KSPSetInitialGuessNonzero (ksp, PETSC_TRUE);
484 const MPI_Comm &mpi_communicator,
488 additional_data (data)
495 PetscErrorCode ierr = KSPSetType (ksp, KSPCGS);
501 ierr = KSPSetInitialGuessNonzero (ksp, PETSC_TRUE);
509 const MPI_Comm &mpi_communicator,
513 additional_data (data)
520 PetscErrorCode ierr = KSPSetType (ksp, KSPTFQMR);
526 ierr = KSPSetInitialGuessNonzero (ksp, PETSC_TRUE);
534 const MPI_Comm &mpi_communicator,
538 additional_data (data)
545 PetscErrorCode ierr = KSPSetType (ksp, KSPTCQMR);
551 ierr = KSPSetInitialGuessNonzero (ksp, PETSC_TRUE);
559 const MPI_Comm &mpi_communicator,
563 additional_data (data)
570 PetscErrorCode ierr = KSPSetType (ksp, KSPCR);
576 ierr = KSPSetInitialGuessNonzero (ksp, PETSC_TRUE);
584 const MPI_Comm &mpi_communicator,
588 additional_data (data)
595 PetscErrorCode ierr = KSPSetType (ksp, KSPLSQR);
601 ierr = KSPSetInitialGuessNonzero (ksp, PETSC_TRUE);
609 const MPI_Comm &mpi_communicator,
613 additional_data (data)
620 PetscErrorCode ierr = KSPSetType (ksp, KSPPREONLY);
636 ierr = KSPSetInitialGuessNonzero (ksp, PETSC_FALSE);
667 PetscErrorCode ierr = KSPSetType (ksp, KSPPREONLY);
682 ierr = KSPSetInitialGuessNonzero (ksp, PETSC_FALSE);
691 #ifdef PETSC_HAVE_MUMPS 702 PetscInt ival=2, icntl=7;
715 if (solver_data.get() == 0)
730 #if DEAL_II_PETSC_VERSION_LT(3, 5, 0) 731 ierr = KSPSetOperators (solver_data->ksp, A, A,
732 DIFFERENT_NONZERO_PATTERN);
734 ierr = KSPSetOperators (solver_data->ksp, A, A);
746 ierr = KSPGetPC (solver_data->ksp, & solver_data->pc);
754 ierr = PCSetType (solver_data->pc, PCCHOLESKY);
756 ierr = PCSetType (solver_data->pc, PCLU);
773 #if DEAL_II_PETSC_VERSION_GTE(3,2,0) 774 ierr = PCFactorSetMatSolverPackage (solver_data->pc, MATSOLVERMUMPS);
776 ierr = PCFactorSetMatSolverPackage (solver_data->pc, MAT_SOLVER_MUMPS);
783 ierr = PCFactorSetUpMatSolverPackage (solver_data->pc);
791 ierr = PCFactorGetMatrix(solver_data->pc, &F);
797 ierr = MatMumpsSetIcntl (F, icntl, ival);
803 ierr = KSPSetOptionsPrefix(solver_data->ksp,
prefix_name.c_str());
810 ierr = KSPSetFromOptions (solver_data->ksp);
818 PetscErrorCode ierr = KSPSolve (solver_data->ksp, b, x);
835 ierr = KSPGetIterationNumber (solver_data->ksp, &its);
837 ierr = KSPGetResidualNorm (solver_data->ksp, &rnorm);
841 #else // PETSC_HAVE_MUMPS 843 ExcMessage (
"Your PETSc installation does not include a copy of " 844 "the MUMPS package necessary for this solver. You will need to configure " 845 "PETSc so that it includes MUMPS, recompile it, and then re-configure " 846 "and recompile deal.II as well."));
857 const PetscInt iteration,
858 const PetscReal residual_norm,
859 KSPConvergedReason *reason,
860 void *solver_control_x)
869 case ::SolverControl::iterate:
870 *reason = KSP_CONVERGED_ITERATING;
873 case ::SolverControl::success:
874 *reason =
static_cast<KSPConvergedReason
>(1);
877 case ::SolverControl::failure:
879 *reason = KSP_DIVERGED_ITS;
881 *reason = KSP_DIVERGED_DTOL;
899 DEAL_II_NAMESPACE_CLOSE
901 #endif // DEAL_II_WITH_PETSC void solve(const MatrixBase &A, VectorBase &x, const VectorBase &b, const PreconditionerBase &preconditioner)
virtual void set_solver_type(KSP &ksp) const =0
virtual State check(const unsigned int step, const double check_value)
const AdditionalData additional_data
virtual void set_solver_type(KSP &ksp) const
AdditionalData(const unsigned int restart_parameter=30, const bool right_preconditioning=false)
SolverRichardson(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
virtual void set_solver_type(KSP &ksp) const
std_cxx11::shared_ptr< SolverData > solver_data
const AdditionalData additional_data
virtual void set_solver_type(KSP &ksp) const
#define AssertThrow(cond, exc)
SparseDirectMUMPS(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
static ::ExceptionBase & ExcMessage(std::string arg1)
void set_prefix(const std::string &prefix)
const MPI_Comm mpi_communicator
Stop iteration, goal reached.
static ::ExceptionBase & ExcPETScError(int arg1) 1
bool right_preconditioning
#define Assert(cond, exc)
SolverBase(SolverControl &cn, const MPI_Comm &mpi_communicator)
static PetscErrorCode convergence_test(KSP ksp, const PetscInt iteration, const PetscReal residual_norm, KSPConvergedReason *reason, void *solver_control)
virtual void set_solver_type(KSP &ksp) const
virtual void set_solver_type(KSP &ksp) const
virtual void set_solver_type(KSP &ksp) const
AdditionalData(const double omega=1)
SolverCR(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
SolverControl & solver_control
double last_value() const
SolverBiCG(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
SolverLSQR(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
SolverChebychev(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
SolverCGS(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
virtual void set_solver_type(KSP &ksp) const
unsigned int last_step() const
PetscErrorCode destroy_krylov_solver(KSP &krylov_solver)
void initialize(const PreconditionerBase &preconditioner)
SolverBicgstab(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
const AdditionalData additional_data
void solve(const MatrixBase &A, VectorBase &x, const VectorBase &b)
virtual void set_solver_type(KSP &ksp) const
virtual void set_solver_type(KSP &ksp) const
SolverTCQMR(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
SolverControl & control() const
virtual void set_solver_type(KSP &ksp) const
unsigned int max_steps() const
SolverTFQMR(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
static ::ExceptionBase & ExcNotImplemented()
void set_symmetric_mode(const bool flag)
virtual void set_solver_type(KSP &ksp) const
virtual void set_solver_type(KSP &ksp) const
static PetscErrorCode convergence_test(KSP ksp, const PetscInt iteration, const PetscReal residual_norm, KSPConvergedReason *reason, void *solver_control)
SolverCG(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
SolverPreOnly(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
const PC & get_pc() const
SolverGMRES(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
unsigned int restart_parameter
virtual void set_solver_type(KSP &ksp) const