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 *
>(PCEISENSTAT));
344 ierr = PCSetFromOptions(
pc);
374 matrix =
static_cast<Mat
>(matrix_);
379 PetscErrorCode ierr = PCSetType(
pc,
const_cast<char *
>(PCICC));
386 ierr = PCSetFromOptions(
pc);
415 matrix =
static_cast<Mat
>(matrix_);
420 PetscErrorCode ierr = PCSetType(
pc,
const_cast<char *
>(PCILU));
427 ierr = PCSetFromOptions(
pc);
438 const bool symmetric_operator,
439 const double strong_threshold,
440 const double max_row_sum,
441 const unsigned int aggressive_coarsening_num_levels,
442 const bool output_details,
446 const unsigned int n_sweeps_coarse,
448 const unsigned int max_iter,
450 : symmetric_operator(symmetric_operator)
451 , strong_threshold(strong_threshold)
452 , max_row_sum(max_row_sum)
453 , aggressive_coarsening_num_levels(aggressive_coarsening_num_levels)
454 , output_details(output_details)
455 , relaxation_type_up(relaxation_type_up)
456 , relaxation_type_down(relaxation_type_down)
457 , relaxation_type_coarse(relaxation_type_coarse)
458 , n_sweeps_coarse(n_sweeps_coarse)
466# ifdef DEAL_II_PETSC_WITH_HYPRE
477 std::string string_type;
479 switch (relaxation_type)
482 string_type =
"Jacobi";
486 string_type =
"sequential-Gauss-Seidel";
490 string_type =
"seqboundary-Gauss-Seidel";
493 string_type =
"SOR/Jacobi";
497 string_type =
"backward-SOR/Jacobi";
501 string_type =
"symmetric-SOR/Jacobi";
505 string_type =
" l1scaled-SOR/Jacobi";
509 string_type =
"Gaussian-elimination";
513 string_type =
"l1-Gauss-Seidel";
517 string_type =
"backward-l1-Gauss-Seidel";
523 string_type =
"Chebyshev";
526 string_type =
"FCF-Jacobi";
530 string_type =
"l1scaled-Jacobi";
533 string_type =
"None";
549 PetscErrorCode ierr = PCCreate(
comm, &
pc);
552# ifdef DEAL_II_PETSC_WITH_HYPRE
557 ExcMessage(
"Your PETSc installation does not include a copy of "
558 "the hypre package necessary for this preconditioner."));
573# ifdef DEAL_II_PETSC_WITH_HYPRE
574 PetscErrorCode ierr = PCSetType(
pc,
const_cast<char *
>(PCHYPRE));
577 ierr = PCHYPRESetType(
pc,
"boomeramg");
621 auto relaxation_type_is_symmetric =
638 ExcMessage(
"Use a symmetric smoother for relaxation_type_up"));
643 ExcMessage(
"Use a symmetric smoother for relaxation_type_down"));
648 ExcMessage(
"Use a symmetric smoother for relaxation_type_coarse"));
669 ierr = PCSetFromOptions(
pc);
673 ExcMessage(
"Your PETSc installation does not include a copy of "
674 "the hypre package necessary for this preconditioner."));
682# ifdef DEAL_II_PETSC_WITH_HYPRE
685 matrix =
static_cast<Mat
>(matrix_);
691 PetscErrorCode ierr = PCSetUp(
pc);
696 (void)additional_data_;
698 ExcMessage(
"Your PETSc installation does not include a copy of "
699 "the hypre package necessary for this preconditioner."));
708 const unsigned int n_levels,
709 const double threshold,
711 const bool output_details)
714 , threshold(threshold)
716 , output_details(output_details)
735 matrix =
static_cast<Mat
>(matrix_);
738# ifdef DEAL_II_PETSC_WITH_HYPRE
741 PetscErrorCode ierr = PCSetType(
pc,
const_cast<char *
>(PCHYPRE));
744 ierr = PCHYPRESetType(
pc,
"parasails");
755 "ParaSails parameter symmetric can only be equal to 0, 1, 2!"));
757 std::stringstream ssStream;
763 ssStream <<
"nonsymmetric";
775 ssStream <<
"nonsymmetric,SPD";
783 "ParaSails parameter symmetric can only be equal to 0, 1, 2!"));
799 ierr = PCSetFromOptions(
pc);
808 ExcMessage(
"Your PETSc installation does not include a copy of "
809 "the hypre package necessary for this preconditioner."));
829 matrix =
static_cast<Mat
>(matrix_);
834 PetscErrorCode ierr = PCSetType(
pc,
const_cast<char *
>(PCNONE));
837 ierr = PCSetFromOptions(
pc);
848 const double zero_pivot,
849 const double damping)
851 , zero_pivot(zero_pivot)
870 matrix =
static_cast<Mat
>(matrix_);
875 PetscErrorCode ierr = PCSetType(
pc,
const_cast<char *
>(PCLU));
888 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
PreconditionEisenstat()=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())
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)
std::string to_string(const T &t)
static ::ExceptionBase & ExcInvalidState()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
@ matrix
Contents is actually a matrix.
@ symmetric
Matrix is symmetric.
void set_option_value(const std::string &name, const std::string &value)
*braid_SplitCommworld & comm
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 double omega=1)
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)