![]() |
Reference documentation for deal.II version 8.1.0
|
#include <trilinos_precondition.h>
Classes | |
struct | AdditionalData |
Public Member Functions | |
void | initialize (const SparseMatrix &matrix, const AdditionalData &additional_data=AdditionalData()) |
void | initialize (const SparseMatrix &matrix, const Teuchos::ParameterList &ml_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 | reinit () |
void | clear () |
size_type | memory_consumption () const |
![]() | |
PreconditionBase () | |
PreconditionBase (const PreconditionBase &) | |
~PreconditionBase () | |
void | clear () |
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 (::parallel::distributed::Vector< double > &dst, const ::parallel::distributed::Vector< double > &src) const |
virtual void | Tvmult (::parallel::distributed::Vector< double > &dst, const ::parallel::distributed::Vector< double > &src) const |
DeclException1 (ExcNonMatchingMaps, std::string,<< "The sparse matrix the preconditioner is based on "<< "uses a map that is not compatible to the one in vector "<< arg1<< ". Check preconditioner and matrix setup.") | |
![]() | |
Subscriptor () | |
Subscriptor (const Subscriptor &) | |
virtual | ~Subscriptor () |
Subscriptor & | operator= (const Subscriptor &) |
void | subscribe (const char *identifier=0) const |
void | unsubscribe (const char *identifier=0) const |
unsigned int | n_subscriptions () const |
void | list_subscribers () const |
DeclException3 (ExcInUse, int, char *, std::string &,<< "Object of class "<< arg2<< " is still used by "<< arg1<< " other objects.\n"<< "(Additional information: "<< arg3<< ")\n"<< "Note the entry in the Frequently Asked Questions of "<< "deal.II (linked to from http://www.dealii.org/) for "<< "more information on what this error means.") | |
DeclException2 (ExcNoSubscriber, char *, char *,<< "No subscriber with identifier \""<< arg2<< "\" did subscribe to this object of class "<< arg1) | |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Private Attributes | |
std_cxx1x::shared_ptr < SparseMatrix > | trilinos_matrix |
Additional Inherited Members | |
![]() | |
typedef ::types::global_dof_index | size_type |
![]() | |
std_cxx1x::shared_ptr < Epetra_Operator > | preconditioner |
Epetra_MpiComm | communicator |
std_cxx1x::shared_ptr< Epetra_Map > | vector_distributor |
This class implements an algebraic multigrid (AMG) preconditioner based on the Trilinos ML implementation, which is a black-box preconditioner that works well for many PDE-based linear problems. What this class does is twofold. When the initialize() function is invoked, a ML preconditioner object is created based on the matrix that we want the preconditioner to be based on. A call of the respective vmult
function does call the respective operation in the Trilinos package, where it is called ApplyInverse
. Use of this class is explained in the step-31 tutorial program.
Since the Trilinos objects we want to use are heavily dependent on Epetra objects, we recommend using this class in conjunction with Trilinos (Epetra) sparse matrices and vectors. There is support for use with matrices of the deal.II::SparseMatrix class and corresponding vectors, too, but this requires generating a copy of the matrix, which is slower and takes (much) more memory. When doing such a copy operation, we can still profit from the fact that some of the entries in the preconditioner matrix are zero and hence can be neglected.
The implementation is able to distinguish between matrices from elliptic problems and convection dominated problems. We use the standard options provided by Trilinos ML for elliptic problems, except that we use a Chebyshev smoother instead of a symmetric Gauss-Seidel smoother. For most elliptic problems, Chebyshev provides a better damping of high frequencies (in the algebraic sense) than Gauss-Seidel (SSOR), and is faster (Chebyshev requires only some matrix-vector products, whereas SSOR requires substitutions which are more expensive). Moreover, Chebyshev is perfectly parallel in the sense that it does not degenerate when used on many processors. SSOR, on the other hand, gets more Jacobi-like on many processors.
For proper functionality of this class we recommend using Trilinos v9.0 and higher. Older versions may have problems with generating the coarse-matrix structure when using matrices with many nonzero entries per row (i.e., matrices stemming from higher order finite element discretizations).
Definition at line 1283 of file trilinos_precondition.h.
void TrilinosWrappers::PreconditionAMG::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.
void TrilinosWrappers::PreconditionAMG::initialize | ( | const SparseMatrix & | matrix, |
const Teuchos::ParameterList & | ml_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 all 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. In particular, users need to follow the ML instructions in case a vector-valued problem ought to be solved.
void TrilinosWrappers::PreconditionAMG::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.
void TrilinosWrappers::PreconditionAMG::reinit | ( | ) |
This function can be used for a faster recalculation of the preconditioner construction when the matrix entries underlying the preconditioner have changed, but the matrix sparsity pattern has remained the same. What this function does is taking the already generated coarsening structure, computing the AMG prolongation and restriction according to a smoothed aggregation strategy and then building the whole multilevel hiearchy. This function can be considerably faster than the initialize function, since the coarsening pattern is usually the most difficult thing to do when setting up the AMG ML preconditioner.
void TrilinosWrappers::PreconditionAMG::clear | ( | ) |
Destroys the preconditioner, leaving an object like just after having called the constructor.
size_type TrilinosWrappers::PreconditionAMG::memory_consumption | ( | ) | const |
Prints an estimate of the memory consumption of this class.
|
private |
A copy of the deal.II matrix into Trilinos format.
Definition at line 1506 of file trilinos_precondition.h.