elsa core¶
Table of Contents
DataDescriptor¶
-
class
elsa
::
DataDescriptor
: public elsa::Cloneable<DataDescriptor>¶ Base class for representing metadata for linearized n-dimensional signal stored in memory.
This class provides an interface for metadata about a signal that is stored in memory. This base class provides other descriptor subclasses with a fundamental interface to access the important parameters (i.e. dimensionality, size and spacing) of the metadata.
- Author
Matthias Wieczorek - initial code
- Author
Tobias Lasser - modularization, modernization
- Author
Maximilian Hornung - various enhancements
- Author
David Frank - inheritance restructuring
Subclassed by elsa::BlockDescriptor, elsa::DetectorDescriptor, elsa::VolumeDescriptor
Public Functions
-
DataDescriptor
() = delete¶ delete default constructor (having no metadata is invalid)
-
~DataDescriptor
() = 0¶ Pure virtual destructor.
-
DataDescriptor
(IndexVector_t numberOfCoefficientsPerDimension)¶ Constructor for DataDescriptor, accepts dimension and size.
- Parameters
[in] numberOfCoefficientsPerDimension
: vector containing the number of coefficients per dimension, (dimension is set implicitly from the size of the vector)
- Exceptions
std::invalid_argument
: if any number of coefficients is non-positive
-
DataDescriptor
(IndexVector_t numberOfCoefficientsPerDimension, RealVector_t spacingPerDimension)¶ Constructor for DataDescriptor, accepts dimension, size and spacing.
- Parameters
[in] numberOfCoefficientsPerDimension
: vector containing the number of coefficients per dimension, (dimension is set implicitly from the size of the vector)[in] spacingPerDimension
: vector containing the spacing per dimension
- Exceptions
std::invalid_argument
: if any number of coefficients is non-positive, or sizes of numberOfCoefficientsPerDimension and spacingPerDimension do not match
-
index_t
getNumberOfDimensions
() const¶ return the number of dimensions
-
index_t
getNumberOfCoefficients
() const¶ return the total number of coefficients
-
IndexVector_t
getNumberOfCoefficientsPerDimension
() const¶ return the number of coefficients per dimension
-
RealVector_t
getSpacingPerDimension
() const¶ return the spacing per dimension
-
RealVector_t
getLocationOfOrigin
() const¶ return the location of the origin (typically the center)
-
index_t
getIndexFromCoordinate
(IndexVector_t coordinate) const¶ computes the linearized index in the data vector from local coordinates
The local coordinates are integers, running from 0 to _numberOfCoefficientsPerDimension[i]-1 for every dimension i = 0,…,_numberOfDimensions. Linearization is assumed to be done in order of the dimensions.
- Return
the index into the linearized data vector
- Parameters
[in] coordinate
: vector containing the local coordinate
-
IndexVector_t
getCoordinateFromIndex
(index_t index) const¶ computes the local coordinates from a linearized index
The local coordinates are integers, running from 0 to _numberOfCoefficientsPerDimension[i]-1 for every dimension i = 0,…,_numberOfDimensions. Linearization is assumed to be done in order of the dimensions.
- Return
the local coordinate corresponding to the index
- Parameters
[in] index
: into the linearized data vector
Protected Functions
-
DataDescriptor
(const DataDescriptor&) = default¶ default copy constructor, hidden from non-derived classes to prevent potential slicing
-
DataDescriptor
(DataDescriptor&&) = default¶ default move constructor, hidden from non-derived classes to prevent potential slicing
-
bool
isEqual
(const DataDescriptor &other) const override¶ implement the polymorphic comparison operation
Protected Attributes
-
index_t
_numberOfDimensions
¶ Number of dimensions.
-
IndexVector_t
_numberOfCoefficientsPerDimension
¶ vector containing the number of coefficients per dimension
-
RealVector_t
_spacingPerDimension
¶ vector containing the spacing per dimension
-
IndexVector_t
_productOfCoefficientsPerDimension
¶ vector containing precomputed partial products of coefficients per dimension for index computations
-
RealVector_t
_locationOfOrigin
¶ vector containing the origin of the described volume (typically the center)
BlockDescriptor¶
-
class
elsa
::
BlockDescriptor
: public elsa::DataDescriptor¶ Abstract class defining the interface of all block descriptors.
A block descriptor provides metadata about a signal that is stored in memory (typically a
DataContainer). This signal can be n-dimensional, and will be stored in memory in a linearized fashion in blocks. The blocks can be used to support various operations (like blocked operators or ordered subsets), however, the blocks have to lie in memory one after the other (i.e. no stride is supported).- Author
Matthias Wieczorek - initial code
- Author
David Frank - rewrite
- Author
Tobias Lasser - rewrite, modularization, modernization
- Author
Nikola Dinev - rework into abstract class
Subclassed by elsa::IdenticalBlocksDescriptor, elsa::PartitionDescriptor, elsa::RandomBlocksDescriptor
Public Functions
-
~BlockDescriptor
() override = default¶ default destructor
-
index_t
getNumberOfBlocks
() const = 0¶ return the number of blocks
-
const DataDescriptor &
getDescriptorOfBlock
(index_t i) const = 0¶ return the DataDescriptor of the i-th block
-
index_t
getOffsetOfBlock
(index_t i) const = 0¶ return the offset to access the data of the i-th block
Protected Functions
-
BlockDescriptor
(DataDescriptor &&base)¶ used by derived classes to initialize the DataDescriptor base
-
BlockDescriptor
(const DataDescriptor &base)¶ used by derived classes to initialize the DataDescriptor base
IdenticalBlocksDescriptor¶
-
class
elsa
::
IdenticalBlocksDescriptor
: public elsa::BlockDescriptor¶ Class representing a series of identical descriptors concatenated along a new dimension (the last dimension of the full descriptor).
The blocks are, essentially, slices (though not necessarily two-dimensional) of the full descriptor along its last dimension. The last dimension of the full descriptor serves solely for the indexing of the different blocks, and will always have a spacing of one and a number of coefficients corresponding to the number of blocks.
- Author
Nikola Dinev
This descriptor should be the preferred choice when dealing with vector fields.
Public Functions
-
IdenticalBlocksDescriptor
(index_t numberOfBlocks, const DataDescriptor &dataDescriptor)¶ Create a new descriptor, replicating the dataDescriptor numberOfBlocks times along a new dimension.
- Parameters
[in] numberOfBlocks
: is the desired number of blocks[in] dataDescriptor
: is the descriptor that will be replicated numberOfBlocks times along a new dimension
- Exceptions
std::invalid_argument
: if numberOfBlocks is non-positive
-
IdenticalBlocksDescriptor
(const IdenticalBlocksDescriptor&) = delete¶ make copy constructor deletion explicit
-
~IdenticalBlocksDescriptor
() override = default¶ default destructor
-
index_t
getNumberOfBlocks
() const override¶ return the number of blocks
-
const DataDescriptor &
getDescriptorOfBlock
(index_t i) const override¶ return the DataDescriptor of the i-th block
-
index_t
getOffsetOfBlock
(index_t i) const override¶ return the offset to access the data of the i-th block
Protected Functions
-
IdenticalBlocksDescriptor *
cloneImpl
() const override¶ implement the polymorphic clone operation
-
bool
isEqual
(const DataDescriptor &other) const override¶ implement the polymorphic comparison operation
Protected Attributes
-
std::unique_ptr<DataDescriptor>
_blockDescriptor
¶ descriptor of a single block
-
index_t
_numberOfBlocks
¶ the total number of identical blocks
Private Functions
-
VolumeDescriptor
initBase
(index_t numberOfBlocks, const DataDescriptor &dataDescriptor)¶ generates the
PartitionDescriptor¶
-
class
elsa
::
PartitionDescriptor
: public elsa::BlockDescriptor¶ Class representing a descriptor obtained after the partition of a normal data descriptor into blocks.
A single block of the
PartitionDescriptor represents a linear segment containing one or more slices of the original descriptor, taken along its last dimension. This also means that the number of coefficients per block can only vary in the last dimension.- Author
Nikola Dinev
The PartitionDescriptor has the same dimension, the same number of coefficients and spacing per dimension as the original.
Public Functions
-
PartitionDescriptor
(const DataDescriptor &dataDescriptor, index_t numberOfBlocks)¶ Construct a PartitionDescriptor by partitioning a given descriptor into blocks of fairly equal sizes.
If the given descriptor has a size of
\( N \) in its last dimension, when dividing it into \( m \) blocks and \( N \) is not evenly divisible by \( m \), the last \( N \bmod m \) blocks will have a size of the last dimension one bigger than that of the others.- Parameters
[in] dataDescriptor
: the descriptor to be partitioned[in] numberOfBlocks
: the number of blocks
- Exceptions
std::invalid_argument
: if numberOfBlocks is less than 2 or greater than the number of coefficients in the last dimension
Note: if the passed in DataDescriptor is a block descriptor, the block information is ignored when generating the new PartitionDescriptor.
-
PartitionDescriptor
(const DataDescriptor &dataDescriptor, IndexVector_t slicesInBlock)¶ Construct a PartitionDescriptor by specifying the number of slices contained in each block.
Note: if the passed in
DataDescriptor is a block descriptor, the block information is ignored when generating the new PartitionDescriptor.- Parameters
[in] dataDescriptor
: the descriptor to be partitioned[in] slicesInBlock
: the number of slices in each block
- Exceptions
std::invalid_argument
: if slicesInBlock does not specify a valid partition scheme for the given descriptor
-
~PartitionDescriptor
() override = default¶ default destructor
-
index_t
getNumberOfBlocks
() const override¶ return the number of blocks
-
const DataDescriptor &
getDescriptorOfBlock
(index_t i) const override¶ return the DataDescriptor of the i-th block
-
index_t
getOffsetOfBlock
(index_t i) const override¶ return the offset to access the data of the i-th block
Protected Functions
-
PartitionDescriptor
(const PartitionDescriptor &other)¶ protected copy constructor; used for cloning
-
PartitionDescriptor *
cloneImpl
() const override¶ implement the polymorphic clone operation
-
bool
isEqual
(const DataDescriptor &other) const override¶ implement the polymorphic comparison operation
Protected Attributes
-
IndexVector_t
_indexMap
¶ maps a block index to the index of the corresponding descriptor in _blockDescriptors
-
std::vector<std::unique_ptr<DataDescriptor>>
_blockDescriptors
¶ vector of unique DataDescriptors describing the individual blocks
-
IndexVector_t
_blockOffsets
¶ vector of the individual block data offsets
Private Functions
-
std::unique_ptr<VolumeDescriptor>
generateDescriptorOfPartition
(index_t numberOfSlices) const¶ generates the descriptor of a partition containing numberOfSlices slices
RandomBlocksDescriptor¶
-
class
elsa
::
RandomBlocksDescriptor
: public elsa::BlockDescriptor¶ Class representing a block descriptor whose different blocks may have completely different descriptors.
There are no restrictions whatsoever imposed on the descriptors of different blocks. Different blocks may even have different number of dimensions.
- Author
Matthias Wieczorek - initial code
- Author
David Frank - rewrite
- Author
Nikola Dinev - various enhancements
- Author
Tobias Lasser - rewrite, modularization, modernization
The full descriptor will always be one-dimensional, and with a spacing of one. The size of it will be the sum of the sizes of all the descriptors, i.e. the sizes returned by DataDescriptor::getNumberOfCoefficients() for each descriptor in the list.
Public Functions
-
RandomBlocksDescriptor
(const std::vector<std::unique_ptr<DataDescriptor>> &blockDescriptors)¶ Construct a RandomBlocksDescriptor from a list of descriptors.
- Parameters
[in] blockDescriptors
: the list of descriptors of each block
- Exceptions
std::invalid_argument
: if the list is empty
-
RandomBlocksDescriptor
(std::vector<std::unique_ptr<DataDescriptor>> &&blockDescriptors)¶ Construct a RandomBlocksDescriptor from a list of descriptors.
- Parameters
[in] blockDescriptors
: the list of descriptors of each block
- Exceptions
std::invalid_argument
: if the list is empty
-
RandomBlocksDescriptor
(const RandomBlocksDescriptor&) = delete¶ make copy constructor deletion explicit
-
~RandomBlocksDescriptor
() override = default¶ default desctructor
-
index_t
getNumberOfBlocks
() const override¶ return the number of blocks
-
const DataDescriptor &
getDescriptorOfBlock
(index_t i) const override¶ return the DataDescriptor of the i-th block
-
index_t
getOffsetOfBlock
(index_t i) const override¶ return the offset to access the data of the i-th block
Protected Functions
-
RandomBlocksDescriptor *
cloneImpl
() const override¶ implement the polymorphic clone operation
-
bool
isEqual
(const DataDescriptor &other) const override¶ implement the polymorphic comparison operation
Protected Attributes
-
std::vector<std::unique_ptr<DataDescriptor>>
_blockDescriptors
¶ vector of DataDescriptors describing the individual blocks
-
IndexVector_t
_blockOffsets
¶ vector of the individual block data offsets
Private Functions
-
IndexVector_t
determineSize
(const std::vector<std::unique_ptr<DataDescriptor>> &blockDescriptors)¶ return the total number of coefficients of the descriptors in the list as an IndexVector
DataContainer¶
-
template<typename
data_t
>
classelsa
::
DataContainer
¶ class representing and storing a linearized n-dimensional signal
forward declaration for predicates and friend test function
This class provides a container for a signal that is stored in memory. This signal can be n-dimensional, and will be stored in memory in a linearized fashion. The information on how this linearization is performed is provided by an associated
DataDescriptor.- Author
Matthias Wieczorek - initial code
- Author
Tobias Lasser - rewrite, modularization, modernization
- Author
David Frank - added DataHandler concept, iterators
- Author
Nikola Dinev - add block support
- Author
Jens Petit - expression templates
- Template Parameters
data_t
: - data type that is stored in the DataContainer, defaulting to real_t.
Public Types
-
using
iterator
= DataContainerIterator<DataContainer<data_t>>¶ iterator for DataContainer (random access and continuous)
-
using
const_iterator
= ConstDataContainerIterator<DataContainer<data_t>>¶ const iterator for DataContainer (random access and continuous)
-
using
const_reverse_iterator
= std::reverse_iterator<const_iterator>¶ alias for const reverse iterator
-
using
value_type
= data_t¶ value_type of the DataContainer elements for iterators
-
using
pointer
= data_t*¶ pointer type of DataContainer elements for iterators
-
using
const_pointer
= const data_t*¶ const pointer type of DataContainer elements for iterators
-
using
reference
= data_t&¶ reference type of DataContainer elements for iterators
-
using
const_reference
= const data_t&¶ const reference type of DataContainer elements for iterators
-
using
difference_type
= std::ptrdiff_t¶ difference type for iterators
Public Functions
-
DataContainer
() = delete¶ delete default constructor (without metadata there can be no valid container)
-
DataContainer
(const DataDescriptor &dataDescriptor, DataHandlerType handlerType = defaultHandlerType)¶ Constructor for empty DataContainer, initializing the data with zeros.
- Parameters
[in] dataDescriptor
: containing the associated metadata[in] handlerType
: the data handler (default: CPU)
-
DataContainer
(const DataDescriptor &dataDescriptor, const Eigen::Matrix<data_t, Eigen::Dynamic, 1> &data, DataHandlerType handlerType = defaultHandlerType)¶ Constructor for DataContainer, initializing it with a DataVector.
- Parameters
[in] dataDescriptor
: containing the associated metadata[in] data
: vector containing the initialization data[in] handlerType
: the data handler (default: CPU)
-
DataContainer
(const DataContainer<data_t> &other)¶ Copy constructor for DataContainer.
- Parameters
[in] other
: DataContainer to copy
-
DataContainer<data_t> &
operator=
(const DataContainer<data_t> &other)¶ copy assignment for DataContainer
Note that a copy assignment with a
DataContainer on a different device (CPU vs GPU) will result in an “infectious” copy which means that afterwards the current container will use the same device as “other”.- Parameters
[in] other
: DataContainer to copy
-
DataContainer
(DataContainer<data_t> &&other) noexcept¶ Move constructor for DataContainer.
The moved-from objects remains in a valid state. However, as preconditions are not fulfilled for any member functions, the object should not be used. After move- or copy- assignment, this is possible again.
- Parameters
[in] other
: DataContainer to move from
-
DataContainer<data_t> &
operator=
(DataContainer<data_t> &&other)¶ Move assignment for DataContainer.
The moved-from objects remains in a valid state. However, as preconditions are not fulfilled for any member functions, the object should not be used. After move- or copy- assignment, this is possible again.
- Parameters
[in] other
: DataContainer to move from
Note that a copy assignment with a DataContainer on a different device (CPU vs GPU) will result in an “infectious” copy which means that afterwards the current container will use the same device as “other”.
-
template<typename
Source
, typename = std::enable_if_t<isExpression<Source>>>
DataContainer<data_t> &operator=
(Source const &source)¶ Expression evaluation assignment for DataContainer.
This evaluates an expression template term into the underlying data member of the
DataHandler in use.- Parameters
[in] source
: expression to evaluate
-
template<typename
Source
, typename = std::enable_if_t<isExpression<Source>>>DataContainer
(Source const &source)¶ Expression constructor.
It creates a new
DataContainer out of an expression. For this the meta information which is saved in the expression is used.- Parameters
[in] source
: expression to evaluate
-
const DataDescriptor &
getDataDescriptor
() const¶ return the current DataDescriptor
-
index_t
getSize
() const¶ return the size of the stored data (i.e. the number of elements in the linearized signal)
-
data_t &
operator[]
(index_t index)¶ return the index-th element of linearized signal (not bounds-checked!)
-
const data_t &
operator[]
(index_t index) const¶ return the index-th element of the linearized signal as read-only (not bounds-checked!)
-
data_t &
operator()
(IndexVector_t coordinate)¶ return an element by n-dimensional coordinate (not bounds-checked!)
-
const data_t &
operator()
(IndexVector_t coordinate) const¶ return an element by n-dimensional coordinate as read-only (not bounds-checked!)
-
template<typename
idx0_t
, typename ...idx_t
, typename = std::enable_if_t<std::is_integral_v<idx0_t> && (... && std::is_integral_v<idx_t>)>>
data_t &operator()
(idx0_t idx0, idx_t... indices)¶ return an element by its coordinates (not bounds-checked!)
-
template<typename
idx0_t
, typename ...idx_t
, typename = std::enable_if_t<std::is_integral_v<idx0_t> && (... && std::is_integral_v<idx_t>)>>
const data_t &operator()
(idx0_t idx0, idx_t... indices) const¶ return an element by its coordinates as read-only (not bounds-checked!)
-
data_t
dot
(const DataContainer<data_t> &other) const¶ return the dot product of this signal with the one from container other
-
template<typename
Source
, typename = std::enable_if_t<isExpression<Source>>>
data_tdot
(const Source &source) const¶ return the dot product of this signal with the one from an expression
-
GetFloatingPointType_t<data_t>
squaredL2Norm
() const¶ return the squared l2 norm of this signal (dot product with itself)
-
GetFloatingPointType_t<data_t>
l1Norm
() const¶ return the l1 norm of this signal (sum of absolute values)
-
GetFloatingPointType_t<data_t>
lInfNorm
() const¶ return the linf norm of this signal (maximum of absolute values)
-
DataContainer<data_t> &
operator+=
(const DataContainer<data_t> &dc)¶ compute in-place element-wise addition of another container
-
template<typename
Source
, typename = std::enable_if_t<isExpression<Source>>>
DataContainer<data_t> &operator+=
(Source const &source)¶ compute in-place element-wise addition with another expression
-
DataContainer<data_t> &
operator-=
(const DataContainer<data_t> &dc)¶ compute in-place element-wise subtraction of another container
-
template<typename
Source
, typename = std::enable_if_t<isExpression<Source>>>
DataContainer<data_t> &operator-=
(Source const &source)¶ compute in-place element-wise subtraction with another expression
-
DataContainer<data_t> &
operator*=
(const DataContainer<data_t> &dc)¶ compute in-place element-wise multiplication with another container
-
template<typename
Source
, typename = std::enable_if_t<isExpression<Source>>>
DataContainer<data_t> &operator*=
(Source const &source)¶ compute in-place element-wise multiplication with another expression
-
DataContainer<data_t> &
operator/=
(const DataContainer<data_t> &dc)¶ compute in-place element-wise division by another container
-
template<typename
Source
, typename = std::enable_if_t<isExpression<Source>>>
DataContainer<data_t> &operator/=
(Source const &source)¶ compute in-place element-wise division with another expression
-
DataContainer<data_t> &
operator+=
(data_t scalar)¶ compute in-place addition of a scalar
-
DataContainer<data_t> &
operator-=
(data_t scalar)¶ compute in-place subtraction of a scalar
-
DataContainer<data_t> &
operator*=
(data_t scalar)¶ compute in-place multiplication with a scalar
-
DataContainer<data_t> &
operator/=
(data_t scalar)¶ compute in-place division by a scalar
-
DataContainer<data_t> &
operator=
(data_t scalar)¶ assign a scalar to the DataContainer
-
bool
operator==
(const DataContainer<data_t> &other) const¶ comparison with another DataContainer
-
bool
operator!=
(const DataContainer<data_t> &other) const¶ comparison with another DataContainer
-
DataContainer<data_t>
getBlock
(index_t i)¶ returns a reference to the i-th block, wrapped in a DataContainer
-
const DataContainer<data_t>
getBlock
(index_t i) const¶ returns a const reference to the i-th block, wrapped in a DataContainer
-
DataContainer<data_t>
viewAs
(const DataDescriptor &dataDescriptor)¶ return a view of this DataContainer with a different descriptor
-
const DataContainer<data_t>
viewAs
(const DataDescriptor &dataDescriptor) const¶ return a const view of this DataContainer with a different descriptor
-
const_iterator
begin
() const¶ returns const iterator to the first element of the container (cannot mutate data)
-
const_iterator
cbegin
() const¶ returns const iterator to the first element of the container (cannot mutate data)
-
const_iterator
end
() const¶ returns const iterator to one past the last element of the container (cannot mutate data)
-
const_iterator
cend
() const¶ returns const iterator to one past the last element of the container (cannot mutate data)
-
reverse_iterator
rbegin
()¶ returns reversed iterator to the last element of the container
-
const_reverse_iterator
rbegin
() const¶ returns const reversed iterator to the last element of the container (cannot mutate data)
-
const_reverse_iterator
crbegin
() const¶ returns const reversed iterator to the last element of the container (cannot mutate data)
-
reverse_iterator
rend
()¶ returns reversed iterator to one past the first element of container
-
const_reverse_iterator
rend
() const¶ returns const reversed iterator to one past the first element of container (cannot mutate data)
-
const_reverse_iterator
crend
() const¶ returns const reversed iterator to one past the first element of container (cannot mutate data)
-
DataHandlerType
getDataHandlerType
() const¶ returns the type of the DataHandler in use
-
DataContainer
loadToGPU
()¶ Factory function which returns GPU based DataContainers.
Note that if this function is called on a container which is already GPU based, it will throw an exception.
- Return
the GPU based DataContainer
-
DataContainer
loadToCPU
()¶ Factory function which returns CPU based DataContainers.
Note that if this function is called on a container which is already CPU based, it will throw an exception.
- Return
the CPU based DataContainer
Private Functions
-
template<typename ...
Args
>
std::unique_ptr<DataHandler<data_t>>createDataHandler
(DataHandlerType handlerType, Args&&... args)¶ factory method to create DataHandlers based on handlerType with perfect forwarding of constructor arguments
-
DataContainer
(const DataDescriptor &dataDescriptor, std::unique_ptr<DataHandler<data_t>> dataHandler, DataHandlerType dataType = defaultHandlerType)¶ private constructor accepting a DataDescriptor and a DataHandler
-
bool
canAssign
(DataHandlerType handlerType)¶ Helper function to indicate if a regular assignment or a clone should be performed.
An assignment operation with a
DataContainer which does not use the same device (CPU / GPU) has to be handled differently. This helper function indicates if a regular assignment should be performed or not.- Return
true if a regular assignment of the pointed to DataHandlers should be done
- Parameters
[in] handlerType
: the member variable of the other container in copy-/move-assignment
Private Members
-
std::unique_ptr<DataDescriptor>
_dataDescriptor
¶ the current DataDescriptor
-
std::unique_ptr<DataHandler<data_t>>
_dataHandler
¶ the current DataHandler
-
DataHandlerType
_dataHandlerType
¶ the current DataHandlerType
LinearOperator¶
-
template<typename
data_t
= real_t>
classelsa
::
LinearOperator
: public elsa::Cloneable<LinearOperator<real_t>>¶ Base class representing a linear operator A. Also implements operator expression functionality.
This class represents a linear operator A, expressed through its apply/applyAdjoint methods, which implement Ax and A^ty for DataContainers x,y of appropriate sizes. Concrete implementations of linear operators will derive from this class and override the applyImpl/applyAdjointImpl methods.
- Author
Matthias Wieczorek - initial code
- Author
Maximilian Hornung - composite rewrite
- Author
Tobias Lasser - rewrite, modularization, modernization
- Template Parameters
data_t
: data type for the domain and range of the operator, defaulting to real_t
LinearOperator also provides functionality to support constructs like the operator expression A^t*B+C, where A,B,C are linear operators. This operator composition is implemented via evaluation trees.
LinearOperator and all its derived classes are expected to be light-weight and easily copyable/cloneable, due to the implementation of evaluation trees. Hence any pre-computations/caching should only be done in a lazy manner (e.g. during the first call of apply), and not in the constructor.
Public Functions
-
LinearOperator
(const DataDescriptor &domainDescriptor, const DataDescriptor &rangeDescriptor)¶ Constructor for the linear operator A, mapping from domain to range.
- Parameters
[in] domainDescriptor
: DataDescriptor describing the domain of the operator[in] rangeDescriptor
: DataDescriptor describing the range of the operator
-
~LinearOperator
() override = default¶ default destructor
-
LinearOperator
(const LinearOperator<data_t> &other)¶ copy construction
-
LinearOperator<data_t> &
operator=
(const LinearOperator<data_t> &other)¶ copy assignment
-
const DataDescriptor &
getDomainDescriptor
() const¶ return the domain DataDescriptor
-
const DataDescriptor &
getRangeDescriptor
() const¶ return the range DataDescriptor
-
DataContainer<data_t>
apply
(const DataContainer<data_t> &x) const¶ apply the operator A to an element in the operator’s domain
Please note: this method uses apply(x, Ax) to perform the actual operation.
- Return
Ax DataContainer containing the application of operator A to data x, i.e. in the range of the operator
- Parameters
[in] x
: input DataContainer (in the domain of the operator)
-
void
apply
(const DataContainer<data_t> &x, DataContainer<data_t> &Ax) const¶ apply the operator A to an element in the operator’s domain
Please note: this method calls the method applyImpl that has to be overridden in derived classes. (Why is this method not virtual itself? Because you cannot have a non-virtual function overloading a virtual one [apply with one vs. two arguments]).
- Parameters
[in] x
: input DataContainer (in the domain of the operator)[out] Ax
: output DataContainer (in the range of the operator)
-
DataContainer<data_t>
applyAdjoint
(const DataContainer<data_t> &y) const¶ apply the adjoint of operator A to an element of the operator’s range
Please note: this method uses applyAdjoint(y, Aty) to perform the actual operation.
- Return
A^ty DataContainer containing the application of A^t to data y, i.e. in the domain of the operator
- Parameters
[in] y
: input DataContainer (in the range of the operator)
-
void
applyAdjoint
(const DataContainer<data_t> &y, DataContainer<data_t> &Aty) const¶ apply the adjoint of operator A to an element of the operator’s range
Please note: this method calls the method applyAdjointImpl that has to be overridden in derived classes. (Why is this method not virtual itself? Because you cannot have a non-virtual function overloading a virtual one [applyAdjoint with one vs. two args]).
- Parameters
[in] y
: input DataContainer (in the range of the operator)[out] Aty
: output DataContainer (in the domain of the operator)
Protected Functions
-
LinearOperator<data_t> *
cloneImpl
() const override¶ implement the polymorphic clone operation
-
bool
isEqual
(const LinearOperator<data_t> &other) const override¶ implement the polymorphic comparison operation
-
void
applyImpl
(const DataContainer<data_t> &x, DataContainer<data_t> &Ax) const¶ the apply method that has to be overridden in derived classes
-
void
applyAdjointImpl
(const DataContainer<data_t> &y, DataContainer<data_t> &Aty) const¶ the applyAdjoint method that has to be overridden in derived classes
Protected Attributes
-
std::unique_ptr<DataDescriptor>
_domainDescriptor
¶ the data descriptor of the domain of the operator
-
std::unique_ptr<DataDescriptor>
_rangeDescriptor
¶ the data descriptor of the range of the operator
Private Types
Private Functions
-
LinearOperator
(const LinearOperator<data_t> &op, bool isAdjoint)¶ constructor to produce an adjoint leaf node
-
LinearOperator
(const LinearOperator<data_t> &lhs, const LinearOperator<data_t> &rhs, CompositeMode mode)¶ constructor to produce a composite (internal node) of the evaluation tree
Private Members
-
std::unique_ptr<LinearOperator<data_t>>
_lhs
= {}¶ pointers to nodes in the evaluation tree
-
bool
_isLeaf
= {false}¶ flag whether this is a leaf-node
-
bool
_isAdjoint
= {false}¶ flag whether this is a leaf-node to implement an adjoint operator
-
bool
_isComposite
= {false}¶ flag whether this is a composite (internal node) of the evaluation tree
-
CompositeMode
_mode
= {CompositeMode::MULT}¶ variable storing the composition mode (+, *)
Friends
-
friend LinearOperator<data_t>
operator+
(const LinearOperator<data_t> &lhs, const LinearOperator<data_t> &rhs)¶ friend operator+ to support composition of LinearOperators (and its derivatives)
-
friend LinearOperator<data_t>
operator*
(const LinearOperator<data_t> &lhs, const LinearOperator<data_t> &rhs)¶ friend operator* to support composition of LinearOperators (and its derivatives)
-
friend LinearOperator<data_t>
adjoint
(const LinearOperator<data_t> &op)¶ friend function to return the adjoint of a LinearOperator (and its derivatives)
-
friend LinearOperator<data_t>
leaf
(const LinearOperator<data_t> &op)¶ friend function to return a leaf node of a LinearOperator (and its derivatives)
Implementation Details¶
Cloneable¶
-
template<typename
Derived
>
classelsa
::
Cloneable
¶ Class implementing polymorphic clones with smart pointers and CRTP, as well as comparison operators.
This class provides a clone method using CRTP to support covariance with smart pointers. For details see
https://www.fluentcpp.com/2017/09/12/how-to-return-a-smart-pointer-and-use-covariance/.- Author
Tobias Lasser
Public Functions
-
Cloneable
() = default¶ default constructor
-
~Cloneable
() = default¶ virtual default destructor
Protected Functions
-
Derived *
cloneImpl
() const = 0¶ actual clone implementation method, abstract to force override in derived classes
-
bool
isEqual
(const Derived &other) const = 0¶ actual comparison method, abstract to force override in derived classes
-
Cloneable
(const Cloneable&) = default¶ default copy constructor, protected to not be publicly available (but available for cloneImpl()
DataHandler¶
-
template<typename
data_t
= real_t>
classelsa
::
DataHandler
: public elsa::Cloneable<DataHandler<real_t>>¶ Base class encapsulating data handling. The data is stored transparently, for example on CPU or GPU.
This abstract base class serves as an interface for data handlers, which encapsulate the actual data being stored e.g. in main memory of the CPU or in various memory types of GPUs. The data itself is treated as a vector, i.e. an array of data_t elements (which usually comes from linearized n-dimensional signals).
- Author
David Frank - initial code
- Author
Tobias Lasser - modularization, modernization
- Author
Nikola Dinev - add block support
Caveat: If data is not stored in main memory (e.g. on GPUs), then some operations may trigger an automatic synchronization of GPU to main memory. Please see the GPU-based handlers’ documentation for details.
Public Types
Public Functions
-
index_t
getSize
() const = 0¶ return the size of the stored data (i.e. number of elements in linearized data vector)
-
data_t &
operator[]
(index_t index) = 0¶ return the index-th element of the data vector (not bounds-checked!)
-
const data_t &
operator[]
(index_t index) const = 0¶ return the index-th element of the data vector as read-only (not bounds-checked!)
-
data_t
dot
(const DataHandler<data_t> &v) const = 0¶ return the dot product of the data vector with vector v
-
GetFloatingPointType_t<data_t>
squaredL2Norm
() const = 0¶ return the squared l2 norm of the data vector (dot product with itself)
-
GetFloatingPointType_t<data_t>
l1Norm
() const = 0¶ return the l1 norm of the data vector (sum of absolute values)
-
GetFloatingPointType_t<data_t>
lInfNorm
() const = 0¶ return the linf norm of the data vector (maximum of absolute values)
-
DataHandler<data_t> &
operator+=
(const DataHandler<data_t> &v) = 0¶ compute in-place element-wise addition of another vector v
-
DataHandler<data_t> &
operator-=
(const DataHandler<data_t> &v) = 0¶ compute in-place element-wise subtraction of another vector v
-
DataHandler<data_t> &
operator*=
(const DataHandler<data_t> &v) = 0¶ compute in-place element-wise multiplication by another vector v
-
DataHandler<data_t> &
operator/=
(const DataHandler<data_t> &v) = 0¶ compute in-place element-wise division by another vector v
-
DataHandler<data_t> &
operator+=
(data_t scalar) = 0¶ compute in-place addition of a scalar
-
DataHandler<data_t> &
operator-=
(data_t scalar) = 0¶ compute in-place subtraction of a scalar
-
DataHandler<data_t> &
operator*=
(data_t scalar) = 0¶ compute in-place multiplication by a scalar
-
DataHandler<data_t> &
operator/=
(data_t scalar) = 0¶ compute in-place division by a scalar
-
DataHandler<data_t> &
operator=
(data_t scalar) = 0¶ assign a scalar to all elements of the data vector
-
DataHandler<data_t> &
operator=
(const DataHandler<data_t> &other)¶ copy assignment operator
-
DataHandler<data_t> &
operator=
(DataHandler<data_t> &&other)¶ move assignment operator
-
std::unique_ptr<DataHandler<data_t>>
getBlock
(index_t startIndex, index_t numberOfElements) = 0¶ return a reference to the sequential block starting at startIndex and containing numberOfElements elements
-
std::unique_ptr<const DataHandler<data_t>>
getBlock
(index_t startIndex, index_t numberOfElements) const = 0¶ return a const reference to the sequential block starting at startIndex and containing numberOfElements elements
Protected Types
-
using
DataVector_t
= Eigen::Matrix<data_t, Eigen::Dynamic, 1>¶ convenience typedef for the Eigen::Matrix data vector
-
using
DataMap_t
= Eigen::Map<DataVector_t>¶ convenience typedef for the Eigen::Map
Protected Functions
-
data_t
slowDotProduct
(const DataHandler<data_t> &v) const¶ slow element-wise dot product fall-back for when DataHandler types do not match
-
void
slowAddition
(const DataHandler<data_t> &v)¶ slow element-wise addition fall-back for when DataHandler types do not match
-
void
slowSubtraction
(const DataHandler<data_t> &v)¶ slow element-wise subtraction fall-back for when DataHandler types do not match
-
void
slowMultiplication
(const DataHandler<data_t> &v)¶ slow element-wise multiplication fall-back for when DataHandler types do not match
-
void
slowDivision
(const DataHandler<data_t> &v)¶ slow element-wise division fall-back for when DataHandler types do not match
-
void
slowAssign
(const DataHandler<data_t> &other)¶ slow element-wise assignment fall-back for when DataHandler types do not match
-
void
assign
(const DataHandler<data_t> &other) = 0¶ derived classes should override this method to implement copy assignment
-
void
assign
(DataHandler<data_t> &&other) = 0¶ derived classes should override this method to implement move assignment
DataHandlerCPU¶
-
template<typename
data_t
>
classelsa
::
DataHandlerCPU
¶ Class representing and owning a vector stored in CPU main memory (using Eigen::Matrix).
forward declaration, allows mutual friending
The class implements copy-on-write. Therefore any non-const functions should call the detach() function first to trigger the copy-on-write mechanism.
- Author
David Frank - main code
- Author
Tobias Lasser - modularization and modernization
- Author
Nikola Dinev - integration of map and copy-on-write concepts
- Template Parameters
data_t
: - data type that is stored, defaulting to real_t.
DataHandlerCPU and DataHandlerMapCPU are mutual friend classes allowing for the vectorization of arithmetic operations with the help of Eigen. A strong bidirectional link exists between the two classes. A Map is associated with the DataHandlerCPU from which it was created for the entirety of its lifetime. If the DataHandlerCPU starts managing a new vector (e.g. through a call to detach()), all associated Maps will also be updated.
Public Functions
-
DataHandlerCPU
() = delete¶ delete default constructor (having no information makes no sense)
-
~DataHandlerCPU
() override¶ default destructor
-
DataHandlerCPU
(index_t size)¶ Constructor initializing an appropriately sized vector.
Note that the values will not be initialized and therefore contain undefined values.
- Parameters
[in] size
: of the vector
- Exceptions
std::invalid_argument
: if the size is non-positive
-
DataHandlerCPU
(DataVector_t const &vector)¶ Constructor initializing a data vector with a given vector.
- Parameters
[in] vector
: that is used for initializing the data
-
DataHandlerCPU
(const DataHandlerCPU<data_t> &other)¶ copy constructor
-
DataHandlerCPU
(DataHandlerCPU<data_t> &&other)¶ move constructor
-
index_t
getSize
() const override¶ return the size of the vector
-
data_t &
operator[]
(index_t index) override¶ return the index-th element of the data vector (not bounds checked!)
-
const data_t &
operator[]
(index_t index) const override¶ return the index-th element of the data vector as read-only (not bound checked!)
-
data_t
dot
(const DataHandler<data_t> &v) const override¶ return the dot product of the data vector with vector v
-
GetFloatingPointType_t<data_t>
squaredL2Norm
() const override¶ return the squared l2 norm of the data vector (dot product with itself)
-
GetFloatingPointType_t<data_t>
l1Norm
() const override¶ return the l1 norm of the data vector (sum of absolute values)
-
GetFloatingPointType_t<data_t>
lInfNorm
() const override¶ return the linf norm of the data vector (maximum of absolute values)
-
DataHandlerCPU<data_t> &
operator=
(const DataHandlerCPU<data_t> &v)¶ copy assign another DataHandlerCPU to this, other types handled in assign()
-
DataHandlerCPU<data_t> &
operator=
(DataHandlerCPU<data_t> &&v)¶ move assign another DataHandlerCPU to this, other types handled in assign()
-
DataHandler<data_t> &
operator+=
(const DataHandler<data_t> &v) override¶ compute in-place element-wise addition of another vector v
-
DataHandler<data_t> &
operator-=
(const DataHandler<data_t> &v) override¶ compute in-place element-wise subtraction of another vector v
-
DataHandler<data_t> &
operator*=
(const DataHandler<data_t> &v) override¶ compute in-place element-wise multiplication by another vector v
-
DataHandler<data_t> &
operator/=
(const DataHandler<data_t> &v) override¶ compute in-place element-wise division by another vector v
-
DataHandler<data_t> &
operator+=
(data_t scalar) override¶ compute in-place addition of a scalar
-
DataHandler<data_t> &
operator-=
(data_t scalar) override¶ compute in-place subtraction of a scalar
-
DataHandler<data_t> &
operator*=
(data_t scalar) override¶ compute in-place multiplication by a scalar
-
DataHandler<data_t> &
operator/=
(data_t scalar) override¶ compute in-place division by a scalar
-
DataHandler<data_t> &
operator=
(data_t scalar) override¶ assign a scalar to all elements of the data vector
-
std::unique_ptr<DataHandler<data_t>>
getBlock
(index_t startIndex, index_t numberOfElements) override¶ return a reference to the sequential block starting at startIndex and containing numberOfElements elements
-
std::unique_ptr<const DataHandler<data_t>>
getBlock
(index_t startIndex, index_t numberOfElements) const override¶ return a const reference to the sequential block starting at startIndex and containing numberOfElements elements
Protected Types
-
using
DataVector_t
= Eigen::Matrix<data_t, Eigen::Dynamic, 1>¶ convenience typedef for the Eigen::Matrix data vector
-
using
DataMap_t
= Eigen::Map<DataVector_t>¶ convenience typedef for the Eigen::Map
Protected Functions
-
DataHandlerCPU<data_t> *
cloneImpl
() const override¶ implement the polymorphic clone operation
-
bool
isEqual
(const DataHandler<data_t> &other) const override¶ implement the polymorphic comparison operation
-
void
assign
(const DataHandler<data_t> &other) override¶ copy the data stored in other
-
void
assign
(DataHandler<data_t> &&other) override¶ move the data stored in other if other is of the same type, otherwise copy the data
Protected Attributes
-
std::shared_ptr<DataVector_t>
_data
¶ the vector storing the data
Private Functions
-
void
detach
()¶ creates the deep copy for the copy-on-write mechanism
-
void
detachWithUninitializedBlock
(index_t startIndex, index_t numberOfElements)¶ same as detach() but leaving an uninitialized block of numberOfElements elements starting at index startIndex
change the vector being handled
change the vector being handled (rvalue version)
Private Members
-
friend DataHandlerMapCPU< data_t >
declare DataHandlerMapCPU as friend, allows the use of Eigen for improved performance
-
friend DataContainer< data_t >
for enabling accessData()
Friends
-
friend long
useCount
(const DataHandlerCPU<data_t> &dh)¶ used for testing only and defined in test file
-
template<bool
GPU
, classOperand
, std::enable_if_t<isDataContainer<Operand>, int>>
friend constexpr autoevaluateOrReturn
(Operand const &operand)¶ friend constexpr function to implement expression templates
friend constexpr function to implement expression templates
Compile time switch to return data in container.
DataHandlerGPU¶
-
template<typename
data_t
>
classelsa
::
DataHandlerGPU
¶ Class representing and owning a vector stored in CPU main memory (using Eigen::Matrix).
forward declaration, allows mutual friending
The class implements copy-on-write. Therefore any non-const functions should call the detach() function first to trigger the copy-on-write mechanism.
- Author
David Frank - main code
- Author
Tobias Lasser - modularization and modernization
- Author
Nikola Dinev - integration of map and copy-on-write concepts
- Template Parameters
data_t
: - data type that is stored, defaulting to real_t.
DataHandlerGPU and DataHandlerMapCPU are mutual friend classes allowing for the vectorization of arithmetic operations with the help of Eigen. A strong bidirectional link exists between the two classes. A Map is associated with the DataHandlerGPU from which it was created for the entirety of its lifetime. If the DataHandlerGPU starts managing a new vector (e.g. through a call to detach()), all associated Maps will also be updated.
Public Functions
-
DataHandlerGPU
() = delete¶ delete default constructor (having no information makes no sense)
-
~DataHandlerGPU
() override¶ default destructor
-
DataHandlerGPU
(index_t size)¶ Constructor initializing an appropriately sized vector with zeros.
- Parameters
[in] size
: of the vector[in] initialize
: - set to false if you do not need initialization with zeros (default: true)
- Exceptions
std::invalid_argument
: if the size is non-positive
-
DataHandlerGPU
(DataVector_t const &vector)¶ Constructor initializing a data vector with a given vector.
- Parameters
[in] vector
: that is used for initializing the data
-
DataHandlerGPU
(quickvec::Vector<data_t> const &vector)¶ Constructor initializing a data vector with a given vector.
- Parameters
[in] vector
: that is used for initializing the data
-
DataHandlerGPU
(const DataHandlerGPU<data_t> &other)¶ copy constructor
-
DataHandlerGPU
(DataHandlerGPU<data_t> &&other) noexcept¶ move constructor
-
index_t
getSize
() const override¶ return the size of the vector
-
data_t &
operator[]
(index_t index) override¶ return the index-th element of the data vector (not bounds checked!)
-
const data_t &
operator[]
(index_t index) const override¶ return the index-th element of the data vector as read-only (not bound checked!)
-
data_t
dot
(const DataHandler<data_t> &v) const override¶ return the dot product of the data vector with vector v
-
GetFloatingPointType_t<data_t>
squaredL2Norm
() const override¶ return the squared l2 norm of the data vector (dot product with itself)
-
GetFloatingPointType_t<data_t>
l1Norm
() const override¶ return the l1 norm of the data vector (sum of absolute values)
-
GetFloatingPointType_t<data_t>
lInfNorm
() const override¶ return the linf norm of the data vector (maximum of absolute values)
-
DataHandlerGPU<data_t> &
operator=
(const DataHandlerGPU<data_t> &v)¶ copy assign another DataHandlerGPU
-
DataHandlerGPU<data_t> &
operator=
(DataHandlerGPU<data_t> &&v)¶ move assign another DataHandlerGPU
-
DataHandler<data_t> &
operator+=
(const DataHandler<data_t> &v) override¶ compute in-place element-wise addition of another vector v
-
DataHandler<data_t> &
operator-=
(const DataHandler<data_t> &v) override¶ compute in-place element-wise subtraction of another vector v
-
DataHandler<data_t> &
operator*=
(const DataHandler<data_t> &v) override¶ compute in-place element-wise multiplication by another vector v
-
DataHandler<data_t> &
operator/=
(const DataHandler<data_t> &v) override¶ compute in-place element-wise division by another vector v
-
DataHandler<data_t> &
operator+=
(data_t scalar) override¶ compute in-place addition of a scalar
-
DataHandler<data_t> &
operator-=
(data_t scalar) override¶ compute in-place subtraction of a scalar
-
DataHandler<data_t> &
operator*=
(data_t scalar) override¶ compute in-place multiplication by a scalar
-
DataHandler<data_t> &
operator/=
(data_t scalar) override¶ compute in-place division by a scalar
-
DataHandler<data_t> &
operator=
(data_t scalar) override¶ assign a scalar to all elements of the data vector
-
std::unique_ptr<DataHandler<data_t>>
getBlock
(index_t startIndex, index_t numberOfElements) override¶ return a reference to the sequential block starting at startIndex and containing numberOfElements elements
-
std::unique_ptr<const DataHandler<data_t>>
getBlock
(index_t startIndex, index_t numberOfElements) const override¶ return a const reference to the sequential block starting at startIndex and containing numberOfElements elements
Protected Types
-
using
DataVector_t
= Eigen::Matrix<data_t, Eigen::Dynamic, 1>¶ convenience typedef for the Eigen::Matrix data vector
-
using
DataMap_t
= Eigen::Map<DataVector_t>¶ convenience typedef for the Eigen::Map
Protected Functions
-
DataHandlerGPU<data_t> *
cloneImpl
() const override¶ implement the polymorphic clone operation
-
bool
isEqual
(const DataHandler<data_t> &other) const override¶ implement the polymorphic comparison operation
-
void
assign
(const DataHandler<data_t> &other) override¶ copy the data stored in other
-
void
assign
(DataHandler<data_t> &&other) override¶ move the data stored in other if other is of the same type, otherwise copy the data
Protected Attributes
Private Functions
-
void
detach
()¶ creates the deep copy for the copy-on-write mechanism
-
void
detachWithUninitializedBlock
(index_t startIndex, index_t numberOfElements)¶ same as detach() but leaving an uninitialized block of numberOfElements elements starting at index startIndex
change the vector being handled
change the vector being handled (rvalue version)
Private Members
-
friend DataContainer< data_t >
for enabling accessData()
-
friend DataHandlerMapGPU< data_t >
declare DataHandlerMapGPU as friend, allows the use of Eigen for improved performance
Friends
-
friend long
useCount
(const DataHandlerGPU<data_t> &dh)¶ used for testing only and defined in test file
-
template<bool
GPU
, classOperand
, std::enable_if_t<isDataContainer<Operand>, int>>
friend constexpr autoevaluateOrReturn
(Operand const &operand)¶ friend constexpr function to implement expression templates
friend constexpr function to implement expression templates
Compile time switch to return data in container.
Expression¶
In elsa, we are using expression templates for fast and efficient mathematical operations while at the same time having intuitive syntax. This technique is also known as lazy-evaluation which means that computations are delayed until the results are actually needed.
Take the example
auto expression = dataContainer1 + dataContainer2;
where the type of expression
is elsa::Expression<operator+, DataContainer<real_t>, DataContainer<real_t>>
.
The expression
contains the work to be done rather than the actual result. Only when assigning to a DataContainer
or constructing a new DataContainer, the expression gets evaluated automatically
DataContainer result = dataContainer1 + dataContainer2;
into result
.
Nesting is easily possible like
auto expression = dataContainer1 + dataContainer2 * expression;
which will store an expression tree as the type information.
Please note also that the operators/member functions available on an Expression
type are different from a DataContainer
.
An expression has a member function eval()
. However, calling eval()
will not return a DataContainer
but depending on the current DataHandler
the underlying raw data computation result.
If single element-wise access is necessary, it is possible to call expression.eval()[index]
. Note that this is computational very expensive as the whole expression gets evaluated at every index individually.
-
template<typename
Callable
, typename ...Operands
>
classelsa
::
Expression
¶ Temporary expression type which enables lazy-evaluation of expression.
Forward declaration for predicates.
- Author
Jens Petit
- Template Parameters
Callable
: - the operation to be performedOperands
: - the objects on which the operation is performed
Public Types
-
using
MetaInfo_t
= std::pair<DataDescriptor const&, DataHandlerType>¶ type which bundles the meta information to create a new DataContainer
Public Functions
-
template<bool
GPU
= false>
autoeval
() const¶ Evaluates the expression.
Private Functions
-
MetaInfo_t
initDescriptor
(Operands const&... args)¶ correctly returns the DataContainer descriptor based on the operands (either expressions or Datacontainers)
-
std::optional<MetaInfo_t>
getMetaInfoFromContainers
()¶ base recursive case if no DataContainer as operand
-
template<class
T
, class ...Ts
>
std::optional<MetaInfo_t>getMetaInfoFromContainers
(T &arg, Ts&... args)¶ recursive traversal of all contained DataContainers
-
std::optional<MetaInfo_t>
getMetaInfoFromExpressions
()¶ base recursive case if no Expression as operand
-
template<class
T
, class ...Ts
>
std::optional<MetaInfo_t>getMetaInfoFromExpressions
(T &arg, Ts&... args)¶ recursive traversal of all contained Expressions
Private Members
-
std::tuple<typename ReferenceOrNot<Operands>::type...>
_args
¶ Contains all operands saved as const references (DataContainers) or copies (Expressions and arithmetic types)
-
const MetaInfo_t
_dataMetaInfo
¶ saves the meta information to create a new DataContainer out of an expression