Reference documentation for deal.II version 8.5.1
|
#include <deal.II/lac/trilinos_solver.h>
Classes | |
struct | AdditionalData |
Public Types | |
enum | SolverName { cg, cgs, gmres, bicgstab, tfqmr } |
Public Member Functions | |
SolverBase (SolverControl &cn) | |
SolverBase (const enum SolverName solver_name, SolverControl &cn) | |
virtual | ~SolverBase () |
void | solve (const SparseMatrix &A, VectorBase &x, const VectorBase &b, const PreconditionBase &preconditioner) |
void | solve (const Epetra_Operator &A, VectorBase &x, const VectorBase &b, const PreconditionBase &preconditioner) |
void | solve (const Epetra_Operator &A, VectorBase &x, const VectorBase &b, const Epetra_Operator &preconditioner) |
void | solve (const Epetra_Operator &A, Epetra_MultiVector &x, const Epetra_MultiVector &b, const PreconditionBase &preconditioner) |
void | solve (const Epetra_Operator &A, Epetra_MultiVector &x, const Epetra_MultiVector &b, const Epetra_Operator &preconditioner) |
void | solve (const SparseMatrix &A, ::Vector< double > &x, const ::Vector< double > &b, const PreconditionBase &preconditioner) |
void | solve (Epetra_Operator &A, ::Vector< double > &x, const ::Vector< double > &b, const PreconditionBase &preconditioner) |
void | solve (const SparseMatrix &A, ::LinearAlgebra::distributed::Vector< double > &x, const ::LinearAlgebra::distributed::Vector< double > &b, const PreconditionBase &preconditioner) |
void | solve (Epetra_Operator &A, ::LinearAlgebra::distributed::Vector< double > &x, const ::LinearAlgebra::distributed::Vector< double > &b, const PreconditionBase &preconditioner) |
SolverControl & | control () const |
Static Public Member Functions | |
static ::ExceptionBase & | ExcTrilinosError (int arg1) |
Protected Attributes | |
SolverControl & | solver_control |
Private Member Functions | |
template<typename Preconditioner > | |
void | do_solve (const Preconditioner &preconditioner) |
template<typename Preconditioner > | |
void | set_preconditioner (AztecOO &solver, const Preconditioner &preconditioner) |
Private Attributes | |
std_cxx11::shared_ptr< Epetra_LinearProblem > | linear_problem |
std_cxx11::shared_ptr< AztecOO_StatusTest > | status_test |
AztecOO | solver |
const AdditionalData | additional_data |
Base class for solver classes using the Trilinos solvers. Since solvers in Trilinos are selected based on flags passed to a generic solver object, basically all the actual solver calls happen in this class, and derived classes simply set the right flags to select one solver or another, or to set certain parameters for individual solvers. For a general discussion on the Trilinos solver package AztecOO, we refer to the AztecOO user guide.
This solver class can also be used as a standalone class, where the respective Krylov method is set via the flag solver_name
. This can be done at runtime (e.g., when parsing the solver from a ParameterList) and is similar to the deal.II class SolverSelector.
Definition at line 67 of file trilinos_solver.h.
Enumeration object that is set in the constructor of the derived classes and tells Trilinos which solver to use. This option can also be set in the user program, so one might use this base class instead of one of the specialized derived classes when the solver should be set at runtime. Currently enabled options are:
Definition at line 78 of file trilinos_solver.h.
TrilinosWrappers::SolverBase::SolverBase | ( | SolverControl & | cn | ) |
Constructor. Takes the solver control object and creates the solver.
Definition at line 50 of file trilinos_solver.cc.
TrilinosWrappers::SolverBase::SolverBase | ( | const enum SolverName | solver_name, |
SolverControl & | cn | ||
) |
Second constructor. This constructor takes an enum object that specifies the solver name and sets the appropriate Krylov method.
Definition at line 58 of file trilinos_solver.cc.
|
virtual |
Destructor.
Definition at line 67 of file trilinos_solver.cc.
void TrilinosWrappers::SolverBase::solve | ( | const SparseMatrix & | A, |
VectorBase & | x, | ||
const VectorBase & | b, | ||
const PreconditionBase & | preconditioner | ||
) |
Solve the linear system Ax=b
. Depending on the information provided by derived classes and the object passed as a preconditioner, one of the linear solvers and preconditioners of Trilinos is chosen.
Definition at line 81 of file trilinos_solver.cc.
void TrilinosWrappers::SolverBase::solve | ( | const Epetra_Operator & | A, |
VectorBase & | x, | ||
const VectorBase & | b, | ||
const PreconditionBase & | preconditioner | ||
) |
Solve the linear system Ax=b
where A
is an operator. This function can be used for matrix free computation. Depending on the information provided by derived classes and the object passed as a preconditioner, one of the linear solvers and preconditioners of Trilinos is chosen.
Definition at line 103 of file trilinos_solver.cc.
void TrilinosWrappers::SolverBase::solve | ( | const Epetra_Operator & | A, |
VectorBase & | x, | ||
const VectorBase & | b, | ||
const Epetra_Operator & | preconditioner | ||
) |
Solve the linear system Ax=b
where both A
and its precondtioner
are an operator. This function can be used when both A
and the preconditioner
are LinearOperators derived from a TrilinosPayload. Depending on the information provided by derived classes and the object passed as a preconditioner, one of the linear solvers and preconditioners of Trilinos is chosen.
Definition at line 125 of file trilinos_solver.cc.
void TrilinosWrappers::SolverBase::solve | ( | const Epetra_Operator & | A, |
Epetra_MultiVector & | x, | ||
const Epetra_MultiVector & | b, | ||
const PreconditionBase & | preconditioner | ||
) |
Solve the linear system Ax=b
where A
is an operator, and the vectors x
and b
are native Trilinos vector types. This function can be used when A
is a LinearOperators derived from a TrilinosPayload. Depending on the information provided by derived classes and the object passed as a preconditioner, one of the linear solvers and preconditioners of Trilinos is chosen.
Definition at line 147 of file trilinos_solver.cc.
void TrilinosWrappers::SolverBase::solve | ( | const Epetra_Operator & | A, |
Epetra_MultiVector & | x, | ||
const Epetra_MultiVector & | b, | ||
const Epetra_Operator & | preconditioner | ||
) |
Solve the linear system Ax=b
where both A
and its precondtioner
are an operator, and the vectors x
and b
are native Trilinos vector types. This function can be used when both A
and the preconditioner
are LinearOperators derived from a TrilinosPayload. Depending on the information provided by derived classes and the object passed as a preconditioner, one of the linear solvers and preconditioners of Trilinos is chosen.
Definition at line 169 of file trilinos_solver.cc.
void TrilinosWrappers::SolverBase::solve | ( | const SparseMatrix & | A, |
::Vector< double > & | x, | ||
const ::Vector< double > & | b, | ||
const PreconditionBase & | preconditioner | ||
) |
Solve the linear system Ax=b
. Depending on the information provided by derived classes and the object passed as a preconditioner, one of the linear solvers and preconditioners of Trilinos is chosen. This class works with matrices according to the TrilinosWrappers format, but can take deal.II vectors as argument. Since deal.II are serial vectors (not distributed), this function does only what you expect in case the matrix is locally owned. Otherwise, an exception will be thrown.
Definition at line 189 of file trilinos_solver.cc.
void TrilinosWrappers::SolverBase::solve | ( | Epetra_Operator & | A, |
::Vector< double > & | x, | ||
const ::Vector< double > & | b, | ||
const PreconditionBase & | preconditioner | ||
) |
Solve the linear system Ax=b
where A
is an operator. This function can be used for matrix free computations. Depending on the information provided by derived classes and the object passed as a preconditioner, one of the linear solvers and preconditioners of Trilinos is chosen. This class works with matrices according to the TrilinosWrappers format, but can take deal.II vectors as argument. Since deal.II are serial vectors (not distributed), this function does only what you expect in case the matrix is locally owned. Otherwise, an exception will be thrown.
Definition at line 222 of file trilinos_solver.cc.
void TrilinosWrappers::SolverBase::solve | ( | const SparseMatrix & | A, |
::LinearAlgebra::distributed::Vector< double > & | x, | ||
const ::LinearAlgebra::distributed::Vector< double > & | b, | ||
const PreconditionBase & | preconditioner | ||
) |
Solve the linear system Ax=b
for deal.II's parallel distributed vectors. Depending on the information provided by derived classes and the object passed as a preconditioner, one of the linear solvers and preconditioners of Trilinos is chosen.
Definition at line 242 of file trilinos_solver.cc.
void TrilinosWrappers::SolverBase::solve | ( | Epetra_Operator & | A, |
::LinearAlgebra::distributed::Vector< double > & | x, | ||
const ::LinearAlgebra::distributed::Vector< double > & | b, | ||
const PreconditionBase & | preconditioner | ||
) |
Solve the linear system Ax=b
where A
is an operator. This function can be used for matrix free computation. Depending on the information provided by derived classes and the object passed as a preconditioner, one of the linear solvers and preconditioners of Trilinos is chosen.
Definition at line 271 of file trilinos_solver.cc.
SolverControl & TrilinosWrappers::SolverBase::control | ( | ) | const |
Access to object that controls convergence.
Definition at line 73 of file trilinos_solver.cc.
|
private |
The solve function is used to set properly the Epetra_LinearProblem, once it is done this function solves the linear problem.
Definition at line 434 of file trilinos_solver.cc.
|
private |
A function that sets the preconditioner that the solver will apply
|
protected |
Reference to the object that controls convergence of the iterative solver. In fact, for these Trilinos wrappers, Trilinos does so itself, but we copy the data from this object before starting the solution process, and copy the data back into it afterwards.
Definition at line 301 of file trilinos_solver.h.
|
private |
A structure that collects the Trilinos sparse matrix, the right hand side vector and the solution vector, which is passed down to the Trilinos solver.
Definition at line 324 of file trilinos_solver.h.
|
private |
A structure that contains a Trilinos object that can query the linear solver and determine whether the convergence criterion have been met.
Definition at line 330 of file trilinos_solver.h.
|
private |
A structure that contains the Trilinos solver and preconditioner objects.
Definition at line 336 of file trilinos_solver.h.
|
private |
Store a copy of the flags for this particular solver.
Definition at line 341 of file trilinos_solver.h.