Reference documentation for deal.II version 9.6.0
|
#include <deal.II/lac/petsc_precondition.h>
Public Types | |
enum class | RelaxationType { Jacobi , sequentialGaussSeidel , seqboundaryGaussSeidel , SORJacobi , backwardSORJacobi , symmetricSORJacobi , l1scaledSORJacobi , GaussianElimination , l1GaussSeidel , backwardl1GaussSeidel , CG , Chebyshev , FCFJacobi , l1scaledJacobi , None } |
Public Member Functions | |
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) | |
Public Attributes | |
bool | symmetric_operator |
double | strong_threshold |
double | max_row_sum |
unsigned int | aggressive_coarsening_num_levels |
bool | output_details |
RelaxationType | relaxation_type_up |
RelaxationType | relaxation_type_down |
RelaxationType | relaxation_type_coarse |
unsigned int | n_sweeps_coarse |
double | tol |
unsigned int | max_iter |
bool | w_cycle |
Standardized data struct to pipe additional flags to the preconditioner.
Definition at line 627 of file petsc_precondition.h.
|
strong |
Defines the available relaxation types for BoomerAMG.
Definition at line 632 of file petsc_precondition.h.
PETScWrappers::PreconditionBoomerAMG::AdditionalData::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 ) |
Constructor. Note that BoomerAMG offers a lot more options to set than what is exposed here.
Definition at line 427 of file petsc_precondition.cc.
bool PETScWrappers::PreconditionBoomerAMG::AdditionalData::symmetric_operator |
Set this flag to true if you have a symmetric system matrix and you want to use a solver which assumes a symmetric preconditioner like CG.
Definition at line 675 of file petsc_precondition.h.
double PETScWrappers::PreconditionBoomerAMG::AdditionalData::strong_threshold |
Threshold of when nodes are considered strongly connected. See HYPRE_BoomerAMGSetStrongThreshold(). Recommended values are 0.25 for 2d and 0.5 for 3d problems, but it is problem dependent.
Definition at line 682 of file petsc_precondition.h.
double PETScWrappers::PreconditionBoomerAMG::AdditionalData::max_row_sum |
If set to a value smaller than 1.0 then diagonally dominant parts of the matrix are treated as having no strongly connected nodes. If the row sum weighted by the diagonal entry is bigger than the given value, it is considered diagonally dominant. This feature is turned of by setting the value to 1.0. This is the default as some matrices can result in having only diagonally dominant entries and thus no multigrid levels are constructed. The default in BoomerAMG for this is 0.9. When you try this, check for a reasonable number of levels created.
Definition at line 695 of file petsc_precondition.h.
unsigned int PETScWrappers::PreconditionBoomerAMG::AdditionalData::aggressive_coarsening_num_levels |
Number of levels of aggressive coarsening. Increasing this value reduces the construction time and memory requirements but may decrease effectiveness.
Definition at line 702 of file petsc_precondition.h.
bool PETScWrappers::PreconditionBoomerAMG::AdditionalData::output_details |
Setting this flag to true produces debug output from HYPRE, when the preconditioner is constructed.
Definition at line 708 of file petsc_precondition.h.
RelaxationType PETScWrappers::PreconditionBoomerAMG::AdditionalData::relaxation_type_up |
Choose relaxation type up.
Definition at line 713 of file petsc_precondition.h.
RelaxationType PETScWrappers::PreconditionBoomerAMG::AdditionalData::relaxation_type_down |
Choose relaxation type down.
Definition at line 718 of file petsc_precondition.h.
RelaxationType PETScWrappers::PreconditionBoomerAMG::AdditionalData::relaxation_type_coarse |
Choose relaxation type coarse.
Definition at line 723 of file petsc_precondition.h.
unsigned int PETScWrappers::PreconditionBoomerAMG::AdditionalData::n_sweeps_coarse |
Choose number of sweeps on coarse grid.
Definition at line 728 of file petsc_precondition.h.
double PETScWrappers::PreconditionBoomerAMG::AdditionalData::tol |
Choose BommerAMG tolerance.
Definition at line 733 of file petsc_precondition.h.
unsigned int PETScWrappers::PreconditionBoomerAMG::AdditionalData::max_iter |
Choose BommerAMG maximum number of cycles.
Definition at line 738 of file petsc_precondition.h.
bool PETScWrappers::PreconditionBoomerAMG::AdditionalData::w_cycle |
Defines whether a w-cycle should be used instead of the standard setting of a v-cycle.
Definition at line 744 of file petsc_precondition.h.