CTL
0.6.1
Computed Tomography Library
|
Specialized sub-class of Matrix<3, 4> to represent a projection matrix. More...
#include <projectionmatrix.h>
Public Types | |
enum | NormalizationMode { NoNormalization = 0, NormalizeAsUnitVector = 1, NormalizeByX = 2, NormalizeByY = 3, NormalizeByChannel = 2, NormalizeByRow = 3 } |
Public Member Functions | |
ProjectionMatrix ()=default | |
ProjectionMatrix (const Matrix< 3, 4 > &other) | |
void | shiftDetectorOrigin (const Matrix< 2, 1 > &translation) |
void | shiftDetectorOrigin (double translationX, double translationY) |
void | changeDetectorResolution (double resamplingFactor) |
void | changeDetectorResolution (double resamplingFactorX, double resamplingFactorY) |
void | normalize () |
ProjectionMatrix | normalized () const |
bool | isNormalized () const |
Matrix< 3, 3 > | M () const |
Matrix< 3, 1 > | p4 () const |
Matrix< 3, 1 > | directionSourceToPixel (const Matrix< 2, 1 > &pixelCoordinates, NormalizationMode normalizationMode=NoNormalization) const |
Matrix< 3, 1 > | directionSourceToPixel (double x, double y, NormalizationMode normalizationMode=NoNormalization) const |
Matrix< 3, 1 > | principalRayDirection () const |
Matrix< 3, 1 > | translationCTS () const |
Matrix< 3, 1 > | sourcePosition () const |
Matrix< 3, 3 > | rotationMatR () const |
Matrix< 3, 3 > | intrinsicMatK () const |
Matrix< 2, 1 > | principalPoint () const |
Matrix< 2, 1 > | focalLength () const |
double | skewCoefficient () const |
double | skewAngle () const |
double | magnificationX (double X, double Y, double Z) const |
double | magnificationX (const Matrix< 3, 1 > &worldCoordinate={ 0.0, 0.0, 0.0 }) const |
double | magnificationY (double X, double Y, double Z) const |
double | magnificationY (const Matrix< 3, 1 > &worldCoordinate={ 0.0, 0.0, 0.0 }) const |
Matrix< 2, 1 > | projectOntoDetector (double X, double Y, double Z) const |
Matrix< 2, 1 > | projectOntoDetector (const Matrix< 3, 1 > &worldCoordinate) const |
![]() | |
Matrix (double fillValue) | |
Matrix (const double(&initArray)[Rows *Cols]) | |
Matrix (double firstElement, Doubles... matrixElements) | |
auto | subMat () const -> Matrix< details::rangeDim(fromRow, toRow), details::rangeDim(fromCol, toCol)> |
auto | subMat () const -> Matrix< details::vecRangeDim< Rows >(from, to), details::vecRangeDim< Cols >(from, to)> |
Matrix< 1, Cols > | row () const |
Matrix< Rows, 1 > | column () const |
void | fill (double fillValue) |
void | normalize () |
Matrix< Rows, Cols > | normalized () const |
double | trace () const |
Matrix< Cols, Rows > | transposed () const |
Matrix< NewRows, NewCols > | reshaped () const |
Matrix< Rows *Cols, 1 > | reshapedAsVector () const |
Matrix< Rows, Cols > | operator- () const |
Matrix< Rows, Cols > | operator- (const Matrix< Rows, Cols > &rhs) const |
Matrix< Rows, Cols > & | operator *= (double scalar) |
Matrix< Rows, Cols > & | operator/= (double scalar) |
Matrix< Rows, Cols > & | operator+= (const Matrix< Rows, Cols > &rhs) |
Matrix< Rows, Cols > & | operator-= (const Matrix< Rows, Cols > &rhs) |
Matrix< Rows, Cols > & | operator%= (const Matrix< Rows, Cols > &rhs) |
Matrix< Rows, Cols > & | operator|= (const Matrix< Rows, Cols > &rhs) |
Matrix< Rows, Cols > | operator * (double scalar) const |
Matrix< Rows, Cols2 > | operator * (const Matrix< Cols, Cols2 > &rhs) const |
Matrix< Rows, Cols > | operator/ (double scalar) const |
Matrix< Rows, Cols > | operator+ (const Matrix< Rows, Cols > &rhs) const |
Matrix< Rows, Cols > | operator% (const Matrix< Rows, Cols > &rhs) const |
Matrix< Rows, Cols > | operator| (const Matrix< Rows, Cols > &rhs) const |
![]() | |
MatrixBase (double fillValue) | |
MatrixBase (const double(&initArray)[Rows *Cols]) | |
template<typename... Doubles> | |
MatrixBase (double firstElement, Doubles... matrixElements) | |
double * | operator[] (uint row) |
const double * | operator[] (uint row) const |
double & | operator() (uint row, uint column) |
double | operator() (uint row, uint column) const |
double & | at (uint row, uint column) noexcept(false) |
double | at (uint row, uint column) const noexcept(false) |
template<uint row, uint column> | |
double & | get () noexcept |
template<uint row, uint column> | |
double | get () const noexcept |
double & | operator() (uint n) |
double | operator() (uint n) const |
double & | at (uint n) |
double | at (uint n) const |
template<uint n> | |
double & | get () noexcept |
template<uint n> | |
double | get () const noexcept |
double * | data () |
const double * | data () const |
const double * | constData () const |
double * | begin () |
const double * | begin () const |
const double * | constBegin () const |
double * | end () |
const double * | end () const |
const double * | constEnd () const |
constexpr size_t | size () const |
std::string | info (const char *lineModifier="") const |
double | norm () const |
bool | operator== (const MatrixBase< Rows, Cols > &rhs) const |
bool | operator!= (const MatrixBase< Rows, Cols > &rhs) const |
Static Public Member Functions | |
static ProjectionMatrix | compose (const Matrix< 3, 3 > &M, const Matrix< 3, 1 > &p4) |
static ProjectionMatrix | compose (const Matrix< 3, 3 > &K, const Matrix< 3, 3 > &R, const Matrix< 3, 1 > &source) |
![]() | |
static Matrix< Rows, Cols > | fromContainer (const Container &vector, size_t NthMat, bool *ok=nullptr) |
Additional Inherited Members | |
![]() | |
static char & | SEPARATOR_CHARACTER_FOR_INFO_STRING = details::SeparatorCharacter<>::c |
global separator character for formatting | |
![]() | |
Matrix< Rows, Cols > | operator * (double scalar, const Matrix< Rows, Cols > &rhs) |
double | dot (const Matrix< N, 1 > &vector1, const Matrix< N, 1 > &vector2) |
double | dot (const Matrix< 1, N > &vector1, const Matrix< 1, N > &vector2) |
double | dot (const Matrix< N, 1 > &vector1, const Matrix< 1, N > &vector2) |
double | dot (const Matrix< 1, N > &vector1, const Matrix< N, 1 > &vector2) |
Matrix< N, N > | diag (const Matrix< N, 1 > &diagElements) |
Matrix< Rows, Cols > | diag (const Matrix< details::minDim(Rows, Cols), 1 > &diagElements) |
auto | diag (const Matrix< Rows, Cols > &m) -> Matrix< details::minDim(Rows, Cols), 1 > |
Matrix< N, N > | eye () |
Matrix< Rows, Cols > | eye () |
Matrix< Rows, Cols1+Cols2 > | horzcat (const Matrix< Rows, Cols1 > &m1, const Matrix< Rows, Cols2 > &m2) |
Matrix< Rows1+Rows2, Cols > | vertcat (const Matrix< Rows1, Cols > &m1, const Matrix< Rows2, Cols > &m2) |
Specialized sub-class of Matrix<3, 4> to represent a projection matrix.
A ProjectionMatrix is used to discribe a projective mapping from 3d to a 2d image. The projection model is a finite (pinhole) camera, which is capable to represent an ideal cone beam geometry. Here, the ProjectionMatrix encodes the whole information that is required for the projective mapping of a view with one X-ray point source and a flat panel detector. However, it is also applicable for a curved or spherical detector, if it is composed of serveral flat panel modules.
Note that a ProjectionMatrix can be multiplied with a scalar value without changing the geometry of the scene. Due to this degree of freedom, it is preferred to standardize the scaling of the projection matrices using the mehtods normalize() or normalized().
The projective mapping works with homogeneous coordinates. A 3d point in homogeneous world coordinates is multiplied by a ProjectionMatrix \(P\in\mathbb{R}^{3\times4}\):
\(\left[\begin{array}{c} \tilde{x}\\ \tilde{y}\\ w \end{array}\right]=P\left[\begin{array}{c} X\\ Y\\ Z\\ 1 \end{array}\right]\)
Converting this vector back to cartesian coordinates, one obtains the 2d image/detector pixel coordinates \(\left[\begin{array}{c} x\\ y \end{array}\right]=\left[\begin{array}{c} \frac{\tilde{x}}{w}\\ \frac{\tilde{y}}{w} \end{array}\right]\).
Usually, the coordinate \(\left[\begin{array}{c} x=0\\ y=0 \end{array}\right]\) denotes the center of the upper left pixel of the image (x-axis points from left to right and the y-axis from top to bottom).
The elements of the ProjectionMatrix are commonly subdivided into two blocks \(M\in\mathbb{R}^{3\times3}\) and \(\mathbf{p}_{4}\in\mathbb{R}^{3\times1}\):
\(P=\left[\begin{array}{cccc} P_{11} & P_{12} & P_{13} & P_{14}\\ P_{21} & P_{22} & P_{23} & P_{24}\\ P_{31} & P_{32} & P_{33} & P_{34} \end{array}\right]=\left[\begin{array}{cc} M & \mathbf{p}_{4}\end{array}\right]\).
Furthermore, a ProjectionMatrix can be composed/decomposed from/in extrinsic and intrinsic parameters:
Here, the (de)composition can be stated as \(P=K\,\left[\begin{array}{cc} R & \mathbf{t}\end{array}\right]\). The relation between translation vector and source position is given by \(\mathbf{t}=-R\mathbf{c}\). The intrisic matrix \(K\) is composed as follows:
\(K=c\left[\begin{array}{ccc} f_{x} & s & p_{x}\\ 0 & f_{y} & p_{y}\\ 0 & 0 & 1 \end{array}\right]\).
The principal point is the intersection of the detector plane and the principal ray, which is the ray that hits perpendicularly the detector.
This enumeration is a switch for the method directionSourceToPixel() that sets the procedure for the normalization of the returned vector.
Type | Type-ID |
---|---|
ProjectionMatrix::NoNormalization | 0 |
ProjectionMatrix::NormalizeAsUnitVector | 1 |
ProjectionMatrix::NormalizeByX | 2 |
ProjectionMatrix::NormalizeByY | 3 |
ProjectionMatrix::NormalizeByChannel | 2 |
ProjectionMatrix::NormalizeByRow | 3 |
|
default |
Default constructor. Constructs a projection matrix without any initialization.
|
inline |
Constructs a projection matrix from an other Matrix<3,4> (base class object).
void CTL::mat::ProjectionMatrix::changeDetectorResolution | ( | double | resamplingFactor | ) |
Increases the number of detector pixels by the same resamplingFactor for both dimensions.
void CTL::mat::ProjectionMatrix::changeDetectorResolution | ( | double | resamplingFactorX, |
double | resamplingFactorY | ||
) |
Increases the number of detector pixels by a resamplingFactorX and resamplingFactorY for each dimensions respectively.
|
static |
Composes a ProjectionMatrix from sub-blocks M and p4: \(P = \left[\begin{array}{cc} M & \mathbf{p}_{4}\end{array}\right]\)
|
static |
Composes a ProjectionMatrix from intrinsic calibration matrix K, extrinsic rotation matrix R und source position \(\mathbf{c}\): \(P = KR\,\left[\begin{array}{cc} I & -\mathbf{c}\end{array}\right]\), where \(I\) is a 3x3 identity matrix.
Matrix< 3, 1 > CTL::mat::ProjectionMatrix::directionSourceToPixel | ( | const Matrix< 2, 1 > & | pixelCoordinates, |
NormalizationMode | normalizationMode = NoNormalization |
||
) | const |
Convenience function that does the same as directionSourceToPixel(double x, double y, NormalizationMode normalizationMode) const.
Matrix< 3, 1 > CTL::mat::ProjectionMatrix::directionSourceToPixel | ( | double | x, |
double | y, | ||
NormalizationMode | normalizationMode = NoNormalization |
||
) | const |
Returns a direction vector of a ray from source to detector pixel (x, y) by calculating \(\textrm{sign}\left(\det M\right)M^{-1}\left[x,y,1\right]^{T}\) using RQ decomposition of \(M\): \(M^{-1}=Q^{T}R^{-1}\). The normalizationMode specifies the length of the returned vector:
NoNormalization
(default): No normalization is performed. The length depend on the scaling of the ProjectionMatrix. NormalizeAsUnitVector
: The length equals to 1 NormalizeByX
: The returned vector points from the source to the detector pixel. This vector is in units of the pixel width in x-/channel-direction. If you multiply with the physical pixel spacing in x-direction, you will obtain a vector with the actual physical length. NormalizeByY
: Same as NormalizeByX, but w.r.t. the y-/row-direction. However, NormalizeByX is recommended, since it is considered to be faster. NormalizeByChannel
: Identical to NormalizeByX. NormalizeByRow
: Identical to NormalizeByY. Matrix< 2, 1 > CTL::mat::ProjectionMatrix::focalLength | ( | ) | const |
Returns the focal length, i.e. the distance between source and detector in units of the pixel spacing in x and y direction.
Matrix< 3, 3 > CTL::mat::ProjectionMatrix::intrinsicMatK | ( | ) | const |
Returns the normalized calibration matrix (intrinsic parameters).
bool CTL::mat::ProjectionMatrix::isNormalized | ( | ) | const |
Return true if the projection matrix is normalized according to the principal ray direction vector \([P_{31}, P_{32}, P_{33}]\) (also referred to as the normal vector of the principal plane), i.e. the norm of this vector is (fuzzy) compared to 1.0
.
Matrix< 3, 3 > CTL::mat::ProjectionMatrix::M | ( | ) | const |
Returns submatrix \(M\), whereby \(P = \left[\begin{array}{cc} M & \mathbf{p}_{4}\end{array}\right]\).
double CTL::mat::ProjectionMatrix::magnificationX | ( | double | X, |
double | Y, | ||
double | Z | ||
) | const |
Same as magnificationX(const Matrix<3, 1>& worldCoordinate).
double CTL::mat::ProjectionMatrix::magnificationX | ( | const Matrix< 3, 1 > & | worldCoordinate = { 0.0, 0.0, 0.0 } | ) | const |
Return the magnification factor w.r.t. to the x dimension of the detector (channel direction). This factor M
describes how strong a point in world coordinates worldCoordinate is magnified when it is projected onto the detector. Precisely, this means that an extent epsilon
mm that is parallel to the detector's x-axis is enlarged to M*epsilon
pixels.
double CTL::mat::ProjectionMatrix::magnificationY | ( | double | X, |
double | Y, | ||
double | Z | ||
) | const |
Same as magnificationY(const Matrix<3, 1>& worldCoordinate).
double CTL::mat::ProjectionMatrix::magnificationY | ( | const Matrix< 3, 1 > & | worldCoordinate = { 0.0, 0.0, 0.0 } | ) | const |
Return the magnification factor w.r.t. to the y dimension of the detector (row direction). This factor M
describes how strong a point in world coordinates worldCoordinate is magnified when it is projected onto the detector. Precisely, this means that an extent epsilon
mm that is parallel to the detector's y-axis is enlarged to M*epsilon
pixels. Note that if there are square pixels and no skew factor, this magnification is identical to the magnification in the x dimension of the detector, s.a. magnificationX().
void CTL::mat::ProjectionMatrix::normalize | ( | ) |
Normalizes the current projection matrix in place.
Note that this function may shadow the Matrix::normalize
function, which normalizes the current matrix using its Frobenius norm (if ENABLE_FROBENIUS_NORM
has been defined).
ProjectionMatrix CTL::mat::ProjectionMatrix::normalized | ( | ) | const |
Returns a normalized ProjectionMatrix by dividing by the norm of the principal ray direction \([P_{31}, P_{32}, P_{33}]\) (also referred to as the normal vector of the principal plane). Moreover, it normalizes the sign, so that the determinant of \(M\) is positive \((P = \left[\begin{array}{cc} M & \mathbf{p}_{4}\end{array}\right])\).
Note that this function may shadow the Matrix::normalized
function, which computes a normalized matrix using its Frobenius norm (if ENABLE_FROBENIUS_NORM
has been defined).
Matrix< 3, 1 > CTL::mat::ProjectionMatrix::p4 | ( | ) | const |
Returns vector \(\mathbf{p}_{4}\), whereby \(P = \left[\begin{array}{cc} M & \mathbf{p}_{4}\end{array}\right]\).
Matrix< 2, 1 > CTL::mat::ProjectionMatrix::principalPoint | ( | ) | const |
Returns the principal point (intrinsic parameters), i.e. the pixel coordinates corresponding to the principal ray (ray that hits perpendicularly the detector).
Matrix< 3, 1 > CTL::mat::ProjectionMatrix::principalRayDirection | ( | ) | const |
Returns unit direction vector of the ray that hits perpendicularly the detector.
Matrix< 2, 1 > CTL::mat::ProjectionMatrix::projectOntoDetector | ( | double | X, |
double | Y, | ||
double | Z | ||
) | const |
Maps a point in (cartesian) world coordinates onto the detector/image plane. This 2-dimensional point (in pixel coordinates) will be returned by the function. The mapping is performed by multiplying the projection matrix (from the left hand side) to the 4d vector [x, y, z, 1] and afterwards the resulting homogeneous 3d vector is normalized back to cartesian coordinates (so that the 3rd element is 1) and only the first two elements are returned.
\(\left[\begin{array}{c} \tilde{x}\\ \tilde{y}\\ w \end{array}\right]=P\left[\begin{array}{c} X\\ Y\\ Z\\ 1 \end{array}\right]\)
Then, the return value is \(\left[\begin{array}{c} x\\ y \end{array}\right]=\left[\begin{array}{c} \frac{\tilde{x}}{w}\\ \frac{\tilde{y}}{w} \end{array}\right]\).
Matrix< 2, 1 > CTL::mat::ProjectionMatrix::projectOntoDetector | ( | const Matrix< 3, 1 > & | worldCoordinate | ) | const |
This function does the same as projectOntoDetector(double X, double Y, double Z).
Matrix< 3, 3 > CTL::mat::ProjectionMatrix::rotationMatR | ( | ) | const |
Returns the rotation matrix (extrinsic parameters).
void CTL::mat::ProjectionMatrix::shiftDetectorOrigin | ( | const Matrix< 2, 1 > & | translation | ) |
Shift image (detector) coordiante origin (upper left image corner), e.g. when cropping the image in a way that the upper left corner is affected. The translation = \([tx,ty]\) of the image frame is performed by a matrix multiplication from the left-hand side:
\(P\rightarrow\left(\begin{array}{ccc} 1 & 0 & -tx\\ 0 & 1 & -ty\\ 0 & 0 & 1 \end{array}\right)P\)
Note that the translation vector is in units of detector pixels (not in mm).
void CTL::mat::ProjectionMatrix::shiftDetectorOrigin | ( | double | translationX, |
double | translationY | ||
) |
Convenience function that does the same as shiftDetectorOrigin(const Matrix<2, 1>& translation).
double CTL::mat::ProjectionMatrix::skewAngle | ( | ) | const |
Returns the skew angle, which refers to the angle between the oblique pixel edge and the (imaginary) edge of a corresponding rectangular pixel.
This is zero if the detector coordinate system is orthogonal.
double CTL::mat::ProjectionMatrix::skewCoefficient | ( | ) | const |
Returns the skew coefficient, which is zero if the detector coordinate system is orthogonal.
Matrix< 3, 1 > CTL::mat::ProjectionMatrix::sourcePosition | ( | ) | const |
Returns the source position (extrinsic parameters).
Matrix< 3, 1 > CTL::mat::ProjectionMatrix::translationCTS | ( | ) | const |
Returns the translation of the CT system (after rotation) \(\mathbf{t}=-R\mathbf{c}\) (with source position \(\mathbf{c}\)). It is the vector \(\mathbf{t}\) in the decomposition \(P=K\,\left[\begin{array}{cc} R & \mathbf{t}\end{array}\right]=KR\,\left[\begin{array}{cc} I & -\mathbf{c}\end{array}\right]\).