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/exceptions.h> 22 # include <deal.II/lac/petsc_compatibility.h> 23 # include <deal.II/lac/petsc_matrix_base.h> 24 # include <deal.II/lac/petsc_vector_base.h> 25 # include <deal.II/lac/petsc_solver.h> 26 # include <petscconf.h> 29 DEAL_II_NAMESPACE_OPEN
35 pc(NULL), matrix(NULL)
55 #if DEAL_II_PETSC_VERSION_LT(3,2,0) 56 PetscErrorCode ierr = PCDestroy(
pc);
58 PetscErrorCode ierr = PCDestroy(&
pc);
72 const PetscErrorCode ierr = PCApply(
pc, src, dst);
88 PetscErrorCode ierr = PetscObjectGetComm(reinterpret_cast<PetscObject>(
matrix),
92 ierr = PCCreate(comm, &
pc);
95 #if DEAL_II_PETSC_VERSION_LT(3, 5, 0) 111 PreconditionerBase::operator Mat ()
const 123 PetscErrorCode ierr = PCCreate(comm, &
pc);
145 PetscErrorCode ierr = PCSetType (
pc, const_cast<char *>(PCJACOBI));
148 ierr = PCSetFromOptions (
pc);
158 matrix =
static_cast<Mat
>(matrix_);
164 PetscErrorCode ierr = PCSetUp (
pc);
175 PetscErrorCode ierr = PCCreate(comm, &
pc);
196 PetscErrorCode ierr = PCSetType (
pc, const_cast<char *>(PCBJACOBI));
199 ierr = PCSetFromOptions (
pc);
210 matrix =
static_cast<Mat
>(matrix_);
216 PetscErrorCode ierr = PCSetUp (
pc);
247 matrix =
static_cast<Mat
>(matrix_);
252 PetscErrorCode ierr = PCSetType (
pc, const_cast<char *>(PCSOR));
259 ierr = PCSetFromOptions (
pc);
293 matrix =
static_cast<Mat
>(matrix_);
298 PetscErrorCode ierr = PCSetType (
pc, const_cast<char *>(PCSOR));
306 ierr = PCSORSetSymmetric (
pc, SOR_SYMMETRIC_SWEEP);
309 ierr = PCSetFromOptions (
pc);
343 matrix =
static_cast<Mat
>(matrix_);
348 PetscErrorCode ierr = PCSetType (
pc, const_cast<char *>(PCEISENSTAT));
355 ierr = PCSetFromOptions (
pc);
390 matrix =
static_cast<Mat
>(matrix_);
395 PetscErrorCode ierr = PCSetType (
pc, const_cast<char *>(PCICC));
402 ierr = PCSetFromOptions (
pc);
436 matrix =
static_cast<Mat
>(matrix_);
441 PetscErrorCode ierr = PCSetType (
pc, const_cast<char *>(PCILU));
448 ierr = PCSetFromOptions (
pc);
460 const double strong_threshold,
461 const double max_row_sum,
462 const unsigned int aggressive_coarsening_num_levels,
463 const bool output_details
466 symmetric_operator(symmetric_operator),
467 strong_threshold(strong_threshold),
468 max_row_sum(max_row_sum),
469 aggressive_coarsening_num_levels(aggressive_coarsening_num_levels),
470 output_details(output_details)
482 PetscErrorCode ierr = PCCreate(comm, &
pc);
485 #ifdef PETSC_HAVE_HYPRE 487 #else // PETSC_HAVE_HYPRE 490 ExcMessage (
"Your PETSc installation does not include a copy of " 491 "the hypre package necessary for this preconditioner."));
505 #ifdef PETSC_HAVE_HYPRE 506 PetscErrorCode ierr = PCSetType (
pc, const_cast<char *>(PCHYPRE));
509 ierr = PCHYPRESetType(
pc,
"boomeramg");
520 std::stringstream ssStream;
530 set_option_value(
"-pc_hypre_boomeramg_relax_type_up",
"symmetric-SOR/Jacobi");
531 set_option_value(
"-pc_hypre_boomeramg_relax_type_down",
"symmetric-SOR/Jacobi");
532 set_option_value(
"-pc_hypre_boomeramg_relax_type_coarse",
"Gaussian-elimination");
538 set_option_value(
"-pc_hypre_boomeramg_relax_type_coarse",
"Gaussian-elimination");
541 ierr = PCSetFromOptions (
pc);
545 ExcMessage (
"Your PETSc installation does not include a copy of " 546 "the hypre package necessary for this preconditioner."));
554 #ifdef PETSC_HAVE_HYPRE 557 matrix =
static_cast<Mat
>(matrix_);
563 PetscErrorCode ierr = PCSetUp (
pc);
566 #else // PETSC_HAVE_HYPRE 568 (void)additional_data_;
570 ExcMessage (
"Your PETSc installation does not include a copy of " 571 "the hypre package necessary for this preconditioner."));
580 const unsigned int n_levels,
581 const double threshold,
583 const bool output_details)
585 symmetric(symmetric),
587 threshold(threshold),
589 output_details(output_details)
610 matrix =
static_cast<Mat
>(matrix_);
613 #ifdef PETSC_HAVE_HYPRE 616 PetscErrorCode ierr = PCSetType (
pc, const_cast<char *>(PCHYPRE));
619 ierr = PCHYPRESetType(
pc,
"parasails");
630 ExcMessage(
"ParaSails parameter symmetric can only be equal to 0, 1, 2!"));
632 std::stringstream ssStream;
638 ssStream <<
"nonsymmetric";
650 ssStream <<
"nonsymmetric,SPD";
656 ExcMessage(
"ParaSails parameter symmetric can only be equal to 0, 1, 2!"));
672 ierr = PCSetFromOptions (
pc);
678 #else // PETSC_HAVE_HYPRE 681 ExcMessage (
"Your PETSc installation does not include a copy of " 682 "the hypre package necessary for this preconditioner."));
706 matrix =
static_cast<Mat
>(matrix_);
711 PetscErrorCode ierr = PCSetType (
pc, const_cast<char *>(PCNONE));
714 ierr = PCSetFromOptions (
pc);
726 const double zero_pivot,
727 const double damping)
730 zero_pivot (zero_pivot),
752 matrix =
static_cast<Mat
>(matrix_);
757 PetscErrorCode ierr = PCSetType (
pc, const_cast<char *>(PCLU));
761 #if DEAL_II_PETSC_VERSION_LT(3,0,1) 771 #if DEAL_II_PETSC_VERSION_LT(3,0,1) 778 ierr = PCSetFromOptions (
pc);
787 DEAL_II_NAMESPACE_CLOSE
789 #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()
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
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
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