|
Reference documentation for deal.II version 9.2.0
|
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\)
\(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\)
\(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\)
\(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Go to the documentation of this file.
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 PreconditionerBase::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)
443 : symmetric_operator(symmetric_operator)
444 , strong_threshold(strong_threshold)
445 , max_row_sum(max_row_sum)
446 , aggressive_coarsening_num_levels(aggressive_coarsening_num_levels)
447 , output_details(output_details)
458 PetscErrorCode ierr = PCCreate(comm, &
pc);
461 # ifdef DEAL_II_PETSC_WITH_HYPRE
463 # else // DEAL_II_PETSC_WITH_HYPRE
466 ExcMessage(
"Your PETSc installation does not include a copy of "
467 "the hypre package necessary for this preconditioner."));
482 # ifdef DEAL_II_PETSC_WITH_HYPRE
483 PetscErrorCode ierr = PCSetType(
pc,
const_cast<char *
>(PCHYPRE));
486 ierr = PCHYPRESetType(
pc,
"boomeramg");
498 std::stringstream ssStream;
509 "symmetric-SOR/Jacobi");
511 "symmetric-SOR/Jacobi");
513 "Gaussian-elimination");
520 "Gaussian-elimination");
523 ierr = PCSetFromOptions(
pc);
527 ExcMessage(
"Your PETSc installation does not include a copy of "
528 "the hypre package necessary for this preconditioner."));
536 # ifdef DEAL_II_PETSC_WITH_HYPRE
539 matrix =
static_cast<Mat
>(matrix_);
545 PetscErrorCode ierr = PCSetUp(
pc);
548 # else // DEAL_II_PETSC_WITH_HYPRE
550 (void)additional_data_;
552 ExcMessage(
"Your PETSc installation does not include a copy of "
553 "the hypre package necessary for this preconditioner."));
562 const unsigned int n_levels,
563 const double threshold,
565 const bool output_details)
568 , threshold(threshold)
570 , output_details(output_details)
589 matrix =
static_cast<Mat
>(matrix_);
592 # ifdef DEAL_II_PETSC_WITH_HYPRE
595 PetscErrorCode ierr = PCSetType(
pc,
const_cast<char *
>(PCHYPRE));
598 ierr = PCHYPRESetType(
pc,
"parasails");
609 "ParaSails parameter symmetric can only be equal to 0, 1, 2!"));
611 std::stringstream ssStream;
617 ssStream <<
"nonsymmetric";
629 ssStream <<
"nonsymmetric,SPD";
637 "ParaSails parameter symmetric can only be equal to 0, 1, 2!"));
653 ierr = PCSetFromOptions(
pc);
659 # else // DEAL_II_PETSC_WITH_HYPRE
662 ExcMessage(
"Your PETSc installation does not include a copy of "
663 "the hypre package necessary for this preconditioner."));
683 matrix =
static_cast<Mat
>(matrix_);
688 PetscErrorCode ierr = PCSetType(
pc,
const_cast<char *
>(PCNONE));
691 ierr = PCSetFromOptions(
pc);
702 const double zero_pivot,
703 const double damping)
705 , zero_pivot(zero_pivot)
724 matrix =
static_cast<Mat
>(matrix_);
729 PetscErrorCode ierr = PCSetType(
pc,
const_cast<char *
>(PCLU));
742 ierr = PCSetFromOptions(
pc);
753 #endif // DEAL_II_WITH_PETSC
AdditionalData(const double omega=1)
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()=default
void vmult(VectorBase &dst, const VectorBase &src) const
AdditionalData additional_data
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())
PreconditionILU()=default
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
AdditionalData(const unsigned int levels=0)
void initialize(const MatrixBase &matrix, const AdditionalData &additional_data=AdditionalData())
void set_option_value(const std::string &name, const std::string &value)
static ::ExceptionBase & ExcInvalidState()
void initialize(const MatrixBase &matrix, const AdditionalData &additional_data=AdditionalData())
PreconditionSOR()=default
AdditionalData(const double omega=1)
std::string to_string(const T &t)
void initialize(const MatrixBase &matrix, const AdditionalData &additional_data=AdditionalData())
PreconditionEisenstat()=default
PreconditionNone()=default
AdditionalData(const double pivoting=1.e-6, const double zero_pivot=1.e-12, const double damping=0.0)
AdditionalData additional_data
PreconditionICC()=default
unsigned int aggressive_coarsening_num_levels
PreconditionSSOR()=default
@ symmetric
Matrix is symmetric.
PreconditionParaSails()=default
static ::ExceptionBase & ExcMessage(std::string arg1)
#define DEAL_II_NAMESPACE_OPEN
AdditionalData additional_data
AdditionalData additional_data
@ matrix
Contents is actually a matrix.
PreconditionJacobi()=default
AdditionalData additional_data
AdditionalData additional_data
void initialize(const MatrixBase &matrix, const AdditionalData &additional_data=AdditionalData())
#define Assert(cond, exc)
AdditionalData additional_data
virtual ~PreconditionerBase()
void Tvmult(VectorBase &dst, const VectorBase &src) const
const PC & get_pc() const
#define DEAL_II_NAMESPACE_CLOSE
AdditionalData(const unsigned int levels=0)
#define AssertThrow(cond, exc)
void initialize(const MatrixBase &matrix, const AdditionalData &additional_data=AdditionalData())
PreconditionBoomerAMG()=default
AdditionalData additional_data