ITK  6.0.0
Insight Toolkit
itkChangeInformationImageFilter.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 itkChangeInformationImageFilter_h
19#define itkChangeInformationImageFilter_h
20
22#include "itkFixedArray.h"
23
24namespace itk
25{
48template <typename TInputImage>
49class ITK_TEMPLATE_EXPORT ChangeInformationImageFilter : public ImageToImageFilter<TInputImage, TInputImage>
50{
51public:
52 ITK_DISALLOW_COPY_AND_MOVE(ChangeInformationImageFilter);
53
59
61 itkNewMacro(Self);
62
63 using InputImageType = TInputImage;
64 using OutputImageType = TInputImage;
65
69
73
75 using OutputImagePixelType = typename OutputImageType::PixelType;
76 using InputImagePixelType = typename InputImageType::PixelType;
77
81 using OutputImageOffsetType = typename OutputImageType::OffsetType;
86 using InputImageOffsetType = typename InputImageType::OffsetType;
88
90 static constexpr unsigned int ImageDimension = InputImageType::ImageDimension;
91
93 using SpacingType = typename InputImageType::SpacingType;
96
98 itkOverrideGetNameOfClassMacro(ChangeInformationImageFilter);
99
102 void
104 {
105 if (image != m_ReferenceImage)
106 {
107 m_ReferenceImage = image;
108 this->ProcessObject::SetNthInput(1, const_cast<InputImageType *>(image));
109 this->Modified();
110 }
111 }
114 itkGetModifiableObjectMacro(ReferenceImage, TInputImage);
115
116 itkSetMacro(UseReferenceImage, bool);
117 itkBooleanMacro(UseReferenceImage);
118 itkGetConstMacro(UseReferenceImage, bool);
119
123 itkSetMacro(OutputSpacing, SpacingType);
124 itkGetConstReferenceMacro(OutputSpacing, SpacingType);
130 itkSetMacro(OutputOrigin, PointType);
131 itkGetConstReferenceMacro(OutputOrigin, PointType);
137 itkSetMacro(OutputDirection, DirectionType);
138 itkGetConstReferenceMacro(OutputDirection, DirectionType);
149 itkSetMacro(OutputOffset, OutputImageOffsetType);
150 itkGetConstReferenceMacro(OutputOffset, OutputImageOffsetType);
151 itkSetVectorMacro(OutputOffset, OutputImageOffsetValueType, ImageDimension);
155 void
157 {
158 this->ChangeSpacingOn();
159 this->ChangeOriginOn();
160 this->ChangeDirectionOn();
161 this->ChangeRegionOn();
162 }
167 void
169 {
170 this->ChangeSpacingOff();
171 this->ChangeOriginOff();
172 this->ChangeDirectionOff();
173 this->ChangeRegionOff();
174 }
183 itkSetMacro(ChangeSpacing, bool);
184 itkBooleanMacro(ChangeSpacing);
185 itkGetConstMacro(ChangeSpacing, bool);
186
193 itkSetMacro(ChangeOrigin, bool);
194 itkBooleanMacro(ChangeOrigin);
195 itkGetConstMacro(ChangeOrigin, bool);
196
203 itkSetMacro(ChangeDirection, bool);
204 itkBooleanMacro(ChangeDirection);
205 itkGetConstMacro(ChangeDirection, bool);
206
209 itkSetMacro(ChangeRegion, bool);
210 itkBooleanMacro(ChangeRegion);
211 itkGetConstMacro(ChangeRegion, bool);
212
216 itkSetMacro(CenterImage, bool);
217 itkBooleanMacro(CenterImage);
218 itkGetConstMacro(CenterImage, bool);
222 void
224
226 void
228
230 void
231 GenerateData() override;
232
233protected:
235 ~ChangeInformationImageFilter() override = default;
236
237 void
238 PrintSelf(std::ostream & os, Indent indent) const override;
239
245 void
246 VerifyInputInformation() const override
247 {}
248
249private:
250 InputImagePointer m_ReferenceImage{ nullptr };
251
252 bool m_CenterImage{ false };
253 bool m_ChangeSpacing{ false };
254 bool m_ChangeOrigin{ false };
255 bool m_ChangeDirection{ false };
256 bool m_ChangeRegion{ false };
257 bool m_UseReferenceImage{ false };
258
259 SpacingType m_OutputSpacing{ MakeFilled<SpacingType>(1.0) };
260 PointType m_OutputOrigin{};
261 DirectionType m_OutputDirection{ DirectionType::GetIdentity() };
262
263 OutputImageOffsetType m_OutputOffset{};
265};
266} // end namespace itk
267
268#ifndef ITK_MANUAL_INSTANTIATION
269# include "itkChangeInformationImageFilter.hxx"
270#endif
271
272#endif
Change the origin, spacing and/or region of an Image.
typename InputImageType::DirectionType InputImageDirectionType
typename OutputImageType::OffsetType OutputImageOffsetType
typename InputImageType::PixelType InputImagePixelType
typename OutputImageType::RegionType OutputImageRegionType
typename OutputImageType::OffsetValueType OutputImageOffsetValueType
void PrintSelf(std::ostream &os, Indent indent) const override
typename InputImageType::SizeType InputImageSizeType
void GenerateInputRequestedRegion() override
typename InputImageType::IndexType InputImageIndexType
void GenerateOutputInformation() override
typename InputImageType::ConstPointer InputImageConstPointer
typename InputImageType::RegionType InputImageRegionType
typename OutputImageType::PixelType OutputImagePixelType
typename OutputImageType::SizeType OutputImageSizeType
typename InputImageType::PointType PointType
typename OutputImageType::IndexType OutputImageIndexType
typename InputImageType::Pointer InputImagePointer
typename InputImageType::OffsetType InputImageOffsetType
typename OutputImageType::DirectionType OutputImageDirectionType
~ChangeInformationImageFilter() override=default
typename InputImageType::DirectionType DirectionType
typename InputImageType::SpacingType SpacingType
Base class for filters that take an image as input and produce an image as output.
Control indentation during Print() invocation.
Definition: itkIndent.h:50
Light weight base class for most itk classes.
virtual void SetNthInput(DataObjectPointerArraySizeType idx, DataObject *input)
SmartPointer< const Self > ConstPointer
SmartPointer< Self > Pointer
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
long OffsetValueType
Definition: itkIntTypes.h:97