16 #include <deal.II/lac/petsc_precondition.h> 18 #ifdef DEAL_II_WITH_PETSC 20 # include <deal.II/base/utilities.h> 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_solver.h> 26 # include <deal.II/lac/petsc_vector_base.h> 28 # include <petscconf.h> 32 DEAL_II_NAMESPACE_OPEN
58 PetscErrorCode ierr = PCDestroy(&
pc);
70 const PetscErrorCode ierr = PCApply(
pc, src, dst);
87 PetscObjectGetComm(reinterpret_cast<PetscObject>(
matrix), &comm);
90 ierr = PCCreate(comm, &
pc);
93 # if DEAL_II_PETSC_VERSION_LT(3, 5, 0) 109 PreconditionerBase::operator Mat()
const 121 PetscErrorCode ierr = PCCreate(comm, &
pc);
140 PetscErrorCode ierr = PCSetType(
pc, const_cast<char *>(PCJACOBI));
143 ierr = PCSetFromOptions(
pc);
153 matrix =
static_cast<Mat
>(matrix_);
159 PetscErrorCode ierr = PCSetUp(
pc);
171 PetscErrorCode ierr = PCCreate(comm, &
pc);
189 PetscErrorCode ierr = PCSetType(
pc, const_cast<char *>(PCBJACOBI));
192 ierr = PCSetFromOptions(
pc);
203 matrix =
static_cast<Mat
>(matrix_);
209 PetscErrorCode ierr = PCSetUp(
pc);
235 matrix =
static_cast<Mat
>(matrix_);
240 PetscErrorCode ierr = PCSetType(
pc, const_cast<char *>(PCSOR));
247 ierr = PCSetFromOptions(
pc);
276 matrix =
static_cast<Mat
>(matrix_);
281 PetscErrorCode ierr = PCSetType(
pc, const_cast<char *>(PCSOR));
289 ierr = PCSORSetSymmetric(
pc, SOR_SYMMETRIC_SWEEP);
292 ierr = PCSetFromOptions(
pc);
322 matrix =
static_cast<Mat
>(matrix_);
327 PetscErrorCode ierr = PCSetType(
pc, const_cast<char *>(PCEISENSTAT));
334 ierr = PCSetFromOptions(
pc);
364 matrix =
static_cast<Mat
>(matrix_);
369 PetscErrorCode ierr = PCSetType(
pc, const_cast<char *>(PCICC));
376 ierr = PCSetFromOptions(
pc);
405 matrix =
static_cast<Mat
>(matrix_);
410 PetscErrorCode ierr = PCSetType(
pc, const_cast<char *>(PCILU));
417 ierr = PCSetFromOptions(
pc);
428 const bool symmetric_operator,
429 const double strong_threshold,
430 const double max_row_sum,
431 const unsigned int aggressive_coarsening_num_levels,
432 const bool output_details)
433 : symmetric_operator(symmetric_operator)
434 , strong_threshold(strong_threshold)
435 , max_row_sum(max_row_sum)
436 , aggressive_coarsening_num_levels(aggressive_coarsening_num_levels)
437 , output_details(output_details)
448 PetscErrorCode ierr = PCCreate(comm, &
pc);
451 # ifdef DEAL_II_PETSC_WITH_HYPRE 453 # else // DEAL_II_PETSC_WITH_HYPRE 456 ExcMessage(
"Your PETSc installation does not include a copy of " 457 "the hypre package necessary for this preconditioner."));
472 # ifdef DEAL_II_PETSC_WITH_HYPRE 473 PetscErrorCode ierr = PCSetType(
pc, const_cast<char *>(PCHYPRE));
476 ierr = PCHYPRESetType(
pc,
"boomeramg");
488 std::stringstream ssStream;
499 "symmetric-SOR/Jacobi");
501 "symmetric-SOR/Jacobi");
503 "Gaussian-elimination");
510 "Gaussian-elimination");
513 ierr = PCSetFromOptions(
pc);
517 ExcMessage(
"Your PETSc installation does not include a copy of " 518 "the hypre package necessary for this preconditioner."));
526 # ifdef DEAL_II_PETSC_WITH_HYPRE 529 matrix =
static_cast<Mat
>(matrix_);
535 PetscErrorCode ierr = PCSetUp(
pc);
538 # else // DEAL_II_PETSC_WITH_HYPRE 540 (void)additional_data_;
542 ExcMessage(
"Your PETSc installation does not include a copy of " 543 "the hypre package necessary for this preconditioner."));
551 const unsigned int symmetric,
552 const unsigned int n_levels,
553 const double threshold,
555 const bool output_details)
556 : symmetric(symmetric)
558 , threshold(threshold)
560 , output_details(output_details)
579 matrix =
static_cast<Mat
>(matrix_);
582 # ifdef DEAL_II_PETSC_WITH_HYPRE 585 PetscErrorCode ierr = PCSetType(
pc, const_cast<char *>(PCHYPRE));
588 ierr = PCHYPRESetType(
pc,
"parasails");
599 "ParaSails parameter symmetric can only be equal to 0, 1, 2!"));
601 std::stringstream ssStream;
607 ssStream <<
"nonsymmetric";
619 ssStream <<
"nonsymmetric,SPD";
627 "ParaSails parameter symmetric can only be equal to 0, 1, 2!"));
643 ierr = PCSetFromOptions(
pc);
649 # else // DEAL_II_PETSC_WITH_HYPRE 652 ExcMessage(
"Your PETSc installation does not include a copy of " 653 "the hypre package necessary for this preconditioner."));
673 matrix =
static_cast<Mat
>(matrix_);
678 PetscErrorCode ierr = PCSetType(
pc, const_cast<char *>(PCNONE));
681 ierr = PCSetFromOptions(
pc);
692 const double zero_pivot,
693 const double damping)
695 , zero_pivot(zero_pivot)
714 matrix =
static_cast<Mat
>(matrix_);
719 PetscErrorCode ierr = PCSetType(
pc, const_cast<char *>(PCLU));
732 ierr = PCSetFromOptions(
pc);
741 DEAL_II_NAMESPACE_CLOSE
743 #endif // DEAL_II_WITH_PETSC PreconditionILU()=default
void vmult(VectorBase &dst, const VectorBase &src) const
AdditionalData(const double pivoting=1.e-6, const double zero_pivot=1.e-12, const double damping=0.0)
AdditionalData additional_data
void initialize(const MatrixBase &matrix, const AdditionalData &additional_data=AdditionalData())
void initialize(const MatrixBase &matrix, const AdditionalData &additional_data=AdditionalData())
AdditionalData(const double omega=1)
void initialize(const MatrixBase &matrix, const AdditionalData &additional_data=AdditionalData())
PreconditionSOR()=default
AdditionalData(const unsigned int levels=0)
PreconditionICC()=default
PreconditionParaSails()=default
#define AssertThrow(cond, exc)
std::string to_string(const number value, const unsigned int digits=numbers::invalid_unsigned_int)
void initialize(const MatrixBase &matrix, const AdditionalData &additional_data=AdditionalData())
void set_option_value(const std::string &name, const std::string &value)
AdditionalData additional_data
PreconditionEisenstat()=default
void initialize(const MatrixBase &matrix, const AdditionalData &additional_data=AdditionalData())
PreconditionSSOR()=default
AdditionalData additional_data
AdditionalData additional_data
PreconditionJacobi()=default
static ::ExceptionBase & ExcInvalidState()
static ::ExceptionBase & ExcMessage(std::string arg1)
AdditionalData(const unsigned int levels=0)
AdditionalData additional_data
#define Assert(cond, exc)
AdditionalData additional_data
PreconditionNone()=default
PreconditionBlockJacobi()=default
AdditionalData additional_data
void initialize(const MatrixBase &matrix, const AdditionalData &additional_data=AdditionalData())
AdditionalData additional_data
unsigned int aggressive_coarsening_num_levels
AdditionalData(const double omega=1)
AdditionalData additional_data
AdditionalData(const unsigned int symmetric=1, const unsigned int n_levels=1, const double threshold=0.1, const double filter=0.05, const bool output_details=false)
void initialize(const MatrixBase &matrix, const AdditionalData &additional_data=AdditionalData())
void initialize(const MatrixBase &matrix, const AdditionalData &additional_data=AdditionalData())
const PC & get_pc() const
AdditionalData(const double omega=1)
virtual ~PreconditionerBase()
AdditionalData(const bool symmetric_operator=false, const double strong_threshold=0.25, const double max_row_sum=0.9, const unsigned int aggressive_coarsening_num_levels=0, const bool output_details=false)
AdditionalData additional_data
PreconditionBoomerAMG()=default