#include <itkOrientImageFilter.h>
Permute axes and then flip images as needed to obtain agreement in coordinateOrientation codes.
This class satisfies a common requirement in medical imaging, which is to properly orient a 3 dimensional image with respect to anatomical features. Due to the wide variety of hardware used to generate 3D images of human anatomy, and the even wider variety of image processing software, it is often necessary to re-orient image volume data.
OrientImageFilter depends on representations of the orientation defined in the AnatomicalOrientation class. The orientation is represented by the following anatomical terms with respect to the subject or patient:
An AnatomicalOrientation object can be constructed unambiguously with the following syntax:
The orientations were previously defined in the itk::SpatialOrientation class. However, the 3 letter code is ambiguous with label and the direction the axis is increasing.
In order to use this filter, you need to supply an input image, the current orientation of the input image (set with SetGivenCoordinateOrientation) and the desired orientation. (set with SetDesiredCoordinateOrientation). You may explicitly set the DesiredOrientation with SetDesiredCoordinateOrientation (if UseImageDirection is "off") or you can use the image's direction cosines to set the DesiredOrientation (if UseImageDirection is "on"). When reading image files that define the coordinate orientation of the image, the current orientation is stored in the MetadataDictionary for the itk::Image object and the Image.Direction direction cosine matrix created from the file.
Definition at line 74 of file itkOrientImageFilter.h.
Static Public Member Functions | |
static Pointer | New () |
![]() | |
static double | GetGlobalDefaultCoordinateTolerance () |
static double | GetGlobalDefaultDirectionTolerance () |
static void | SetGlobalDefaultCoordinateTolerance (double) |
static void | SetGlobalDefaultDirectionTolerance (double) |
![]() | |
static bool | GetGlobalWarningDisplay () |
static void | GlobalWarningDisplayOff () |
static void | GlobalWarningDisplayOn () |
static Pointer | New () |
static void | SetGlobalWarningDisplay (bool val) |
![]() | |
static void | BreakOnError () |
static Pointer | New () |
Static Public Attributes | |
static constexpr unsigned int | InputImageDimension = TInputImage::ImageDimension |
static constexpr unsigned int | OutputImageDimension = TOutputImage::ImageDimension |
![]() | |
static constexpr unsigned int | InputImageDimension = TInputImage::ImageDimension |
static constexpr unsigned int | OutputImageDimension = TOutputImage::ImageDimension |
![]() | |
static constexpr unsigned int | OutputImageDimension = TOutputImage::ImageDimension |
Private Attributes | |
CoordinateOrientationCode | m_DesiredCoordinateOrientation |
FlipAxesArrayType | m_FlipAxes { false } |
CoordinateOrientationCode | m_GivenCoordinateOrientation |
PermuteOrderArrayType | m_PermuteOrder {} |
bool | m_UseImageDirection { false } |
Additional Inherited Members | |
![]() | |
using | InputToOutputRegionCopierType |
using | OutputToInputRegionCopierType |
![]() | |
static const ImageRegionSplitterBase * | GetGlobalDefaultSplitter () |
static ITK_THREAD_RETURN_FUNCTION_CALL_CONVENTION | ThreaderCallback (void *arg) |
![]() | |
template<typename TSourceObject> | |
static void | MakeRequiredOutputs (TSourceObject &sourceObject, const DataObjectPointerArraySizeType numberOfRequiredOutputs) |
static constexpr float | progressFixedToFloat (uint32_t fixed) |
static uint32_t | progressFloatToFixed (float f) |
![]() | |
bool | m_DynamicMultiThreading { true } |
![]() | |
TimeStamp | m_OutputInformationMTime {} |
bool | m_Updating {} |
![]() | |
std::atomic< int > | m_ReferenceCount {} |
using itk::OrientImageFilter< TInputImage, TOutputImage >::ConstPointer = SmartPointer<const Self> |
Definition at line 83 of file itkOrientImageFilter.h.
using itk::OrientImageFilter< TInputImage, TOutputImage >::CoordinateOrientationCode = AnatomicalOrientation |
Definition at line 96 of file itkOrientImageFilter.h.
using itk::OrientImageFilter< TInputImage, TOutputImage >::FlipAxesArrayType = typename FlipperType::FlipAxesArrayType |
Definition at line 104 of file itkOrientImageFilter.h.
using itk::OrientImageFilter< TInputImage, TOutputImage >::FlipperType = FlipImageFilter<TInputImage> |
Axes flipper type.
Definition at line 103 of file itkOrientImageFilter.h.
using itk::OrientImageFilter< TInputImage, TOutputImage >::InputImageConstPointer = typename InputImageType::ConstPointer |
Definition at line 88 of file itkOrientImageFilter.h.
using itk::OrientImageFilter< TInputImage, TOutputImage >::InputImagePixelType = typename InputImageType::PixelType |
Definition at line 90 of file itkOrientImageFilter.h.
using itk::OrientImageFilter< TInputImage, TOutputImage >::InputImagePointer = typename InputImageType::Pointer |
Definition at line 87 of file itkOrientImageFilter.h.
using itk::OrientImageFilter< TInputImage, TOutputImage >::InputImageRegionType = typename InputImageType::RegionType |
Definition at line 89 of file itkOrientImageFilter.h.
using itk::OrientImageFilter< TInputImage, TOutputImage >::InputImageType = TInputImage |
Some convenient type alias.
Definition at line 86 of file itkOrientImageFilter.h.
using itk::OrientImageFilter< TInputImage, TOutputImage >::OutputImageConstPointer = typename OutputImageType::ConstPointer |
Definition at line 93 of file itkOrientImageFilter.h.
using itk::OrientImageFilter< TInputImage, TOutputImage >::OutputImagePixelType = typename OutputImageType::PixelType |
Definition at line 95 of file itkOrientImageFilter.h.
using itk::OrientImageFilter< TInputImage, TOutputImage >::OutputImagePointer = typename OutputImageType::Pointer |
Definition at line 92 of file itkOrientImageFilter.h.
using itk::OrientImageFilter< TInputImage, TOutputImage >::OutputImageRegionType = typename OutputImageType::RegionType |
Definition at line 94 of file itkOrientImageFilter.h.
using itk::OrientImageFilter< TInputImage, TOutputImage >::OutputImageType = TOutputImage |
Definition at line 91 of file itkOrientImageFilter.h.
using itk::OrientImageFilter< TInputImage, TOutputImage >::PermuteOrderArrayType = typename PermuterType::PermuteOrderArrayType |
Definition at line 100 of file itkOrientImageFilter.h.
using itk::OrientImageFilter< TInputImage, TOutputImage >::PermuterType = PermuteAxesImageFilter<TInputImage> |
Axes permuter type.
Definition at line 99 of file itkOrientImageFilter.h.
using itk::OrientImageFilter< TInputImage, TOutputImage >::Pointer = SmartPointer<Self> |
Definition at line 82 of file itkOrientImageFilter.h.
using itk::OrientImageFilter< TInputImage, TOutputImage >::Self = OrientImageFilter |
Standard class type aliases.
Definition at line 80 of file itkOrientImageFilter.h.
using itk::OrientImageFilter< TInputImage, TOutputImage >::Superclass = ImageToImageFilter<TInputImage, TOutputImage> |
Definition at line 81 of file itkOrientImageFilter.h.
|
protected |
Referenced by GetNameOfClass().
|
overrideprotecteddefault |
|
virtual |
Create an object from an instance, potentially deferring to a factory. This method allows you to create an instance of an object that is exactly the same type as the referring object. This is useful in cases where an object has been cast back to a base class.
Reimplemented from itk::Object.
|
protected |
|
overrideprotectedvirtual |
OrientImageFilter will produce the entire output.
Reimplemented from itk::ProcessObject.
|
overrideprotectedvirtual |
Single-threaded version of GenerateData. This filter delegates to PermuteAxesImageFilter and FlipImageFilter.
Reimplemented from itk::ImageSource< TOutputImage >.
|
overrideprotectedvirtual |
OrientImageFilter needs the entire input be available. Thus, it needs to provide an implementation of GenerateInputRequestedRegion().
Reimplemented from itk::ImageToImageFilter< TInputImage, TOutputImage >.
|
overridevirtual |
OrientImageFilter produces an image which is a different dimensionality than its input image, in general. As such, OrientImageFilter needs to provide an implementation for GenerateOutputInformation() in order to inform the pipeline execution model. The original documentation of this method is below.
Reimplemented from itk::ProcessObject.
|
virtual |
|
virtual |
Get flip axes.
|
virtual |
Set/Get the orientation codes to define the coordinate transform.
|
overridevirtual |
Reimplemented from itk::ImageToImageFilter< TInputImage, TOutputImage >.
References OrientImageFilter().
|
virtual |
Get axes permute order.
|
virtual |
Controls how the GivenCoordinateOrientation is determined. If set to On, the direction cosines determine the coordinate orientation. If set to Off, the user must use the SetGivenCoordinateOrientation method to establish the orientation.
For compatibility with the original API, the default value is Off.
|
protected |
Returns true if flipping is required. Returns false otherwise.
|
protected |
Returns true if a permute is required. Returns false otherwise.
|
static |
Standard New method.
|
overrideprotectedvirtual |
Methods invoked by Print() to print information about the object including superclasses. Typically not called by the user (use Print() instead) but used in the hierarchical print process to combine the output of several classes.
Reimplemented from itk::ImageToImageFilter< TInputImage, TOutputImage >.
|
inline |
Definition at line 134 of file itkOrientImageFilter.h.
References SetDesiredCoordinateOrientation().
void itk::OrientImageFilter< TInputImage, TOutputImage >::SetDesiredCoordinateOrientation | ( | CoordinateOrientationCode | newCode | ) |
|
inline |
Convenience methods to set desired slice orientation These methods allow a limited selection of slice orientations without having to specify the SpatialOrientation.
SetDesiredCoordinateOrientationToAxial is equivalent to AnatomicalOrientation::PositiveEnum::LPS.
SetDesiredCoordinateOrientationToCoronal is equivalent to AnatomicalOrientation::PositiveEnum::LIP.
SetDesiredCoordinateOrientationToSagittal is equivalent to AnatomicalOrientation::PositiveEnum::PIR.
Definition at line 171 of file itkOrientImageFilter.h.
References itk::AnatomicalOrientation::AnteriorToPosterior, itk::AnatomicalOrientation::InferiorToSuperior, itk::AnatomicalOrientation::RightToLeft, and SetDesiredCoordinateOrientation().
|
inline |
Convenience methods to set desired slice orientation These methods allow a limited selection of slice orientations without having to specify the SpatialOrientation.
SetDesiredCoordinateOrientationToAxial is equivalent to AnatomicalOrientation::PositiveEnum::LPS.
SetDesiredCoordinateOrientationToCoronal is equivalent to AnatomicalOrientation::PositiveEnum::LIP.
SetDesiredCoordinateOrientationToSagittal is equivalent to AnatomicalOrientation::PositiveEnum::PIR.
Definition at line 179 of file itkOrientImageFilter.h.
References itk::AnatomicalOrientation::AnteriorToPosterior, itk::AnatomicalOrientation::RightToLeft, SetDesiredCoordinateOrientation(), and itk::AnatomicalOrientation::SuperiorToInferior.
|
inline |
Convenience methods to set desired slice orientation These methods allow a limited selection of slice orientations without having to specify the SpatialOrientation.
SetDesiredCoordinateOrientationToAxial is equivalent to AnatomicalOrientation::PositiveEnum::LPS.
SetDesiredCoordinateOrientationToCoronal is equivalent to AnatomicalOrientation::PositiveEnum::LIP.
SetDesiredCoordinateOrientationToSagittal is equivalent to AnatomicalOrientation::PositiveEnum::PIR.
Definition at line 187 of file itkOrientImageFilter.h.
References itk::AnatomicalOrientation::AnteriorToPosterior, itk::AnatomicalOrientation::LeftToRight, SetDesiredCoordinateOrientation(), and itk::AnatomicalOrientation::SuperiorToInferior.
|
inline |
Definition at line 124 of file itkOrientImageFilter.h.
References SetGivenCoordinateOrientation().
void itk::OrientImageFilter< TInputImage, TOutputImage >::SetGivenCoordinateOrientation | ( | CoordinateOrientationCode | newCode | ) |
Set/Get the orientation codes to define the coordinate transform.
Referenced by SetGivenCoordinateDirection().
|
virtual |
Controls how the GivenCoordinateOrientation is determined. If set to On, the direction cosines determine the coordinate orientation. If set to Off, the user must use the SetGivenCoordinateOrientation method to establish the orientation.
For compatibility with the original API, the default value is Off.
|
virtual |
Controls how the GivenCoordinateOrientation is determined. If set to On, the direction cosines determine the coordinate orientation. If set to Off, the user must use the SetGivenCoordinateOrientation method to establish the orientation.
For compatibility with the original API, the default value is Off.
|
virtual |
Controls how the GivenCoordinateOrientation is determined. If set to On, the direction cosines determine the coordinate orientation. If set to Off, the user must use the SetGivenCoordinateOrientation method to establish the orientation.
For compatibility with the original API, the default value is Off.
|
overrideprotectedvirtual |
Verifies that the process object has been configured correctly, that all required inputs are set, and needed parameters are set appropriately. If not valid an exceptions will be thrown.
This method is called before UpdateOutputInformation() is propagated to the inputs.
The ProcessObject's implementation verifies that the m_NumberOfRequiredInputs are set and not null.
Reimplemented from itk::ProcessObject.
|
staticconstexpr |
ImageDimension constants
Definition at line 107 of file itkOrientImageFilter.h.
|
private |
Definition at line 250 of file itkOrientImageFilter.h.
|
private |
Definition at line 258 of file itkOrientImageFilter.h.
|
private |
Definition at line 247 of file itkOrientImageFilter.h.
|
private |
Definition at line 257 of file itkOrientImageFilter.h.
|
private |
Definition at line 255 of file itkOrientImageFilter.h.
|
staticconstexpr |
Definition at line 108 of file itkOrientImageFilter.h.