16 #include <deal.II/lac/petsc_precondition.h> 18 #ifdef DEAL_II_WITH_PETSC 20 # include <deal.II/base/utilities.h> 21 # include <deal.II/lac/petsc_compatibility.h> 22 # include <deal.II/lac/petsc_matrix_base.h> 23 # include <deal.II/lac/petsc_vector_base.h> 24 # include <deal.II/lac/petsc_solver.h> 25 # include <petscconf.h> 28 DEAL_II_NAMESPACE_OPEN
34 pc(NULL), matrix(NULL)
42 #if DEAL_II_PETSC_VERSION_LT(3,2,0) 43 int ierr = PCDestroy(
pc);
45 int ierr = PCDestroy(&
pc);
56 AssertThrow (
pc != NULL, StandardExceptions::ExcInvalidState ());
59 ierr = PCApply(
pc, src, dst);
69 AssertThrow (
pc == NULL, StandardExceptions::ExcInvalidState ());
76 ierr = PetscObjectGetComm(reinterpret_cast<PetscObject>(
matrix), &comm);
79 ierr = PCCreate(comm, &
pc);
82 #if DEAL_II_PETSC_VERSION_LT(3, 5, 0) 98 PreconditionerBase::operator Mat ()
const 108 additional_data = additional_data_;
110 int ierr = PCCreate(comm, &
pc);
124 initialize(matrix, additional_data);
131 ierr = PCSetType (
pc, const_cast<char *>(PCJACOBI));
134 ierr = PCSetFromOptions (
pc);
142 matrix =
static_cast<Mat
>(matrix_);
143 additional_data = additional_data_;
148 int ierr = PCSetUp (
pc);
157 additional_data = additional_data_;
159 int ierr = PCCreate(comm, &
pc);
174 initialize(matrix, additional_data);
181 ierr = PCSetType (
pc, const_cast<char *>(PCBJACOBI));
184 ierr = PCSetFromOptions (
pc);
193 matrix =
static_cast<Mat
>(matrix_);
194 additional_data = additional_data_;
199 int ierr = PCSetUp (
pc);
228 matrix =
static_cast<Mat
>(matrix_);
234 ierr = PCSetType (
pc, const_cast<char *>(PCSOR));
241 ierr = PCSetFromOptions (
pc);
273 matrix =
static_cast<Mat
>(matrix_);
279 ierr = PCSetType (
pc, const_cast<char *>(PCSOR));
287 ierr = PCSORSetSymmetric (
pc, SOR_SYMMETRIC_SWEEP);
290 ierr = PCSetFromOptions (
pc);
322 matrix =
static_cast<Mat
>(matrix_);
328 ierr = PCSetType (
pc, const_cast<char *>(PCEISENSTAT));
335 ierr = PCSetFromOptions (
pc);
368 matrix =
static_cast<Mat
>(matrix_);
374 ierr = PCSetType (
pc, const_cast<char *>(PCICC));
381 ierr = PCSetFromOptions (
pc);
413 matrix =
static_cast<Mat
>(matrix_);
419 ierr = PCSetType (
pc, const_cast<char *>(PCILU));
426 ierr = PCSetFromOptions (
pc);
438 const double strong_threshold,
439 const double max_row_sum,
440 const unsigned int aggressive_coarsening_num_levels,
441 const bool output_details
444 symmetric_operator(symmetric_operator),
445 strong_threshold(strong_threshold),
446 max_row_sum(max_row_sum),
447 aggressive_coarsening_num_levels(aggressive_coarsening_num_levels),
448 output_details(output_details)
460 int ierr = PCCreate(comm, &
pc);
463 #ifdef PETSC_HAVE_HYPRE 465 #else // PETSC_HAVE_HYPRE 468 ExcMessage (
"Your PETSc installation does not include a copy of " 469 "the hypre package necessary for this preconditioner."));
484 ierr = PCSetType (
pc, const_cast<char *>(PCHYPRE));
487 ierr = PCHYPRESetType(
pc,
"boomeramg");
498 std::stringstream ssStream;
508 set_option_value(
"-pc_hypre_boomeramg_relax_type_up",
"symmetric-SOR/Jacobi");
509 set_option_value(
"-pc_hypre_boomeramg_relax_type_down",
"symmetric-SOR/Jacobi");
510 set_option_value(
"-pc_hypre_boomeramg_relax_type_coarse",
"Gaussian-elimination");
516 set_option_value(
"-pc_hypre_boomeramg_relax_type_coarse",
"Gaussian-elimination");
519 ierr = PCSetFromOptions (
pc);
527 matrix =
static_cast<Mat
>(matrix_);
530 #ifdef PETSC_HAVE_HYPRE 534 int ierr = PCSetUp (
pc);
537 #else // PETSC_HAVE_HYPRE 540 ExcMessage (
"Your PETSc installation does not include a copy of " 541 "the hypre package necessary for this preconditioner."));
550 const unsigned int n_levels,
551 const double threshold,
555 symmetric(symmetric),
557 threshold(threshold),
559 output_details(output_details)
578 matrix =
static_cast<Mat
>(matrix_);
581 #ifdef PETSC_HAVE_HYPRE 585 ierr = PCSetType (
pc, const_cast<char *>(PCHYPRE));
588 ierr = PCHYPRESetType(
pc,
"parasails");
599 ExcMessage(
"ParaSails parameter symmetric can only be equal to 0, 1, 2!"));
601 std::stringstream ssStream;
607 ssStream <<
"nonsymmetric";
619 ssStream <<
"nonsymmetric,SPD";
625 ExcMessage(
"ParaSails parameter symmetric can only be equal to 0, 1, 2!"));
641 ierr = PCSetFromOptions (
pc);
647 #else // PETSC_HAVE_HYPRE 650 ExcMessage (
"Your PETSc installation does not include a copy of " 651 "the hypre package necessary for this preconditioner."));
673 matrix =
static_cast<Mat
>(matrix_);
679 ierr = PCSetType (
pc, const_cast<char *>(PCNONE));
682 ierr = PCSetFromOptions (
pc);
694 const double zero_pivot,
695 const double damping)
698 zero_pivot (zero_pivot),
718 matrix =
static_cast<Mat
>(matrix_);
724 ierr = PCSetType (
pc, const_cast<char *>(PCLU));
728 #if DEAL_II_PETSC_VERSION_LT(3,0,1) 738 #if DEAL_II_PETSC_VERSION_LT(3,0,1) 745 ierr = PCSetFromOptions (
pc);
754 DEAL_II_NAMESPACE_CLOSE
756 #endif // DEAL_II_WITH_PETSC
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)
PreconditionBlockJacobi()
::ExceptionBase & ExcMessage(std::string arg1)
void initialize(const MatrixBase &matrix, const AdditionalData &additional_data=AdditionalData())
AdditionalData(const unsigned int levels=0)
#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
void initialize(const MatrixBase &matrix, const AdditionalData &additional_data=AdditionalData())
AdditionalData additional_data
AdditionalData additional_data
AdditionalData(const unsigned int levels=0)
AdditionalData additional_data
#define Assert(cond, exc)
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