Math and statistics » Tensor operators module

Operators specific to tensor images.

Contents

Functions

void dip::AllTensorElements(dip::Image const& in, dip::Image& out)
Determines if all tensor elements are non-zero, producing a binary scalar image.
void dip::Angle(dip::Image const& in, dip::Image& out)
Computes the angle of the vector at each pixel in image in. more...
void dip::AnyTensorElement(dip::Image const& in, dip::Image& out)
Determines if any tensor element is non-zero, producing a binary scalar image.
void dip::CartesianToPolar(dip::Image const& in, dip::Image& out)
Converts the vector at each pixel in image in from Cartesian coordinates to polar (or spherical) coordinates. more...
void dip::ConjugateTranspose(dip::Image const& in, dip::Image& out)
Computes the conjugate transpose of the tensor image in.
void dip::CrossProduct(dip::Image const& lhs, dip::Image const& rhs, dip::Image& out)
Computes the cross product (inner product) of two vector images. more...
void dip::Determinant(dip::Image const& in, dip::Image& out)
Computes the determinant of the square matrix at each pixel in image in.
void dip::DotProduct(dip::Image const& lhs, dip::Image const& rhs, dip::Image& out)
Computes the dot product (inner product) of two vector images.
void dip::EigenDecomposition(dip::Image const& in, dip::Image& out, dip::Image& eigenvectors, dip::String const& method = S::PRECISE)
Computes the eigenvalues and eigenvectors of the square matrix at each pixel in image in. more...
void dip::Eigenvalues(dip::Image const& in, dip::Image& out, dip::String const& method = S::PRECISE)
Computes the eigenvalues of the square matrix at each pixel in image in. more...
void dip::GeometricMeanTensorElement(dip::Image const& in, dip::Image& out)
Computes the geometric mean tensor element value at each pixel, producing a scalar image.
void dip::Identity(dip::Image const& in, dip::Image& out)
Creates an image whose pixels are identity matrices. more...
void dip::Inverse(dip::Image const& in, dip::Image& out)
Computes the inverse of the square matrix at each pixel in image in. more...
void dip::LargestEigenvalue(dip::Image const& in, dip::Image& out, dip::String const& method = S::PRECISE)
Finds the largest eigenvalue of the square matrix at each pixel in image in. more...
void dip::LargestEigenvector(dip::Image const& in, dip::Image& out)
Finds the largest eigenvector of the symmetric matrix at each pixel in image in. more...
void dip::MaximumAbsTensorElement(dip::Image const& in, dip::Image& out)
Takes the maximum absolute tensor element at each pixel, producing a scalar image. For float and complex images only.
void dip::MaximumTensorElement(dip::Image const& in, dip::Image& out)
Takes the maximum tensor element at each pixel, producing a scalar image.
void dip::MeanTensorElement(dip::Image const& in, dip::Image& out)
Computes the mean tensor element value at each pixel, producing a scalar image.
void dip::MinimumAbsTensorElement(dip::Image const& in, dip::Image& out)
Takes the minimum absolute tensor element at each pixel, producing a scalar image. For float and complex images only.
void dip::MinimumTensorElement(dip::Image const& in, dip::Image& out)
Takes the minimum tensor element at each pixel, producing a scalar image.
void dip::Norm(dip::Image const& in, dip::Image& out)
Computes the norm of the vector at each pixel in image in.
void dip::Orientation(dip::Image const& in, dip::Image& out)
Computes the orientation of the vector at each pixel in image in. more...
void dip::PolarToCartesian(dip::Image const& in, dip::Image& out)
Converts the vector at each pixel in image in from polar (or spherical) coordinates to Cartesian coordinates. more...
void dip::ProductTensorElements(dip::Image const& in, dip::Image& out)
Multiplies all tensor elements, producing a scalar image.
void dip::PseudoInverse(dip::Image const& in, dip::Image& out, dip::dfloat tolerance = 1e-7)
Computes the pseudo-inverse of the matrix at each pixel in image in. more...
void dip::Rank(dip::Image const& in, dip::Image& out)
Computes the rank of the square matrix at each pixel in image in. The output is DT_UINT8, under the assumption that we won’t have tensor images with a rank higher than 255.
void dip::SingularValueDecomposition(dip::Image const& A, dip::Image& U, dip::Image& S, dip::Image& V)
Computes the “thin” singular value decomposition of the matrix at each pixel in image in. more...
void dip::SingularValues(dip::Image const& in, dip::Image& out)
Computes the “thin” singular value decomposition of the matrix at each pixel in image in. more...
void dip::SmallestEigenvalue(dip::Image const& in, dip::Image& out, dip::String const& method = S::PRECISE)
Finds the smallest eigenvalue of the square matrix at each pixel in image in. more...
void dip::SmallestEigenvector(dip::Image const& in, dip::Image& out)
Finds the smallest eigenvector of the symmetric matrix at each pixel in image in. more...
void dip::SortTensorElements(dip::Image& out)
Sorts the tensor elements within each pixel from largest to smallest. Works in-place. out must be real-valued.
void dip::SortTensorElementsByMagnitude(dip::Image& out)
Sorts the tensor elements within each pixel by magnitude from largest to smallest. Works in-place. out must be of a floating point or complex type.
void dip::SquareNorm(dip::Image const& in, dip::Image& out)
Computes the square of the norm of the vector at each pixel in image in.
void dip::SumTensorElements(dip::Image const& in, dip::Image& out)
Adds all tensor elements, producing a scalar image.
void dip::Trace(dip::Image const& in, dip::Image& out)
Computes the trace of the square matrix at each pixel in image in.
auto dip::Transpose(dip::Image const& in) -> dip::Image
Transposes the tensor image, the data are not copied.

Function documentation

dip::Image dip::Transpose( dip::Image const& in)

Transposes the tensor image, the data are not copied.

void dip::ConjugateTranspose( dip::Image const& in, dip::Image& out)

Computes the conjugate transpose of the tensor image in.

void dip::DotProduct( dip::Image const& lhs, dip::Image const& rhs, dip::Image& out)

Computes the dot product (inner product) of two vector images.

void dip::CrossProduct( dip::Image const& lhs, dip::Image const& rhs, dip::Image& out)

Computes the cross product (inner product) of two vector images.

Input image tensors must be 2-vectors or 3-vectors. For 3-vectors, the cross product is as commonly defined in 3D. For 2-vectors, we define the cross product as the z-component of the cross product of the 3D vectors obtained by adding a 0 z-component to the inputs. That is, it is the area of the parallelogram formed by the two 2D vectors.

void dip::Norm( dip::Image const& in, dip::Image& out)

Computes the norm of the vector at each pixel in image in.

void dip::SquareNorm( dip::Image const& in, dip::Image& out)

Computes the square of the norm of the vector at each pixel in image in.

void dip::Angle( dip::Image const& in, dip::Image& out)

Computes the angle of the vector at each pixel in image in.

in must be a 2-vector or a 3-vector. For a 2-vector, out is a scalar image representing phi, the angle from the x-axis. For a 3-vector, out has 2 tensor components, corresponding to phi and theta. phi, as in the 2D case, is the angle from the x-axis within the x-y plane (azimuth). theta is the angle from the z-axis (inclination). See dip::CartesianToPolar for more details. This function yields the same output as dip::CartesianToPolar, but without the first tensor component.

void dip::Orientation( dip::Image const& in, dip::Image& out)

Computes the orientation of the vector at each pixel in image in.

Orientation is defined as the angle mapped to the half-circle or half-sphere with positive x-coordinate. That is, in 2D it is an angle in the range (-π/2, π/2), and in 3D the phi component is mapped to that same range. See dip::Angle for more information.

void dip::CartesianToPolar( dip::Image const& in, dip::Image& out)

Converts the vector at each pixel in image in from Cartesian coordinates to polar (or spherical) coordinates.

in must be a 2-vector or a 3-vector. out is a same-size vector containing r and phi in the 2D case, and r, phi and theta in the 3D case. phi is the angle to the x-axis within the x-y plane (azimuth). theta is the angle from the z-axis (inclination).

That is, in 2D:

in[ 0 ] == out[ 0 ] * Cos( out[ 1 ] );
in[ 1 ] == out[ 0 ] * Sin( out[ 1 ] );

and in 3D:

in[ 0 ] == out[ 0 ] * Cos( out[ 1 ] ) * Sin( out[ 2 ] );
in[ 1 ] == out[ 0 ] * Sin( out[ 1 ] ) * Sin( out[ 2 ] );
in[ 2 ] == out[ 0 ] * Cos( out[ 2 ] );

void dip::PolarToCartesian( dip::Image const& in, dip::Image& out)

Converts the vector at each pixel in image in from polar (or spherical) coordinates to Cartesian coordinates.

in must be a 2-vector or a 3-vector. See dip::CartesianToPolar for a description of the polar coordinates used.

void dip::Determinant( dip::Image const& in, dip::Image& out)

Computes the determinant of the square matrix at each pixel in image in.

void dip::Trace( dip::Image const& in, dip::Image& out)

Computes the trace of the square matrix at each pixel in image in.

void dip::Rank( dip::Image const& in, dip::Image& out)

Computes the rank of the square matrix at each pixel in image in. The output is DT_UINT8, under the assumption that we won’t have tensor images with a rank higher than 255.

void dip::Eigenvalues( dip::Image const& in, dip::Image& out, dip::String const& method = S::PRECISE)

Computes the eigenvalues of the square matrix at each pixel in image in.

out is a vector image containing the eigenvalues. If in is symmetric and real-valued, then out is real-valued, otherwise, out is complex-valued. The eigenvalues are sorted by magnitude, in descending order.

method is either "precise" or "fast". The precise method uses an iterative QR decomposition. The fast method uses a closed-form algorithm that is much faster, but potentially less accurate. This faster algorithm is only used for real-valued, 2x2 or 3x3 symmetric tensor images.

void dip::LargestEigenvalue( dip::Image const& in, dip::Image& out, dip::String const& method = S::PRECISE)

Finds the largest eigenvalue of the square matrix at each pixel in image in.

Computes the eigenvalues in the same way as dip::Eigenvalues, but outputs only the eigenvalue with the largest magnitude. See the linked function’s documentation for a description of the method parameter.

void dip::SmallestEigenvalue( dip::Image const& in, dip::Image& out, dip::String const& method = S::PRECISE)

Finds the smallest eigenvalue of the square matrix at each pixel in image in.

Computes the eigenvalues in the same way as dip::Eigenvalues, but outputs only the eigenvalue with the smallest magnitude. See the linked function’s documentation for a description of the method parameter.

void dip::EigenDecomposition( dip::Image const& in, dip::Image& out, dip::Image& eigenvectors, dip::String const& method = S::PRECISE)

Computes the eigenvalues and eigenvectors of the square matrix at each pixel in image in.

The decomposition is such that in * eigenvectors == eigenvectors * out. eigenvectors is almost always invertible, in which case one can write in == eigenvectors * out * Inverse( eigenvectors ).

out is a diagonal matrix image containing the eigenvalues. If in is symmetric and real-valued, then out is real-valued, otherwise, out is complex-valued. The eigenvalues are sorted by magnitude, in descending order.

The eigenvectors are the columns of eigenvectors. It has the same data type as out.

method is either "precise" or "fast". The precise method uses an iterative QR decomposition. The fast method uses a closed-form algorithm that is much faster, but potentially less accurate. This faster algorithm is only used for real-valued, 2x2 or 3x3 symmetric tensor images.

void dip::LargestEigenvector( dip::Image const& in, dip::Image& out)

Finds the largest eigenvector of the symmetric matrix at each pixel in image in.

Computes the eigen decomposition in the same way as dip::EigenDecomposition, but outputs only the eigenvector that corresponds to the eigenvalue with largest magnitude. Always uses the precise algorithm.

in must be symmetric and real-valued.

void dip::SmallestEigenvector( dip::Image const& in, dip::Image& out)

Finds the smallest eigenvector of the symmetric matrix at each pixel in image in.

Computes the eigen decomposition in the same way as dip::EigenDecomposition, but outputs only the eigenvector that corresponds to the eigenvalue with smallest magnitude.

in must be symmetric and real-valued.

void dip::Inverse( dip::Image const& in, dip::Image& out)

Computes the inverse of the square matrix at each pixel in image in.

The result is undetermined if the matrix is not invertible.

void dip::PseudoInverse( dip::Image const& in, dip::Image& out, dip::dfloat tolerance = 1e-7)

Computes the pseudo-inverse of the matrix at each pixel in image in.

Computes the Moore-Penrose pseudo-inverse using tolerance. Singular values smaller than tolerance * max(rows,cols) * p, with p the largest singular value, will be set to zero in the inverse.

void dip::SingularValues( dip::Image const& in, dip::Image& out)

Computes the “thin” singular value decomposition of the matrix at each pixel in image in.

For an input image in with a tensor size of NxP, and with M the smaller of N and P, out is a vector image with M elements, corresponding to the singular values, sorted in decreasing order.

Use dip::SingularValueDecomposition if you need the full decomposition.

This function uses the two-sided Jacobi SVD decomposition algorithm. This is efficient for small matrices only.

void dip::SingularValueDecomposition( dip::Image const& A, dip::Image& U, dip::Image& S, dip::Image& V)

Computes the “thin” singular value decomposition of the matrix at each pixel in image in.

For an input image A with a tensor size of NxP, and with M the smaller of N and P, S is a square diagonal MxM matrix, U is a NxM matrix, and V is a PxM matrix. These matrices satisfy the relation \(A = USV^*\) .

The (diagonal) elements of S are the singular values, sorted in decreasing order. You can use dip::SingularValues if you are not interested in computing U and V.

This function uses the two-sided Jacobi SVD decomposition algorithm. This is efficient for small matrices only.

void dip::Identity( dip::Image const& in, dip::Image& out)

Creates an image whose pixels are identity matrices.

out will have the same sizes as in, and with a tensor representation of a diagonal matrix with a size concordant to that of the tensor representation of in. For example, for an N-vector image, the resulting output matrix image will be NxN. out will be of type dip::DT_SFLOAT.

void dip::SumTensorElements( dip::Image const& in, dip::Image& out)

Adds all tensor elements, producing a scalar image.

void dip::ProductTensorElements( dip::Image const& in, dip::Image& out)

Multiplies all tensor elements, producing a scalar image.

void dip::AllTensorElements( dip::Image const& in, dip::Image& out)

Determines if all tensor elements are non-zero, producing a binary scalar image.

void dip::AnyTensorElement( dip::Image const& in, dip::Image& out)

Determines if any tensor element is non-zero, producing a binary scalar image.

void dip::MaximumTensorElement( dip::Image const& in, dip::Image& out)

Takes the maximum tensor element at each pixel, producing a scalar image.

void dip::MaximumAbsTensorElement( dip::Image const& in, dip::Image& out)

Takes the maximum absolute tensor element at each pixel, producing a scalar image. For float and complex images only.

void dip::MinimumTensorElement( dip::Image const& in, dip::Image& out)

Takes the minimum tensor element at each pixel, producing a scalar image.

void dip::MinimumAbsTensorElement( dip::Image const& in, dip::Image& out)

Takes the minimum absolute tensor element at each pixel, producing a scalar image. For float and complex images only.

void dip::MeanTensorElement( dip::Image const& in, dip::Image& out)

Computes the mean tensor element value at each pixel, producing a scalar image.

void dip::GeometricMeanTensorElement( dip::Image const& in, dip::Image& out)

Computes the geometric mean tensor element value at each pixel, producing a scalar image.

void dip::SortTensorElements( dip::Image& out)

Sorts the tensor elements within each pixel from largest to smallest. Works in-place. out must be real-valued.

void dip::SortTensorElementsByMagnitude( dip::Image& out)

Sorts the tensor elements within each pixel by magnitude from largest to smallest. Works in-place. out must be of a floating point or complex type.