Reference documentation for deal.II version 9.0.0
|
#include <deal.II/lac/linear_operator.h>
Public Member Functions | |
LinearOperator (const Payload &payload=Payload()) | |
LinearOperator (const LinearOperator< Range, Domain, Payload > &)=default | |
template<typename Op , typename = typename std::enable_if<!std::is_base_of<LinearOperator<Range, Domain, Payload>, Op>::value>::type> | |
LinearOperator (const Op &op) | |
LinearOperator< Range, Domain, Payload > & | operator= (const LinearOperator< Range, Domain, Payload > &)=default |
template<typename Op , typename = typename std::enable_if<!std::is_base_of<LinearOperator<Range, Domain, Payload>, Op>::value>::type> | |
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 Functions | |
(Note that these are not member functions.) | |
Indirectly applying constraints to LinearOperator | |
template<typename Range , typename Domain > | |
LinearOperator< Range, Domain > | distribute_constraints_linear_operator (const ConstraintMatrix &constraint_matrix, const LinearOperator< Range, Domain > &exemplar) |
template<typename Range , typename Domain > | |
LinearOperator< Range, Domain > | project_to_constrained_linear_operator (const ConstraintMatrix &constraint_matrix, const LinearOperator< Range, Domain > &exemplar) |
template<typename Range , typename Domain > | |
LinearOperator< Range, Domain > | constrained_linear_operator (const ConstraintMatrix &constraint_matrix, const LinearOperator< Range, Domain > &linop) |
template<typename Range , typename Domain > | |
PackagedOperation< Range > | constrained_right_hand_side (const ConstraintMatrix &constraint_matrix, const LinearOperator< Range, Domain > &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) |
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::TrilinosPayload > | linear_operator (const TrilinosWrappers::SparseMatrix &operator_exemplar, const Matrix &matrix) |
template<typename Range , typename Domain = Range> | |
LinearOperator< Range, Domain, TrilinosWrappers::internal::LinearOperatorImplementation::TrilinosPayload > | linear_operator (const TrilinosWrappers::SparseMatrix &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) |
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:
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
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:
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).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.
Definition at line 44 of file linear_operator.h.
|
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 163 of file linear_operator.h.
|
default |
Default copy constructor.
|
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 218 of file linear_operator.h.
|
default |
Default copy assignment operator.
|
inline |
Templated copy assignment operator for an object op
for which the conversion function linear_operator
is defined.
Definition at line 235 of file linear_operator.h.
|
inline |
Addition with a LinearOperator second_op
with the same Domain
and Range
.
Definition at line 293 of file linear_operator.h.
|
inline |
Subtraction with a LinearOperator second_op
with the same Domain
and Range
.
Definition at line 304 of file linear_operator.h.
|
inline |
Composition of the LinearOperator with an endomorphism second_op
of the Domain
space.
Definition at line 315 of file linear_operator.h.
|
inline |
Scalar multiplication of the LinearOperator with number
from the right.
Definition at line 326 of file linear_operator.h.
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 245 of file linear_operator.h.
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 251 of file linear_operator.h.
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 257 of file linear_operator.h.
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 263 of file linear_operator.h.
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 272 of file linear_operator.h.
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 281 of file linear_operator.h.
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 337 of file linear_operator.h.