RTK
2.6.0
Reconstruction Toolkit
|
#include <rtkThreeDCircularProjectionGeometry.h>
Public Member Functions | |
void | AddProjection (const double sid, const double sdd, const double gantryAngle, const double projOffsetX=0., const double projOffsetY=0., const double outOfPlaneAngle=0., const double inPlaneAngle=0., const double sourceOffsetX=0., const double sourceOffsetY=0.) |
bool | AddProjection (const PointType &sourcePosition, const PointType &detectorPosition, const VectorType &detectorRowVector, const VectorType &detectorColumnVector) |
bool | AddProjection (const HomogeneousProjectionMatrixType &pMat) |
virtual void | AddProjectionInRadians (const double sid, const double sdd, const double gantryAngle, const double projOffsetX=0., const double projOffsetY=0., const double outOfPlaneAngle=0., const double inPlaneAngle=0., const double sourceOffsetX=0., const double sourceOffsetY=0.) |
void | Clear () override |
virtual ::itk::LightObject::Pointer | CreateAnother () const |
const std::vector< double > | GetAngularGaps (const std::vector< double > &angles) |
const std::vector< double > | GetAngularGapsWithNext (const std::vector< double > &angles) const |
virtual double | GetFixAnglesTolerance () const |
const std::vector< Superclass::MatrixType > & | GetMagnificationMatrices () const |
Superclass::MatrixType | GetMagnificationMatrices (const unsigned int i) const |
const ThreeDHomogeneousMatrixType | GetProjectionCoordinatesToDetectorSystemMatrix (const unsigned int i) const |
const ThreeDHomogeneousMatrixType | GetProjectionCoordinatesToFixedSystemMatrix (const unsigned int i) const |
const std::vector< ThreeDHomogeneousMatrixType > & | GetRotationMatrices () const |
ThreeDHomogeneousMatrixType | GetRotationMatrix (const unsigned int i) const |
const std::multimap< double, unsigned int > | GetSortedAngles (const std::vector< double > &angles) const |
const std::vector< double > & | GetSourceAngles () const |
const HomogeneousVectorType | GetSourcePosition (const unsigned int i) const |
const std::vector< ThreeDHomogeneousMatrixType > & | GetSourceTranslationMatrices () const |
ThreeDHomogeneousMatrixType | GetSourceTranslationMatrices (const unsigned int i) const |
const std::vector< double > | GetTiltAngles () const |
const std::map< double, unsigned int > | GetUniqueSortedAngles (const std::vector< double > &angles) const |
void | SetCollimationOfLastProjection (const double uinf, const double usup, const double vinf, const double vsup) |
virtual void | SetFixAnglesTolerance (double _arg) |
double | ToUntiltedCoordinateAtIsocenter (const unsigned int noProj, const double tiltedCoord) const |
const std::vector< double > & | GetGantryAngles () const |
const std::vector< double > & | GetOutOfPlaneAngles () const |
const std::vector< double > & | GetInPlaneAngles () const |
const std::vector< double > & | GetSourceToIsocenterDistances () const |
const std::vector< double > & | GetSourceOffsetsX () const |
const std::vector< double > & | GetSourceOffsetsY () const |
const std::vector< double > & | GetSourceToDetectorDistances () const |
const std::vector< double > & | GetProjectionOffsetsX () const |
const std::vector< double > & | GetProjectionOffsetsY () const |
const std::vector< TwoDHomogeneousMatrixType > & | GetProjectionTranslationMatrices () const |
TwoDHomogeneousMatrixType | GetProjectionTranslationMatrix (const unsigned int i) const |
const std::vector< double > & | GetCollimationUInf () const |
const std::vector< double > & | GetCollimationUSup () const |
const std::vector< double > & | GetCollimationVInf () const |
const std::vector< double > & | GetCollimationVSup () const |
virtual double | GetRadiusCylindricalDetector () const |
virtual void | SetRadiusCylindricalDetector (double _arg) |
virtual double | GetVerifyAnglesTolerance () const |
virtual void | SetVerifyAnglesTolerance (double _arg) |
Public Member Functions inherited from rtk::ProjectionGeometry< 3 > | |
virtual ::itk::LightObject::Pointer | CreateAnother () const |
const std::vector< MatrixType > & | GetMatrices () const |
MatrixType | GetMatrix (const unsigned int i) const |
Static Public Member Functions | |
static Superclass::MatrixType | ComputeProjectionMagnificationMatrix (double sdd, double sid) |
static ThreeDHomogeneousMatrixType | ComputeRotationHomogeneousMatrix (double angleX, double angleY, double angleZ) |
static double | ConvertAngleBetween0And2PIRadians (const double a) |
static double | ConvertAngleBetween0And360Degrees (const double a) |
static double | ConvertAngleBetweenMinusAndPlusPIRadians (const double a) |
static Pointer | New () |
static TwoDHomogeneousMatrixType | ComputeTranslationHomogeneousMatrix (double transX, double transY) |
static ThreeDHomogeneousMatrixType | ComputeTranslationHomogeneousMatrix (double transX, double transY, double transZ) |
Static Public Member Functions inherited from rtk::ProjectionGeometry< 3 > | |
static Pointer | New () |
Protected Member Functions | |
virtual void | AddMagnificationMatrix (const Superclass::MatrixType &m) |
virtual void | AddProjectionTranslationMatrix (const TwoDHomogeneousMatrixType &m) |
virtual void | AddRotationMatrix (const ThreeDHomogeneousMatrixType &m) |
virtual void | AddSourceTranslationMatrix (const ThreeDHomogeneousMatrixType &m) |
itk::LightObject::Pointer | InternalClone () const override |
ThreeDCircularProjectionGeometry () | |
~ThreeDCircularProjectionGeometry () override=default | |
bool | VerifyAngles (const double outOfPlaneAngleRAD, const double gantryAngleRAD, const double inPlaneAngleRAD, const Matrix3x3Type &referenceMatrix) const |
bool | FixAngles (double &outOfPlaneAngleRAD, double &gantryAngleRAD, double &inPlaneAngleRAD, const Matrix3x3Type &referenceMatrix) const |
Protected Member Functions inherited from rtk::ProjectionGeometry< 3 > | |
void | PrintSelf (std::ostream &os, itk::Indent indent) const override |
ProjectionGeometry ()=default | |
~ProjectionGeometry () override=default | |
virtual void | AddMatrix (const MatrixType &m) |
Protected Attributes | |
std::vector< double > | m_CollimationUInf |
std::vector< double > | m_CollimationUSup |
std::vector< double > | m_CollimationVInf |
std::vector< double > | m_CollimationVSup |
double | m_FixAnglesTolerance { 1e-6 } |
std::vector< double > | m_GantryAngles |
std::vector< double > | m_InPlaneAngles |
std::vector< Superclass::MatrixType > | m_MagnificationMatrices |
std::vector< double > | m_OutOfPlaneAngles |
std::vector< double > | m_ProjectionOffsetsX |
std::vector< double > | m_ProjectionOffsetsY |
std::vector< TwoDHomogeneousMatrixType > | m_ProjectionTranslationMatrices |
double | m_RadiusCylindricalDetector { 0. } |
std::vector< ThreeDHomogeneousMatrixType > | m_RotationMatrices |
std::vector< double > | m_SourceAngles |
std::vector< double > | m_SourceOffsetsX |
std::vector< double > | m_SourceOffsetsY |
std::vector< double > | m_SourceToDetectorDistances |
std::vector< double > | m_SourceToIsocenterDistances |
std::vector< ThreeDHomogeneousMatrixType > | m_SourceTranslationMatrices |
double | m_VerifyAnglesTolerance { 1e-4 } |
Projection geometry for a source and a 2-D flat panel.
The source and the detector rotate around a circle paremeterized with the SourceToDetectorDistance and the SourceToIsocenterDistance. The position of each projection along this circle is parameterized by the GantryAngle. The detector can be shifted in plane with the ProjectionOffsetsX and the ProjectionOffsetsY. It can be also rotated with InPlaneAngles and OutOfPlaneAngles. All angles are in radians except for the function AddProjection that takes angles in degrees. The source can be shifted in plane with the SourceOffsetsX and the SourceOffsetsY.
If SDD equals 0., then one is dealing with a parallel geometry.
More information is provided in RTK 3D circular projection geometry.
Definition at line 50 of file rtkThreeDCircularProjectionGeometry.h.
Definition at line 58 of file rtkThreeDCircularProjectionGeometry.h.
using rtk::ThreeDCircularProjectionGeometry::HomogeneousProjectionMatrixType = Superclass::MatrixType |
Definition at line 66 of file rtkThreeDCircularProjectionGeometry.h.
using rtk::ThreeDCircularProjectionGeometry::HomogeneousVectorType = itk::Vector<double, 4> |
Definition at line 61 of file rtkThreeDCircularProjectionGeometry.h.
using rtk::ThreeDCircularProjectionGeometry::Matrix3x3Type = itk::Matrix<double, 3, 3> |
Definition at line 65 of file rtkThreeDCircularProjectionGeometry.h.
Definition at line 57 of file rtkThreeDCircularProjectionGeometry.h.
using rtk::ThreeDCircularProjectionGeometry::PointType = itk::Point<double, 3> |
Definition at line 64 of file rtkThreeDCircularProjectionGeometry.h.
Definition at line 55 of file rtkThreeDCircularProjectionGeometry.h.
Definition at line 56 of file rtkThreeDCircularProjectionGeometry.h.
using rtk::ThreeDCircularProjectionGeometry::ThreeDHomogeneousMatrixType = itk::Matrix<double, 4, 4> |
Definition at line 63 of file rtkThreeDCircularProjectionGeometry.h.
using rtk::ThreeDCircularProjectionGeometry::TwoDHomogeneousMatrixType = itk::Matrix<double, 3, 3> |
Definition at line 62 of file rtkThreeDCircularProjectionGeometry.h.
using rtk::ThreeDCircularProjectionGeometry::VectorType = itk::Vector<double, 3> |
Definition at line 60 of file rtkThreeDCircularProjectionGeometry.h.
|
protected |
|
overrideprotecteddefault |
|
inlineprotectedvirtual |
Definition at line 403 of file rtkThreeDCircularProjectionGeometry.h.
void rtk::ThreeDCircularProjectionGeometry::AddProjection | ( | const double | sid, |
const double | sdd, | ||
const double | gantryAngle, | ||
const double | projOffsetX = 0. , |
||
const double | projOffsetY = 0. , |
||
const double | outOfPlaneAngle = 0. , |
||
const double | inPlaneAngle = 0. , |
||
const double | sourceOffsetX = 0. , |
||
const double | sourceOffsetY = 0. |
||
) |
Add projection to geometry. One projection is defined with the rotation angle in degrees and the in-plane translation of the detector in physical units (e.g. mm). The rotation axis is assumed to be (0,1,0).
bool rtk::ThreeDCircularProjectionGeometry::AddProjection | ( | const PointType & | sourcePosition, |
const PointType & | detectorPosition, | ||
const VectorType & | detectorRowVector, | ||
const VectorType & | detectorColumnVector | ||
) |
Add a REG23-based geometry set to the RTK projections list.
sourcePosition | absolute position of the point source S in WCS |
detectorPosition | absolute position of the detector origin R in WCS |
detectorRowVector | absolute direction vector indicating the orientation of the detector's rows r (sometimes referred to as v1) |
detectorColumnVector | absolute direction vector indicating the orientation of the detector's columns c (sometimes referred to as v2) |
bool rtk::ThreeDCircularProjectionGeometry::AddProjection | ( | const HomogeneousProjectionMatrixType & | pMat | ) |
Add projection from a projection matrix. A projection matrix is defined up to a scaling factor. The function here Assumes that the input matrix pMat is normalized such that pMat*(x,y,z,1)'=(u,v,1)'. This code assumes that the SourceToDetectorDistance is positive. If the projection is for parallel geometry, then the SourceToIsocenterDistance is set to 1000 mm.
Column-1 | Column-2 | Column-3 | Column-4 |
---|---|---|---|
R11 | R12 | R13 | -ProjectionOffsetX |
R21 | R22 | R23 | -ProjectionOffsetY |
0 | 0 | 0 | 1 |
where
Column-1 | Column-2 | Column-3 |
---|---|---|
R11 | R12 | R13 |
R21 | R22 | R23 |
0 | 0 | 0 |
is a cropped 3x3 rotation matrix.
|
virtual |
Idem with angles in radians.
|
inlineprotectedvirtual |
Definition at line 385 of file rtkThreeDCircularProjectionGeometry.h.
|
inlineprotectedvirtual |
Definition at line 391 of file rtkThreeDCircularProjectionGeometry.h.
|
inlineprotectedvirtual |
Definition at line 397 of file rtkThreeDCircularProjectionGeometry.h.
|
overridevirtual |
Empty the geometry object.
Reimplemented from rtk::ProjectionGeometry< 3 >.
|
static |
Compute the magnification matrix from 3D to 2D given a source to detector and to isocenter distance.
|
static |
Compute rotation matrix in homogeneous coordinates from 3 angles in degrees. The convention is the default in itk, i.e. ZXY of Euler angles.
|
static |
Compute translation matrix in homogeneous coordinates from translation parameters.
|
static |
Compute translation matrix in homogeneous coordinates from translation parameters.
|
static |
This function wraps an angle value between 0 and 2*PI radians.
Referenced by CheckGeometries().
|
static |
This function wraps an angle value between 0 and 360 degrees.
|
static |
This function wraps an angle value between -PI and PI radians.
|
virtual |
Reimplemented from itk::Object.
|
protected |
Try to fix Euler angles, which were found incorrect, to match the specified reference matrix.
[out] | outOfPlaneAngleRAD | out-of-plane angle of the detector in radians; if this method returns TRUE, this angle can be safely considered |
[out] | gantryAngleRAD | gantry angle of the detector in radians; if this method returns TRUE, this angle can be safely considered |
[out] | inPlaneAngleRAD | in-plane angle of the detector in radians; if this method returns TRUE, this angle can be safely considered |
referenceMatrix | reference matrix which reflects detector orientation in WCS |
const std::vector<double> rtk::ThreeDCircularProjectionGeometry::GetAngularGaps | ( | const std::vector< double > & | angles | ) |
Get for each projection half the angular distance between the previous and the next projection in radians.
const std::vector<double> rtk::ThreeDCircularProjectionGeometry::GetAngularGapsWithNext | ( | const std::vector< double > & | angles | ) | const |
Get for each projection the angular gaps with next projection in radians.
|
inline |
Get the vector containing the collimation jaw parameters.
Definition at line 300 of file rtkThreeDCircularProjectionGeometry.h.
Referenced by CheckGeometries().
|
inline |
Get the vector containing the collimation jaw parameters.
Definition at line 305 of file rtkThreeDCircularProjectionGeometry.h.
Referenced by CheckGeometries().
|
inline |
Get the vector containing the collimation jaw parameters.
Definition at line 310 of file rtkThreeDCircularProjectionGeometry.h.
Referenced by CheckGeometries().
|
inline |
Get the vector containing the collimation jaw parameters.
Definition at line 315 of file rtkThreeDCircularProjectionGeometry.h.
Referenced by CheckGeometries().
|
virtual |
|
inline |
Get the vector of geometry parameters (one per projection). Angles are in radians.
Definition at line 155 of file rtkThreeDCircularProjectionGeometry.h.
Referenced by CheckGeometries().
|
inline |
Get the vector of geometry parameters (one per projection). Angles are in radians.
Definition at line 165 of file rtkThreeDCircularProjectionGeometry.h.
Referenced by CheckGeometries().
|
inline |
Definition at line 288 of file rtkThreeDCircularProjectionGeometry.h.
|
inline |
Definition at line 293 of file rtkThreeDCircularProjectionGeometry.h.
|
inline |
Get the vector of geometry parameters (one per projection). Angles are in radians.
Definition at line 160 of file rtkThreeDCircularProjectionGeometry.h.
Referenced by CheckGeometries().
const ThreeDHomogeneousMatrixType rtk::ThreeDCircularProjectionGeometry::GetProjectionCoordinatesToDetectorSystemMatrix | ( | const unsigned int | i | ) | const |
Compute the ith matrix to convert projection coordinates to coordinates in the detector coordinate system (u,v,u^v). Note that the matrix is square but the third element of the projection coordinates is ignored because projection coordinates are 2D. This is meant to manipulate more easily stack of projection images.
const ThreeDHomogeneousMatrixType rtk::ThreeDCircularProjectionGeometry::GetProjectionCoordinatesToFixedSystemMatrix | ( | const unsigned int | i | ) | const |
Compute the ith matrix to convert projection coordinates to coordinates in the fixed coordinate system. Note that the matrix is square but the third element of the projection coordinates is ignored because projection coordinates are 2D. This is meant to manipulate more easily stack of projection images.
|
inline |
Get the vector of geometry parameters (one per projection). Angles are in radians.
Definition at line 190 of file rtkThreeDCircularProjectionGeometry.h.
Referenced by CheckGeometries().
|
inline |
Get the vector of geometry parameters (one per projection). Angles are in radians.
Definition at line 195 of file rtkThreeDCircularProjectionGeometry.h.
Referenced by CheckGeometries().
|
inline |
Get the vector containing the sub matrices used to compute the main projection matrix.
Definition at line 254 of file rtkThreeDCircularProjectionGeometry.h.
|
inline |
Get the vector containing the sub matrices used to compute the main projection matrix.
Definition at line 259 of file rtkThreeDCircularProjectionGeometry.h.
|
virtual |
Accessor for the radius of curved detector. The default is 0 and it means a flat detector.
Referenced by CheckGeometries().
|
inline |
Definition at line 266 of file rtkThreeDCircularProjectionGeometry.h.
|
inline |
Definition at line 271 of file rtkThreeDCircularProjectionGeometry.h.
const std::multimap<double, unsigned int> rtk::ThreeDCircularProjectionGeometry::GetSortedAngles | ( | const std::vector< double > & | angles | ) | const |
Get a multimap containing all sorted angles in radians and corresponding index.
|
inline |
Get a vector containing the source angles in radians. The source angle is defined as the angle between the z-axis and the isocenter-source line projected on the central plane.
Definition at line 205 of file rtkThreeDCircularProjectionGeometry.h.
|
inline |
Get the vector of geometry parameters (one per projection). Angles are in radians.
Definition at line 175 of file rtkThreeDCircularProjectionGeometry.h.
Referenced by CheckGeometries().
|
inline |
Get the vector of geometry parameters (one per projection). Angles are in radians.
Definition at line 180 of file rtkThreeDCircularProjectionGeometry.h.
Referenced by CheckGeometries().
const HomogeneousVectorType rtk::ThreeDCircularProjectionGeometry::GetSourcePosition | ( | const unsigned int | i | ) | const |
Get the source position for the ith projection in the fixed reference system and in homogeneous coordinates.
|
inline |
Get the vector of geometry parameters (one per projection). Angles are in radians.
Definition at line 185 of file rtkThreeDCircularProjectionGeometry.h.
Referenced by CheckGeometries().
|
inline |
Get the vector of geometry parameters (one per projection). Angles are in radians.
Definition at line 170 of file rtkThreeDCircularProjectionGeometry.h.
Referenced by CheckGeometries().
|
inline |
Definition at line 277 of file rtkThreeDCircularProjectionGeometry.h.
|
inline |
Definition at line 282 of file rtkThreeDCircularProjectionGeometry.h.
const std::vector<double> rtk::ThreeDCircularProjectionGeometry::GetTiltAngles | ( | ) | const |
Get a vector containing the tilt angles in radians. The tilt angle is defined as the difference between -GantryAngle and the SourceAngle.
const std::map<double, unsigned int> rtk::ThreeDCircularProjectionGeometry::GetUniqueSortedAngles | ( | const std::vector< double > & | angles | ) | const |
Get a map containing unique sorted angles in radians and corresponding index.
|
virtual |
Accessor for the internal tolerance on angles.
|
overrideprotectedvirtual |
Clone the geometry object in a new one.
Reimplemented from itk::LightObject.
|
static |
Method for creation through the object factory.
void rtk::ThreeDCircularProjectionGeometry::SetCollimationOfLastProjection | ( | const double | uinf, |
const double | usup, | ||
const double | vinf, | ||
const double | vsup | ||
) |
Set the collimation of the latest added projection (to be called after AddProjection).
|
virtual |
|
virtual |
Accessor for the radius of curved detector. The default is 0 and it means a flat detector.
|
virtual |
Accessor for the internal tolerance on angles.
double rtk::ThreeDCircularProjectionGeometry::ToUntiltedCoordinateAtIsocenter | ( | const unsigned int | noProj, |
const double | tiltedCoord | ||
) | const |
Changes the coordinate on the projection image to the coordinate on a virtual detector that is perpendicular to the source to isocenter line and positioned at the isocenter. It is assumed that OutOfPlaneAngle=0 and InPlaneAngle=0.
|
protected |
Verify that the specified Euler angles in ZXY result in a rotation matrix which corresponds to the specified detector orientation. Rationale for this utility method is that in some situations numerical instabilities (e.g. if gantry=+90deg,in-plane=-90deg or vice versa, "invalid" angles may be computed using the standard ITK Euler transform) may occur.
outOfPlaneAngleRAD | out-of-plane angle of the detector in radians |
gantryAngleRAD | gantry angle of the detector in radians |
inPlaneAngleRAD | in-plane angle of the detector in radians |
referenceMatrix | reference matrix which reflects detector orientation in WCS |
|
protected |
Parameters of the collimation jaws. The collimation position is with respect to the distance of the m_RotationCenter along
Definition at line 484 of file rtkThreeDCircularProjectionGeometry.h.
|
protected |
Definition at line 485 of file rtkThreeDCircularProjectionGeometry.h.
|
protected |
Definition at line 486 of file rtkThreeDCircularProjectionGeometry.h.
|
protected |
Definition at line 487 of file rtkThreeDCircularProjectionGeometry.h.
|
protected |
Definition at line 497 of file rtkThreeDCircularProjectionGeometry.h.
|
protected |
Circular geometry parameters per projection (angles in degrees between 0 and 360).
Definition at line 463 of file rtkThreeDCircularProjectionGeometry.h.
|
protected |
Definition at line 465 of file rtkThreeDCircularProjectionGeometry.h.
|
protected |
Definition at line 491 of file rtkThreeDCircularProjectionGeometry.h.
|
protected |
Definition at line 464 of file rtkThreeDCircularProjectionGeometry.h.
|
protected |
Definition at line 471 of file rtkThreeDCircularProjectionGeometry.h.
|
protected |
Definition at line 472 of file rtkThreeDCircularProjectionGeometry.h.
|
protected |
Matrices to change coordiate systems.
Definition at line 490 of file rtkThreeDCircularProjectionGeometry.h.
|
protected |
Radius of curved detector. The default is 0 and it means a flat detector.
Definition at line 475 of file rtkThreeDCircularProjectionGeometry.h.
|
protected |
Definition at line 492 of file rtkThreeDCircularProjectionGeometry.h.
|
protected |
Definition at line 466 of file rtkThreeDCircularProjectionGeometry.h.
|
protected |
Definition at line 468 of file rtkThreeDCircularProjectionGeometry.h.
|
protected |
Definition at line 469 of file rtkThreeDCircularProjectionGeometry.h.
|
protected |
Definition at line 470 of file rtkThreeDCircularProjectionGeometry.h.
|
protected |
Definition at line 467 of file rtkThreeDCircularProjectionGeometry.h.
|
protected |
Definition at line 493 of file rtkThreeDCircularProjectionGeometry.h.
|
protected |
Internal tolerance parameters.
Definition at line 496 of file rtkThreeDCircularProjectionGeometry.h.