Reference documentation for deal.II version 9.5.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
Classes | Public Types | Public Member Functions | Static Public Member Functions | Public Attributes | Protected Attributes | Private Member Functions | Private Attributes | List of all members
TrilinosWrappers::SolverBase Class Reference

#include <deal.II/lac/trilinos_solver.h>

Inheritance diagram for TrilinosWrappers::SolverBase:
[legend]

Classes

struct  AdditionalData
 

Public Types

enum  SolverName {
  cg , cgs , gmres , bicgstab ,
  tfqmr
}
 

Public Member Functions

 SolverBase (SolverControl &cn, const AdditionalData &data=AdditionalData())
 
 SolverBase (const enum SolverName solver_name, SolverControl &cn, const AdditionalData &data=AdditionalData())
 
virtual ~SolverBase ()=default
 
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)
 
SolverControlcontrol () const
 
template<>
void set_preconditioner (AztecOO &solver, const PreconditionBase &preconditioner)
 
template<>
void set_preconditioner (AztecOO &solver, const Epetra_Operator &preconditioner)
 

Static Public Member Functions

static ::ExceptionBaseExcTrilinosError (int arg1)
 

Public Attributes

enum TrilinosWrappers::SolverBase::SolverName solver_name
 

Protected Attributes

SolverControlsolver_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::unique_ptr< Epetra_LinearProblem > linear_problem
 
std::unique_ptr< AztecOO_StatusTest > status_test
 
AztecOO solver
 
const AdditionalData additional_data
 

Detailed Description

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 78 of file trilinos_solver.h.

Member Enumeration Documentation

◆ 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 88 of file trilinos_solver.h.

Constructor & Destructor Documentation

◆ SolverBase() [1/2]

SolverBase< VectorType >::SolverBase ( SolverControl cn,
const AdditionalData data = AdditionalData() 
)

Constructor. Takes the solver control object and creates the solver.

Definition at line 49 of file trilinos_solver.cc.

◆ SolverBase() [2/2]

SolverBase< VectorType >::SolverBase ( const enum SolverName  solver_name,
SolverControl cn,
const AdditionalData data = AdditionalData() 
)

Second constructor. This constructor takes an enum object that specifies the solver name and sets the appropriate Krylov method.

Definition at line 57 of file trilinos_solver.cc.

◆ ~SolverBase()

virtual TrilinosWrappers::SolverBase::~SolverBase ( )
virtualdefault

Destructor.

Member Function Documentation

◆ solve() [1/9]

void SolverBase< VectorType >::solve ( const SparseMatrix A,
MPI::Vector x,
const MPI::Vector 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 76 of file trilinos_solver.cc.

◆ solve() [2/9]

void SolverBase< VectorType >::solve ( const Epetra_Operator A,
MPI::Vector x,
const MPI::Vector 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 96 of file trilinos_solver.cc.

◆ solve() [3/9]

void SolverBase< VectorType >::solve ( const Epetra_Operator A,
MPI::Vector x,
const MPI::Vector b,
const Epetra_Operator preconditioner 
)

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 117 of file trilinos_solver.cc.

◆ solve() [4/9]

void SolverBase< VectorType >::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 138 of file trilinos_solver.cc.

◆ solve() [5/9]

void SolverBase< VectorType >::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 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 159 of file trilinos_solver.cc.

◆ solve() [6/9]

void SolverBase< VectorType >::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 178 of file trilinos_solver.cc.

◆ solve() [7/9]

void SolverBase< VectorType >::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 208 of file trilinos_solver.cc.

◆ solve() [8/9]

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.

◆ solve() [9/9]

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.

◆ control()

SolverControl & SolverBase< VectorType >::control ( ) const

Access to object that controls convergence.

Definition at line 68 of file trilinos_solver.cc.

◆ do_solve()

template<typename Preconditioner >
template void SolverBase< VectorType >::do_solve ( const Preconditioner &  preconditioner)
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 415 of file trilinos_solver.cc.

◆ set_preconditioner() [1/3]

template<typename Preconditioner >
void TrilinosWrappers::SolverBase::set_preconditioner ( AztecOO &  solver,
const Preconditioner &  preconditioner 
)
private

A function that sets the preconditioner that the solver will apply

◆ set_preconditioner() [2/3]

template<>
void SolverBase::set_preconditioner ( AztecOO &  solver,
const PreconditionBase preconditioner 
)

Definition at line 565 of file trilinos_solver.cc.

◆ set_preconditioner() [3/3]

template<>
void SolverBase::set_preconditioner ( AztecOO &  solver,
const Epetra_Operator preconditioner 
)

Definition at line 583 of file trilinos_solver.cc.

Member Data Documentation

◆ solver_name

enum TrilinosWrappers::SolverBase::SolverName TrilinosWrappers::SolverBase::solver_name

◆ solver_control

SolverControl& TrilinosWrappers::SolverBase::solver_control
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 312 of file trilinos_solver.h.

◆ linear_problem

std::unique_ptr<Epetra_LinearProblem> TrilinosWrappers::SolverBase::linear_problem
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 335 of file trilinos_solver.h.

◆ status_test

std::unique_ptr<AztecOO_StatusTest> TrilinosWrappers::SolverBase::status_test
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 341 of file trilinos_solver.h.

◆ solver

AztecOO TrilinosWrappers::SolverBase::solver
private

A structure that contains the Trilinos solver and preconditioner objects.

Definition at line 347 of file trilinos_solver.h.

◆ additional_data

const AdditionalData TrilinosWrappers::SolverBase::additional_data
private

Store a copy of the flags for this particular solver.

Definition at line 352 of file trilinos_solver.h.


The documentation for this class was generated from the following files: