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
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;
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
122 inline void
124 {
125 SetGivenCoordinateOrientation(AnatomicalOrientation(GivenDirection));
126 }
127
128 itkGetEnumMacro(DesiredCoordinateOrientation, CoordinateOrientationCode);
129 void
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
201
202#ifdef ITK_USE_CONCEPT_CHECKING
203 // Begin concept checking
207 // End concept checking
208#endif
209
210protected:
212 ~OrientImageFilter() override = default;
213 void
214 PrintSelf(std::ostream & os, Indent indent) const override;
215
219 void
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
232 const CoordinateOrientationCode moving_orient);
233
235 bool
237
239 bool
241
244 void
245 GenerateData() override;
246
247private:
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
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.
Base class for all process objects that output image data.
typename OutputImageType::PixelType OutputImagePixelType
typename OutputImageType::RegionType OutputImageRegionType
TOutputImage OutputImageType
typename OutputImageType::Pointer OutputImagePointer
Base class for filters that take an image as input and produce an image as output.
typename InputImageType::ConstPointer InputImageConstPointer
typename InputImageType::Pointer InputImagePointer
typename InputImageType::PixelType InputImagePixelType
typename InputImageType::RegionType InputImageRegionType
Control indentation during Print() invocation.
Definition: itkIndent.h:50
Permute axes and then flip images as needed to obtain agreement in coordinateOrientation codes.
void VerifyPreconditions() const override
Verifies that the process object has been configured correctly, that all required inputs are set,...
void GenerateInputRequestedRegion() override
void EnlargeOutputRequestedRegion(DataObject *) override
void DeterminePermutationsAndFlips(const CoordinateOrientationCode fixed_orient, const CoordinateOrientationCode moving_orient)
void SetGivenCoordinateOrientation(CoordinateOrientationCode newCode)
typename FlipperType::FlipAxesArrayType FlipAxesArrayType
void SetDesiredCoordinateOrientation(CoordinateOrientationCode newCode)
void GenerateOutputInformation() override
void SetDesiredCoordinateDirection(const typename TOutputImage::DirectionType &DesiredDirection)
void PrintSelf(std::ostream &os, Indent indent) const override
typename OutputImageType::ConstPointer OutputImageConstPointer
void SetGivenCoordinateDirection(const typename TInputImage::DirectionType &GivenDirection)
typename PermuterType::PermuteOrderArrayType PermuteOrderArrayType
void GenerateData() override
~OrientImageFilter() override=default
Permutes the image axes according to a user specified order.
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
SmartPointer< const Self > ConstPointer
SmartPointer< Self > Pointer
#define itkConceptMacro(name, concept)
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....