ITK  5.2.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  * http://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 {
140 template <typename TInputImage, typename TOutputImage>
141 class ITK_TEMPLATE_EXPORT OrientImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
142 {
143 public:
144  ITK_DISALLOW_COPY_AND_MOVE(OrientImageFilter);
145 
151 
153  using InputImageType = TInputImage;
154  using InputImagePointer = typename InputImageType::Pointer;
155  using InputImageConstPointer = typename InputImageType::ConstPointer;
157  using InputImagePixelType = typename InputImageType::PixelType;
158  using OutputImageType = TOutputImage;
159  using OutputImagePointer = typename OutputImageType::Pointer;
160  using OutputImageConstPointer = typename OutputImageType::ConstPointer;
162  using OutputImagePixelType = typename OutputImageType::PixelType;
164 
168 
172 
174  static constexpr unsigned int InputImageDimension = TInputImage::ImageDimension;
175  static constexpr unsigned int OutputImageDimension = TOutputImage::ImageDimension;
176 
178  itkNewMacro(Self);
179 
181  itkTypeMacro(OrientImageFilter, ImageToImageFilter);
182 
184  itkGetEnumMacro(GivenCoordinateOrientation, CoordinateOrientationCode);
185  void
186  SetGivenCoordinateOrientation(CoordinateOrientationCode newCode);
188 
189  inline void
191  {
192  SetGivenCoordinateOrientation(itk::SpatialOrientationAdapter().FromDirectionCosines(GivenDirection));
193  }
194 
195  itkGetEnumMacro(DesiredCoordinateOrientation, CoordinateOrientationCode);
196  void
197  SetDesiredCoordinateOrientation(CoordinateOrientationCode newCode);
198 
199  inline void
201  {
202  SetDesiredCoordinateOrientation(itk::SpatialOrientationAdapter().FromDirectionCosines(DesiredDirection));
203  }
204 
211  itkBooleanMacro(UseImageDirection);
212  itkGetConstMacro(UseImageDirection, bool);
213  itkSetMacro(UseImageDirection, bool);
215 
217  itkGetConstReferenceMacro(PermuteOrder, PermuteOrderArrayType);
218 
220  itkGetConstReferenceMacro(FlipAxes, FlipAxesArrayType);
221 
235  void
237  {
238  this->SetDesiredCoordinateOrientation(SpatialOrientation::ITK_COORDINATE_ORIENTATION_RAI);
239  }
240 
241  void
243  {
244  this->SetDesiredCoordinateOrientation(SpatialOrientation::ITK_COORDINATE_ORIENTATION_RSA);
245  }
246 
247  void
249  {
250  this->SetDesiredCoordinateOrientation(SpatialOrientation::ITK_COORDINATE_ORIENTATION_ASL);
251  }
252 
260  void
261  GenerateOutputInformation() override;
262 
263 #ifdef ITK_USE_CONCEPT_CHECKING
264  // Begin concept checking
268  // End concept checking
269 #endif
270 
271 protected:
273  ~OrientImageFilter() override = default;
274  void
275  PrintSelf(std::ostream & os, Indent indent) const override;
276 
280  void
281  GenerateInputRequestedRegion() override;
282 
284  void
285  EnlargeOutputRequestedRegion(DataObject * itkNotUsed(output)) override;
286 
287 /*** Member functions used by GenerateData: */
288  void
289  DeterminePermutationsAndFlips(const SpatialOrientation::ValidCoordinateOrientationFlags fixed_orient,
291 
293  bool
294  NeedToPermute();
295 
297  bool
298  NeedToFlip();
299 
302  void
303  GenerateData() override;
304 
305 private:
306  std::string
307  GetMajorAxisFromPatientRelativeDirectionCosine(double x, double y, double z);
308 
311  bool m_UseImageDirection{ false };
312 
315 
316  std::map<std::string, CoordinateOrientationCode> m_StringToCode;
317  std::map<CoordinateOrientationCode, std::string> m_CodeToString;
318 }; // end of class
319 } // end namespace itk
320 
321 #ifndef ITK_MANUAL_INSTANTIATION
322 # include "itkOrientImageFilter.hxx"
323 #endif
324 
325 #endif
itkSpatialOrientationAdapter.h
itk::OrientImageFilter::SetDesiredCoordinateOrientationToCoronal
void SetDesiredCoordinateOrientationToCoronal()
Definition: itkOrientImageFilter.h:242
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:141
itk::OrientImageFilter::m_StringToCode
std::map< std::string, CoordinateOrientationCode > m_StringToCode
Definition: itkOrientImageFilter.h:316
itkFlipImageFilter.h
itk::OrientImageFilter::m_CodeToString
std::map< CoordinateOrientationCode, std::string > m_CodeToString
Definition: itkOrientImageFilter.h:317
itk::OrientImageFilter::m_PermuteOrder
PermuteOrderArrayType m_PermuteOrder
Definition: itkOrientImageFilter.h:313
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::OrientImageFilter::OutputImageConstPointer
typename OutputImageType::ConstPointer OutputImageConstPointer
Definition: itkOrientImageFilter.h:160
itk::SpatialOrientationAdapter
Converts SpatialOrientation flags 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::OrientImageFilter::PermuteOrderArrayType
typename PermuterType::PermuteOrderArrayType PermuteOrderArrayType
Definition: itkOrientImageFilter.h:167
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itk::SpatialOrientation::ITK_COORDINATE_ORIENTATION_RSA
Definition: itkSpatialOrientation.h:105
itk::ImageToImageFilter::InputImageType
TInputImage InputImageType
Definition: itkImageToImageFilter.h:129
itk::OrientImageFilter::SetGivenCoordinateDirection
void SetGivenCoordinateDirection(const typename TInputImage::DirectionType &GivenDirection)
Definition: itkOrientImageFilter.h:190
itk::FixedArray< unsigned int, Self::ImageDimension >
itk::OrientImageFilter::SetDesiredCoordinateOrientationToAxial
void SetDesiredCoordinateOrientationToAxial()
Definition: itkOrientImageFilter.h:236
itk::SpatialOrientation::ITK_COORDINATE_ORIENTATION_RIP
Definition: itkSpatialOrientation.h:87
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: 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:171
itk::SpatialOrientation::ValidCoordinateOrientationFlags
ValidCoordinateOrientationFlags
Definition: itkSpatialOrientation.h:84
itk::ProcessObject
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Definition: itkProcessObject.h:138
itk::SpatialOrientation::ITK_COORDINATE_ORIENTATION_ASL
Definition: itkSpatialOrientation.h:233
itk::ImageSource::OutputImagePixelType
typename OutputImageType::PixelType OutputImagePixelType
Definition: itkImageSource.h:93
itk::Concept::Convertible
Definition: itkConceptChecking.h:216
itk::OrientImageFilter::m_FlipAxes
FlipAxesArrayType m_FlipAxes
Definition: itkOrientImageFilter.h:314
itkPermuteAxesImageFilter.h
itk::SpatialOrientation::ITK_COORDINATE_ORIENTATION_RAI
Definition: itkSpatialOrientation.h:143
itk::OrientImageFilter::SetDesiredCoordinateOrientationToSagittal
void SetDesiredCoordinateOrientationToSagittal()
Definition: itkOrientImageFilter.h:248
itk::ImageToImageFilter::InputImageRegionType
typename InputImageType::RegionType InputImageRegionType
Definition: itkImageToImageFilter.h:132
itk::OrientImageFilter::SetDesiredCoordinateDirection
void SetDesiredCoordinateDirection(const typename TOutputImage::DirectionType &DesiredDirection)
Definition: itkOrientImageFilter.h:200
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