deal.II version GIT relicensing-2165-gc91f007519 2024-11-20 01:40:00+00:00
|
#include <deal.II/lac/trilinos_precondition.h>
Public Member Functions | |
AdditionalData (const bool elliptic=true, const bool higher_order_elements=false, const unsigned int n_cycles=1, const bool w_cycle=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") | |
void | set_parameters (Teuchos::ParameterList ¶meter_list, std::unique_ptr< Epetra_MultiVector > &distributed_constant_modes, const Epetra_RowMatrix &matrix) const |
void | set_parameters (Teuchos::ParameterList ¶meter_list, std::unique_ptr< Epetra_MultiVector > &distributed_constant_modes, const SparseMatrix &matrix) const |
void | set_operator_null_space (Teuchos::ParameterList ¶meter_list, std::unique_ptr< Epetra_MultiVector > &distributed_constant_modes, const Epetra_RowMatrix &matrix) const |
void | set_operator_null_space (Teuchos::ParameterList ¶meter_list, std::unique_ptr< Epetra_MultiVector > &distributed_constant_modes, const SparseMatrix &matrix) const |
Public Attributes | |
bool | elliptic |
bool | higher_order_elements |
unsigned int | n_cycles |
bool | w_cycle |
double | aggregation_threshold |
std::vector< std::vector< bool > > | constant_modes |
std::vector< std::vector< double > > | constant_modes_values |
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 ML implementation. A structure of the current type are passed to the constructor of PreconditionAMG.
Definition at line 1329 of file trilinos_precondition.h.
TrilinosWrappers::PreconditionAMG::AdditionalData::AdditionalData | ( | const bool | elliptic = true , |
const bool | higher_order_elements = false , |
||
const unsigned int | n_cycles = 1 , |
||
const bool | w_cycle = 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.
Making use of the DoFTools::extract_constant_modes() function, the constant_modes
vector can be initialized for a given field in the following manner:
Definition at line 38 of file trilinos_precondition_ml.cc.
void TrilinosWrappers::PreconditionAMG::AdditionalData::set_parameters | ( | Teuchos::ParameterList & | parameter_list, |
std::unique_ptr< Epetra_MultiVector > & | distributed_constant_modes, | ||
const Epetra_RowMatrix & | matrix | ||
) | const |
Fill in a parameter_list
that can be used to initialize the AMG preconditioner.
The matrix
is used in conjunction with the constant_modes
to configure the null space settings for the preconditioner. The distributed_constant_modes
are initialized by this function, and must remain in scope until PreconditionAMG::initialize() has been called.
See the documentation for the Trilinos ML package for details on what options are available for modification.
Definition at line 66 of file trilinos_precondition_ml.cc.
void TrilinosWrappers::PreconditionAMG::AdditionalData::set_parameters | ( | Teuchos::ParameterList & | parameter_list, |
std::unique_ptr< Epetra_MultiVector > & | distributed_constant_modes, | ||
const SparseMatrix & | matrix | ||
) | const |
Fill in a parameter list that can be used to initialize the AMG preconditioner.
Definition at line 200 of file trilinos_precondition_ml.cc.
void TrilinosWrappers::PreconditionAMG::AdditionalData::set_operator_null_space | ( | Teuchos::ParameterList & | parameter_list, |
std::unique_ptr< Epetra_MultiVector > & | distributed_constant_modes, | ||
const Epetra_RowMatrix & | matrix | ||
) | const |
Configure the null space setting in the parameter_list
for the input matrix
based on the state of the constant_modes
variable.
Definition at line 126 of file trilinos_precondition_ml.cc.
void TrilinosWrappers::PreconditionAMG::AdditionalData::set_operator_null_space | ( | Teuchos::ParameterList & | parameter_list, |
std::unique_ptr< Epetra_MultiVector > & | distributed_constant_modes, | ||
const SparseMatrix & | matrix | ||
) | const |
Configure the null space setting in the parameter_list
for the input matrix
based on the state of the constant_modes
variable.
Definition at line 213 of file trilinos_precondition_ml.cc.
bool TrilinosWrappers::PreconditionAMG::AdditionalData::elliptic |
Determines whether the AMG preconditioner should be optimized for elliptic problems (ML option smoothed aggregation SA, using a Chebyshev smoother) or for non-elliptic problems (ML option non-symmetric smoothed aggregation NSSA, smoother is SSOR with underrelaxation).
Definition at line 1445 of file trilinos_precondition.h.
bool TrilinosWrappers::PreconditionAMG::AdditionalData::higher_order_elements |
Determines whether the matrix that the preconditioner is built upon is generated from linear or higher-order elements.
Definition at line 1451 of file trilinos_precondition.h.
unsigned int TrilinosWrappers::PreconditionAMG::AdditionalData::n_cycles |
Defines how many multigrid cycles should be performed by the preconditioner.
Definition at line 1457 of file trilinos_precondition.h.
bool TrilinosWrappers::PreconditionAMG::AdditionalData::w_cycle |
Defines whether a w-cycle should be used instead of the standard setting of a v-cycle.
Definition at line 1463 of file trilinos_precondition.h.
double TrilinosWrappers::PreconditionAMG::AdditionalData::aggregation_threshold |
This threshold tells the AMG setup how the coarsening should be performed. In the AMG used by ML, 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 1474 of file trilinos_precondition.h.
std::vector<std::vector<bool> > TrilinosWrappers::PreconditionAMG::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, and default value is OK) or on a vector-valued equation. For vector-valued problem with n_components
components, the provided constant_modes
should fulfill the following requirements:
constant_modes.size() == n_components
constant_modes[*].size()
needs to be equal to either the total number of degrees of freedom associated with the linear system (constant_modes[*].size() == n_dofs
), or the number of locally owned degrees of freedom (constant_modes[*].size() == n_locally_owned_dofs
). In parallel computations, the latter is the more appropriate choice since one does not want to store vectors of global size on individual processes. constant_modes[ic][id] == true
if DoF id
is a degree of freedom that is part of vector component ic
. We obtain the constant_modes
fulfilling the above requirements with the function DoFTools::extract_constant_modes.
Definition at line 1500 of file trilinos_precondition.h.
std::vector<std::vector<double> > TrilinosWrappers::PreconditionAMG::AdditionalData::constant_modes_values |
Same as above, but with values instead of Booleans. This is useful if you want to specify rotational modes in addition to the translational modes. See also: DoFTools::extract_rigid_body_modes and DoFTools::extract_level_rigid_body_modes.
Definition at line 1509 of file trilinos_precondition.h.
unsigned int TrilinosWrappers::PreconditionAMG::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 1521 of file trilinos_precondition.h.
unsigned int TrilinosWrappers::PreconditionAMG::AdditionalData::smoother_overlap |
Determines the overlap in the SSOR/Chebyshev error smoother when run in parallel.
Definition at line 1527 of file trilinos_precondition.h.
bool TrilinosWrappers::PreconditionAMG::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 1534 of file trilinos_precondition.h.
const char* TrilinosWrappers::PreconditionAMG::AdditionalData::smoother_type |
Determines which smoother to use for the AMG cycle. Possibilities for smoother_type are the following:
Definition at line 1570 of file trilinos_precondition.h.
const char* TrilinosWrappers::PreconditionAMG::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 1576 of file trilinos_precondition.h.