Reference documentation for deal.II version GIT relicensing-136-gb80d0be4af 2024-03-18 08:20:02+00:00
\(\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
Public Member Functions | Public Attributes | List of all members
LinearOperator< Range, Domain, Payload > Class Template Reference

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

Inheritance diagram for LinearOperator< Range, Domain, Payload >:
Inheritance graph
[legend]

Public Member Functions

 LinearOperator (const Payload &payload=Payload())
 
 LinearOperator (const LinearOperator< Range, Domain, Payload > &)=default
 
template<typename Op , typename = std::enable_if_t< !std::is_base_of_v<LinearOperator<Range, Domain, Payload>, Op>>>
 LinearOperator (const Op &op)
 
LinearOperator< Range, Domain, Payload > & operator= (const LinearOperator< Range, Domain, Payload > &)=default
 
template<typename Op , typename = std::enable_if_t< !std::is_base_of_v<LinearOperator<Range, Domain, Payload>, Op>>>
LinearOperator< Range, Domain, Payload > & operator= (const Op &op)
 

Public Attributes

std::function< void(Range &v, const Domain &u)> vmult
 
std::function< void(Range &v, const Domain &u)> vmult_add
 
std::function< void(Domain &v, const Range &u)> Tvmult
 
std::function< void(Domain &v, const Range &u)> Tvmult_add
 
std::function< void(Range &v, bool omit_zeroing_entries)> reinit_range_vector
 
std::function< void(Domain &v, bool omit_zeroing_entries)> reinit_domain_vector
 

Related Symbols

(Note that these are not member symbols.)

Manipulation of a BlockLinearOperator
template<typename Range = BlockVector<double>, typename Domain = Range, typename BlockPayload = internal::BlockLinearOperatorImplementation::EmptyBlockPayload<>>
LinearOperator< Domain, Range, typename BlockPayload::BlockType > block_forward_substitution (const BlockLinearOperator< Range, Domain, BlockPayload > &block_operator, const BlockLinearOperator< Domain, Range, BlockPayload > &diagonal_inverse)
 
Indirectly applying constraints to LinearOperator
template<typename Range , typename Domain , typename Payload >
LinearOperator< Range, Domain, Payload > distribute_constraints_linear_operator (const AffineConstraints< typename Range::value_type > &constraints, const LinearOperator< Range, Domain, Payload > &exemplar)
 
template<typename Range , typename Domain , typename Payload >
LinearOperator< Range, Domain, Payload > project_to_constrained_linear_operator (const AffineConstraints< typename Range::value_type > &constraints, const LinearOperator< Range, Domain, Payload > &exemplar)
 
template<typename Range , typename Domain , typename Payload >
LinearOperator< Range, Domain, Payload > constrained_linear_operator (const AffineConstraints< typename Range::value_type > &constraints, const LinearOperator< Range, Domain, Payload > &linop)
 
template<typename Range , typename Domain , typename Payload >
PackagedOperation< Range > constrained_right_hand_side (const AffineConstraints< typename Range::value_type > &constraints, const LinearOperator< Range, Domain, Payload > &linop, const Range &right_hand_side)
 
Vector space operations
template<typename Range , typename Domain , typename Payload >
LinearOperator< Range, Domain, Payload > operator+ (const LinearOperator< Range, Domain, Payload > &first_op, const LinearOperator< Range, Domain, Payload > &second_op)
 
template<typename Range , typename Domain , typename Payload >
LinearOperator< Range, Domain, Payload > operator- (const LinearOperator< Range, Domain, Payload > &first_op, const LinearOperator< Range, Domain, Payload > &second_op)
 
template<typename Range , typename Domain , typename Payload >
LinearOperator< Range, Domain, Payload > operator* (typename Range::value_type number, const LinearOperator< Range, Domain, Payload > &op)
 
template<typename Range , typename Domain , typename Payload >
LinearOperator< Range, Domain, Payload > operator* (const LinearOperator< Range, Domain, Payload > &op, typename Domain::value_type number)
 
Composition and manipulation of a LinearOperator
template<typename Range , typename Intermediate , typename Domain , typename Payload >
LinearOperator< Range, Domain, Payload > operator* (const LinearOperator< Range, Intermediate, Payload > &first_op, const LinearOperator< Intermediate, Domain, Payload > &second_op)
 
template<typename Range , typename Domain , typename Payload >
LinearOperator< Domain, Range, Payload > transpose_operator (const LinearOperator< Range, Domain, Payload > &op)
 
template<typename Payload , typename Solver , typename Preconditioner , typename Range = typename Solver::vector_type, typename Domain = Range>
LinearOperator< Domain, Range, Payload > inverse_operator (const LinearOperator< Range, Domain, Payload > &op, Solver &solver, const Preconditioner &preconditioner)
 
template<typename Payload , typename Solver , typename Range = typename Solver::vector_type, typename Domain = Range>
LinearOperator< Domain, Range, Payload > inverse_operator (const LinearOperator< Range, Domain, Payload > &op, Solver &solver, const LinearOperator< Range, Domain, Payload > &preconditioner)
 
template<typename Payload , typename Solver , typename Range = typename Solver::vector_type, typename Domain = Range>
LinearOperator< Domain, Range, Payload > inverse_operator (const LinearOperator< Range, Domain, Payload > &op, Solver &solver)
 
template<typename Payload , typename Solver , typename Range = typename Solver::vector_type, typename Domain = Range>
LinearOperator< Domain, Range, Payload > inverse_operator (const LinearOperator< Range, Domain, Payload > &op, Solver &solver, const PreconditionIdentity &)
 
Creation of a LinearOperator
template<typename Range , typename Payload = internal::LinearOperatorImplementation::EmptyPayload>
LinearOperator< Range, Range, Payload > identity_operator (const std::function< void(Range &, bool)> &reinit_vector)
 
template<typename Range , typename Domain , typename Payload >
LinearOperator< Range, Domain, Payload > identity_operator (const LinearOperator< Range, Domain, Payload > &op)
 
template<typename Range , typename Domain , typename Payload >
LinearOperator< Range, Domain, Payload > null_operator (const LinearOperator< Range, Domain, Payload > &op)
 
template<typename Range , typename Payload = internal::LinearOperatorImplementation::EmptyPayload>
LinearOperator< Range, Range, Payload > mean_value_filter (const std::function< void(Range &, bool)> &reinit_vector)
 
template<typename Range , typename Domain , typename Payload >
LinearOperator< Range, Domain, Payload > mean_value_filter (const LinearOperator< Range, Domain, Payload > &op)
 
template<typename Range , typename Domain , typename Payload , typename Matrix >
LinearOperator< Range, Domain, Payload > linear_operator (const Matrix &matrix)
 
template<typename Range , typename Domain , typename Payload , typename OperatorExemplar , typename Matrix >
LinearOperator< Range, Domain, Payload > linear_operator (const OperatorExemplar &operator_exemplar, const Matrix &matrix)
 
template<typename Range , typename Domain , typename Payload , typename Matrix >
LinearOperator< Range, Domain, Payload > linear_operator (const LinearOperator< Range, Domain, Payload > &operator_exemplar, const Matrix &matrix)
 
Creation of a LinearOperator related to the Schur Complement
template<typename Range_1 , typename Domain_1 , typename Range_2 , typename Domain_2 , typename Payload >
LinearOperator< Range_2, Domain_2, Payload > schur_complement (const LinearOperator< Domain_1, Range_1, Payload > &A_inv, const LinearOperator< Range_1, Domain_2, Payload > &B, const LinearOperator< Range_2, Domain_1, Payload > &C, const LinearOperator< Range_2, Domain_2, Payload > &D)
 
Creation of a LinearOperator
template<typename Range , typename Domain = Range, typename Matrix >
LinearOperator< Range, Domain, TrilinosWrappers::internal::LinearOperatorImplementation::TrilinosPayloadlinear_operator (const TrilinosWrappers::SparseMatrix &operator_exemplar, const Matrix &matrix)
 
template<typename Range , typename Domain = Range>
LinearOperator< Range, Domain, TrilinosWrappers::internal::LinearOperatorImplementation::TrilinosPayloadlinear_operator (const TrilinosWrappers::SparseMatrix &matrix)
 
template<typename Range , typename Domain , typename Matrix >
LinearOperator< Range, Domain, TrilinosWrappers::internal::LinearOperatorImplementation::TrilinosPayloadlinear_operator (const LinearOperator< Range, Domain, TrilinosWrappers::internal::LinearOperatorImplementation::TrilinosPayload > &operator_exemplar, const Matrix &matrix)
 

In-place vector space operations

bool is_null_operator
 
LinearOperator< Range, Domain, Payload > & operator+= (const LinearOperator< Range, Domain, Payload > &second_op)
 
LinearOperator< Range, Domain, Payload > & operator-= (const LinearOperator< Range, Domain, Payload > &second_op)
 
LinearOperator< Range, Domain, Payload > & operator*= (const LinearOperator< Domain, Domain, Payload > &second_op)
 
LinearOperator< Range, Domain, Payload > operator*= (typename Domain::value_type number)
 

Detailed Description

template<typename Range, typename Domain, typename Payload>
class LinearOperator< Range, Domain, Payload >

A class to store the abstract concept of a linear operator.

The class essentially consists of std::function objects that store the knowledge of how to apply the linear operator by implementing the abstract Matrix interface:

std::function<void(Range &, const Domain &)> vmult;
std::function<void(Range &, const Domain &)> vmult_add;
std::function<void(Domain &, const Range &)> Tvmult;
std::function<void(Domain &, const Range &)> Tvmult_add;
std::function< void(Range &v, const Domain &u)> vmult_add
std::function< void(Domain &v, const Range &u)> Tvmult
std::function< void(Range &v, const Domain &u)> vmult
std::function< void(Domain &v, const Range &u)> Tvmult_add

But, in contrast to a usual matrix object, the domain and range of the linear operator are also bound to the LinearOperator class on the type level. Because of this, LinearOperator<Range, Domain> has two additional function objects

std::function<void(Range &, bool)> reinit_range_vector;
std::function<void(Domain &, bool)> reinit_domain_vector;
std::function< void(Domain &v, bool omit_zeroing_entries)> reinit_domain_vector
std::function< void(Range &v, bool omit_zeroing_entries)> reinit_range_vector

that store the knowledge how to initialize (resize + internal data structures) an arbitrary vector of the Range and Domain space.

The primary purpose of this class is to provide syntactic sugar for complex matrix-vector operations and free the user from having to create, set up and handle intermediate storage locations by hand.

As an example consider the operation \((A+k\,B)\,C\), where \(A\), \(B\) and \(C\) denote (possible different) matrices. In order to construct a LinearOperator op that stores the knowledge of this operation, one can write:

const double k = ...;
// Setup and assembly of matrices
const auto op_a = linear_operator(A);
const auto op_b = linear_operator(B);
const auto op_c = linear_operator(C);
const auto op = (op_a + k * op_b) * op_c;
LinearOperator< Range, Domain, Payload > linear_operator(const Matrix &matrix)
Note
This class makes heavy use of std::function objects and lambda functions. This flexibility comes with a run-time penalty. Only use this object to encapsulate matrix object of medium to large size (as a rule of thumb, sparse matrices with a size \(1000\times1000\), or larger).
In order to use Trilinos or PETSc sparse matrices and preconditioners in conjunction with the LinearOperator class, it is necessary to extend the functionality of the LinearOperator class by means of an additional Payload.

For example: LinearOperator instances representing matrix inverses usually require calling some linear solver. These solvers may not have interfaces to the LinearOperator (which, for example, may represent a composite operation). The TrilinosWrappers::internal::LinearOperatorImplementation::TrilinosPayload therefore provides an interface extension to the LinearOperator so that it can be passed to the solver and used by the solver as if it were a Trilinos operator. This implies that all of the necessary functionality of the specific Trilinos operator has been overloaded within the Payload class. This includes operator-vector multiplication and inverse operator-vector multiplication, where the operator can be either a TrilinosWrappers::SparseMatrix or a TrilinosWrappers::PreconditionBase and the vector is a native Trilinos vector.

Another case where payloads provide a crucial supplement to the LinearOperator class are when composite operations are constructed (via operator overloading). In this instance, it is again necessary to provide an interface that produces the result of this composite operation that is compatible with Trilinos operator used by Trilinos solvers.

Note
Many use cases of LinearOperator lead to intermediate expressions requiring a PackagedOperation. In order to include all necessary header files in one go consider using

In order to use the full LinearOperator and PackagedOperation

Note
To ensure that the correct payload is provided, wrapper functions for linear operators have been provided within the respective TrilinosWrappers (and, in the future, PETScWrappers) namespaces.

Examples of use

The step-20 tutorial program has a detailed usage example of the LinearOperator class.

Instrumenting operations

It is sometimes useful to know when functions are called, or to inject additional operations. In such cases, what one wants is to replace, for example, the vmult object of this class with one that does the additional operations and then calls what was originally supposed to happen. This can be done with commands such as the following:

auto A_inv = inverse_operator(A, solver_A, preconditioner_A);
A_inv.vmult = [base_vmult = A_inv.vmult](Vector<double> &dst,
const Vector<double> &src) {
std::cout << "Calling A_inv.vmult()" << std::endl;
base_vmult(dst, src);
};
LinearOperator< Domain, Range, Payload > inverse_operator(const LinearOperator< Range, Domain, Payload > &op, Solver &solver, const Preconditioner &preconditioner)

Here, we replace A_inv.vmult with a lambda function that first captures the previous value of A_inv.vmult and stores it in the base_vmult object. The newly installed A_inv.vmult function then first outputs some status information, and then calls the original functionality.

This approach works for all of the other function objects mentioned above as well.

Definition at line 200 of file linear_operator.h.

Constructor & Destructor Documentation

◆ LinearOperator() [1/3]

template<typename Range , typename Domain , typename Payload >
LinearOperator< Range, Domain, Payload >::LinearOperator ( const Payload &  payload = Payload())
inline

Create an empty LinearOperator object. When a payload is passed to this constructor, the resulting operator is constructed with a functional payload. In either case, this constructor yields an object that can not actually be used for any linear operator operations, and will throw an exception upon invocation.

Definition at line 211 of file linear_operator.h.

◆ LinearOperator() [2/3]

template<typename Range , typename Domain , typename Payload >
LinearOperator< Range, Domain, Payload >::LinearOperator ( const LinearOperator< Range, Domain, Payload > &  )
default

Default copy constructor.

◆ LinearOperator() [3/3]

template<typename Range , typename Domain , typename Payload >
template<typename Op , typename = std::enable_if_t< !std::is_base_of_v<LinearOperator<Range, Domain, Payload>, Op>>>
LinearOperator< Range, Domain, Payload >::LinearOperator ( const Op &  op)
inline

Templated copy constructor that creates a LinearOperator object from an object op for which the conversion function linear_operator is defined.

Definition at line 265 of file linear_operator.h.

Member Function Documentation

◆ operator=() [1/2]

template<typename Range , typename Domain , typename Payload >
LinearOperator< Range, Domain, Payload > & LinearOperator< Range, Domain, Payload >::operator= ( const LinearOperator< Range, Domain, Payload > &  )
default

Default copy assignment operator.

◆ operator=() [2/2]

template<typename Range , typename Domain , typename Payload >
template<typename Op , typename = std::enable_if_t< !std::is_base_of_v<LinearOperator<Range, Domain, Payload>, Op>>>
LinearOperator< Range, Domain, Payload > & LinearOperator< Range, Domain, Payload >::operator= ( const Op &  op)
inline

Templated copy assignment operator for an object op for which the conversion function linear_operator is defined.

Definition at line 284 of file linear_operator.h.

◆ operator+=()

template<typename Range , typename Domain , typename Payload >
LinearOperator< Range, Domain, Payload > & LinearOperator< Range, Domain, Payload >::operator+= ( const LinearOperator< Range, Domain, Payload > &  second_op)
inline

Addition with a LinearOperator second_op with the same Domain and Range.

Definition at line 343 of file linear_operator.h.

◆ operator-=()

template<typename Range , typename Domain , typename Payload >
LinearOperator< Range, Domain, Payload > & LinearOperator< Range, Domain, Payload >::operator-= ( const LinearOperator< Range, Domain, Payload > &  second_op)
inline

Subtraction with a LinearOperator second_op with the same Domain and Range.

Definition at line 354 of file linear_operator.h.

◆ operator*=() [1/2]

template<typename Range , typename Domain , typename Payload >
LinearOperator< Range, Domain, Payload > & LinearOperator< Range, Domain, Payload >::operator*= ( const LinearOperator< Domain, Domain, Payload > &  second_op)
inline

Composition of the LinearOperator with an endomorphism second_op of the Domain space.

Definition at line 365 of file linear_operator.h.

◆ operator*=() [2/2]

template<typename Range , typename Domain , typename Payload >
LinearOperator< Range, Domain, Payload > LinearOperator< Range, Domain, Payload >::operator*= ( typename Domain::value_type  number)
inline

Scalar multiplication of the LinearOperator with number from the right.

Definition at line 376 of file linear_operator.h.

Member Data Documentation

◆ vmult

template<typename Range , typename Domain , typename Payload >
std::function<void(Range &v, const Domain &u)> LinearOperator< Range, Domain, Payload >::vmult

Application of the LinearOperator object to a vector u of the Domain space giving a vector v of the Range space.

Definition at line 294 of file linear_operator.h.

◆ vmult_add

template<typename Range , typename Domain , typename Payload >
std::function<void(Range &v, const Domain &u)> LinearOperator< Range, Domain, Payload >::vmult_add

Application of the LinearOperator object to a vector u of the Domain space. The result is added to the vector v.

Definition at line 300 of file linear_operator.h.

◆ Tvmult

template<typename Range , typename Domain , typename Payload >
std::function<void(Domain &v, const Range &u)> LinearOperator< Range, Domain, Payload >::Tvmult

Application of the transpose LinearOperator object to a vector u of the Range space giving a vector v of the Domain space.

Definition at line 306 of file linear_operator.h.

◆ Tvmult_add

template<typename Range , typename Domain , typename Payload >
std::function<void(Domain &v, const Range &u)> LinearOperator< Range, Domain, Payload >::Tvmult_add

Application of the transpose LinearOperator object to a vector u of the Range space.The result is added to the vector v.

Definition at line 312 of file linear_operator.h.

◆ reinit_range_vector

template<typename Range , typename Domain , typename Payload >
std::function<void(Range &v, bool omit_zeroing_entries)> LinearOperator< Range, Domain, Payload >::reinit_range_vector

Initializes a vector v of the Range space to be directly usable as the destination parameter in an application of vmult. Similar to the reinit functions of the vector classes, the boolean determines whether a fast initialization is done, i.e., if it is set to false the content of the vector is set to 0.

Definition at line 321 of file linear_operator.h.

◆ reinit_domain_vector

template<typename Range , typename Domain , typename Payload >
std::function<void(Domain &v, bool omit_zeroing_entries)> LinearOperator< Range, Domain, Payload >::reinit_domain_vector

Initializes a vector of the Domain space to be directly usable as the source parameter in an application of vmult. Similar to the reinit functions of the vector classes, the boolean determines whether a fast initialization is done, i.e., if it is set to false the content of the vector is set to 0.

Definition at line 331 of file linear_operator.h.

◆ is_null_operator

template<typename Range , typename Domain , typename Payload >
bool LinearOperator< Range, Domain, Payload >::is_null_operator

This bool is used to determine whether a linear operator is a null operator. In this case the class is able to optimize some operations like multiplication or addition.

Definition at line 387 of file linear_operator.h.


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