Reference documentation for deal.II version 9.1.1
|
#include <deal.II/lac/trilinos_sparse_matrix.h>
Inherits Epetra_Operator.
Public Types | |
using | VectorType = Epetra_MultiVector |
using | Range = VectorType |
using | Domain = VectorType |
Public Member Functions | |
Constructors / destructor | |
TrilinosPayload () | |
TrilinosPayload (const TrilinosWrappers::SparseMatrix &matrix_exemplar, const TrilinosWrappers::SparseMatrix &matrix) | |
TrilinosPayload (const TrilinosWrappers::SparseMatrix &matrix_exemplar, const TrilinosWrappers::PreconditionBase &preconditioner) | |
TrilinosPayload (const TrilinosWrappers::PreconditionBase &preconditioner_exemplar, const TrilinosWrappers::PreconditionBase &preconditioner) | |
TrilinosPayload (const TrilinosPayload &payload) | |
TrilinosPayload (const TrilinosPayload &first_op, const TrilinosPayload &second_op) | |
virtual | ~TrilinosPayload () override=default |
TrilinosPayload | identity_payload () const |
TrilinosPayload | null_payload () const |
TrilinosPayload | transpose_payload () const |
template<typename Solver , typename Preconditioner > | |
std::enable_if< std::is_base_of< TrilinosWrappers::SolverBase, Solver >::value &&std::is_base_of< TrilinosWrappers::PreconditionBase, Preconditioner >::value, TrilinosPayload >::type | inverse_payload (Solver &, const Preconditioner &) const |
template<typename Solver , typename Preconditioner > | |
std::enable_if< !(std::is_base_of< TrilinosWrappers::SolverBase, Solver >::value &&std::is_base_of< TrilinosWrappers::PreconditionBase, Preconditioner >::value), TrilinosPayload >::type | inverse_payload (Solver &, const Preconditioner &) const |
Core Epetra_Operator functionality | |
virtual bool | UseTranspose () const override |
virtual int | SetUseTranspose (bool UseTranspose) override |
virtual int | Apply (const VectorType &X, VectorType &Y) const override |
virtual int | ApplyInverse (const VectorType &Y, VectorType &X) const override |
Additional Epetra_Operator functionality | |
virtual const char * | Label () const override |
virtual const Epetra_Comm & | Comm () const override |
virtual const Epetra_Map & | OperatorDomainMap () const override |
virtual const Epetra_Map & | OperatorRangeMap () const override |
Private Member Functions | |
virtual bool | HasNormInf () const override |
virtual double | NormInf () const override |
Private Attributes | |
bool | use_transpose |
Epetra_MpiComm | communicator |
Epetra_Map | domain_map |
Epetra_Map | range_map |
LinearOperator functionality | |
std::function< void(VectorType &, const VectorType &)> | vmult |
std::function< void(VectorType &, const VectorType &)> | Tvmult |
std::function< void(VectorType &, const VectorType &)> | inv_vmult |
std::function< void(VectorType &, const VectorType &)> | inv_Tvmult |
IndexSet | locally_owned_domain_indices () const |
IndexSet | locally_owned_range_indices () const |
MPI_Comm | get_mpi_communicator () const |
void | transpose () |
This is an extension class to LinearOperators for Trilinos sparse matrix and preconditioner types. It provides the interface to performing basic operations (vmult
and Tvmult
) on Trilinos vector types. It fulfills the requirements necessary for wrapping a Trilinos solver, which calls Epetra_Operator functions, as a LinearOperator.
vmult
and Tvmult
functions. This object is not thread-safe when the transpose flag is set on it or the Trilinos object to which it refers. See the docuemtation for the TrilinosWrappers::internal::LinearOperatorImplementation::TrilinosPayload::SetUseTranspose() function for further details.Definition at line 2243 of file trilinos_sparse_matrix.h.
using TrilinosWrappers::internal::LinearOperatorImplementation::TrilinosPayload::VectorType = Epetra_MultiVector |
Definition for the internally supported vector type.
Definition at line 2249 of file trilinos_sparse_matrix.h.
Definition for the vector type for the domain space of the operator.
Definition at line 2254 of file trilinos_sparse_matrix.h.
using TrilinosWrappers::internal::LinearOperatorImplementation::TrilinosPayload::Domain = VectorType |
Definition for the vector type for the range space of the operator.
Definition at line 2259 of file trilinos_sparse_matrix.h.
TrilinosWrappers::internal::LinearOperatorImplementation::TrilinosPayload::TrilinosPayload | ( | ) |
Default constructor
Definition at line 2509 of file trilinos_sparse_matrix.cc.
TrilinosWrappers::internal::LinearOperatorImplementation::TrilinosPayload::TrilinosPayload | ( | const TrilinosWrappers::SparseMatrix & | matrix_exemplar, |
const TrilinosWrappers::SparseMatrix & | matrix | ||
) |
Constructor for a sparse matrix based on an exemplary matrix
Definition at line 2548 of file trilinos_sparse_matrix.cc.
TrilinosWrappers::internal::LinearOperatorImplementation::TrilinosPayload::TrilinosPayload | ( | const TrilinosWrappers::SparseMatrix & | matrix_exemplar, |
const TrilinosWrappers::PreconditionBase & | preconditioner | ||
) |
Constructor for a preconditioner based on an exemplary matrix
Definition at line 2633 of file trilinos_sparse_matrix.cc.
TrilinosWrappers::internal::LinearOperatorImplementation::TrilinosPayload::TrilinosPayload | ( | const TrilinosWrappers::PreconditionBase & | preconditioner_exemplar, |
const TrilinosWrappers::PreconditionBase & | preconditioner | ||
) |
Constructor for a preconditioner based on an exemplary preconditioner
Definition at line 2753 of file trilinos_sparse_matrix.cc.
TrilinosWrappers::internal::LinearOperatorImplementation::TrilinosPayload::TrilinosPayload | ( | const TrilinosPayload & | payload | ) |
Default copy constructor
Definition at line 2872 of file trilinos_sparse_matrix.cc.
TrilinosWrappers::internal::LinearOperatorImplementation::TrilinosPayload::TrilinosPayload | ( | const TrilinosPayload & | first_op, |
const TrilinosPayload & | second_op | ||
) |
Composite copy constructor
This is required for PackagedOperations as it sets up the domain and range maps, and composite vmult
and Tvmult
operations based on the combined operation of both operations
Definition at line 2887 of file trilinos_sparse_matrix.cc.
|
overridevirtualdefault |
Destructor
TrilinosPayload TrilinosWrappers::internal::LinearOperatorImplementation::TrilinosPayload::identity_payload | ( | ) | const |
Return a payload configured for identity operations
Definition at line 2900 of file trilinos_sparse_matrix.cc.
TrilinosPayload TrilinosWrappers::internal::LinearOperatorImplementation::TrilinosPayload::null_payload | ( | ) | const |
Return a payload configured for null operations
Definition at line 2926 of file trilinos_sparse_matrix.cc.
TrilinosPayload TrilinosWrappers::internal::LinearOperatorImplementation::TrilinosPayload::transpose_payload | ( | ) | const |
Return a payload configured for transpose operations
Definition at line 2966 of file trilinos_sparse_matrix.cc.
std::enable_if< std::is_base_of<TrilinosWrappers::SolverBase, Solver>::value && std::is_base_of<TrilinosWrappers::PreconditionBase, Preconditioner>::value, TrilinosPayload>::type TrilinosWrappers::internal::LinearOperatorImplementation::TrilinosPayload::inverse_payload | ( | Solver & | , |
const Preconditioner & | |||
) | const |
Return a payload configured for inverse operations
Invoking this factory function will configure two additional functions, namely inv_vmult
and inv_Tvmult
, both of which wrap inverse operations. The vmult
and Tvmult
operations retain the standard definitions inherited from op
.
Solver
derives from TrilinosWrappers::SolverBase, andPreconditioner
derives from TrilinosWrappers::PreconditionBase. std::enable_if< !(std::is_base_of<TrilinosWrappers::SolverBase, Solver>::value && std::is_base_of<TrilinosWrappers::PreconditionBase, Preconditioner>::value), TrilinosPayload>::type TrilinosWrappers::internal::LinearOperatorImplementation::TrilinosPayload::inverse_payload | ( | Solver & | , |
const Preconditioner & | |||
) | const |
Return a payload configured for inverse operations
Invoking this factory function will configure two additional functions, namely inv_vmult
and inv_Tvmult
, both of which are disabled because the Solver
or Preconditioner
are not compatible with Epetra_MultiVector. The vmult
and Tvmult
operations retain the standard definitions inherited from op
.
Solver
does not derive from TrilinosWrappers::SolverBase, andPreconditioner
does not derive from TrilinosWrappers::PreconditionBase. IndexSet TrilinosWrappers::internal::LinearOperatorImplementation::TrilinosPayload::locally_owned_domain_indices | ( | ) | const |
Return an IndexSet that defines the partitioning of the domain space of this matrix, i.e., the partitioning of the vectors this matrix has to be multiplied with / operate on.
Definition at line 2976 of file trilinos_sparse_matrix.cc.
IndexSet TrilinosWrappers::internal::LinearOperatorImplementation::TrilinosPayload::locally_owned_range_indices | ( | ) | const |
Return an IndexSet that defines the partitioning of the range space of this matrix, i.e., the partitioning of the vectors that result from matrix-vector products.
Definition at line 2984 of file trilinos_sparse_matrix.cc.
MPI_Comm TrilinosWrappers::internal::LinearOperatorImplementation::TrilinosPayload::get_mpi_communicator | ( | ) | const |
Return the MPI communicator object in use with this Payload.
Definition at line 2992 of file trilinos_sparse_matrix.cc.
void TrilinosWrappers::internal::LinearOperatorImplementation::TrilinosPayload::transpose | ( | ) |
Sets an internal flag so that all operations performed by the matrix, i.e., multiplications, are done in transposed order.
Definition at line 3004 of file trilinos_sparse_matrix.cc.
|
overridevirtual |
Return the status of the transpose flag for this operator
This overloads the same function from the Trilinos class Epetra_Operator.
Definition at line 3012 of file trilinos_sparse_matrix.cc.
|
overridevirtual |
Sets an internal flag so that all operations performed by the matrix, i.e., multiplications, are done in transposed order.
This overloads the same function from the Trilinos class Epetra_Operator.
Definition at line 3020 of file trilinos_sparse_matrix.cc.
|
overridevirtual |
Apply the vmult operation on a vector X
(of internally defined type VectorType) and store the result in the vector Y
.
This overloads the same function from the Trilinos class Epetra_Operator.
Definition at line 3035 of file trilinos_sparse_matrix.cc.
|
overridevirtual |
Apply the vmult inverse operation on a vector X
(of internally defined type VectorType) and store the result in the vector Y
.
In practise, this function is only called from a Trilinos solver if the wrapped object is to act as a preconditioner.
This overloads the same function from the Trilinos class Epetra_Operator.
Definition at line 3046 of file trilinos_sparse_matrix.cc.
|
overridevirtual |
Return a label to describe this class.
This overloads the same function from the Trilinos class Epetra_Operator.
Definition at line 3057 of file trilinos_sparse_matrix.cc.
|
overridevirtual |
Return a reference to the underlying MPI communicator for this object.
This overloads the same function from the Trilinos class Epetra_Operator.
Definition at line 3065 of file trilinos_sparse_matrix.cc.
|
overridevirtual |
Return the partitioning of the domain space of this matrix, i.e., the partitioning of the vectors this matrix has to be multiplied with.
This overloads the same function from the Trilinos class Epetra_Operator.
Definition at line 3073 of file trilinos_sparse_matrix.cc.
|
overridevirtual |
Return the partitioning of the range space of this matrix, i.e., the partitioning of the vectors that are result from matrix-vector products.
This overloads the same function from the Trilinos class Epetra_Operator.
Definition at line 3081 of file trilinos_sparse_matrix.cc.
|
overrideprivatevirtual |
Return a flag that describes whether this operator can return the computation of the infinity norm. Since in general this is not the case, this always returns a negetive result.
This overloads the same function from the Trilinos class Epetra_Operator.
Definition at line 3089 of file trilinos_sparse_matrix.cc.
|
overrideprivatevirtual |
Return the infinity norm of this operator. Throws an error since, in general, we cannot compute this value.
This overloads the same function from the Trilinos class Epetra_Operator.
Definition at line 3097 of file trilinos_sparse_matrix.cc.
std::function<void(VectorType &, const VectorType &)> TrilinosWrappers::internal::LinearOperatorImplementation::TrilinosPayload::vmult |
The standard matrix-vector operation to be performed by the payload when Apply is called.
Definition at line 2427 of file trilinos_sparse_matrix.h.
std::function<void(VectorType &, const VectorType &)> TrilinosWrappers::internal::LinearOperatorImplementation::TrilinosPayload::Tvmult |
The standard transpose matrix-vector operation to be performed by the payload when Apply is called.
Definition at line 2436 of file trilinos_sparse_matrix.h.
std::function<void(VectorType &, const VectorType &)> TrilinosWrappers::internal::LinearOperatorImplementation::TrilinosPayload::inv_vmult |
The inverse matrix-vector operation to be performed by the payload when ApplyInverse is called.
Definition at line 2446 of file trilinos_sparse_matrix.h.
std::function<void(VectorType &, const VectorType &)> TrilinosWrappers::internal::LinearOperatorImplementation::TrilinosPayload::inv_Tvmult |
The inverse transpose matrix-vector operation to be performed by the payload when ApplyInverse is called.
Definition at line 2456 of file trilinos_sparse_matrix.h.
|
private |
A flag recording whether the operator is to perform standard matrix-vector multiplication, or the transpose operation.
Definition at line 2579 of file trilinos_sparse_matrix.h.
|
private |
Internal communication pattern in case the matrix needs to be copied from deal.II format.
Definition at line 2586 of file trilinos_sparse_matrix.h.
|
private |
Epetra_Map that sets the partitioning of the domain space of this operator.
Definition at line 2595 of file trilinos_sparse_matrix.h.
|
private |
Epetra_Map that sets the partitioning of the range space of this operator.
Definition at line 2601 of file trilinos_sparse_matrix.h.