17 #include <deal.II/base/logstream.h> 18 #include <deal.II/base/std_cxx14/memory.h> 20 #include <deal.II/lac/petsc_solver.h> 22 #ifdef DEAL_II_WITH_PETSC 24 # include <deal.II/lac/exceptions.h> 25 # include <deal.II/lac/petsc_compatibility.h> 26 # include <deal.II/lac/petsc_matrix_base.h> 27 # include <deal.II/lac/petsc_precondition.h> 28 # include <deal.II/lac/petsc_vector_base.h> 30 # include <petscversion.h> 34 DEAL_II_NAMESPACE_OPEN
90 const Mat B = preconditioner;
92 ExcMessage(
"PETSc preconditioner should have an" 93 "associated matrix set to be used in solver."));
97 # if DEAL_II_PETSC_VERSION_LT(3, 5, 0) 103 SAME_PRECONDITIONER);
105 ierr = KSPSetOperators(
solver_data->ksp, A, preconditioner);
122 PetscErrorCode ierr =
176 const PetscInt iteration,
177 const PetscReal residual_norm,
178 KSPConvergedReason *reason,
179 void * solver_control_x)
189 case ::SolverControl::iterate:
190 *reason = KSP_CONVERGED_ITERATING;
193 case ::SolverControl::success:
194 *reason =
static_cast<KSPConvergedReason
>(1);
197 case ::SolverControl::failure:
199 *reason = KSP_DIVERGED_ITS;
201 *reason = KSP_DIVERGED_DTOL;
217 solver_data = std_cxx14::make_unique<SolverData>();
269 PetscErrorCode ierr = KSPSetType(ksp, KSPRICHARDSON);
279 ierr = KSPSetInitialGuessNonzero(ksp, PETSC_TRUE);
295 ierr = KSPSetTolerances(ksp,
307 const MPI_Comm & mpi_communicator,
310 , additional_data(data)
317 PetscErrorCode ierr = KSPSetType(ksp, KSPCHEBYSHEV);
323 ierr = KSPSetInitialGuessNonzero(ksp, PETSC_TRUE);
331 const MPI_Comm & mpi_communicator,
334 , additional_data(data)
341 PetscErrorCode ierr = KSPSetType(ksp, KSPCG);
347 ierr = KSPSetInitialGuessNonzero(ksp, PETSC_TRUE);
355 const MPI_Comm & mpi_communicator,
358 , additional_data(data)
365 PetscErrorCode ierr = KSPSetType(ksp, KSPBICG);
371 ierr = KSPSetInitialGuessNonzero(ksp, PETSC_TRUE);
379 const unsigned int restart_parameter,
380 const bool right_preconditioning)
381 : restart_parameter(restart_parameter)
382 , right_preconditioning(right_preconditioning)
398 PetscErrorCode ierr = KSPSetType(ksp, KSPGMRES);
421 int (*fun_ptr)(KSP, int);
422 ierr = PetscObjectQueryFunction(reinterpret_cast<PetscObject>(ksp),
423 "KSPGMRESSetRestart_C",
424 reinterpret_cast<void (**)()
>(&fun_ptr));
434 ierr = KSPSetPCSide(ksp, PC_RIGHT);
441 ierr = KSPSetInitialGuessNonzero(ksp, PETSC_TRUE);
449 const MPI_Comm & mpi_communicator,
452 , additional_data(data)
459 PetscErrorCode ierr = KSPSetType(ksp, KSPBCGS);
465 ierr = KSPSetInitialGuessNonzero(ksp, PETSC_TRUE);
473 const MPI_Comm & mpi_communicator,
476 , additional_data(data)
483 PetscErrorCode ierr = KSPSetType(ksp, KSPCGS);
489 ierr = KSPSetInitialGuessNonzero(ksp, PETSC_TRUE);
497 const MPI_Comm & mpi_communicator,
500 , additional_data(data)
507 PetscErrorCode ierr = KSPSetType(ksp, KSPTFQMR);
513 ierr = KSPSetInitialGuessNonzero(ksp, PETSC_TRUE);
521 const MPI_Comm & mpi_communicator,
524 , additional_data(data)
531 PetscErrorCode ierr = KSPSetType(ksp, KSPTCQMR);
537 ierr = KSPSetInitialGuessNonzero(ksp, PETSC_TRUE);
545 const MPI_Comm & mpi_communicator,
548 , additional_data(data)
555 PetscErrorCode ierr = KSPSetType(ksp, KSPCR);
561 ierr = KSPSetInitialGuessNonzero(ksp, PETSC_TRUE);
569 const MPI_Comm & mpi_communicator,
572 , additional_data(data)
579 PetscErrorCode ierr = KSPSetType(ksp, KSPLSQR);
585 ierr = KSPSetInitialGuessNonzero(ksp, PETSC_TRUE);
593 const MPI_Comm & mpi_communicator,
596 , additional_data(data)
603 PetscErrorCode ierr = KSPSetType(ksp, KSPPREONLY);
619 ierr = KSPSetInitialGuessNonzero(ksp, PETSC_FALSE);
651 PetscErrorCode ierr = KSPSetType(ksp, KSPPREONLY);
666 ierr = KSPSetInitialGuessNonzero(ksp, PETSC_FALSE);
675 # ifdef DEAL_II_PETSC_WITH_MUMPS 686 PetscInt ival = 2, icntl = 7;
699 if (solver_data ==
nullptr)
701 solver_data = std_cxx14::make_unique<SolverDataMUMPS>();
714 # if DEAL_II_PETSC_VERSION_LT(3, 5, 0) 716 KSPSetOperators(solver_data->ksp, A, A, DIFFERENT_NONZERO_PATTERN);
718 ierr = KSPSetOperators(solver_data->ksp, A, A);
730 ierr = KSPGetPC(solver_data->ksp, &solver_data->pc);
738 ierr = PCSetType(solver_data->pc, PCCHOLESKY);
740 ierr = PCSetType(solver_data->pc, PCLU);
747 ierr = KSPSetConvergenceTest(solver_data->ksp,
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);
806 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 834 "Your PETSc installation does not include a copy of " 835 "the MUMPS package necessary for this solver. You will need to configure " 836 "PETSc so that it includes MUMPS, recompile it, and then re-configure " 837 "and recompile deal.II as well."));
848 const PetscInt iteration,
849 const PetscReal residual_norm,
850 KSPConvergedReason *reason,
851 void * solver_control_x)
861 case ::SolverControl::iterate:
862 *reason = KSP_CONVERGED_ITERATING;
865 case ::SolverControl::success:
866 *reason =
static_cast<KSPConvergedReason
>(1);
869 case ::SolverControl::failure:
871 *reason = KSP_DIVERGED_ITS;
873 *reason = KSP_DIVERGED_DTOL;
891 DEAL_II_NAMESPACE_CLOSE
893 #endif // DEAL_II_WITH_PETSC virtual void set_solver_type(KSP &ksp) const override
virtual void set_solver_type(KSP &ksp) const =0
virtual void set_solver_type(KSP &ksp) const override
virtual State check(const unsigned int step, const double check_value)
virtual void set_solver_type(KSP &ksp) const override
const AdditionalData additional_data
void set_prefix(const std::string &prefix)
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())
const AdditionalData additional_data
void initialize(const PreconditionerBase &preconditioner)
#define AssertThrow(cond, exc)
virtual void set_solver_type(KSP &ksp) const override
SparseDirectMUMPS(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
static ::ExceptionBase & ExcMessage(std::string arg1)
const MPI_Comm mpi_communicator
Stop iteration, goal reached.
std::unique_ptr< SolverData > solver_data
virtual void set_solver_type(KSP &ksp) const override
bool right_preconditioning
#define Assert(cond, exc)
void solve(const MatrixBase &A, VectorBase &x, const VectorBase &b, const PreconditionerBase &preconditioner)
virtual void set_solver_type(KSP &ksp) const override
virtual void set_solver_type(KSP &ksp) const override
SolverBase(SolverControl &cn, const MPI_Comm &mpi_communicator)
virtual void set_solver_type(KSP &ksp) const override
virtual void set_solver_type(KSP &ksp) const override
AdditionalData(const double omega=1)
virtual void set_solver_type(KSP &ksp) const override
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())
unsigned int last_step() const
PetscErrorCode destroy_krylov_solver(KSP &krylov_solver)
SolverBicgstab(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
const AdditionalData additional_data
SolverControl & control() const
virtual void set_solver_type(KSP &ksp) const override
void solve(const MatrixBase &A, VectorBase &x, const VectorBase &b)
SolverTCQMR(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
unsigned int max_steps() const
SolverTFQMR(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
static PetscErrorCode convergence_test(KSP ksp, const PetscInt iteration, const PetscReal residual_norm, KSPConvergedReason *reason, void *solver_control)
static ::ExceptionBase & ExcNotImplemented()
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)
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 override
virtual void set_solver_type(KSP &ksp) const override