elsa core¶
Table of Contents
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
Tobias Lasser - rewrite, modularization, modernization
David Frank - added DataHandler concept, iterators, integrated unified memory
Nikola Dinev - add block support
Jens Petit - expression templates
Jonas Jelten - various enhancements, fft, complex handling, pretty formatting
- Template Parameters
data_t
: data type that is stored in the DataContainer, defaulting to real_t.
Public Types
-
using
iterator
= typename ContiguousStorage<data_t>::iterator¶ iterator for DataContainer (random access and continuous)
-
using
const_iterator
= typename ContiguousStorage<data_t>::const_iterator¶ const iterator for DataContainer (random access and continuous)
-
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
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)¶ Constructor for empty DataContainer, no initialisation is performed, but the underlying space is allocated.
- Parameters
[in] dataDescriptor
: containing the associated metadata
-
DataContainer
(const DataDescriptor &dataDescriptor, const Eigen::Matrix<data_t, Eigen::Dynamic, 1> &data)¶ Constructor for DataContainer, initializing it with a DataVector.
- Parameters
[in] dataDescriptor
: containing the associated metadata[in] data
: vector containing the initialization data
-
DataContainer
(const DataDescriptor &dataDescriptor, const ContiguousStorage<data_t> &storage)¶ constructor accepting a DataDescriptor and a DataHandler
-
DataContainer
(const DataDescriptor &dataDescriptor, ContiguousStorageView<data_t> storage)¶ constructor accepting a DataDescriptor and a DataHandler
-
template<mr::StorageType
tag
>DataContainer
(const DataDescriptor &dataDescriptor, NdViewTagged<data_t, tag> &view)¶ constructor accepting a DataDescriptor and an NdView. Dimensions in the DataDescriptor must match the dimensions of the NdView.
-
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) noexcept¶ 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”.
-
const DataDescriptor &
getDataDescriptor
() const¶ return the current DataDescriptor
-
bool
isOwning
() const¶ return true, if the current DataContainer is owning its memory
-
bool
isView
() const¶ return true, if the current DataContainer is a view, i.e. is not owning its memory
-
index_t
getSize
() const¶ return the size of the stored data (i.e. the number of elements in the linearized signal)
-
index_t
getNumberOfBlocks
() const¶ get the number of blocks the signal is created from. If the descriptor is not of type
BlockDescriptor
1 is returned;
-
reference
operator[]
(index_t index)¶ return the index-th element of linearized signal (not bounds-checked!)
-
const_reference
operator[]
(index_t index) const¶ return the index-th element of the linearized signal as read-only (not bounds-checked!)
-
reference
operator()
(const IndexVector_t &coordinate)¶ return an element by n-dimensional coordinate (not bounds-checked!) The indexing follows the convention
(x1, x2, ..., xn)
. Specifically, in 2D(x, y)
and 3D(x, y, z)
.
-
const_reference
operator()
(const IndexVector_t &coordinate) const¶ return an element by n-dimensional coordinate as read-only (not bounds-checked!) The indexing follows the convention
(x1, x2, ..., xn)
. Specifically, in 2D(x, y)
and 3D(x, y, z)
.
-
template<typename
idx0_t
, typename ...idx_t
, typename = std::enable_if_t<std::is_integral_v<idx0_t> && (... && std::is_integral_v<idx_t>)>>
referenceoperator()
(idx0_t idx0, idx_t... indices)¶ return an element by its coordinates (not bounds-checked!) The indexing follows the convention
(x1, x2, ..., xn)
. Specifically, in 2D(x, y)
and 3D(x, y, z)
.
-
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_referenceoperator()
(idx0_t idx0, idx_t... indices) const¶ return an element by its coordinates as read-only (not bounds-checked!) The indexing follows the convention
(x1, x2, ..., xn)
. Specifically, in 2D(x, y)
and 3D(x, y, z)
.
-
data_t
dot
(const DataContainer<data_t> &other) const¶ return the dot product of this signal with the one from container other
-
GetFloatingPointType_t<data_t>
squaredL2Norm
() const¶ return the squared l2 norm of this signal (dot product with itself)
-
GetFloatingPointType_t<data_t>
l2Norm
() const¶ return the l2 norm of this signal (square root of dot product with itself)
-
DataContainer<GetFloatingPointType_t<data_t>>
pL2Norm
() const¶ return the pointwise l2 norm of this signal. Pointwise norms assume blocked DataContainers, and then the norm is taken column-wise. The pointwise l2 norm is equivalent to the following NumPy code snippet:
np.sqrt(np.sum(x ** 2, axis=0))
, assuming a signal, where the first axis indexes each block of the data.
-
GetFloatingPointType_t<data_t>
l1Norm
() const¶ return the l1 norm of this signal (sum of absolute values)
-
DataContainer<GetFloatingPointType_t<data_t>>
pL1Norm
() const¶ return the pointwise l1 norm of this signal. Pointwise norms assume blocked DataContainers, and then the norm is taken column-wise.
The pointwise l1 norm is equivalent to the following NumPy code snippet:
np.sum(np.abs(x), axis=0)
, assuming a signal, where the first axis indexes each block of the data.
-
GetFloatingPointType_t<data_t>
lInfNorm
() const¶ return the linf norm of this signal (maximum of absolute values)
-
data_t
l21SmoothMixedNorm
(data_t epsilon) const¶ return the mixed L21 norm of this signal with smoothing parameter epsilon
-
void
fft
(FFTNorm norm, FFTPolicy policy = FFTPolicy::AUTO)¶ convert to the fourier transformed signal.
- See
FFTPolicy
- Parameters
policy
: controls the choice of implementation. If the requested policy cannot be applied, a runtime error is generated.
-
void
ifft
(FFTNorm norm, FFTPolicy policy = FFTPolicy::AUTO)¶ convert to the inverse fourier transformed signal.
- See
FFTPolicy
- Parameters
policy
: controls the choice of implementation. If the requested policy cannot be applied, a runtime error is generated.
-
void
assign
(const DataContainer<data_t> &other)¶ copy the values from other DataContainer to this DataContainer
-
DataContainer<data_t> &
zero
() &¶ Set all values of the DataContainer to zero.
-
DataContainer<data_t>
zero
() &&¶ Set all values of the DataContainer to zero.
-
DataContainer<data_t> &
one
() &¶ Set all values of the DataContainer to one.
-
DataContainer<data_t>
one
() &&¶ Set all values of the DataContainer to one.
-
DataContainer<data_t> &
fill
(SelfType_t<data_t> value) &¶ Set all values of the DataContainer to the given value.
-
DataContainer<data_t>
fill
(SelfType_t<data_t> value) &&¶ Set all values of the DataContainer to the given value.
-
DataContainer<add_complex_t<data_t>>
asComplex
() const¶ if the datacontainer is already complex, return itself.
-
DataContainer<data_t> &
operator+=
(const DataContainer<data_t> &dc)¶ compute in-place element-wise addition of another container
-
DataContainer<data_t> &
operator-=
(const DataContainer<data_t> &dc)¶ compute in-place element-wise subtraction of another container
-
DataContainer<data_t> &
operator*=
(const DataContainer<data_t> &dc)¶ compute in-place element-wise multiplication with another container
-
DataContainer<data_t> &
operator/=
(const DataContainer<data_t> &dc)¶ compute in-place element-wise division by another container
-
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 DataContainer<data_t>
slice
(index_t i) const¶ Slice the container in the last dimension.
Access a portion of the container via a slice. The slicing always is in the last dimension. So for a 3D volume, the slice would be an sliced in the z direction and would be a part of the x-y plane.
A slice is always the same dimension as the original DataContainer, but with a thickness of 1 in the last dimension (i.e. the coefficient of the last dimension is 1)
-
DataContainer<data_t>
slice
(index_t i)¶ Slice the container in the last dimension, non-const overload.
This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
-
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)
-
void
format
(std::ostream &os, format_config cfg = {}) const¶ write a pretty-formatted string representation to stream
Private Members
-
std::unique_ptr<DataDescriptor>
_dataDescriptor
¶ the current DataDescriptor
-
std::variant<ContiguousStorage<data_t>, ContiguousStorageView<data_t>>
storage_
¶ the current DataHandler
Transformations¶
-
template<typename
data_t
>
DataContainer<data_t>elsa
::
exp
(const DataContainer<data_t> &dc)¶ Compute a coefficient wise exponential for each element of the
DataContainer
-
template<typename
data_t
>
DataContainer<data_t>elsa
::
log
(const DataContainer<data_t> &dc)¶ Compute a coefficient wise log for each element of the
DataContainer
-
template<typename
data_t
>
DataContainer<data_t>elsa
::
square
(const DataContainer<data_t> &dc)¶ Compute a coefficient wise square for each element of the
DataContainer
-
template<typename
data_t
>
DataContainer<data_t>elsa
::
sqrt
(const DataContainer<data_t> &dc)¶ Compute a coefficient wise square root for each element of the
DataContainer
-
template<class
data_t
>
DataContainer<data_t>elsa
::
bessel_log_0
(const DataContainer<data_t> &dc)¶ Compute an element-wise log of modified bessel function of the first kind of order 0 for each element of the
DataContainer
-
template<class
data_t
>
DataContainer<data_t>elsa
::
bessel_1_0
(const DataContainer<data_t> &dc)¶ Compute an element-wise modified bessel function of the first kind of order 1 divided by that of the order 0 for each element of the
DataContainer
-
template<typename
data_t
>
DataContainer<data_t>elsa
::
minimum
(const DataContainer<data_t> &dc, SelfType_t<data_t> scalar)¶ Compute a coefficient wise minimum with a scalar. For each element in
x_i
the givenDataContainer
, computemin(x_i, scalar)
-
template<typename
data_t
>
DataContainer<data_t>elsa
::
maximum
(const DataContainer<data_t> &dc, SelfType_t<data_t> scalar)¶ Compute a coefficient wise maximum with a scalar. For each element in
x_i
the givenDataContainer
, computemax(x_i, scalar)
-
template<typename
data_t
>
DataContainer<value_type_of_t<data_t>>elsa
::
cwiseAbs
(const DataContainer<data_t> &dc)¶ Compute the absolute value for each of the coefficients of the given DataContainer.
-
template<typename
xdata_t
, typenameydata_t
>
DataContainer<value_type_of_t<std::common_type_t<xdata_t, ydata_t>>>elsa
::
cwiseMin
(const DataContainer<xdata_t> &lhs, const DataContainer<ydata_t> &rhs)¶
-
template<typename
xdata_t
, typenameydata_t
>
DataContainer<value_type_of_t<std::common_type_t<xdata_t, ydata_t>>>elsa
::
cwiseMax
(const DataContainer<xdata_t> &lhs, const DataContainer<ydata_t> &rhs)¶
-
template<typename
data_t
>
DataContainer<value_type_of_t<data_t>>elsa
::
sign
(const DataContainer<data_t> &dc)¶ compute the sign of each entry of the input DataContainer. The function is defined as:
\[ \operatorname{sign}(x_i) = \begin{cases} 1 & \text{if } x_i > 0 \\ -1 & \text{if } x_i < 0 \\ 0 & \text{elsa} \end{cases} \quad \forall i \]For complex numbers, the definition is givens as:\[ \operatorname{sign}(x_i) = \begin{cases} 1 & \text{if } \mathrm{Re}(x_i) > 0 \\ -1 & \text{if } \mathrm{Re}(x_i) < 0 \\ sign(\mathrm{Im}(x_i)) & \text{elsa} \end{cases} \quad \forall i \]
-
template<typename
data_t
>
DataContainer<value_type_of_t<data_t>>elsa
::
real
(const DataContainer<data_t> &dc)¶ Return the real part a complex DataContainer. If the input DataContainer is real, it is assumed that the imaginary part is zero and a copy is returned.
-
template<typename
data_t
>
DataContainer<value_type_of_t<data_t>>elsa
::
imag
(const DataContainer<data_t> &dc)¶ Return the imaginary part a complex DataContainer. If the input DataContainer is real, it is assumed that the imaginary part is zero and zero initialized DataContainer of the same size and DataDescriptor is returned.
-
template<typename
data_t
>
DataContainer<data_t>elsa
::
clip
(const DataContainer<data_t> &dc, data_t min, data_t max)¶ clip the container values outside of the interval, to the interval edges
-
template<typename
data_t
>
DataContainer<add_complex_t<data_t>>elsa
::
asComplex
(const DataContainer<data_t> &dc)¶ Return a DataContainer with complex data type. If the input vector is real, all coefficients will be converted to complex values with a zero imaginary part. If the input is complex already, a copy of the container is returned.
Creation¶
-
template<class
data_t
>
DataContainer<data_t>elsa
::
zeros
(const DataDescriptor &desc)¶ Create a DataContainer filled with zeros and the given DataDescriptor.
-
template<class
data_t
>
DataContainer<data_t>elsa
::
zeroslike
(const DataContainer<data_t> &dc)¶ Create a DataContainer filled with zeros and the DataDescriptor of the given DataContainer.
-
template<class
data_t
>
DataContainer<data_t>elsa
::
ones
(const DataDescriptor &desc)¶ Create a DataContainer filled with ones values and the given DataDescriptor.
-
template<class
data_t
>
DataContainer<data_t>elsa
::
oneslike
(const DataContainer<data_t> &dc)¶ Create a DataContainer filled with ones and the DataDescriptor of the given DataContainer.
-
template<class
data_t
>
DataContainer<data_t>elsa
::
full
(const DataDescriptor &desc, SelfType_t<data_t> value)¶ Create a DataContainer filled with the given value and DataDescriptor.
-
template<class
data_t
>
DataContainer<data_t>elsa
::
fulllike
(const DataContainer<data_t> &dc, SelfType_t<data_t> value)¶ Create a DataContainer filled with the given value and the DataDescriptor of the given DataContainer
-
template<class
data_t
>
DataContainer<data_t>elsa
::
empty
(const DataDescriptor &desc)¶ Create an uninitialized DataContainer, the caller is responsible to fill DataContainer before use
-
template<class
data_t
>
DataContainer<data_t>elsa
::
emptylike
(const DataContainer<data_t> &dc)¶ Create an uninitialized DataContainer with the DataDescriptor of the given DataContainer. The caller is responsible to fill DataContainer before use
Others¶
-
template<class
data_t
>
DataContainer<data_t>elsa
::
lincomb
(SelfType_t<data_t> a, const DataContainer<data_t> &x, SelfType_t<data_t> b, const DataContainer<data_t> &y)¶ Compute the linear combination of $a * x + b * y$.
This function can be used as a memory efficient version for the computation of the above expression, as for such an expression (without expression template) multiple copies need to be created and allocated.
The function throws, if x and y do not have the same data descriptor
-
template<class
data_t
>
voidelsa
::
lincomb
(SelfType_t<data_t> a, const DataContainer<data_t> &x, SelfType_t<data_t> b, const DataContainer<data_t> &y, DataContainer<data_t> &out)¶ Compute the linear combination of $a * x + b * y$, and write it to the output variable.
This function can be used as a memory efficient version for the computation of the above expression, as for such an expression (without expression template) multiple copies need to be created and allocated.
The function throws, if x, y and out do not have the same data descriptor
-
template<class
data_t
>
DataContainer<data_t>elsa
::
materialize
(const DataContainer<data_t> &x)¶ Return an owning DataContainer, if given an non-owning one, the data is copied to a new owning buffer.
-
template<typename
data_t
>
DataContainer<data_t>elsa
::
concatenate
(const DataContainer<data_t> &dc1, const DataContainer<data_t> &dc2)¶ Concatenate two DataContainers to one (requires copying of both)
-
template<typename
data_t
>
DataContainer<data_t>elsa
::
fftShift2D
(const DataContainer<data_t> &dc)¶ Perform the FFT shift operation to the provided signal. Refer to https://numpy.org/doc/stable/reference/generated/numpy.fft.fftshift.html for further details.
-
template<typename
data_t
>
DataContainer<data_t>elsa
::
ifftShift2D
(const DataContainer<data_t> &dc)¶ Perform the IFFT shift operation to the provided signal. Refer to https://numpy.org/doc/stable/reference/generated/numpy.fft.ifftshift.html for further details.
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.
Subclassed by elsa::axdt::SphericalHarmonicsTransform< data_t >, elsa::AXDTOperator< data_t >
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
-
LinearOperator
(data_t scalar, const LinearOperator<data_t> &op)¶ 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>
operator*
(data_t scalar, const LinearOperator<data_t> &op)¶ friend operator* to support composition of a scalar and a LinearOperator
-
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)
Descriptors¶
Descriptors are storing metda data for the DataContainer. They give meaning to the actual data signal. Whether or not the data is the volume, measurements or some blocked data format depends on the specific type of descriptor.
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
InvalidArgumentError
: 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
InvalidArgumentError
: if any number of coefficients is non-positive, or sizes of numberOfCoefficientsPerDimension and spacingPerDimension do not match
-
IndexVector_t
getNumberOfCoefficientsPerDimension
() const¶ return the number of coefficients per dimension
-
IndexVector_t
getProductOfCoefficientsPerDimension
() const¶ return the product 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
(const elsa::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
-
template<class
data_t
>
DataContainer<data_t>element
() const¶ create a DataContainer for the given DataDescriptor. By default an uninitialized element is returned, and the caller is responsible to properly initialize the element
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
-
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)
VolumeDescriptor¶
-
class
elsa
::
VolumeDescriptor
: public elsa::DataDescriptor¶ Class representing metadata for linearized n-dimensional signal stored in memory.
This class 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.- Author
Matthias Wieczorek - initial code
- Author
Tobias Lasser - modularization, modernization
- Author
Maximilian Hornung - various enhancements
- Author
David Frank - inheritance restructuring
Public Functions
-
VolumeDescriptor
() = delete¶ delete default constructor (having no metadata is invalid)
-
~VolumeDescriptor
() override = default¶ default destructor
-
VolumeDescriptor
(IndexVector_t numberOfCoefficientsPerDimension)¶ Constructor for DataDescriptor, accepts vector for coefficients per dimensions.
- Parameters
[in] numberOfCoefficientsPerDimension
: vector containing the number of coefficients per dimension, (dimension is set implicitly from the size of the vector)
- Exceptions
InvalidArgumentError
: if any number of coefficients is non-positive
-
VolumeDescriptor
(std::initializer_list<index_t> numberOfCoefficientsPerDimension)¶ Constructs VolumeDescriptor from initializer list for the coefficients per dimensions.
- Parameters
[in] numberOfCoefficientsPerDimension
: initializer list containing the number of coefficients per dimension (dimension is set implicitly from the size of the list)
- Exceptions
InvalidArgumentError
: if any number of coefficients is non-positive
-
VolumeDescriptor
(IndexVector_t numberOfCoefficientsPerDimension, RealVector_t spacingPerDimension)¶ Constructor for DataDescriptor, accepts vectors for coefficients 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
InvalidArgumentError
: if any number of coefficients or spacing is non-positive, or sizes of numberOfCoefficientsPerDimension and spacingPerDimension do not match
-
VolumeDescriptor
(std::initializer_list<index_t> numberOfCoefficientsPerDimension, std::initializer_list<real_t> spacingPerDimension)¶ Constructs VolumeDescriptor from two initializer lists for coefficients and spacing.
- Parameters
[in] numberOfCoefficientsPerDimension
: initializer list containing the number of coefficients per dimension (dimension is set implicitly from the size of the list)[in] spacingPerDimension
: initializer list containing the spacing per dimension
- Exceptions
InvalidArgumentError
: if any number of coefficients or spacing is non-positive, or sizes of numberOfCoefficientsPerDimension and spacingPerDimension do not match
Private Functions
-
VolumeDescriptor *
cloneImpl
() const override¶ implement the polymorphic clone operation
-
bool
isEqual
(const DataDescriptor &other) const override¶ implement the polymorphic comparison operation
DetectorDescriptor¶
-
class
elsa
::
DetectorDescriptor
: public elsa::DataDescriptor¶ Class representing metadata for linearized n-dimensional signal stored in memory. It is a base class for different type signals captured by some kind of detectors (i.e. a planar detector, curved or some other shaped detector). It additionally stores information about the different poses of a trajectory.
Subclassed by elsa::CurvedDetectorDescriptor, elsa::PlanarDetectorDescriptor, elsa::XGIDetectorDescriptor
Public Functions
-
DetectorDescriptor
() = delete¶ There is not default signal.
-
~DetectorDescriptor
() override = default¶ Default destructor.
-
DetectorDescriptor
(const IndexVector_t &numOfCoeffsPerDim, const std::vector<Geometry> &geometryList)¶ Construct a DetectorDescriptor with a number of coefficients for each dimension and a list of geometry poses.
-
DetectorDescriptor
(const IndexVector_t &numOfCoeffsPerDim, const RealVector_t &spacingPerDim, const std::vector<Geometry> &geometryList)¶ Construct a DetectorDescriptor with a number of coefficients and spacing for each dimension and a list of geometry poses.
-
RealRay_t
computeRayFromDetectorCoord
(const index_t detectorIndex) const¶ Overload of computeRayToDetector with a single detectorIndex. Compute the pose and coord index using getCoordinateFromIndex and call overload.
-
RealRay_t
computeRayFromDetectorCoord
(const IndexVector_t coord) const¶ Overload of computeRayToDetector with a single coord vector. This vector encodes the pose index and the detector coordinate. So for a 1D detector, this will be 2D and the second dimension, will reference the pose index.
-
RealRay_t
computeRayFromDetectorCoord
(const RealVector_t &detectorCoord, const index_t poseIndex) const = 0¶ Compute a ray from the source from a pose to the given detector coordinate.
- Parameters
[in] detectorCoord
: Vector of size dim - 1, specifying the coordinate the ray should hit[in] poseIndex
: index into geometryList array, which pose to use for ray computation
-
std::pair<RealVector_t, real_t>
projectAndScaleVoxelOnDetector
(const RealVector_t &voxelCoord, const index_t poseIndex) const¶ Computes the projection of the center of a voxel to the detector and its scaling !!be aware that this function is not optimized, as it uses dynamic sized arrays!!
- Return
std::pair<RealVector_t, real_t> detector coordinate and scaling on detector
- Parameters
[in] voxelCoord
: coordinate of the voxel to be projected in volume coordinates[in] poseIndex
: index into geometryList array, which pose to use for projection
-
std::unique_ptr<DetectorDescriptor>
cloneWithGeometry
(std::vector<Geometry> geometries) const = 0¶ Create another DetectorDescriptor of same type sharing all parameters except geometries Useful to split datasets for stochastic methods
Protected Functions
-
bool
isEqual
(const DataDescriptor &other) const override¶ implement the polymorphic comparison operation
-
-
class
elsa
::
PlanarDetectorDescriptor
: public elsa::DetectorDescriptor¶ Class representing metadata for lineraized n-dimensional signal stored in memory. It specifically describes signals, which were captured by a planar detector and stores additional information such as different poses.
- Author
David Frank - initial code
Public Functions
-
PlanarDetectorDescriptor
(const IndexVector_t &numOfCoeffsPerDim, const RealVector_t &spacingPerDim, const std::vector<Geometry> &geometryList)¶ Construct a PlanarDetectorDescriptor with given number of coefficients and spacing per dimension and a list of geometry poses in the trajectory.
-
PlanarDetectorDescriptor
(const IndexVector_t &numOfCoeffsPerDim, const std::vector<Geometry> &geometryList)¶ Construct a PlanatDetectorDescriptor with given number of coefficients per dimension and a list of geometry poses in the trajectory.
-
RealRay_t
computeRayFromDetectorCoord
(const RealVector_t &detectorCoord, const index_t poseIndex) const override¶ Override function to compute rays for a planar detector.
-
std::unique_ptr<DetectorDescriptor>
cloneWithGeometry
(std::vector<Geometry> geometries) const override¶ Create another DetectorDescriptor of same type sharing all parameters except geometries Useful to split datasets for stochastic methods
-
RealRay_t
computeRayFromDetectorCoord
(const index_t detectorIndex) const¶ Overload of computeRayToDetector with a single detectorIndex. Compute the pose and coord index using getCoordinateFromIndex and call overload.
-
RealRay_t
computeRayFromDetectorCoord
(const IndexVector_t coord) const¶ Overload of computeRayToDetector with a single coord vector. This vector encodes the pose index and the detector coordinate. So for a 1D detector, this will be 2D and the second dimension, will reference the pose index.
-
RealRay_t
computeRayFromDetectorCoord
(const RealVector_t &detectorCoord, const index_t poseIndex) const = 0 Compute a ray from the source from a pose to the given detector coordinate.
- Parameters
[in] detectorCoord
: Vector of size dim - 1, specifying the coordinate the ray should hit[in] poseIndex
: index into geometryList array, which pose to use for ray computation
Private Functions
-
bool
isEqual
(const DataDescriptor &other) const override¶ implement the polymorphic comparison operation
-
class
elsa
::
CurvedDetectorDescriptor
: public elsa::DetectorDescriptor¶ Class representing a curved detector surface. It uses a virtual PlanarDetectorDescriptor in the background by mapping the coordinates of the curved detector detector to coordinates of the planar one.
The idea behind the current approach is to reduce the problem of computing intersections between rays and a curved detector to the planar detector case. We project coordinates of the curved detector onto coordinates on a virtual flat detector behind. Conceptually, the flat detector has the same amount of pixels but they become increasingly narrower towards the endpoints of the flat detector. The coordinates for the principal ray are the same for the flat and the curved detector.
- Author
Julia Spindler, Robert Imschweiler - adapt PlanarDetectorDescriptor for CurvedDetectorDescriptor
Public Functions
-
CurvedDetectorDescriptor
(const IndexVector_t &numOfCoeffsPerDim, const RealVector_t &spacingPerDim, const std::vector<Geometry> &geometryList, const geometry::Radian angle, const real_t s2d)¶ Construct a CurvedDetectorDescriptor with given number of coefficients and spacing per dimension, a list of geometry poses in the trajectory, an angle in radians, and the length from the source to the detector.
The information needed to model a curved detector includes the information required for a planar detector. Additionally, we need a parameter describing the curvature. Currently, this is implemented by providing the fanout angle (in radians), which we use internally to compute the radius of the curved detector. Furthermore, we need a parameter for the distance of the source to the detector (so, the sum of the distances “source to center” and “center to detector”).
-
CurvedDetectorDescriptor
(const IndexVector_t &numOfCoeffsPerDim, const std::vector<Geometry> &geometryList, const geometry::Radian angle, const real_t s2d)¶ Construct a CurvedDetectorDescriptor with given number of coefficients per dimension, a list of geometry poses in the trajectory, an angle in radians, and the length from the source to the detector.
The information needed to model a curved detector includes the information required for a planar detector. Additionally, we need a parameter describing the curvature. Currently, this is implemented by providing the fanout angle (in radians), which we use internally to compute the radius of the curved detector. Furthermore, we need a parameter for the distance of the source to the detector (so, the sum of the distances “source to center” and “center to detector”).
-
RealRay_t
computeRayFromDetectorCoord
(const RealVector_t &detectorCoord, const index_t poseIndex) const override¶ The ray computation of the curved detector descriptor is a wrapper around the ray computation of the planar detector descriptor. The wrapper is responsible for mapping the user-provided coordinates of the curved detector to the corresponding coordinates of the virtual planar detector descriptor. This overhead is encapsulated in the mapCurvedCoordToPlanarCoord method, which “extends” the rays hitting the curved detector to the virtual flat detector behind. Most importantly, the CurvedDetectorDescriptor class overrides the computeRayFromDetectorCoord function. This function receives the parameter detectorCoord which is treated as the index of a pixel on the curved detector. If the x-coordinate of detectorCoord is of the form x.5, i.e. references the pixel center, we can use the coordinates that were precomputed in the constructor. This should usually be the case. Otherwise, the planar detector coordinate needs to be dynamically computed by mapCurvedCoordToPlanarCoord. Finally, we receive a coordinate in the detector space of the flat detector which we pass to the parent implementation of computeRayFromDetectorCoord.
-
const std::vector<RealVector_t> &
getPlanarCoords
() const¶ Return the coordinates of the planar detector descriptor which operates in background.
-
std::unique_ptr<DetectorDescriptor>
cloneWithGeometry
(std::vector<Geometry> geometries) const override¶ Create another DetectorDescriptor of same type sharing all parameters except geometries Useful to split datasets for stochastic methods
-
RealRay_t
computeRayFromDetectorCoord
(const index_t detectorIndex) const¶ Overload of computeRayToDetector with a single detectorIndex. Compute the pose and coord index using getCoordinateFromIndex and call overload.
-
RealRay_t
computeRayFromDetectorCoord
(const IndexVector_t coord) const¶ Overload of computeRayToDetector with a single coord vector. This vector encodes the pose index and the detector coordinate. So for a 1D detector, this will be 2D and the second dimension, will reference the pose index.
-
RealRay_t
computeRayFromDetectorCoord
(const RealVector_t &detectorCoord, const index_t poseIndex) const = 0 Compute a ray from the source from a pose to the given detector coordinate.
- Parameters
[in] detectorCoord
: Vector of size dim - 1, specifying the coordinate the ray should hit[in] poseIndex
: index into geometryList array, which pose to use for ray computation
Private Functions
-
bool
isEqual
(const DataDescriptor &other) const override¶ actual comparison method, abstract to force override in derived classes
-
void
setup
()¶ setup function that is called in the constructor. Precomputes the coordinates of the pixel midpoints on the hypothetical planar detector
-
RealVector_t
mapCurvedCoordToPlanarCoord
(const RealVector_t &coord) const¶ Map a given coordinate of the curved detector to the corresponding coordinate of the virtual flat detector.
- Return
real_t
- Parameters
coord
:
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
-
const DataDescriptor &
getDescriptorOfBlock
(index_t i) const = 0¶ return the DataDescriptor 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.
Subclassed by elsa::SphericalCoefficientsDescriptor
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
InvalidArgumentError
: if numberOfBlocks is non-positive
-
IdenticalBlocksDescriptor
(const IdenticalBlocksDescriptor&) = delete¶ make copy constructor deletion explicit
-
~IdenticalBlocksDescriptor
() override = default¶ default destructor
-
const DataDescriptor &
getDescriptorOfBlock
(index_t i) const override¶ return the DataDescriptor 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
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
InvalidArgumentError
: 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
InvalidArgumentError
: if slicesInBlock does not specify a valid partition scheme for the given descriptor
-
~PartitionDescriptor
() override = default¶ default destructor
-
const DataDescriptor &
getDescriptorOfBlock
(index_t i) const override¶ return the DataDescriptor 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
InvalidArgumentError
: 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
InvalidArgumentError
: if the list is empty
-
RandomBlocksDescriptor
(const RandomBlocksDescriptor&) = delete¶ make copy constructor deletion explicit
-
~RandomBlocksDescriptor
() override = default¶ default desctructor
-
const DataDescriptor &
getDescriptorOfBlock
(index_t i) const override¶ return the DataDescriptor 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
Core Type Declarations¶
-
namespace
elsa
¶ Typedefs
-
using
real_t
= float¶ global type for real numbers
-
using
index_t
= std::ptrdiff_t¶ global type for indices
-
using
RealVector_t
= Eigen::Matrix<real_t, Eigen::Dynamic, 1>¶ global type for vectors of real numbers
-
using
ComplexVector_t
= Eigen::Matrix<complex_t, Eigen::Dynamic, 1>¶ global type for vectors of complex numbers
-
using
BooleanVector_t
= Eigen::Matrix<bool, Eigen::Dynamic, 1>¶ global type for vectors of booleans
-
template<typename
data_t
>
usingVector_t
= Eigen::Matrix<data_t, Eigen::Dynamic, 1>¶ global type for vectors of data_t
-
template<int
dim
>
usingIndexArray_t
= Eigen::Array<index_t, dim, 1>¶ global type for arrays of index_t of size dim
-
template<int
dim
>
usingRealArray_t
= Eigen::Array<real_t, dim, 1>¶ global type for arrays of real_t of size dim
-
template<int
dim
>
usingBooleanArray_t
= Eigen::Array<bool, dim, 1>¶ global type for arrays of bol of size dim
-
template<class
data_t
>
usingMatrix_t
= Eigen::Matrix<data_t, Eigen::Dynamic, Eigen::Dynamic>¶ global type for matrices of real numbers
-
template<typename
data_t
>
usingRay_t
= Eigen::ParametrizedLine<data_t, Eigen::Dynamic>¶ global type alias for rays
-
template<typename
T
>
usingGetFloatingPointType_t
= typename GetFloatingPointType<T>::type¶ helper typedef to facilitate usage
-
template<class
T
>
usingRemoveCvRef_t
= typename RemoveCvRef<T>::type¶ Helper to make type available.
Enums
Variables
-
template<typename
T
>
constexpr autopi
= static_cast<T>(3.14159265358979323846)¶ template global constexpr for the number pi
-
template<typename
T
>
constexpr boolisComplex
= std::is_same<RemoveCvRef_t<T>, complex<float>>::value || std::is_same<RemoveCvRef_t<T>, complex<double>>::value¶ Predicate to check if of complex type.
-
template<typename
T
>
structGetFloatingPointType
¶ - #include <elsaDefines.h>
base case for deducing floating point type of std::complex
-
template<typename
T
>
structGetFloatingPointType
<complex<T>>¶ - #include <elsaDefines.h>
partial specialization to derive correct floating point type
-
template<typename
T
>
structRemoveCvRef
¶ - #include <elsaDefines.h>
Remove cv qualifiers as well as reference of given type.
-
template<class
T
>
structSelfType
¶ - #include <elsaDefines.h>
With C++20 this can be replaced by std::type_identity.
-
namespace
axdt
¶ Typedefs
Enums
-
enum
Symmetry
¶ Symmetry influences the coefficients of the underlying spherical harmonics even -> odd degrees will be deemed to have zero coefficients; regular -> all non-zero (no simplification) Purely odd functions make no physical sense and can hence be ignored
Values:
-
enumerator
even
¶
-
enumerator
regular
¶
-
enumerator
-
enum
-
using
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()
Utilities¶
-
template<typename
T
>
classelsa
::
Badge
¶ Make public interfaces only callable by certain types.
If type T need insight into certain function of another class U, but if the it doesn’t make sense to make it publicly available to all classes, you could make T a friend of U. But then access to the complete private implementation is granted.
This Badge types is some middle way. The specific member functions needed by T are still public, but class T has to show its badge to get access to it.
This is specifically useful, if the private function is a private constructor of U and used to call
std::make_unique
(such that only T can create objects of class U with this specific constructor). Making T a friend of U, doesn’t help, asstd::make_unique
still can’t access the private implementation, but with the badge it works.This is a small example of the Badge:
class Bar; // Forward declaration class Foo { public: void foo(Bar bar) { // This is a legal call bar.internal_needed_by_foo({}); } }; class Bar { public: void internal_needed_by_foo(Badge<Foo>) };
Then
internal_needed_by_foo
is only callable fromFoo
, but no other class.References: https://awesomekling.github.io/Serenity-C++-patterns-The-Badge/
- Author
David Frank - initial code
- Template Parameters
T
:
Private Functions
-
Badge
() = default¶ Private constructor accessible by T.
Private Members
-
friend T
Make T a friend of
Badge
-
namespace
elsa
Functions
-
template<typename
Container
>
constexpr autocalculateMeanStddev
(Container v) -> std::pair<typename Container::value_type, typename Container::value_type>¶ Calculate mean and standard deviation of a container.
- Return
a pair of mean and standard deviation (of type
Conatiner::value_type
)- Parameters
v
: Any container, such asstd::vector
- Template Parameters
Container
: Container type of argument (e.g.std::vector
)
-
template<typename
data_t
= real_t>
constexpr autoconfidenceInterval95
(std::size_t n, data_t mean, data_t stddev) -> std::pair<real_t, real_t>¶ Compute the 95% confidence interval for a given number of samples
n
and the mean and standard deviation of the measurements.Compute it as $mean - c(n) * stddev, mean + c(n) * stddev$, where $c(n)$ is the n-th entry in the two tails T distribution table. For $n > 30$, it is assumed that $n = 1.96$.
- Return
pair of lower and upper bound of 95% confidence interval
- Parameters
n
: Number of of samplesmean
: mean of samplesstddev
: standard deviation of samples
-
template<typename