17 #include <deal.II/base/logstream.h> 18 #include <deal.II/base/std_cxx14/memory.h> 19 #include <deal.II/lac/petsc_solver.h> 21 #ifdef DEAL_II_WITH_PETSC 23 # include <deal.II/lac/exceptions.h> 24 # include <deal.II/lac/petsc_compatibility.h> 25 # include <deal.II/lac/petsc_matrix_base.h> 26 # include <deal.II/lac/petsc_vector_base.h> 27 # include <deal.II/lac/petsc_precondition.h> 30 #include <petscversion.h> 32 DEAL_II_NAMESPACE_OPEN
85 const Mat B = preconditioner;
87 ExcMessage(
"PETSc preconditioner should have an" 88 "associated matrix set to be used in solver."));
92 #if DEAL_II_PETSC_VERSION_LT(3, 5, 0) 95 ierr = KSPSetOperators (
solver_data->ksp, A, preconditioner,
98 ierr = KSPSetOperators (
solver_data->ksp, A, preconditioner);
166 const PetscInt iteration,
167 const PetscReal residual_norm,
168 KSPConvergedReason *reason,
169 void *solver_control_x)
178 case ::SolverControl::iterate:
179 *reason = KSP_CONVERGED_ITERATING;
182 case ::SolverControl::success:
183 *reason =
static_cast<KSPConvergedReason
>(1);
186 case ::SolverControl::failure:
188 *reason = KSP_DIVERGED_ITS;
190 *reason = KSP_DIVERGED_DTOL;
206 solver_data = std_cxx14::make_unique<SolverData>();
260 PetscErrorCode ierr = KSPSetType (ksp, KSPRICHARDSON);
270 ierr = KSPSetInitialGuessNonzero (ksp, PETSC_TRUE);
295 const MPI_Comm &mpi_communicator,
299 additional_data (data)
306 PetscErrorCode ierr = KSPSetType (ksp, KSPCHEBYSHEV);
312 ierr = KSPSetInitialGuessNonzero (ksp, PETSC_TRUE);
320 const MPI_Comm &mpi_communicator,
324 additional_data (data)
331 PetscErrorCode ierr = KSPSetType (ksp, KSPCG);
337 ierr = KSPSetInitialGuessNonzero (ksp, PETSC_TRUE);
345 const MPI_Comm &mpi_communicator,
349 additional_data (data)
356 PetscErrorCode ierr = KSPSetType (ksp, KSPBICG);
362 ierr = KSPSetInitialGuessNonzero (ksp, PETSC_TRUE);
371 const bool right_preconditioning)
373 restart_parameter (restart_parameter),
374 right_preconditioning (right_preconditioning)
391 PetscErrorCode ierr = KSPSetType (ksp, KSPGMRES);
414 int (*fun_ptr)(KSP,int);
415 ierr = PetscObjectQueryFunction((PetscObject)(ksp),
416 "KSPGMRESSetRestart_C",
417 (
void (* *)())&fun_ptr);
427 ierr = KSPSetPCSide(ksp, PC_RIGHT);
434 ierr = KSPSetInitialGuessNonzero (ksp, PETSC_TRUE);
442 const MPI_Comm &mpi_communicator,
446 additional_data (data)
453 PetscErrorCode ierr = KSPSetType (ksp, KSPBCGS);
459 ierr = KSPSetInitialGuessNonzero (ksp, PETSC_TRUE);
467 const MPI_Comm &mpi_communicator,
471 additional_data (data)
478 PetscErrorCode ierr = KSPSetType (ksp, KSPCGS);
484 ierr = KSPSetInitialGuessNonzero (ksp, PETSC_TRUE);
492 const MPI_Comm &mpi_communicator,
496 additional_data (data)
503 PetscErrorCode ierr = KSPSetType (ksp, KSPTFQMR);
509 ierr = KSPSetInitialGuessNonzero (ksp, PETSC_TRUE);
517 const MPI_Comm &mpi_communicator,
521 additional_data (data)
528 PetscErrorCode ierr = KSPSetType (ksp, KSPTCQMR);
534 ierr = KSPSetInitialGuessNonzero (ksp, PETSC_TRUE);
542 const MPI_Comm &mpi_communicator,
546 additional_data (data)
553 PetscErrorCode ierr = KSPSetType (ksp, KSPCR);
559 ierr = KSPSetInitialGuessNonzero (ksp, PETSC_TRUE);
567 const MPI_Comm &mpi_communicator,
571 additional_data (data)
578 PetscErrorCode ierr = KSPSetType (ksp, KSPLSQR);
584 ierr = KSPSetInitialGuessNonzero (ksp, PETSC_TRUE);
592 const MPI_Comm &mpi_communicator,
596 additional_data (data)
603 PetscErrorCode ierr = KSPSetType (ksp, KSPPREONLY);
619 ierr = KSPSetInitialGuessNonzero (ksp, PETSC_FALSE);
652 PetscErrorCode ierr = KSPSetType (ksp, KSPPREONLY);
667 ierr = KSPSetInitialGuessNonzero (ksp, PETSC_FALSE);
676 #ifdef DEAL_II_PETSC_WITH_MUMPS 687 PetscInt ival=2, icntl=7;
700 if (solver_data ==
nullptr)
702 solver_data = std_cxx14::make_unique<SolverDataMUMPS >();
715 #if DEAL_II_PETSC_VERSION_LT(3, 5, 0) 716 ierr = KSPSetOperators (solver_data->ksp, A, A,
717 DIFFERENT_NONZERO_PATTERN);
719 ierr = KSPSetOperators (solver_data->ksp, A, A);
731 ierr = KSPGetPC (solver_data->ksp, & solver_data->pc);
739 ierr = PCSetType (solver_data->pc, PCCHOLESKY);
741 ierr = PCSetType (solver_data->pc, PCLU);
758 #if DEAL_II_PETSC_VERSION_LT(3, 9, 0) 759 ierr = PCFactorSetMatSolverPackage (solver_data->pc, MATSOLVERMUMPS);
761 ierr = PCFactorSetMatSolverType (solver_data->pc, MATSOLVERMUMPS);
768 #if DEAL_II_PETSC_VERSION_LT(3, 9, 0) 769 ierr = PCFactorSetUpMatSolverPackage (solver_data->pc);
771 ierr = PCFactorSetUpMatSolverType (solver_data->pc);
780 ierr = PCFactorGetMatrix(solver_data->pc, &F);
786 ierr = MatMumpsSetIcntl (F, icntl, ival);
792 ierr = KSPSetOptionsPrefix(solver_data->ksp,
prefix_name.c_str());
799 ierr = KSPSetFromOptions (solver_data->ksp);
807 PetscErrorCode ierr = KSPSolve (solver_data->ksp, b, x);
824 ierr = KSPGetIterationNumber (solver_data->ksp, &its);
826 ierr = KSPGetResidualNorm (solver_data->ksp, &rnorm);
830 #else // DEAL_II_PETSC_WITH_MUMPS 832 ExcMessage (
"Your PETSc installation does not include a copy of " 833 "the MUMPS package necessary for this solver. You will need to configure " 834 "PETSc so that it includes MUMPS, recompile it, and then re-configure " 835 "and recompile deal.II as well."));
846 const PetscInt iteration,
847 const PetscReal residual_norm,
848 KSPConvergedReason *reason,
849 void *solver_control_x)
858 case ::SolverControl::iterate:
859 *reason = KSP_CONVERGED_ITERATING;
862 case ::SolverControl::success:
863 *reason =
static_cast<KSPConvergedReason
>(1);
866 case ::SolverControl::failure:
868 *reason = KSP_DIVERGED_ITS;
870 *reason = KSP_DIVERGED_DTOL;
888 DEAL_II_NAMESPACE_CLOSE
890 #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
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.
std::unique_ptr< SolverData > solver_data
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