Reference documentation for deal.II version 8.5.1
|
#include <deal.II/lac/trilinos_precondition.h>
Classes | |
struct | AdditionalData |
Public Member Functions | |
PreconditionAMGMueLu () | |
~PreconditionAMGMueLu () | |
void | initialize (const SparseMatrix &matrix, const AdditionalData &additional_data=AdditionalData()) |
void | initialize (const Epetra_CrsMatrix &matrix, const AdditionalData &additional_data=AdditionalData()) |
void | initialize (const SparseMatrix &matrix, Teuchos::ParameterList &muelu_parameters) |
void | initialize (const Epetra_CrsMatrix &matrix, Teuchos::ParameterList &muelu_parameters) |
template<typename number > | |
void | initialize (const ::SparseMatrix< number > &deal_ii_sparse_matrix, const AdditionalData &additional_data=AdditionalData(), const double drop_tolerance=1e-13, const ::SparsityPattern *use_this_sparsity=0) |
void | clear () |
size_type | memory_consumption () const |
Public Member Functions inherited from TrilinosWrappers::PreconditionBase | |
PreconditionBase () | |
PreconditionBase (const PreconditionBase &) | |
~PreconditionBase () | |
void | clear () |
MPI_Comm | get_mpi_communicator () const |
void | transpose () |
virtual void | vmult (VectorBase &dst, const VectorBase &src) const |
virtual void | Tvmult (VectorBase &dst, const VectorBase &src) const |
virtual void | vmult (::Vector< double > &dst, const ::Vector< double > &src) const |
virtual void | Tvmult (::Vector< double > &dst, const ::Vector< double > &src) const |
virtual void | vmult (::LinearAlgebra::distributed::Vector< double > &dst, const ::LinearAlgebra::distributed::Vector< double > &src) const |
virtual void | Tvmult (::LinearAlgebra::distributed::Vector< double > &dst, const ::LinearAlgebra::distributed::Vector< double > &src) const |
Epetra_Operator & | trilinos_operator () const |
IndexSet | locally_owned_domain_indices () const |
IndexSet | locally_owned_range_indices () const |
Public Member Functions inherited from Subscriptor | |
Subscriptor () | |
Subscriptor (const Subscriptor &) | |
Subscriptor (Subscriptor &&) | |
virtual | ~Subscriptor () |
Subscriptor & | operator= (const Subscriptor &) |
Subscriptor & | operator= (Subscriptor &&) |
void | subscribe (const char *identifier=0) const |
void | unsubscribe (const char *identifier=0) const |
unsigned int | n_subscriptions () const |
void | list_subscribers () const |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Private Attributes | |
std_cxx11::shared_ptr< SparseMatrix > | trilinos_matrix |
Additional Inherited Members | |
Public Types inherited from TrilinosWrappers::PreconditionBase | |
typedef ::types::global_dof_index | size_type |
Static Public Member Functions inherited from TrilinosWrappers::PreconditionBase | |
static ::ExceptionBase & | ExcNonMatchingMaps (std::string arg1) |
Static Public Member Functions inherited from Subscriptor | |
static ::ExceptionBase & | ExcInUse (int arg1, char *arg2, std::string &arg3) |
static ::ExceptionBase & | ExcNoSubscriber (char *arg1, char *arg2) |
Protected Attributes inherited from TrilinosWrappers::PreconditionBase | |
std_cxx11::shared_ptr< Epetra_Operator > | preconditioner |
Epetra_MpiComm | communicator |
std_cxx11::shared_ptr< Epetra_Map > | vector_distributor |
This class implements an algebraic multigrid (AMG) preconditioner based on the Trilinos MueLu implementation, which is a black-box preconditioner that works well for many PDE-based linear problems. The interface of PreconditionerAMGMueLu is the same as the interface of PreconditionerAMG except for the higher_order_elements parameter which does not exist in PreconditionerAMGMueLu.
Definition at line 1630 of file trilinos_precondition.h.
TrilinosWrappers::PreconditionAMGMueLu::PreconditionAMGMueLu | ( | ) |
Constructor.
Definition at line 104 of file trilinos_precondition_muelu.cc.
TrilinosWrappers::PreconditionAMGMueLu::~PreconditionAMGMueLu | ( | ) |
Destructor.
Definition at line 113 of file trilinos_precondition_muelu.cc.
void TrilinosWrappers::PreconditionAMGMueLu::initialize | ( | const SparseMatrix & | matrix, |
const AdditionalData & | additional_data = AdditionalData() |
||
) |
Let Trilinos compute a multilevel hierarchy for the solution of a linear system with the given matrix. The function uses the matrix format specified in TrilinosWrappers::SparseMatrix.
Definition at line 122 of file trilinos_precondition_muelu.cc.
void TrilinosWrappers::PreconditionAMGMueLu::initialize | ( | const Epetra_CrsMatrix & | matrix, |
const AdditionalData & | additional_data = AdditionalData() |
||
) |
Let Trilinos compute a multilevel hierarchy for the solution of a linear system with the given matrix. As opposed to the other initialize function above, this function uses an object of type Epetra_CrsMatrixCrs.
Definition at line 131 of file trilinos_precondition_muelu.cc.
void TrilinosWrappers::PreconditionAMGMueLu::initialize | ( | const SparseMatrix & | matrix, |
Teuchos::ParameterList & | muelu_parameters | ||
) |
Let Trilinos compute a multilevel hierarchy for the solution of a linear system with the given matrix. The function uses the matrix format specified in TrilinosWrappers::SparseMatrix.
This function is similar to the one above, but allows the user to set most of the options of the Trilinos ML preconditioner. In order to find out about all the options for ML, we refer to the ML user's guide. Not all ML options have a corresponding MueLu option.
Definition at line 231 of file trilinos_precondition_muelu.cc.
void TrilinosWrappers::PreconditionAMGMueLu::initialize | ( | const Epetra_CrsMatrix & | matrix, |
Teuchos::ParameterList & | muelu_parameters | ||
) |
Let Trilinos compute a multilevel hierarchy for the solution of a linear system with the given matrix. As opposed to the other initialize function above, this function uses an object of type Epetra_CrsMatrix.
Definition at line 240 of file trilinos_precondition_muelu.cc.
void TrilinosWrappers::PreconditionAMGMueLu::initialize | ( | const ::SparseMatrix< number > & | deal_ii_sparse_matrix, |
const AdditionalData & | additional_data = AdditionalData() , |
||
const double | drop_tolerance = 1e-13 , |
||
const ::SparsityPattern * | use_this_sparsity = 0 |
||
) |
Let Trilinos compute a multilevel hierarchy for the solution of a linear system with the given matrix. This function takes a deal.ii matrix and copies the content into a Trilinos matrix, so the function can be considered rather inefficient.
Definition at line 278 of file trilinos_precondition_muelu.cc.
void TrilinosWrappers::PreconditionAMGMueLu::clear | ( | ) |
Destroys the preconditioner, leaving an object like just after having called the constructor.
Definition at line 305 of file trilinos_precondition_muelu.cc.
PreconditionAMGMueLu::size_type TrilinosWrappers::PreconditionAMGMueLu::memory_consumption | ( | ) | const |
Prints an estimate of the memory consumption of this class.
Definition at line 314 of file trilinos_precondition_muelu.cc.
|
private |
A copy of the deal.II matrix into Trilinos format.
Definition at line 1842 of file trilinos_precondition.h.