ITK  6.0.0
Insight Toolkit
itkOrientImageFilter.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright NumFOCUS
4  *
5  * Licensed under the Apache License, Version 2.0 (the "License");
6  * you may not use this file except in compliance with the License.
7  * You may obtain a copy of the License at
8  *
9  * https://www.apache.org/licenses/LICENSE-2.0.txt
10  *
11  * Unless required by applicable law or agreed to in writing, software
12  * distributed under the License is distributed on an "AS IS" BASIS,
13  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14  * See the License for the specific language governing permissions and
15  * limitations under the License.
16  *
17  *=========================================================================*/
18 #ifndef itkOrientImageFilter_h
19 #define itkOrientImageFilter_h
20 
22 #include "itkFlipImageFilter.h"
24 #include <map>
25 #include <string>
26 
27 namespace itk
28 {
73 template <typename TInputImage, typename TOutputImage>
74 class ITK_TEMPLATE_EXPORT OrientImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
75 {
76 public:
77  ITK_DISALLOW_COPY_AND_MOVE(OrientImageFilter);
78 
84 
86  using InputImageType = TInputImage;
90  using InputImagePixelType = typename InputImageType::PixelType;
91  using OutputImageType = TOutputImage;
95  using OutputImagePixelType = typename OutputImageType::PixelType;
97 
101 
105 
107  static constexpr unsigned int InputImageDimension = TInputImage::ImageDimension;
108  static constexpr unsigned int OutputImageDimension = TOutputImage::ImageDimension;
109 
111  itkNewMacro(Self);
112 
114  itkOverrideGetNameOfClassMacro(OrientImageFilter);
115 
117  itkGetEnumMacro(GivenCoordinateOrientation, CoordinateOrientationCode);
118  void
119  SetGivenCoordinateOrientation(CoordinateOrientationCode newCode);
122  inline void
124  {
125  SetGivenCoordinateOrientation(AnatomicalOrientation(GivenDirection));
126  }
127 
128  itkGetEnumMacro(DesiredCoordinateOrientation, CoordinateOrientationCode);
129  void
130  SetDesiredCoordinateOrientation(CoordinateOrientationCode newCode);
131 
132  inline void
134  {
135  SetDesiredCoordinateOrientation(AnatomicalOrientation(DesiredDirection));
136  }
137 
146  itkBooleanMacro(UseImageDirection);
147  itkGetConstMacro(UseImageDirection, bool);
148  itkSetMacro(UseImageDirection, bool);
152  itkGetConstReferenceMacro(PermuteOrder, PermuteOrderArrayType);
153 
155  itkGetConstReferenceMacro(FlipAxes, FlipAxesArrayType);
156 
167  void
169  {
170  this->SetDesiredCoordinateOrientation({ AnatomicalOrientation::CoordinateEnum::RightToLeft,
173  }
176  void
178  {
179  this->SetDesiredCoordinateOrientation({ AnatomicalOrientation::CoordinateEnum::RightToLeft,
182  }
183 
184  void
186  {
187  this->SetDesiredCoordinateOrientation({ AnatomicalOrientation::CoordinateEnum::AnteriorToPosterior,
190  }
191 
199  void
200  GenerateOutputInformation() override;
201 
202 #ifdef ITK_USE_CONCEPT_CHECKING
203  // Begin concept checking
207  // End concept checking
208 #endif
209 
210 protected:
212  ~OrientImageFilter() override = default;
213  void
214  PrintSelf(std::ostream & os, Indent indent) const override;
215 
219  void
220  GenerateInputRequestedRegion() override;
221 
223  void
224  EnlargeOutputRequestedRegion(DataObject * itkNotUsed(output)) override;
225 
226  void
227  VerifyPreconditions() const override;
228 
229  /*** Member functions used by GenerateData: */
230  void
231  DeterminePermutationsAndFlips(const CoordinateOrientationCode fixed_orient,
232  const CoordinateOrientationCode moving_orient);
233 
235  bool
236  NeedToPermute();
237 
239  bool
240  NeedToFlip();
241 
244  void
245  GenerateData() override;
246 
247 private:
251  CoordinateOrientationCode m_DesiredCoordinateOrientation{
255  };
256  bool m_UseImageDirection{ false };
257 
258  PermuteOrderArrayType m_PermuteOrder{};
259  FlipAxesArrayType m_FlipAxes{ false };
260 
261 }; // end of class
262 } // end namespace itk
263 
264 #ifndef ITK_MANUAL_INSTANTIATION
265 # include "itkOrientImageFilter.hxx"
266 #endif
267 
268 #endif
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
ConstPointer
SmartPointer< const Self > ConstPointer
Definition: itkAddImageFilter.h:94
itk::OrientImageFilter::SetDesiredCoordinateOrientationToCoronal
void SetDesiredCoordinateOrientationToCoronal()
Definition: itkOrientImageFilter.h:177
itk::GTest::TypedefsAndConstructors::Dimension2::DirectionType
ImageBaseType::DirectionType DirectionType
Definition: itkGTestTypedefsAndConstructors.h:52
itk::ImageSource::OutputImagePointer
typename OutputImageType::Pointer OutputImagePointer
Definition: itkImageSource.h:91
itk::AnatomicalOrientation::CoordinateEnum::SuperiorToInferior
to foot
itk::OrientImageFilter
Permute axes and then flip images as needed to obtain agreement in coordinateOrientation codes.
Definition: itkOrientImageFilter.h:74
itk::AnatomicalOrientation::CoordinateEnum::LeftToRight
itkFlipImageFilter.h
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::OrientImageFilter::OutputImageConstPointer
typename OutputImageType::ConstPointer OutputImageConstPointer
Definition: itkOrientImageFilter.h:93
itk::AnatomicalOrientation::CoordinateEnum::PosteriorToAnterior
to front - 0b0100
itk::Concept::SameDimension
Definition: itkConceptChecking.h:697
itk::ImageToImageFilter::InputImagePixelType
typename InputImageType::PixelType InputImagePixelType
Definition: itkImageToImageFilter.h:133
itk::PermuteAxesImageFilter
Permutes the image axes according to a user specified order.
Definition: itkPermuteAxesImageFilter.h:51
itk::ImageToImageFilter
Base class for filters that take an image as input and produce an image as output.
Definition: itkImageToImageFilter.h:108
itk::ImageSource
Base class for all process objects that output image data.
Definition: itkImageSource.h:67
itk::ImageToImageFilter::InputImagePointer
typename InputImageType::Pointer InputImagePointer
Definition: itkImageToImageFilter.h:130
itk::OrientImageFilter::PermuteOrderArrayType
typename PermuterType::PermuteOrderArrayType PermuteOrderArrayType
Definition: itkOrientImageFilter.h:100
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itkAnatomicalOrientation.h
itk::ImageToImageFilter::InputImageType
TInputImage InputImageType
Definition: itkImageToImageFilter.h:129
itk::AnatomicalOrientation::CoordinateEnum::AnteriorToPosterior
to back
itk::OrientImageFilter::SetGivenCoordinateDirection
void SetGivenCoordinateDirection(const typename TInputImage::DirectionType &GivenDirection)
Definition: itkOrientImageFilter.h:123
itk::FixedArray< unsigned int, Self::ImageDimension >
itk::OrientImageFilter::SetDesiredCoordinateOrientationToAxial
void SetDesiredCoordinateOrientationToAxial()
Definition: itkOrientImageFilter.h:168
itk::ImageSource::OutputImageRegionType
typename OutputImageType::RegionType OutputImageRegionType
Definition: itkImageSource.h:92
itkConceptMacro
#define itkConceptMacro(name, concept)
Definition: itkConceptChecking.h:65
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnatomicalOrientation.h:29
itk::FlipImageFilter
Flips an image across user specified axes.
Definition: itkFlipImageFilter.h:53
itk::OrientImageFilter::FlipAxesArrayType
typename FlipperType::FlipAxesArrayType FlipAxesArrayType
Definition: itkOrientImageFilter.h:104
itk::ProcessObject
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Definition: itkProcessObject.h:139
itk::AnatomicalOrientation::CoordinateEnum::RightToLeft
0b0010
itk::ImageSource::OutputImagePixelType
typename OutputImageType::PixelType OutputImagePixelType
Definition: itkImageSource.h:93
itk::Concept::Convertible
Definition: itkConceptChecking.h:217
itkPermuteAxesImageFilter.h
itk::AnatomicalOrientation::CoordinateEnum::InferiorToSuperior
to head - 0b1000
itk::AnatomicalOrientation
Representations of anatomical orientations and methods to convert between conventions.
Definition: itkAnatomicalOrientation.h:53
itk::OrientImageFilter::SetDesiredCoordinateOrientationToSagittal
void SetDesiredCoordinateOrientationToSagittal()
Definition: itkOrientImageFilter.h:185
itk::ImageToImageFilter::InputImageRegionType
typename InputImageType::RegionType InputImageRegionType
Definition: itkImageToImageFilter.h:132
itk::OrientImageFilter::SetDesiredCoordinateDirection
void SetDesiredCoordinateDirection(const typename TOutputImage::DirectionType &DesiredDirection)
Definition: itkOrientImageFilter.h:133
itk::ImageToImageFilter::InputImageConstPointer
typename InputImageType::ConstPointer InputImageConstPointer
Definition: itkImageToImageFilter.h:131
itk::ImageSource::OutputImageType
TOutputImage OutputImageType
Definition: itkImageSource.h:90
itk::DataObject
Base class for all data objects in ITK.
Definition: itkDataObject.h:293