ITK  5.2.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  * http://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 
23 namespace itk
24 {
42 template <typename TPixel, unsigned int VImageDimension = 2>
43 class ITK_TEMPLATE_EXPORT ImportImageFilter : public ImageSource<Image<TPixel, VImageDimension>>
44 {
45 public:
46  ITK_DISALLOW_COPY_AND_MOVE(ImportImageFilter);
47 
54 
60 
62  itkNewMacro(Self);
63 
65  itkTypeMacro(ImportImageFilter, ImageSource);
66 
69 
72 
76 
78  using OutputImagePixelType = TPixel;
79 
81  TPixel *
82  GetImportPointer();
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  }
108 
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);
125 
128  itkSetMacro(Origin, OriginType);
129  itkGetConstReferenceMacro(Origin, OriginType);
130  itkSetVectorMacro(Origin, const float, VImageDimension);
132 
134 
137  virtual void
138  SetDirection(const DirectionType & direction);
139 
142  itkGetConstReferenceMacro(Direction, DirectionType);
143 
144 protected:
146  ~ImportImageFilter() override = default;
147  void
148  PrintSelf(std::ostream & os, Indent indent) const override;
149 
152  void
153  GenerateData() override;
154 
158  void
159  GenerateOutputInformation() override;
160 
168  void
169  EnlargeOutputRequestedRegion(DataObject * output) override;
170 
171 private:
176 
179 };
180 } // end namespace itk
181 
182 #ifndef ITK_MANUAL_INSTANTIATION
183 # include "itkImportImageFilter.hxx"
184 #endif
185 
186 #endif
itk::ImportImageFilter::m_Region
RegionType m_Region
Definition: itkImportImageFilter.h:172
itk::ImportImageFilter::m_Origin
OriginType m_Origin
Definition: itkImportImageFilter.h:174
itk::Index< VImageDimension >
itk::ImportImageFilter::OutputImagePixelType
TPixel OutputImagePixelType
Definition: itkImportImageFilter.h:78
itkImageSource.h
itk::Size
Represent a n-dimensional size (bounds) of a n-dimensional image.
Definition: itkSize.h:69
itk::ImageRegion
An image region represents a structured region of data.
Definition: itkImageRegion.h:69
itk::GTest::TypedefsAndConstructors::Dimension2::PointType
ImageBaseType::PointType PointType
Definition: itkGTestTypedefsAndConstructors.h:51
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::ImportImageFilter::SetRegion
void SetRegion(const RegionType &region)
Definition: itkImportImageFilter.h:99
itk::ImageSource
Base class for all process objects that output image data.
Definition: itkImageSource.h:67
itk::ImportImageFilter::GetRegion
const RegionType & GetRegion() const
Definition: itkImportImageFilter.h:114
itk::ImportImageFilter::OutputImagePointer
typename OutputImageType::Pointer OutputImagePointer
Definition: itkImportImageFilter.h:50
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itk::ImportImageFilter::SpacingType
typename OutputImageType::SpacingType SpacingType
Definition: itkImportImageFilter.h:51
itk::ImportImageContainer
Defines an itk::Image front-end to a standard C-array.
Definition: itkImportImageContainer.h:45
itk::ImportImageFilter::OriginType
typename OutputImageType::PointType OriginType
Definition: itkImportImageFilter.h:52
itk::ImportImageFilter
Import data from a standard C array into an itk::Image.
Definition: itkImportImageFilter.h:43
itk::Matrix< SpacePrecisionType, VImageDimension, VImageDimension >
itk::ImportImageFilter::m_Spacing
SpacingType m_Spacing
Definition: itkImportImageFilter.h:173
itk::Image::SpacingType
typename Superclass::SpacingType SpacingType
Definition: itkImage.h:154
itk::ImportImageFilter::m_ImportImageContainer
ImportImageContainerType::Pointer m_ImportImageContainer
Definition: itkImportImageFilter.h:177
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itk::ImportImageFilter::m_Direction
DirectionType m_Direction
Definition: itkImportImageFilter.h:175
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:86
itk::ImportImageFilter::m_Size
SizeValueType m_Size
Definition: itkImportImageFilter.h:178
itk::SizeValueType
unsigned long SizeValueType
Definition: itkIntTypes.h:83
itk::DataObject
Base class for all data objects in ITK.
Definition: itkDataObject.h:293