#include <deal.II/lac/trilinos_solver.h>
|
| SolverGMRES (SolverControl &cn, const AdditionalData &data=AdditionalData()) |
|
void | solve (const SparseMatrix &A, MPI::Vector &x, const MPI::Vector &b, const PreconditionBase &preconditioner) |
|
void | solve (const Epetra_Operator &A, MPI::Vector &x, const MPI::Vector &b, const PreconditionBase &preconditioner) |
|
void | solve (const Epetra_Operator &A, MPI::Vector &x, const MPI::Vector &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 |
|
template<> |
void | set_preconditioner (AztecOO &solver, const PreconditionBase &preconditioner) |
|
template<> |
void | set_preconditioner (AztecOO &solver, const Epetra_Operator &preconditioner) |
|
template<> |
void | set_preconditioner (AztecOO &solver, const PreconditionBase &preconditioner) |
|
template<> |
void | set_preconditioner (AztecOO &solver, const Epetra_Operator &preconditioner) |
|
|
template<typename Preconditioner > |
void | do_solve (const Preconditioner &preconditioner) |
|
template<typename Preconditioner > |
void | set_preconditioner (AztecOO &solver, const Preconditioner &preconditioner) |
|
An implementation of the solver interface using the Trilinos GMRES solver.
Definition at line 413 of file trilinos_solver.h.
◆ SolverName
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:
Enumerator |
---|
cg | Use the conjugate gradient (CG) algorithm.
|
cgs | Use the conjugate gradient squared (CGS) algorithm.
|
gmres | Use the generalized minimum residual (GMRES) algorithm.
|
bicgstab | Use the biconjugate gradient stabilized (BICGStab) algorithm.
|
tfqmr | Use the transpose-free quasi-minimal residual (TFQMR) method.
|
Definition at line 89 of file trilinos_solver.h.
◆ SolverGMRES()
Constructor. In contrast to deal.II's own solvers, there is no need to give a vector memory object.
The last argument takes a structure with additional, solver dependent flags for tuning.
Definition at line 602 of file trilinos_solver.cc.
◆ solve() [1/9]
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 75 of file trilinos_solver.cc.
◆ solve() [2/9]
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 95 of file trilinos_solver.cc.
◆ solve() [3/9]
Solve the linear system Ax=b
where both A
and its preconditioner
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 116 of file trilinos_solver.cc.
◆ solve() [4/9]
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 137 of file trilinos_solver.cc.
◆ solve() [5/9]
Solve the linear system Ax=b
where both A
and its preconditioner
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 158 of file trilinos_solver.cc.
◆ solve() [6/9]
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 177 of file trilinos_solver.cc.
◆ solve() [7/9]
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 207 of file trilinos_solver.cc.
◆ solve() [8/9]
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 227 of file trilinos_solver.cc.
◆ solve() [9/9]
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 255 of file trilinos_solver.cc.
◆ control()
◆ do_solve()
template<typename Preconditioner >
template void SolverBase< VectorType >::do_solve |
( |
const Preconditioner & | preconditioner | ) |
|
|
privateinherited |
The solve function is used to set properly the Epetra_LinearProblem, once it is done this function solves the linear problem.
Definition at line 416 of file trilinos_solver.cc.
◆ set_preconditioner() [1/5]
template<typename Preconditioner >
void TrilinosWrappers::SolverBase::set_preconditioner |
( |
AztecOO & | solver, |
|
|
const Preconditioner & | preconditioner ) |
|
privateinherited |
A function that sets the preconditioner that the solver will apply
◆ set_preconditioner() [2/5]
template<>
void TrilinosWrappers::SolverBase::set_preconditioner |
( |
AztecOO & | solver, |
|
|
const PreconditionBase & | preconditioner ) |
|
inherited |
◆ set_preconditioner() [3/5]
template<>
void TrilinosWrappers::SolverBase::set_preconditioner |
( |
AztecOO & | solver, |
|
|
const Epetra_Operator & | preconditioner ) |
|
inherited |
◆ set_preconditioner() [4/5]
template<>
void TrilinosWrappers::SolverBase::set_preconditioner |
( |
AztecOO & | solver, |
|
|
const PreconditionBase & | preconditioner ) |
|
inherited |
◆ set_preconditioner() [5/5]
template<>
void TrilinosWrappers::SolverBase::set_preconditioner |
( |
AztecOO & | solver, |
|
|
const Epetra_Operator & | preconditioner ) |
|
inherited |
◆ solver_name
◆ solver_control
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 313 of file trilinos_solver.h.
◆ linear_problem
std::unique_ptr<Epetra_LinearProblem> TrilinosWrappers::SolverBase::linear_problem |
|
privateinherited |
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 336 of file trilinos_solver.h.
◆ status_test
std::unique_ptr<AztecOO_StatusTest> TrilinosWrappers::SolverBase::status_test |
|
privateinherited |
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 342 of file trilinos_solver.h.
◆ solver
AztecOO TrilinosWrappers::SolverBase::solver |
|
privateinherited |
A structure that contains the Trilinos solver and preconditioner objects.
Definition at line 348 of file trilinos_solver.h.
◆ additional_data
Store a copy of the flags for this particular solver.
Definition at line 353 of file trilinos_solver.h.
The documentation for this class was generated from the following files: