ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
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
27namespace itk
28{
73template <typename TInputImage, typename TOutputImage>
74class ITK_TEMPLATE_EXPORT OrientImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
75{
76public:
77 ITK_DISALLOW_COPY_AND_MOVE(OrientImageFilter);
78
84
86 using InputImageType = TInputImage;
87 using InputImagePointer = typename InputImageType::Pointer;
88 using InputImageConstPointer = typename InputImageType::ConstPointer;
89 using InputImageRegionType = typename InputImageType::RegionType;
90 using InputImagePixelType = typename InputImageType::PixelType;
91 using OutputImageType = TOutputImage;
92 using OutputImagePointer = typename OutputImageType::Pointer;
93 using OutputImageConstPointer = typename OutputImageType::ConstPointer;
94 using OutputImageRegionType = typename OutputImageType::RegionType;
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
121
122 inline void
123 SetGivenCoordinateDirection(const typename TInputImage::DirectionType & GivenDirection)
124 {
126 }
127
128 itkGetEnumMacro(DesiredCoordinateOrientation, CoordinateOrientationCode);
129 void
131
132 inline void
133 SetDesiredCoordinateDirection(const typename TOutputImage::DirectionType & DesiredDirection)
134 {
136 }
137
146 itkBooleanMacro(UseImageDirection);
147 itkGetConstMacro(UseImageDirection, bool);
148 itkSetMacro(UseImageDirection, bool);
150
152 itkGetConstReferenceMacro(PermuteOrder, PermuteOrderArrayType);
153
155 itkGetConstReferenceMacro(FlipAxes, FlipAxesArrayType);
156
167 void
174
175
176 void
183
184 void
191
199 void
201
205
206protected:
208 ~OrientImageFilter() override = default;
209 void
210 PrintSelf(std::ostream & os, Indent indent) const override;
211
215 void
217
219 void
220 EnlargeOutputRequestedRegion(DataObject * itkNotUsed(output)) override;
221
222 void
223 VerifyPreconditions() const override;
224
225 /*** Member functions used by GenerateData: */
226 void
228 const CoordinateOrientationCode moving_orient);
229
231 bool
233
235 bool
237
240 void
241 GenerateData() override;
242
243private:
252 bool m_UseImageDirection{ false };
253
256
257}; // end of class
258} // end namespace itk
259
260#ifndef ITK_MANUAL_INSTANTIATION
261# include "itkOrientImageFilter.hxx"
262#endif
263
264#endif
Representations of anatomical orientations and methods to convert between conventions.
Base class for all data objects in ITK.
Flips an image across user specified axes.
FixedArray< bool, Self::ImageDimension > FlipAxesArrayType
Control indentation during Print() invocation.
Definition itkIndent.h:50
void VerifyPreconditions() const override
Verifies that the process object has been configured correctly, that all required inputs are set,...
static constexpr unsigned int InputImageDimension
void GenerateInputRequestedRegion() override
typename InputImageType::PixelType InputImagePixelType
typename InputImageType::ConstPointer InputImageConstPointer
void DeterminePermutationsAndFlips(const CoordinateOrientationCode fixed_orient, const CoordinateOrientationCode moving_orient)
void SetGivenCoordinateOrientation(CoordinateOrientationCode newCode)
PermuteOrderArrayType m_PermuteOrder
void EnlargeOutputRequestedRegion(DataObject *output) override
typename FlipperType::FlipAxesArrayType FlipAxesArrayType
typename InputImageType::Pointer InputImagePointer
void SetDesiredCoordinateOrientation(CoordinateOrientationCode newCode)
void GenerateOutputInformation() override
PermuteAxesImageFilter< TInputImage > PermuterType
typename OutputImageType::PixelType OutputImagePixelType
void SetDesiredCoordinateDirection(const typename TOutputImage::DirectionType &DesiredDirection)
CoordinateOrientationCode m_GivenCoordinateOrientation
typename InputImageType::RegionType InputImageRegionType
typename OutputImageType::Pointer OutputImagePointer
void PrintSelf(std::ostream &os, Indent indent) const override
ImageToImageFilter< TInputImage, TOutputImage > Superclass
SmartPointer< const Self > ConstPointer
typename OutputImageType::ConstPointer OutputImageConstPointer
void SetGivenCoordinateDirection(const typename TInputImage::DirectionType &GivenDirection)
typename PermuterType::PermuteOrderArrayType PermuteOrderArrayType
AnatomicalOrientation CoordinateOrientationCode
FlipImageFilter< TInputImage > FlipperType
static constexpr unsigned int OutputImageDimension
CoordinateOrientationCode m_DesiredCoordinateOrientation
SmartPointer< Self > Pointer
typename OutputImageType::RegionType OutputImageRegionType
void GenerateData() override
~OrientImageFilter() override=default
Permutes the image axes according to a user specified order.
FixedArray< unsigned int, Self::ImageDimension > PermuteOrderArrayType
Implements transparent reference counting.
#define itkConceptMacro(name, concept)
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....