Reference documentation for deal.II version 8.5.1
Public Types | Private Member Functions | Private Attributes | List of all members
TrilinosWrappers::internal::LinearOperator::TrilinosPayload Class Reference

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

Inherits Epetra_Operator.

Public Types

typedef Epetra_MultiVector VectorType
 
typedef VectorType Range
 
typedef VectorType Domain
 

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 ()
 
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
 
virtual int SetUseTranspose (bool UseTranspose)
 
virtual int Apply (const VectorType &X, VectorType &Y) const
 
virtual int ApplyInverse (const VectorType &Y, VectorType &X) const
 
Additional Epetra_Operator functionality
virtual const char * Label () const
 
virtual const Epetra_Comm & Comm () const
 
virtual const Epetra_Map & OperatorDomainMap () const
 
virtual const Epetra_Map & OperatorRangeMap () const
 

Private Member Functions

virtual bool HasNormInf () const
 
virtual double NormInf () const
 

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 ()
 

Detailed Description

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.

Note
The TrilinosWrappers::SparseMatrix or TrilinosWrappers::PreconditionBase that this payload wraps is passed by reference to the 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::LinearOperator::TrilinosPayload::SetUseTranspose() function for further details.
Author
Jean-Paul Pelteret, 2016

Definition at line 2133 of file trilinos_sparse_matrix.h.

Member Typedef Documentation

◆ VectorType

Definition for the internally supported vector type.

Definition at line 2141 of file trilinos_sparse_matrix.h.

◆ Range

Definition for the vector type for the domain space of the operator.

Definition at line 2146 of file trilinos_sparse_matrix.h.

◆ Domain

Definition for the vector type for the range space of the operator.

Definition at line 2151 of file trilinos_sparse_matrix.h.

Constructor & Destructor Documentation

◆ TrilinosPayload() [1/6]

TrilinosWrappers::internal::LinearOperator::TrilinosPayload::TrilinosPayload ( )

Default constructor

Note
By design, the resulting object is inoperable since there is insufficient information with which to construct the domain and range maps.

Definition at line 2470 of file trilinos_sparse_matrix.cc.

◆ TrilinosPayload() [2/6]

TrilinosWrappers::internal::LinearOperator::TrilinosPayload::TrilinosPayload ( const TrilinosWrappers::SparseMatrix matrix_exemplar,
const TrilinosWrappers::SparseMatrix matrix 
)

Constructor for a sparse matrix based on an exemplary matrix

Definition at line 2512 of file trilinos_sparse_matrix.cc.

◆ TrilinosPayload() [3/6]

TrilinosWrappers::internal::LinearOperator::TrilinosPayload::TrilinosPayload ( const TrilinosWrappers::SparseMatrix matrix_exemplar,
const TrilinosWrappers::PreconditionBase preconditioner 
)

Constructor for a preconditioner based on an exemplary matrix

Definition at line 2576 of file trilinos_sparse_matrix.cc.

◆ TrilinosPayload() [4/6]

TrilinosWrappers::internal::LinearOperator::TrilinosPayload::TrilinosPayload ( const TrilinosWrappers::PreconditionBase preconditioner_exemplar,
const TrilinosWrappers::PreconditionBase preconditioner 
)

Constructor for a preconditioner based on an exemplary preconditioner

Definition at line 2659 of file trilinos_sparse_matrix.cc.

◆ TrilinosPayload() [5/6]

TrilinosWrappers::internal::LinearOperator::TrilinosPayload::TrilinosPayload ( const TrilinosPayload payload)

Default copy constructor

Definition at line 2742 of file trilinos_sparse_matrix.cc.

◆ TrilinosPayload() [6/6]

TrilinosWrappers::internal::LinearOperator::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 2757 of file trilinos_sparse_matrix.cc.

◆ ~TrilinosPayload()

TrilinosWrappers::internal::LinearOperator::TrilinosPayload::~TrilinosPayload ( )
virtual

Destructor

Definition at line 2768 of file trilinos_sparse_matrix.cc.

Member Function Documentation

◆ identity_payload()

TrilinosPayload TrilinosWrappers::internal::LinearOperator::TrilinosPayload::identity_payload ( ) const

Returns a payload configured for identity operations

Definition at line 2774 of file trilinos_sparse_matrix.cc.

◆ null_payload()

TrilinosPayload TrilinosWrappers::internal::LinearOperator::TrilinosPayload::null_payload ( ) const

Returns a payload configured for null operations

Definition at line 2808 of file trilinos_sparse_matrix.cc.

◆ transpose_payload()

TrilinosPayload TrilinosWrappers::internal::LinearOperator::TrilinosPayload::transpose_payload ( ) const

Returns a payload configured for transpose operations

Definition at line 2854 of file trilinos_sparse_matrix.cc.

◆ inverse_payload() [1/2]

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 TrilinosWrappers::internal::LinearOperator::TrilinosPayload::inverse_payload ( Solver ,
const Preconditioner &   
) const

Returns 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.

Note
This function is enabled only if the solver and preconditioner derive from the respective TrilinosWrappers base classes. The C++ compiler will therefore only consider this function if the following criterion are satisfied:
  1. the Solver derives from TrilinosWrappers::SolverBase, and
  2. the Preconditioner derives from TrilinosWrappers::PreconditionBase.

◆ inverse_payload() [2/2]

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 TrilinosWrappers::internal::LinearOperator::TrilinosPayload::inverse_payload ( Solver ,
const Preconditioner &   
) const

Returns 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.

Note
The C++ compiler will only consider this function if the following criterion are satisfied:
  1. the Solver does not derive from TrilinosWrappers::SolverBase, and
  2. the Preconditioner does not derive from TrilinosWrappers::PreconditionBase.

◆ locally_owned_domain_indices()

IndexSet TrilinosWrappers::internal::LinearOperator::TrilinosPayload::locally_owned_domain_indices ( ) const

Returns 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 2864 of file trilinos_sparse_matrix.cc.

◆ locally_owned_range_indices()

IndexSet TrilinosWrappers::internal::LinearOperator::TrilinosPayload::locally_owned_range_indices ( ) const

Returns 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 2872 of file trilinos_sparse_matrix.cc.

◆ get_mpi_communicator()

MPI_Comm TrilinosWrappers::internal::LinearOperator::TrilinosPayload::get_mpi_communicator ( ) const

Return the MPI communicator object in use with this Payload.

Definition at line 2880 of file trilinos_sparse_matrix.cc.

◆ transpose()

void TrilinosWrappers::internal::LinearOperator::TrilinosPayload::transpose ( )

Sets an internal flag so that all operations performed by the matrix, i.e., multiplications, are done in transposed order.

Note
This does not reshape the matrix to transposed form directly, so care should be taken when using this flag.

Definition at line 2892 of file trilinos_sparse_matrix.cc.

◆ UseTranspose()

bool TrilinosWrappers::internal::LinearOperator::TrilinosPayload::UseTranspose ( ) const
virtual

Return the status of the transpose flag for this operator

This overloads the same function from the Trilinos class Epetra_Operator.

Definition at line 2900 of file trilinos_sparse_matrix.cc.

◆ SetUseTranspose()

int TrilinosWrappers::internal::LinearOperator::TrilinosPayload::SetUseTranspose ( bool  UseTranspose)
virtual

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.

Note
This does not reshape the matrix to transposed form directly, so care should be taken when using this flag. When the flag is set to true (either here or directly on the underlying Trilinos object itself), this object is no longer thread-safe. In essense, it is not possible ensure that the transposed state of the LinearOperator and the underlying Trilinos object remain synchronized throughout all operations that may occur on different threads simultaneously.

Definition at line 2908 of file trilinos_sparse_matrix.cc.

◆ Apply()

int TrilinosWrappers::internal::LinearOperator::TrilinosPayload::Apply ( const VectorType X,
VectorType Y 
) const
virtual

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.

Note
The intended operation depends on the status of the internal transpose flag. If this flag is set to true, the result will be the equivalent of performing a Tvmult operation.

Definition at line 2923 of file trilinos_sparse_matrix.cc.

◆ ApplyInverse()

int TrilinosWrappers::internal::LinearOperator::TrilinosPayload::ApplyInverse ( const VectorType Y,
VectorType X 
) const
virtual

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.

Note
This function will only be operable if the payload has been initalised with an InverseOperator, or is a wrapper to a preconditioner. If not, then using this function will lead to an error being thrown.
The intended operation depends on the status of the internal transpose flag. If this flag is set to true, the result will be the equivalent of performing a Tvmult operation.

Definition at line 2935 of file trilinos_sparse_matrix.cc.

◆ Label()

const char * TrilinosWrappers::internal::LinearOperator::TrilinosPayload::Label ( ) const
virtual

Returns a label to describe this class.

This overloads the same function from the Trilinos class Epetra_Operator.

Definition at line 2947 of file trilinos_sparse_matrix.cc.

◆ Comm()

const Epetra_Comm & TrilinosWrappers::internal::LinearOperator::TrilinosPayload::Comm ( ) const
virtual

Returns a reference to the underlying MPI communicator for this object.

This overloads the same function from the Trilinos class Epetra_Operator.

Definition at line 2955 of file trilinos_sparse_matrix.cc.

◆ OperatorDomainMap()

const Epetra_Map & TrilinosWrappers::internal::LinearOperator::TrilinosPayload::OperatorDomainMap ( ) const
virtual

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 2963 of file trilinos_sparse_matrix.cc.

◆ OperatorRangeMap()

const Epetra_Map & TrilinosWrappers::internal::LinearOperator::TrilinosPayload::OperatorRangeMap ( ) const
virtual

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 2971 of file trilinos_sparse_matrix.cc.

◆ HasNormInf()

bool TrilinosWrappers::internal::LinearOperator::TrilinosPayload::HasNormInf ( ) const
privatevirtual

Returns 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 2979 of file trilinos_sparse_matrix.cc.

◆ NormInf()

double TrilinosWrappers::internal::LinearOperator::TrilinosPayload::NormInf ( ) const
privatevirtual

Returns 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 2987 of file trilinos_sparse_matrix.cc.

Member Data Documentation

◆ vmult

std::function<void(VectorType &, const VectorType &)> TrilinosWrappers::internal::LinearOperator::TrilinosPayload::vmult

The standard matrix-vector operation to be performed by the payload when Apply is called.

Note
This is not called by a LinearOperator, but rather by Trilinos functions that expect this to mimic the action of the LinearOperator.

Definition at line 2311 of file trilinos_sparse_matrix.h.

◆ Tvmult

std::function<void(VectorType &, const VectorType &)> TrilinosWrappers::internal::LinearOperator::TrilinosPayload::Tvmult

The standard transpose matrix-vector operation to be performed by the payload when Apply is called.

Note
This is not called by a LinearOperator, but rather by Trilinos functions that expect this to mimic the action of the LinearOperator.

Definition at line 2320 of file trilinos_sparse_matrix.h.

◆ inv_vmult

std::function<void(VectorType &, const VectorType &)> TrilinosWrappers::internal::LinearOperator::TrilinosPayload::inv_vmult

The inverse matrix-vector operation to be performed by the payload when ApplyInverse is called.

Note
This is not called by a LinearOperator, but rather by Trilinos functions that expect this to mimic the action of the InverseOperator.

Definition at line 2329 of file trilinos_sparse_matrix.h.

◆ inv_Tvmult

std::function<void(VectorType &, const VectorType &)> TrilinosWrappers::internal::LinearOperator::TrilinosPayload::inv_Tvmult

The inverse transpose matrix-vector operation to be performed by the payload when ApplyInverse is called.

Note
This is not called by a LinearOperator, but rather by Trilinos functions that expect this to mimic the action of the InverseOperator.

Definition at line 2338 of file trilinos_sparse_matrix.h.

◆ use_transpose

bool TrilinosWrappers::internal::LinearOperator::TrilinosPayload::use_transpose
private

A flag recording whether the operator is to perform standard matrix-vector multiplication, or the transpose operation.

Definition at line 2463 of file trilinos_sparse_matrix.h.

◆ communicator

Epetra_MpiComm TrilinosWrappers::internal::LinearOperator::TrilinosPayload::communicator
private

Internal communication pattern in case the matrix needs to be copied from deal.II format.

Definition at line 2470 of file trilinos_sparse_matrix.h.

◆ domain_map

Epetra_Map TrilinosWrappers::internal::LinearOperator::TrilinosPayload::domain_map
private

Epetra_Map that sets the partitioning of the domain space of this operator.

Definition at line 2479 of file trilinos_sparse_matrix.h.

◆ range_map

Epetra_Map TrilinosWrappers::internal::LinearOperator::TrilinosPayload::range_map
private

Epetra_Map that sets the partitioning of the range space of this operator.

Definition at line 2485 of file trilinos_sparse_matrix.h.


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