Reference documentation for deal.II version 9.1.1
|
#include <deal.II/lac/trilinos_solver.h>
Classes | |
struct | AdditionalData |
Public Member Functions | |
SolverDirect (SolverControl &cn, const AdditionalData &data=AdditionalData()) | |
virtual | ~SolverDirect ()=default |
void | initialize (const SparseMatrix &A) |
void | solve (MPI::Vector &x, const MPI::Vector &b) |
void | solve (::LinearAlgebra::distributed::Vector< double > &x, const ::LinearAlgebra::distributed::Vector< double > &b) |
void | solve (const SparseMatrix &A, MPI::Vector &x, const MPI::Vector &b) |
void | solve (const SparseMatrix &A, ::Vector< double > &x, const ::Vector< double > &b) |
void | solve (const SparseMatrix &A, ::LinearAlgebra::distributed::Vector< double > &x, const ::LinearAlgebra::distributed::Vector< double > &b) |
SolverControl & | control () const |
Static Public Member Functions | |
static ::ExceptionBase & | ExcTrilinosError (int arg1) |
Private Member Functions | |
void | do_solve () |
Private Attributes | |
SolverControl & | solver_control |
std::unique_ptr< Epetra_LinearProblem > | linear_problem |
std::unique_ptr< Amesos_BaseSolver > | solver |
const AdditionalData | additional_data |
An implementation of Trilinos direct solvers (using the Amesos package). The data field AdditionalData::solver_type can be used to specify the type of solver. It allows the use of built-in solvers Amesos_Klu as well as third-party solvers Amesos_Superludist or Amesos_Mumps.
For instructions on how to install Trilinos for use with direct solvers other than KLU, see the link to the Trilinos installation instructions linked to from the deal.II ReadMe file.
Definition at line 568 of file trilinos_solver.h.
TrilinosWrappers::SolverDirect::SolverDirect | ( | SolverControl & | cn, |
const AdditionalData & | data = AdditionalData() |
||
) |
Constructor. Takes the solver control object and creates the solver.
Definition at line 684 of file trilinos_solver.cc.
|
virtualdefault |
Destructor.
void TrilinosWrappers::SolverDirect::initialize | ( | const SparseMatrix & | A | ) |
Initializes the direct solver for the matrix A
and creates a factorization for it with the package chosen from the additional data structure. Note that there is no need for a preconditioner here and solve() is not called.
Definition at line 700 of file trilinos_solver.cc.
void TrilinosWrappers::SolverDirect::solve | ( | MPI::Vector & | x, |
const MPI::Vector & | b | ||
) |
Solve the linear system Ax=b
based on the package set in initialize(). Note the matrix is not refactorized during this call.
Definition at line 745 of file trilinos_solver.cc.
void TrilinosWrappers::SolverDirect::solve | ( | ::LinearAlgebra::distributed::Vector< double > & | x, |
const ::LinearAlgebra::distributed::Vector< double > & | b | ||
) |
Solve the linear system Ax=b
based on the package set in initialize() for deal.II's own parallel vectors. Note the matrix is not refactorized during this call.
Definition at line 773 of file trilinos_solver.cc.
void TrilinosWrappers::SolverDirect::solve | ( | const SparseMatrix & | A, |
MPI::Vector & | x, | ||
const MPI::Vector & | b | ||
) |
Solve the linear system Ax=b
. Creates a factorization of the matrix with the package chosen from the additional data structure and performs the solve. Note that there is no need for a preconditioner here.
Definition at line 859 of file trilinos_solver.cc.
void TrilinosWrappers::SolverDirect::solve | ( | const SparseMatrix & | A, |
::Vector< double > & | x, | ||
const ::Vector< double > & | b | ||
) |
Solve the linear system Ax=b
. This class works with Trilinos matrices, but takes deal.II serial vectors as argument. Since these vectors are not distributed, this function does only what you expect in case the matrix is serial (i.e., locally owned). Otherwise, an exception will be thrown.
Definition at line 876 of file trilinos_solver.cc.
void TrilinosWrappers::SolverDirect::solve | ( | const SparseMatrix & | A, |
::LinearAlgebra::distributed::Vector< double > & | x, | ||
const ::LinearAlgebra::distributed::Vector< double > & | b | ||
) |
Solve the linear system Ax=b
for deal.II's own parallel vectors. Creates a factorization of the matrix with the package chosen from the additional data structure and performs the solve. Note that there is no need for a preconditioner here.
Definition at line 902 of file trilinos_solver.cc.
SolverControl & TrilinosWrappers::SolverDirect::control | ( | ) | const |
Access to object that controls convergence.
Definition at line 692 of file trilinos_solver.cc.
|
private |
Actually performs the operations for solving the linear system, including the factorization and forward and backward substitution.
Definition at line 808 of file trilinos_solver.cc.
|
private |
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 707 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 714 of file trilinos_solver.h.
|
private |
A structure that contains the Trilinos solver and preconditioner objects.
Definition at line 720 of file trilinos_solver.h.
|
private |
Store a copy of the flags for this particular solver.
Definition at line 725 of file trilinos_solver.h.