ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
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 }
107
112 const RegionType &
113 GetRegion() const
114 {
115 return m_Region;
116 }
117
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);
133
136 virtual void
137 SetDirection(const DirectionType & direction);
138
141 itkGetConstReferenceMacro(Direction, DirectionType);
142
143protected:
144 ImportImageFilter() = default;
145 ~ImportImageFilter() override = default;
146 void
147 PrintSelf(std::ostream & os, Indent indent) const override;
148
151 void
152 GenerateData() override;
153
157 void
159
167 void
169
170private:
175
178};
179} // end namespace itk
180
181#ifndef ITK_MANUAL_INSTANTIATION
182# include "itkImportImageFilter.hxx"
183#endif
184
185#endif
Base class for all data objects in ITK.
An image region represents a structured region of data.
Templated n-dimensional image class.
Definition itkImage.h:89
Point< PointValueType, VImageDimension > PointType
Vector< SpacingValueType, VImageDimension > SpacingType
Defines an itk::Image front-end to a standard C-array.
typename OutputImageType::Pointer OutputImagePointer
ImportImageContainerType::Pointer m_ImportImageContainer
~ImportImageFilter() override=default
Index< VImageDimension > IndexType
virtual void SetDirection(const DirectionType &direction)
const RegionType & GetRegion() const
Size< VImageDimension > SizeType
void SetImportPointer(TPixel *ptr, SizeValueType num, bool LetImageContainerManageMemory)
Matrix< SpacePrecisionType, VImageDimension, VImageDimension > DirectionType
SmartPointer< const Self > ConstPointer
void PrintSelf(std::ostream &os, Indent indent) const override
TPixel * GetImportPointer()
void GenerateOutputInformation() override
ImportImageContainer< SizeValueType, TPixel > ImportImageContainerType
void EnlargeOutputRequestedRegion(DataObject *output) override
ImageSource< OutputImageType > Superclass
void GenerateData() override
ImageRegion< VImageDimension > RegionType
Image< TPixel, VImageDimension > OutputImageType
void SetRegion(const RegionType &region)
SmartPointer< Self > Pointer
typename OutputImageType::SpacingType SpacingType
typename OutputImageType::PointType OriginType
Control indentation during Print() invocation.
Definition itkIndent.h:50
A templated class holding a M x N size Matrix.
Definition itkMatrix.h:53
virtual void Modified() const
Implements transparent reference counting.
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
constexpr TContainer MakeFilled(typename TContainer::const_reference value)
unsigned long SizeValueType
Definition itkIntTypes.h:86
Represent a n-dimensional index in a n-dimensional image.
Definition itkIndex.h:69
Represent a n-dimensional size (bounds) of a n-dimensional image.
Definition itkSize.h:70