ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkImportImageFilter.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright Insight Software Consortium
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:
44  public ImageSource< Image< TPixel, VImageDimension > >
45 {
46 public:
47  ITK_DISALLOW_COPY_AND_ASSIGN(ImportImageFilter);
48 
55 
61 
63  itkNewMacro(Self);
64 
66  itkTypeMacro(ImportImageFilter, ImageSource);
67 
70 
73 
77 
79  using OutputImagePixelType = TPixel;
80 
82  TPixel * GetImportPointer();
83 
91  void SetImportPointer(TPixel *ptr, SizeValueType num,
92  bool LetImageContainerManageMemory);
93 
98  void SetRegion(const RegionType & region)
99  { if ( m_Region != region ) { m_Region = region; this->Modified(); } }
100 
105  const RegionType & GetRegion() const
106  { return m_Region; }
107 
110  itkSetMacro(Spacing, SpacingType);
111  itkGetConstReferenceMacro(Spacing, SpacingType);
112  itkSetVectorMacro(Spacing, const float, VImageDimension);
114 
117  itkSetMacro(Origin, OriginType);
118  itkGetConstReferenceMacro(Origin, OriginType);
119  itkSetVectorMacro(Origin, const float, VImageDimension);
121 
123 
126  virtual void SetDirection(const DirectionType & direction);
127 
130  itkGetConstReferenceMacro(Direction, DirectionType);
131 
132 protected:
134  ~ImportImageFilter() override = default;
135  void PrintSelf(std::ostream & os, Indent indent) const override;
136 
139  void GenerateData() override;
140 
144  void GenerateOutputInformation() override;
145 
153  void EnlargeOutputRequestedRegion(DataObject *output) override;
154 
155 private:
160 
163 };
164 } // end namespace itk
165 
166 #ifndef ITK_MANUAL_INSTANTIATION
167 #include "itkImportImageFilter.hxx"
168 #endif
169 
170 #endif
typename OutputImageType::SpacingType SpacingType
typename OutputImageType::PointType OriginType
unsigned long SizeValueType
Definition: itkIntTypes.h:83
An image region represents a structured region of data.
Base class for all process objects that output image data.
typename OutputImageType::Pointer OutputImagePointer
Import data from a standard C array into an itk::Image.
Represent a n-dimensional size (bounds) of a n-dimensional image.
Definition: itkSize.h:68
typename Superclass::SpacingType SpacingType
Definition: itkImage.h:143
const RegionType & GetRegion() const
ImportImageContainerType::Pointer m_ImportImageContainer
Control indentation during Print() invocation.
Definition: itkIndent.h:49
void SetRegion(const RegionType &region)
Base class for all data objects in ITK.
Templated n-dimensional image class.
Definition: itkImage.h:75
Defines an itk::Image front-end to a standard C-array.