ITK  5.4.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 {
141 template <typename TInputImage, typename TOutputImage>
142 class ITK_TEMPLATE_EXPORT OrientImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
143 {
144 public:
145  ITK_DISALLOW_COPY_AND_MOVE(OrientImageFilter);
146 
152 
154  using InputImageType = TInputImage;
158  using InputImagePixelType = typename InputImageType::PixelType;
159  using OutputImageType = TOutputImage;
163  using OutputImagePixelType = typename OutputImageType::PixelType;
165 
169 
173 
175  static constexpr unsigned int InputImageDimension = TInputImage::ImageDimension;
176  static constexpr unsigned int OutputImageDimension = TOutputImage::ImageDimension;
177 
179  itkNewMacro(Self);
180 
182  itkOverrideGetNameOfClassMacro(OrientImageFilter);
183 
185  itkGetEnumMacro(GivenCoordinateOrientation, CoordinateOrientationCode);
186  void
187  SetGivenCoordinateOrientation(CoordinateOrientationCode newCode);
190  inline void
192  {
193  SetGivenCoordinateOrientation(itk::SpatialOrientationAdapter().FromDirectionCosines(GivenDirection));
194  }
195 
196  itkGetEnumMacro(DesiredCoordinateOrientation, CoordinateOrientationCode);
197  void
198  SetDesiredCoordinateOrientation(CoordinateOrientationCode newCode);
199 
200  inline void
202  {
203  SetDesiredCoordinateOrientation(itk::SpatialOrientationAdapter().FromDirectionCosines(DesiredDirection));
204  }
205 
212  itkBooleanMacro(UseImageDirection);
213  itkGetConstMacro(UseImageDirection, bool);
214  itkSetMacro(UseImageDirection, bool);
218  itkGetConstReferenceMacro(PermuteOrder, PermuteOrderArrayType);
219 
221  itkGetConstReferenceMacro(FlipAxes, FlipAxesArrayType);
222 
236  void
238  {
239  this->SetDesiredCoordinateOrientation(
241  }
242 
243  void
245  {
246  this->SetDesiredCoordinateOrientation(
248  }
249 
250  void
252  {
253  this->SetDesiredCoordinateOrientation(
255  }
256 
264  void
265  GenerateOutputInformation() override;
266 
267 #ifdef ITK_USE_CONCEPT_CHECKING
268  // Begin concept checking
272  // End concept checking
273 #endif
274 
275 protected:
277  ~OrientImageFilter() override = default;
278  void
279  PrintSelf(std::ostream & os, Indent indent) const override;
280 
284  void
285  GenerateInputRequestedRegion() override;
286 
288  void
289  EnlargeOutputRequestedRegion(DataObject * itkNotUsed(output)) override;
290 
291  /*** Member functions used by GenerateData: */
292  void
293  DeterminePermutationsAndFlips(const SpatialOrientationEnums::ValidCoordinateOrientations fixed_orient,
295 
297  bool
298  NeedToPermute();
299 
301  bool
302  NeedToFlip();
303 
306  void
307  GenerateData() override;
308 
309 private:
310  std::string
311  GetMajorAxisFromPatientRelativeDirectionCosine(double x, double y, double z);
312 
313  CoordinateOrientationCode m_GivenCoordinateOrientation{
315  };
316  CoordinateOrientationCode m_DesiredCoordinateOrientation{
318  };
319  bool m_UseImageDirection{ false };
320 
321  PermuteOrderArrayType m_PermuteOrder{};
322  FlipAxesArrayType m_FlipAxes{};
323 
324  std::map<std::string, CoordinateOrientationCode> m_StringToCode{};
325  std::map<CoordinateOrientationCode, std::string> m_CodeToString{};
326 }; // end of class
327 } // end namespace itk
328 
329 #ifndef ITK_MANUAL_INSTANTIATION
330 # include "itkOrientImageFilter.hxx"
331 #endif
332 
333 #endif
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
ConstPointer
SmartPointer< const Self > ConstPointer
Definition: itkAddImageFilter.h:94
itkSpatialOrientationAdapter.h
itk::OrientImageFilter::SetDesiredCoordinateOrientationToCoronal
void SetDesiredCoordinateOrientationToCoronal()
Definition: itkOrientImageFilter.h:244
itk::GTest::TypedefsAndConstructors::Dimension2::DirectionType
ImageBaseType::DirectionType DirectionType
Definition: itkGTestTypedefsAndConstructors.h:52
itk::ImageSource::OutputImagePointer
typename OutputImageType::Pointer OutputImagePointer
Definition: itkImageSource.h:91
itk::OrientImageFilter
Permute axes and then flip images as needed to obtain agreement in coordinateOrientation codes.
Definition: itkOrientImageFilter.h:142
itk::SpatialOrientationEnums::ValidCoordinateOrientations::ITK_COORDINATE_ORIENTATION_RIP
itkFlipImageFilter.h
itk::SpatialOrientationEnums::ValidCoordinateOrientations
ValidCoordinateOrientations
Definition: itkSpatialOrientation.h:111
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::OrientImageFilter::OutputImageConstPointer
typename OutputImageType::ConstPointer OutputImageConstPointer
Definition: itkOrientImageFilter.h:161
itk::SpatialOrientationAdapter
Converts SpatialOrientationEnums to/from direction cosines.
Definition: itkSpatialOrientationAdapter.h:81
itk::Concept::SameDimension
Definition: itkConceptChecking.h:694
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::SpatialOrientationEnums::ValidCoordinateOrientations::ITK_COORDINATE_ORIENTATION_ASL
itk::OrientImageFilter::PermuteOrderArrayType
typename PermuterType::PermuteOrderArrayType PermuteOrderArrayType
Definition: itkOrientImageFilter.h:168
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itk::ImageToImageFilter::InputImageType
TInputImage InputImageType
Definition: itkImageToImageFilter.h:129
itk::OrientImageFilter::SetGivenCoordinateDirection
void SetGivenCoordinateDirection(const typename TInputImage::DirectionType &GivenDirection)
Definition: itkOrientImageFilter.h:191
itk::FixedArray< unsigned int, Self::ImageDimension >
itk::OrientImageFilter::SetDesiredCoordinateOrientationToAxial
void SetDesiredCoordinateOrientationToAxial()
Definition: itkOrientImageFilter.h:237
itk::ImageSource::OutputImageRegionType
typename OutputImageType::RegionType OutputImageRegionType
Definition: itkImageSource.h:92
itk::SpatialOrientationEnums::ValidCoordinateOrientations::ITK_COORDINATE_ORIENTATION_RAI
itkConceptMacro
#define itkConceptMacro(name, concept)
Definition: itkConceptChecking.h:65
itk::SpatialOrientationEnums::ValidCoordinateOrientations::ITK_COORDINATE_ORIENTATION_RSA
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itk::FlipImageFilter
Flips an image across user specified axes.
Definition: itkFlipImageFilter.h:53
itk::OrientImageFilter::FlipAxesArrayType
typename FlipperType::FlipAxesArrayType FlipAxesArrayType
Definition: itkOrientImageFilter.h:172
itk::ProcessObject
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Definition: itkProcessObject.h:139
itk::ImageSource::OutputImagePixelType
typename OutputImageType::PixelType OutputImagePixelType
Definition: itkImageSource.h:93
itk::Concept::Convertible
Definition: itkConceptChecking.h:216
itkPermuteAxesImageFilter.h
itk::OrientImageFilter::SetDesiredCoordinateOrientationToSagittal
void SetDesiredCoordinateOrientationToSagittal()
Definition: itkOrientImageFilter.h:251
itk::ImageToImageFilter::InputImageRegionType
typename InputImageType::RegionType InputImageRegionType
Definition: itkImageToImageFilter.h:132
itk::OrientImageFilter::SetDesiredCoordinateDirection
void SetDesiredCoordinateDirection(const typename TOutputImage::DirectionType &DesiredDirection)
Definition: itkOrientImageFilter.h:201
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