Reference documentation for deal.II version 9.1.1
|
#include <deal.II/lac/trilinos_precondition.h>
Public Member Functions | |
AdditionalData (const bool elliptic=true, const unsigned int n_cycles=1, const bool w_cyle=false, const double aggregation_threshold=1e-4, const std::vector< std::vector< bool >> &constant_modes=std::vector< std::vector< bool >>(0), const unsigned int smoother_sweeps=2, const unsigned int smoother_overlap=0, const bool output_details=false, const char *smoother_type="Chebyshev", const char *coarse_type="Amesos-KLU") | |
Public Attributes | |
bool | elliptic |
unsigned int | n_cycles |
bool | w_cycle |
double | aggregation_threshold |
std::vector< std::vector< bool > > | constant_modes |
unsigned int | smoother_sweeps |
unsigned int | smoother_overlap |
bool | output_details |
const char * | smoother_type |
const char * | coarse_type |
A data structure that is used to control details of how the algebraic multigrid is set up. The flags detailed in here are then passed to the Trilinos MueLu implementation. A structure of the current type are passed to the constructor of PreconditionAMGMueLu.
Definition at line 1743 of file trilinos_precondition.h.
TrilinosWrappers::PreconditionAMGMueLu::AdditionalData::AdditionalData | ( | const bool | elliptic = true , |
const unsigned int | n_cycles = 1 , |
||
const bool | w_cyle = false , |
||
const double | aggregation_threshold = 1e-4 , |
||
const std::vector< std::vector< bool >> & | constant_modes = std::vector<std::vector<bool>>(0) , |
||
const unsigned int | smoother_sweeps = 2 , |
||
const unsigned int | smoother_overlap = 0 , |
||
const bool | output_details = false , |
||
const char * | smoother_type = "Chebyshev" , |
||
const char * | coarse_type = "Amesos-KLU" |
||
) |
Constructor. By default, we pretend to work on elliptic problems with linear finite elements on a scalar equation.
Definition at line 31 of file trilinos_precondition_muelu.cc.
bool TrilinosWrappers::PreconditionAMGMueLu::AdditionalData::elliptic |
Determines whether the AMG preconditioner should be optimized for elliptic problems (MueLu option smoothed aggregation SA, using a Chebyshev smoother) or for non-elliptic problems (MueLu option non- symmetric smoothed aggregation NSSA, smoother is SSOR with underrelaxation).
Definition at line 1768 of file trilinos_precondition.h.
unsigned int TrilinosWrappers::PreconditionAMGMueLu::AdditionalData::n_cycles |
Defines how many multigrid cycles should be performed by the preconditioner.
Definition at line 1774 of file trilinos_precondition.h.
bool TrilinosWrappers::PreconditionAMGMueLu::AdditionalData::w_cycle |
Defines whether a w-cycle should be used instead of the standard setting of a v-cycle.
Definition at line 1780 of file trilinos_precondition.h.
double TrilinosWrappers::PreconditionAMGMueLu::AdditionalData::aggregation_threshold |
This threshold tells the AMG setup how the coarsening should be performed. In the AMG used by MueLu, all points that strongly couple with the tentative coarse-level point form one aggregate. The term strong coupling is controlled by the variable aggregation_threshold
, meaning that all elements that are not smaller than aggregation_threshold
times the diagonal element do couple strongly.
Definition at line 1791 of file trilinos_precondition.h.
std::vector<std::vector<bool> > TrilinosWrappers::PreconditionAMGMueLu::AdditionalData::constant_modes |
Specifies the constant modes (near null space) of the matrix. This parameter tells AMG whether we work on a scalar equation (where the near null space only consists of ones) or on a vector-valued equation.
Definition at line 1799 of file trilinos_precondition.h.
unsigned int TrilinosWrappers::PreconditionAMGMueLu::AdditionalData::smoother_sweeps |
Determines how many sweeps of the smoother should be performed. When the flag elliptic
is set to true
, i.e., for elliptic or almost elliptic problems, the polynomial degree of the Chebyshev smoother is set to smoother_sweeps
. The term sweeps refers to the number of matrix-vector products performed in the Chebyshev case. In the non-elliptic case, smoother_sweeps
sets the number of SSOR relaxation sweeps for post-smoothing to be performed.
Definition at line 1811 of file trilinos_precondition.h.
unsigned int TrilinosWrappers::PreconditionAMGMueLu::AdditionalData::smoother_overlap |
Determines the overlap in the SSOR/Chebyshev error smoother when run in parallel.
Definition at line 1817 of file trilinos_precondition.h.
bool TrilinosWrappers::PreconditionAMGMueLu::AdditionalData::output_details |
If this flag is set to true
, then internal information from the ML preconditioner is printed to screen. This can be useful when debugging the preconditioner.
Definition at line 1824 of file trilinos_precondition.h.
const char* TrilinosWrappers::PreconditionAMGMueLu::AdditionalData::smoother_type |
Determines which smoother to use for the AMG cycle. Possibilities for smoother_type are the following:
Definition at line 1860 of file trilinos_precondition.h.
const char* TrilinosWrappers::PreconditionAMGMueLu::AdditionalData::coarse_type |
Determines which solver to use on the coarsest level. The same settings as for the smoother type are possible.
Definition at line 1866 of file trilinos_precondition.h.