18#ifdef DEAL_II_WITH_PETSC
28# include <petscconf.h>
58 PetscErrorCode ierr = PCDestroy(&
pc);
70 const PetscErrorCode ierr = PCApply(
pc, src, dst);
80 const PetscErrorCode ierr = PCApplyTranspose(
pc, src, dst);
97 PetscObjectGetComm(
reinterpret_cast<PetscObject
>(
matrix), &
comm);
100 ierr = PCCreate(
comm, &
pc);
103# if DEAL_II_PETSC_VERSION_LT(3, 5, 0)
119 PreconditionBase::operator Mat()
const
131 PetscErrorCode ierr = PCCreate(
comm, &
pc);
150 PetscErrorCode ierr = PCSetType(
pc,
const_cast<char *
>(PCJACOBI));
153 ierr = PCSetFromOptions(
pc);
163 matrix =
static_cast<Mat
>(matrix_);
169 PetscErrorCode ierr = PCSetUp(
pc);
181 PetscErrorCode ierr = PCCreate(
comm, &
pc);
199 PetscErrorCode ierr = PCSetType(
pc,
const_cast<char *
>(PCBJACOBI));
202 ierr = PCSetFromOptions(
pc);
213 matrix =
static_cast<Mat
>(matrix_);
219 PetscErrorCode ierr = PCSetUp(
pc);
245 matrix =
static_cast<Mat
>(matrix_);
250 PetscErrorCode ierr = PCSetType(
pc,
const_cast<char *
>(PCSOR));
257 ierr = PCSetFromOptions(
pc);
286 matrix =
static_cast<Mat
>(matrix_);
291 PetscErrorCode ierr = PCSetType(
pc,
const_cast<char *
>(PCSOR));
299 ierr = PCSORSetSymmetric(
pc, SOR_SYMMETRIC_SWEEP);
302 ierr = PCSetFromOptions(
pc);
332 matrix =
static_cast<Mat
>(matrix_);
337 PetscErrorCode ierr = PCSetType(
pc,
const_cast<char *
>(PCICC));
344 ierr = PCSetFromOptions(
pc);
373 matrix =
static_cast<Mat
>(matrix_);
378 PetscErrorCode ierr = PCSetType(
pc,
const_cast<char *
>(PCILU));
385 ierr = PCSetFromOptions(
pc);
396 const bool symmetric_operator,
397 const double strong_threshold,
398 const double max_row_sum,
399 const unsigned int aggressive_coarsening_num_levels,
400 const bool output_details,
404 const unsigned int n_sweeps_coarse,
406 const unsigned int max_iter,
408 : symmetric_operator(symmetric_operator)
409 , strong_threshold(strong_threshold)
410 , max_row_sum(max_row_sum)
411 , aggressive_coarsening_num_levels(aggressive_coarsening_num_levels)
412 , output_details(output_details)
413 , relaxation_type_up(relaxation_type_up)
414 , relaxation_type_down(relaxation_type_down)
415 , relaxation_type_coarse(relaxation_type_coarse)
416 , n_sweeps_coarse(n_sweeps_coarse)
424# ifdef DEAL_II_PETSC_WITH_HYPRE
435 std::string string_type;
437 switch (relaxation_type)
440 string_type =
"Jacobi";
444 string_type =
"sequential-Gauss-Seidel";
448 string_type =
"seqboundary-Gauss-Seidel";
451 string_type =
"SOR/Jacobi";
455 string_type =
"backward-SOR/Jacobi";
459 string_type =
"symmetric-SOR/Jacobi";
463 string_type =
" l1scaled-SOR/Jacobi";
467 string_type =
"Gaussian-elimination";
471 string_type =
"l1-Gauss-Seidel";
475 string_type =
"backward-l1-Gauss-Seidel";
481 string_type =
"Chebyshev";
484 string_type =
"FCF-Jacobi";
488 string_type =
"l1scaled-Jacobi";
491 string_type =
"None";
507 PetscErrorCode ierr = PCCreate(
comm, &
pc);
510# ifdef DEAL_II_PETSC_WITH_HYPRE
515 ExcMessage(
"Your PETSc installation does not include a copy of "
516 "the hypre package necessary for this preconditioner."));
531# ifdef DEAL_II_PETSC_WITH_HYPRE
532 PetscErrorCode ierr = PCSetType(
pc,
const_cast<char *
>(PCHYPRE));
535 ierr = PCHYPRESetType(
pc,
"boomeramg");
579 auto relaxation_type_is_symmetric =
596 ExcMessage(
"Use a symmetric smoother for relaxation_type_up"));
601 ExcMessage(
"Use a symmetric smoother for relaxation_type_down"));
606 ExcMessage(
"Use a symmetric smoother for relaxation_type_coarse"));
627 ierr = PCSetFromOptions(
pc);
631 ExcMessage(
"Your PETSc installation does not include a copy of "
632 "the hypre package necessary for this preconditioner."));
640# ifdef DEAL_II_PETSC_WITH_HYPRE
643 matrix =
static_cast<Mat
>(matrix_);
649 PetscErrorCode ierr = PCSetUp(
pc);
654 (void)additional_data_;
656 ExcMessage(
"Your PETSc installation does not include a copy of "
657 "the hypre package necessary for this preconditioner."));
665 const unsigned int symmetric,
666 const unsigned int n_levels,
667 const double threshold,
669 const bool output_details)
670 : symmetric(symmetric)
672 , threshold(threshold)
674 , output_details(output_details)
693 matrix =
static_cast<Mat
>(matrix_);
696# ifdef DEAL_II_PETSC_WITH_HYPRE
699 PetscErrorCode ierr = PCSetType(
pc,
const_cast<char *
>(PCHYPRE));
702 ierr = PCHYPRESetType(
pc,
"parasails");
713 "ParaSails parameter symmetric can only be equal to 0, 1, 2!"));
715 std::stringstream ssStream;
721 ssStream <<
"nonsymmetric";
733 ssStream <<
"nonsymmetric,SPD";
741 "ParaSails parameter symmetric can only be equal to 0, 1, 2!"));
757 ierr = PCSetFromOptions(
pc);
766 ExcMessage(
"Your PETSc installation does not include a copy of "
767 "the hypre package necessary for this preconditioner."));
787 matrix =
static_cast<Mat
>(matrix_);
792 PetscErrorCode ierr = PCSetType(
pc,
const_cast<char *
>(PCNONE));
795 ierr = PCSetFromOptions(
pc);
806 const double zero_pivot,
807 const double damping)
809 , zero_pivot(zero_pivot)
828 matrix =
static_cast<Mat
>(matrix_);
833 PetscErrorCode ierr = PCSetType(
pc,
const_cast<char *
>(PCLU));
846 ierr = PCSetFromOptions(
pc);
void Tvmult(VectorBase &dst, const VectorBase &src) const
const PC & get_pc() const
virtual ~PreconditionBase()
void vmult(VectorBase &dst, const VectorBase &src) const
PreconditionBlockJacobi()=default
AdditionalData additional_data
PreconditionBoomerAMG()=default
AdditionalData additional_data
AdditionalData additional_data
void initialize(const MatrixBase &matrix, const AdditionalData &additional_data=AdditionalData())
PreconditionICC()=default
AdditionalData additional_data
PreconditionILU()=default
void initialize(const MatrixBase &matrix, const AdditionalData &additional_data=AdditionalData())
PreconditionJacobi()=default
AdditionalData additional_data
void initialize(const MatrixBase &matrix, const AdditionalData &additional_data=AdditionalData())
AdditionalData additional_data
void initialize(const MatrixBase &matrix, const AdditionalData &additional_data=AdditionalData())
PreconditionNone()=default
void initialize(const MatrixBase &matrix, const AdditionalData &additional_data=AdditionalData())
PreconditionParaSails()=default
AdditionalData additional_data
AdditionalData additional_data
PreconditionSOR()=default
void initialize(const MatrixBase &matrix, const AdditionalData &additional_data=AdditionalData())
PreconditionSSOR()=default
AdditionalData additional_data
void initialize(const MatrixBase &matrix, const AdditionalData &additional_data=AdditionalData())
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
static ::ExceptionBase & ExcInvalidState()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
void set_option_value(const std::string &name, const std::string &value)
RelaxationType relaxation_type_up
unsigned int n_sweeps_coarse
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, const RelaxationType relaxation_type_up=RelaxationType::SORJacobi, const RelaxationType relaxation_type_down=RelaxationType::SORJacobi, const RelaxationType relaxation_type_coarse=RelaxationType::GaussianElimination, const unsigned int n_sweeps_coarse=1, const double tol=0.0, const unsigned int max_iter=1, const bool w_cycle=false)
unsigned int aggressive_coarsening_num_levels
RelaxationType relaxation_type_down
RelaxationType relaxation_type_coarse
AdditionalData(const unsigned int levels=0)
AdditionalData(const unsigned int levels=0)
AdditionalData(const double pivoting=1.e-6, const double zero_pivot=1.e-12, const double damping=0.0)
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)
AdditionalData(const double omega=1)
AdditionalData(const double omega=1)