ITK  6.0.0
Insight Toolkit
itkImportImageFilter.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 itkImportImageFilter_h
19#define itkImportImageFilter_h
20
21#include "itkImageSource.h"
22
23namespace itk
24{
42template <typename TPixel, unsigned int VImageDimension = 2>
43class ITK_TEMPLATE_EXPORT ImportImageFilter : public ImageSource<Image<TPixel, VImageDimension>>
44{
45public:
46 ITK_DISALLOW_COPY_AND_MOVE(ImportImageFilter);
47
54
60
62 itkNewMacro(Self);
63
65 itkOverrideGetNameOfClassMacro(ImportImageFilter);
66
69
72
76
78 using OutputImagePixelType = TPixel;
79
81 TPixel *
83
91 void
92 SetImportPointer(TPixel * ptr, SizeValueType num, bool LetImageContainerManageMemory);
93
98 void
99 SetRegion(const RegionType & region)
100 {
101 if (m_Region != region)
102 {
103 m_Region = region;
104 this->Modified();
105 }
106 }
113 const RegionType &
114 GetRegion() const
115 {
116 return m_Region;
117 }
118
121 itkSetMacro(Spacing, SpacingType);
122 itkGetConstReferenceMacro(Spacing, SpacingType);
123 itkSetVectorMacro(Spacing, const float, VImageDimension);
128 itkSetMacro(Origin, OriginType);
129 itkGetConstReferenceMacro(Origin, OriginType);
130 itkSetVectorMacro(Origin, const float, VImageDimension);
134
137 virtual void
138 SetDirection(const DirectionType & direction);
139
142 itkGetConstReferenceMacro(Direction, DirectionType);
143
144protected:
145 ImportImageFilter() = default;
146 ~ImportImageFilter() override = default;
147 void
148 PrintSelf(std::ostream & os, Indent indent) const override;
149
152 void
153 GenerateData() override;
154
158 void
160
168 void
170
171private:
172 RegionType m_Region{};
173 SpacingType m_Spacing{ MakeFilled<SpacingType>(1.0) };
174 OriginType m_Origin{};
175 DirectionType m_Direction{ DirectionType::GetIdentity() };
176
177 typename ImportImageContainerType::Pointer m_ImportImageContainer{};
179};
180} // end namespace itk
181
182#ifndef ITK_MANUAL_INSTANTIATION
183# include "itkImportImageFilter.hxx"
184#endif
185
186#endif
Base class for all data objects in ITK.
An image region represents a structured region of data.
Base class for all process objects that output image data.
Templated n-dimensional image class.
Definition: itkImage.h:89
Defines an itk::Image front-end to a standard C-array.
Import data from a standard C array into an itk::Image.
typename OutputImageType::Pointer OutputImagePointer
~ImportImageFilter() override=default
virtual void SetDirection(const DirectionType &direction)
const RegionType & GetRegion() const
void SetImportPointer(TPixel *ptr, SizeValueType num, bool LetImageContainerManageMemory)
void PrintSelf(std::ostream &os, Indent indent) const override
TPixel * GetImportPointer()
void GenerateOutputInformation() override
void EnlargeOutputRequestedRegion(DataObject *output) override
void GenerateData() override
void SetRegion(const RegionType &region)
typename OutputImageType::SpacingType SpacingType
typename OutputImageType::PointType OriginType
Control indentation during Print() invocation.
Definition: itkIndent.h:50
SmartPointer< Self > Pointer
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
unsigned long SizeValueType
Definition: itkIntTypes.h:86