ITK  4.1.0
Insight Segmentation and Registration Toolkit
Public Types | Public Member Functions | Static Public Member Functions | Protected Member Functions | Private Member Functions | Private Attributes
itk::OrientImageFilter< TInputImage, TOutputImage > Class Template Reference

#include <itkOrientImageFilter.h>

+ Inheritance diagram for itk::OrientImageFilter< TInputImage, TOutputImage >:
+ Collaboration diagram for itk::OrientImageFilter< TInputImage, TOutputImage >:

List of all members.

Public Types

typedef SmartPointer< const SelfConstPointer
typedef
SpatialOrientation::ValidCoordinateOrientationFlags 
CoordinateOrientationCode
typedef
FlipperType::FlipAxesArrayType 
FlipAxesArrayType
typedef FlipImageFilter
< TInputImage > 
FlipperType
typedef
InputImageType::ConstPointer 
InputImageConstPointer
typedef InputImageType::PixelType InputImagePixelType
typedef InputImageType::Pointer InputImagePointer
typedef InputImageType::RegionType InputImageRegionType
typedef TInputImage InputImageType
typedef
OutputImageType::ConstPointer 
OutputImageConstPointer
typedef OutputImageType::PixelType OutputImagePixelType
typedef OutputImageType::Pointer OutputImagePointer
typedef OutputImageType::RegionType OutputImageRegionType
typedef TOutputImage OutputImageType
typedef
PermuterType::PermuteOrderArrayType 
PermuteOrderArrayType
typedef PermuteAxesImageFilter
< TInputImage > 
PermuterType
typedef SmartPointer< SelfPointer
typedef OrientImageFilter Self
typedef ImageToImageFilter
< TInputImage, TOutputImage > 
Superclass

Public Member Functions

virtual ::itk::LightObject::Pointer CreateAnother (void) const
virtual void GenerateOutputInformation ()
virtual CoordinateOrientationCode GetDesiredCoordinateOrientation () const
virtual const FlipAxesArrayTypeGetFlipAxes ()
virtual const char * GetNameOfClass () const
virtual const
PermuteOrderArrayType
GetPermuteOrder ()
void SetDesiredCoordinateDirection (const typename TOutputImage::DirectionType &DesiredDirection)
void SetDesiredCoordinateOrientation (CoordinateOrientationCode newCode)
void SetDesiredCoordinateOrientationToAxial ()
void SetDesiredCoordinateOrientationToCoronal ()
void SetDesiredCoordinateOrientationToSagittal ()
void SetGivenCoordinateDirection (const typename TInputImage::DirectionType &GivenDirection)
 typedef (Concept::Convertible< InputImagePixelType, OutputImagePixelType >) InputConvertibleToOutput
 typedef (Concept::SameDimension< itkGetStaticConstMacro(InputImageDimension), itkGetStaticConstMacro(OutputImageDimension) >) SameDimension
 typedef (Concept::SameDimension< itkGetStaticConstMacro(InputImageDimension), 3 >) DimensionShouldBe3
virtual CoordinateOrientationCode GetGivenCoordinateOrientation () const
void SetGivenCoordinateOrientation (CoordinateOrientationCode newCode)
virtual void UseImageDirectionOn ()
virtual void UseImageDirectionOff ()
virtual bool GetUseImageDirection () const
virtual void SetUseImageDirection (bool _arg)

Static Public Member Functions

static Pointer New ()

Static Public Attributes

static const unsigned int InputImageDimension = TInputImage::ImageDimension
static const unsigned int OutputImageDimension = TOutputImage::ImageDimension

Protected Member Functions

void DeterminePermutationsAndFlips (const SpatialOrientation::ValidCoordinateOrientationFlags fixed_orient, const SpatialOrientation::ValidCoordinateOrientationFlags moving_orient)
void EnlargeOutputRequestedRegion (DataObject *)
void GenerateData ()
void GenerateInputRequestedRegion ()
bool NeedToFlip ()
bool NeedToPermute ()
 OrientImageFilter ()
 ~OrientImageFilter ()
void PrintSelf (std::ostream &os, Indent indent) const

Private Member Functions

std::string GetMajorAxisFromPatientRelativeDirectionCosine (double x, double y, double z)
void operator= (const Self &)
 OrientImageFilter (const Self &)

Private Attributes

std::map
< CoordinateOrientationCode,
std::string > 
m_CodeToString
CoordinateOrientationCode m_DesiredCoordinateOrientation
FlipAxesArrayType m_FlipAxes
CoordinateOrientationCode m_GivenCoordinateOrientation
PermuteOrderArrayType m_PermuteOrder
std::map< std::string,
CoordinateOrientationCode
m_StringToCode
bool m_UseImageDirection

Detailed Description

template<class TInputImage, class TOutputImage>
class itk::OrientImageFilter< TInputImage, TOutputImage >

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 a set of constants that describe all possible labeled according to the following scheme: Directions are labeled in terms of following pairs:

The initials of these directions are used in a 3 letter code in the enumerated type itk::SpatialOrientation::ValidCoordinateOrientationFlags. The initials are given fastest moving index first, second fastest second, third fastest third. Examples:

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.

As an example, if you wished to keep all images within your program in the orientation corresponding to the Analyze file format's 'CORONAL' orientation you could do something like the following

 // DEPRECATED -- using metadata for orientation is no longer supported
 //
 #include "itkAnalyzeImageIO.h"
 #include "itkMetaDataObject.h"
 #include "itkImage.h"
 #include "itkOrientImageFilter.h"
 typedef itk::Image<unsigned char,3> ImageType;
 typedef itk::ImageFileReader< TstImageType > ImageReaderType;
 ImageType::Pointer ReadAnalyzeFile(const char *path)
 {
   itk::AnalyzeImageIO::Pointer io = itk::AnalyzeImageIO::New();
   ImageReaderType::Pointer fileReader = ImageReaderType::New();
   fileReader->SetImageIO(io);
   fileReader->SetFileName(path);
   fileReader->Update();
   ImageType::Pointer rval = fileReader->GetOutput();

 Deprecated -- use direction cosines
  itk::SpatialOrientation::ValidCoordinateOrientationFlags fileOrientation;
  itk::ExposeMetaData<itk::SpatialOrientation::ValidCoordinateOrientationFlags>
    (rval->GetMetaDataDictionary(),itk::ITK_CoordinateOrientation,fileOrientation);
   itk::OrientImageFilter<ImageType,ImageType>::Pointer orienter =
     itk::OrientImageFilter<ImageType,ImageType>::New();
   orienter->SetGivenCoordinateOrientation(fileOrientation); // deprecated

   orienter->SetDesiredCoordinateOrientation(itk::SpatialOrientation::ITK_COORDINATE_ORIENTATION_RIP);
   orienter->SetInput(rval);
   orienter->Update();
   rval = orienter->GetOutput();
   return rval;
 }

Or, using the direction cosines of the image,

 #include "itkAnalyzeImageIO.h"
 #include "itkImage.h"
 #include "itkOrientImageFilter.h"
 typedef itk::Image<unsigned char,3> ImageType;
 typedef itk::ImageFileReader< TstImageType > ImageReaderType;
 ImageType::Pointer ReadAnalyzeFile(const char *path)
 {
   itk::AnalyzeImageIO::Pointer io = itk::AnalyzeImageIO::New();
   ImageReaderType::Pointer fileReader = ImageReaderType::New();
   fileReader->SetImageIO(io);
   fileReader->SetFileName(path);
   fileReader->Update();
   ImageType::Pointer rval = fileReader->GetOutput();

   itk::OrientImageFilter<ImageType,ImageType>::Pointer orienter =
     itk::OrientImageFilter<ImageType,ImageType>::New();
   orienter->UseImageDirectionOn();
   orienter->SetDesiredCoordinateOrientation(itk::SpatialOrientation::ITK_COORDINATE_ORIENTATION_RIP);
   orienter->SetInput(rval);
   orienter->Update();
   rval = orienter->GetOutput();
   return rval;
 }

Definition at line 140 of file itkOrientImageFilter.h.


Member Typedef Documentation

template<class TInputImage , class TOutputImage >
typedef SmartPointer< const Self > itk::OrientImageFilter< TInputImage, TOutputImage >::ConstPointer

Reimplemented from itk::ImageToImageFilter< TInputImage, TOutputImage >.

Definition at line 148 of file itkOrientImageFilter.h.

template<class TInputImage , class TOutputImage >
typedef SpatialOrientation::ValidCoordinateOrientationFlags itk::OrientImageFilter< TInputImage, TOutputImage >::CoordinateOrientationCode

Definition at line 162 of file itkOrientImageFilter.h.

template<class TInputImage , class TOutputImage >
typedef FlipperType::FlipAxesArrayType itk::OrientImageFilter< TInputImage, TOutputImage >::FlipAxesArrayType

Definition at line 170 of file itkOrientImageFilter.h.

template<class TInputImage , class TOutputImage >
typedef FlipImageFilter< TInputImage > itk::OrientImageFilter< TInputImage, TOutputImage >::FlipperType

Axes flipper type.

Definition at line 169 of file itkOrientImageFilter.h.

template<class TInputImage , class TOutputImage >
typedef InputImageType::ConstPointer itk::OrientImageFilter< TInputImage, TOutputImage >::InputImageConstPointer

Reimplemented from itk::ImageToImageFilter< TInputImage, TOutputImage >.

Definition at line 153 of file itkOrientImageFilter.h.

template<class TInputImage , class TOutputImage >
typedef InputImageType::PixelType itk::OrientImageFilter< TInputImage, TOutputImage >::InputImagePixelType

Reimplemented from itk::ImageToImageFilter< TInputImage, TOutputImage >.

Definition at line 155 of file itkOrientImageFilter.h.

template<class TInputImage , class TOutputImage >
typedef InputImageType::Pointer itk::OrientImageFilter< TInputImage, TOutputImage >::InputImagePointer

Reimplemented from itk::ImageToImageFilter< TInputImage, TOutputImage >.

Definition at line 152 of file itkOrientImageFilter.h.

template<class TInputImage , class TOutputImage >
typedef InputImageType::RegionType itk::OrientImageFilter< TInputImage, TOutputImage >::InputImageRegionType

Reimplemented from itk::ImageToImageFilter< TInputImage, TOutputImage >.

Definition at line 154 of file itkOrientImageFilter.h.

template<class TInputImage , class TOutputImage >
typedef TInputImage itk::OrientImageFilter< TInputImage, TOutputImage >::InputImageType

Some convenient typedefs.

Reimplemented from itk::ImageToImageFilter< TInputImage, TOutputImage >.

Definition at line 151 of file itkOrientImageFilter.h.

template<class TInputImage , class TOutputImage >
typedef OutputImageType::ConstPointer itk::OrientImageFilter< TInputImage, TOutputImage >::OutputImageConstPointer

Definition at line 158 of file itkOrientImageFilter.h.

template<class TInputImage , class TOutputImage >
typedef OutputImageType::PixelType itk::OrientImageFilter< TInputImage, TOutputImage >::OutputImagePixelType

Reimplemented from itk::ImageToImageFilter< TInputImage, TOutputImage >.

Definition at line 160 of file itkOrientImageFilter.h.

template<class TInputImage , class TOutputImage >
typedef OutputImageType::Pointer itk::OrientImageFilter< TInputImage, TOutputImage >::OutputImagePointer

Reimplemented from itk::ImageSource< TOutputImage >.

Definition at line 157 of file itkOrientImageFilter.h.

template<class TInputImage , class TOutputImage >
typedef OutputImageType::RegionType itk::OrientImageFilter< TInputImage, TOutputImage >::OutputImageRegionType

Superclass typedefs.

Reimplemented from itk::ImageToImageFilter< TInputImage, TOutputImage >.

Definition at line 159 of file itkOrientImageFilter.h.

template<class TInputImage , class TOutputImage >
typedef TOutputImage itk::OrientImageFilter< TInputImage, TOutputImage >::OutputImageType

Some convenient typedefs.

Reimplemented from itk::ImageSource< TOutputImage >.

Definition at line 156 of file itkOrientImageFilter.h.

template<class TInputImage , class TOutputImage >
typedef PermuterType::PermuteOrderArrayType itk::OrientImageFilter< TInputImage, TOutputImage >::PermuteOrderArrayType

Definition at line 166 of file itkOrientImageFilter.h.

template<class TInputImage , class TOutputImage >
typedef PermuteAxesImageFilter< TInputImage > itk::OrientImageFilter< TInputImage, TOutputImage >::PermuterType

Axes permuter type.

Definition at line 165 of file itkOrientImageFilter.h.

template<class TInputImage , class TOutputImage >
typedef SmartPointer< Self > itk::OrientImageFilter< TInputImage, TOutputImage >::Pointer

Reimplemented from itk::ImageToImageFilter< TInputImage, TOutputImage >.

Definition at line 147 of file itkOrientImageFilter.h.

template<class TInputImage , class TOutputImage >
typedef OrientImageFilter itk::OrientImageFilter< TInputImage, TOutputImage >::Self

Standard class typedefs.

Reimplemented from itk::ImageToImageFilter< TInputImage, TOutputImage >.

Definition at line 145 of file itkOrientImageFilter.h.

template<class TInputImage , class TOutputImage >
typedef ImageToImageFilter< TInputImage, TOutputImage > itk::OrientImageFilter< TInputImage, TOutputImage >::Superclass

Reimplemented from itk::ImageToImageFilter< TInputImage, TOutputImage >.

Definition at line 146 of file itkOrientImageFilter.h.


Constructor & Destructor Documentation

template<class TInputImage , class TOutputImage >
itk::OrientImageFilter< TInputImage, TOutputImage >::OrientImageFilter ( ) [protected]

End concept checking

template<class TInputImage , class TOutputImage >
itk::OrientImageFilter< TInputImage, TOutputImage >::~OrientImageFilter ( ) [inline, protected]

End concept checking

Definition at line 274 of file itkOrientImageFilter.h.

template<class TInputImage , class TOutputImage >
itk::OrientImageFilter< TInputImage, TOutputImage >::OrientImageFilter ( const Self ) [private]

Member Function Documentation

template<class TInputImage , class TOutputImage >
virtual::itk::LightObject::Pointer itk::OrientImageFilter< TInputImage, TOutputImage >::CreateAnother ( void  ) const [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.

template<class TInputImage , class TOutputImage >
void itk::OrientImageFilter< TInputImage, TOutputImage >::DeterminePermutationsAndFlips ( const SpatialOrientation::ValidCoordinateOrientationFlags  fixed_orient,
const SpatialOrientation::ValidCoordinateOrientationFlags  moving_orient 
) [protected]
template<class TInputImage , class TOutputImage >
void itk::OrientImageFilter< TInputImage, TOutputImage >::EnlargeOutputRequestedRegion ( DataObject ) [protected, virtual]

OrientImageFilter will produce the entire output.

Reimplemented from itk::ProcessObject.

template<class TInputImage , class TOutputImage >
void itk::OrientImageFilter< TInputImage, TOutputImage >::GenerateData ( ) [protected, virtual]

Single-threaded version of GenerateData. This filter delegates to PermuteAxesImageFilter and FlipImageFilter.

Reimplemented from itk::ImageSource< TOutputImage >.

template<class TInputImage , class TOutputImage >
void itk::OrientImageFilter< TInputImage, TOutputImage >::GenerateInputRequestedRegion ( ) [protected, virtual]

OrientImageFilter needs the entire input be available. Thus, it needs to provide an implementation of GenerateInputRequestedRegion().

Reimplemented from itk::ImageToImageFilter< TInputImage, TOutputImage >.

template<class TInputImage , class TOutputImage >
virtual void itk::OrientImageFilter< TInputImage, TOutputImage >::GenerateOutputInformation ( ) [virtual]

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.

See also:
ProcessObject::GenerateOutputInformaton()

Reimplemented from itk::ProcessObject.

template<class TInputImage , class TOutputImage >
virtual CoordinateOrientationCode itk::OrientImageFilter< TInputImage, TOutputImage >::GetDesiredCoordinateOrientation ( ) const [virtual]
template<class TInputImage , class TOutputImage >
virtual const FlipAxesArrayType& itk::OrientImageFilter< TInputImage, TOutputImage >::GetFlipAxes ( ) [virtual]

Get flip axes.

template<class TInputImage , class TOutputImage >
virtual CoordinateOrientationCode itk::OrientImageFilter< TInputImage, TOutputImage >::GetGivenCoordinateOrientation ( ) const [virtual]

Set/Get the orientation codes to define the coordinate transform.

template<class TInputImage , class TOutputImage >
std::string itk::OrientImageFilter< TInputImage, TOutputImage >::GetMajorAxisFromPatientRelativeDirectionCosine ( double  x,
double  y,
double  z 
) [private]
template<class TInputImage , class TOutputImage >
virtual const char* itk::OrientImageFilter< TInputImage, TOutputImage >::GetNameOfClass ( ) const [virtual]

Runtime information support.

Reimplemented from itk::ImageToImageFilter< TInputImage, TOutputImage >.

template<class TInputImage , class TOutputImage >
virtual const PermuteOrderArrayType& itk::OrientImageFilter< TInputImage, TOutputImage >::GetPermuteOrder ( ) [virtual]

Get axes permute order.

template<class TInputImage , class TOutputImage >
virtual bool itk::OrientImageFilter< TInputImage, TOutputImage >::GetUseImageDirection ( ) const [virtual]

Controls how the GivenCoordinateOrientation is determined. If "on", the direction cosines determine the coordinate orientation. If "off", the user must use the SetGivenCoordinateOrientation method to establish the orientation. For compatbility with the original API, the default if "off".

template<class TInputImage , class TOutputImage >
bool itk::OrientImageFilter< TInputImage, TOutputImage >::NeedToFlip ( ) [protected]
template<class TInputImage , class TOutputImage >
bool itk::OrientImageFilter< TInputImage, TOutputImage >::NeedToPermute ( ) [protected]
template<class TInputImage , class TOutputImage >
static Pointer itk::OrientImageFilter< TInputImage, TOutputImage >::New ( ) [static]

Standard New method.

Reimplemented from itk::Object.

template<class TInputImage , class TOutputImage >
void itk::OrientImageFilter< TInputImage, TOutputImage >::operator= ( const Self ) [private]

PushBackInput(), PushFronInput() in the public section force the input to be the type expected by an ImageToImageFilter. However, these methods end of "hiding" the versions from the superclass (ProcessObject) whose arguments are DataObjects. Here, we re-expose the versions from ProcessObject to avoid warnings about hiding methods from the superclass.

Reimplemented from itk::ImageToImageFilter< TInputImage, TOutputImage >.

template<class TInputImage , class TOutputImage >
void itk::OrientImageFilter< TInputImage, TOutputImage >::PrintSelf ( std::ostream &  os,
Indent  indent 
) const [protected, virtual]

End concept checking

Reimplemented from itk::ImageToImageFilter< TInputImage, TOutputImage >.

template<class TInputImage , class TOutputImage >
void itk::OrientImageFilter< TInputImage, TOutputImage >::SetDesiredCoordinateDirection ( const typename TOutputImage::DirectionType &  DesiredDirection) [inline]

Definition at line 199 of file itkOrientImageFilter.h.

template<class TInputImage , class TOutputImage >
void itk::OrientImageFilter< TInputImage, TOutputImage >::SetDesiredCoordinateOrientation ( CoordinateOrientationCode  newCode)
template<class TInputImage , class TOutputImage >
void itk::OrientImageFilter< TInputImage, TOutputImage >::SetDesiredCoordinateOrientationToAxial ( ) [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 SetDesiredCoordinateOrientation (ITK_COORDINATE_ORIENTATION_RAI).

SetDesiredCoordinateOrientationToCoronal is equivalent to SetDesiredCoordinateOrientation (ITK_COORDINATE_ORIENTATION_RSA).

SetDesiredCoordinateOrientationToSagittal is equivalent to SetDesiredCoordinateOrientation (ITK_COORDINATE_ORIENTATION_ASL).

Definition at line 236 of file itkOrientImageFilter.h.

References itk::SpatialOrientation::ITK_COORDINATE_ORIENTATION_RAI.

template<class TInputImage , class TOutputImage >
void itk::OrientImageFilter< TInputImage, TOutputImage >::SetDesiredCoordinateOrientationToCoronal ( ) [inline]
template<class TInputImage , class TOutputImage >
void itk::OrientImageFilter< TInputImage, TOutputImage >::SetDesiredCoordinateOrientationToSagittal ( ) [inline]
template<class TInputImage , class TOutputImage >
void itk::OrientImageFilter< TInputImage, TOutputImage >::SetGivenCoordinateDirection ( const typename TInputImage::DirectionType &  GivenDirection) [inline]

Definition at line 190 of file itkOrientImageFilter.h.

template<class TInputImage , class TOutputImage >
void itk::OrientImageFilter< TInputImage, TOutputImage >::SetGivenCoordinateOrientation ( CoordinateOrientationCode  newCode)

Set/Get the orientation codes to define the coordinate transform.

template<class TInputImage , class TOutputImage >
virtual void itk::OrientImageFilter< TInputImage, TOutputImage >::SetUseImageDirection ( bool  _arg) [virtual]

Controls how the GivenCoordinateOrientation is determined. If "on", the direction cosines determine the coordinate orientation. If "off", the user must use the SetGivenCoordinateOrientation method to establish the orientation. For compatbility with the original API, the default if "off".

template<class TInputImage , class TOutputImage >
itk::OrientImageFilter< TInputImage, TOutputImage >::typedef ( Concept::Convertible< InputImagePixelType, OutputImagePixelType )

Begin concept checking This class requires InputConvertibleToOutput in the form of ( Concept::Convertible< InputImagePixelType, OutputImagePixelType > )

template<class TInputImage , class TOutputImage >
itk::OrientImageFilter< TInputImage, TOutputImage >::typedef ( Concept::SameDimension< itkGetStaticConstMacro(InputImageDimension), itkGetStaticConstMacro(OutputImageDimension) >  )

This class requires SameDimension in the form of ( Concept::SameDimension< itkGetStaticConstMacro(InputImageDimension), itkGetStaticConstMacro(OutputImageDimension) > )

template<class TInputImage , class TOutputImage >
itk::OrientImageFilter< TInputImage, TOutputImage >::typedef ( Concept::SameDimension< itkGetStaticConstMacro(InputImageDimension), 3 >  )

This class requires DimensionShouldBe3 in the form of ( Concept::SameDimension< itkGetStaticConstMacro(InputImageDimension), 3 > )

template<class TInputImage , class TOutputImage >
virtual void itk::OrientImageFilter< TInputImage, TOutputImage >::UseImageDirectionOff ( ) [virtual]

Controls how the GivenCoordinateOrientation is determined. If "on", the direction cosines determine the coordinate orientation. If "off", the user must use the SetGivenCoordinateOrientation method to establish the orientation. For compatbility with the original API, the default if "off".

template<class TInputImage , class TOutputImage >
virtual void itk::OrientImageFilter< TInputImage, TOutputImage >::UseImageDirectionOn ( ) [virtual]

Controls how the GivenCoordinateOrientation is determined. If "on", the direction cosines determine the coordinate orientation. If "off", the user must use the SetGivenCoordinateOrientation method to establish the orientation. For compatbility with the original API, the default if "off".


Member Data Documentation

template<class TInputImage , class TOutputImage >
const unsigned int itk::OrientImageFilter< TInputImage, TOutputImage >::InputImageDimension = TInputImage::ImageDimension [static]

ImageDimension constants

Reimplemented from itk::ImageToImageFilter< TInputImage, TOutputImage >.

Definition at line 174 of file itkOrientImageFilter.h.

template<class TInputImage , class TOutputImage >
std::map< CoordinateOrientationCode, std::string > itk::OrientImageFilter< TInputImage, TOutputImage >::m_CodeToString [private]

Definition at line 312 of file itkOrientImageFilter.h.

template<class TInputImage , class TOutputImage >
CoordinateOrientationCode itk::OrientImageFilter< TInputImage, TOutputImage >::m_DesiredCoordinateOrientation [private]

Definition at line 305 of file itkOrientImageFilter.h.

template<class TInputImage , class TOutputImage >
FlipAxesArrayType itk::OrientImageFilter< TInputImage, TOutputImage >::m_FlipAxes [private]

Definition at line 309 of file itkOrientImageFilter.h.

template<class TInputImage , class TOutputImage >
CoordinateOrientationCode itk::OrientImageFilter< TInputImage, TOutputImage >::m_GivenCoordinateOrientation [private]

Definition at line 304 of file itkOrientImageFilter.h.

template<class TInputImage , class TOutputImage >
PermuteOrderArrayType itk::OrientImageFilter< TInputImage, TOutputImage >::m_PermuteOrder [private]

Definition at line 308 of file itkOrientImageFilter.h.

template<class TInputImage , class TOutputImage >
std::map< std::string, CoordinateOrientationCode > itk::OrientImageFilter< TInputImage, TOutputImage >::m_StringToCode [private]

Definition at line 311 of file itkOrientImageFilter.h.

template<class TInputImage , class TOutputImage >
bool itk::OrientImageFilter< TInputImage, TOutputImage >::m_UseImageDirection [private]

Definition at line 306 of file itkOrientImageFilter.h.

template<class TInputImage , class TOutputImage >
const unsigned int itk::OrientImageFilter< TInputImage, TOutputImage >::OutputImageDimension = TOutputImage::ImageDimension [static]

ImageDimension constants

Reimplemented from itk::ImageToImageFilter< TInputImage, TOutputImage >.

Definition at line 176 of file itkOrientImageFilter.h.


The documentation for this class was generated from the following file: