ITK  4.13.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:
53 
59 
61  itkNewMacro(Self);
62 
65 
68 
71 
75 
77  typedef TPixel OutputImagePixelType;
78 
80  TPixel * GetImportPointer();
81 
89  void SetImportPointer(TPixel *ptr, SizeValueType num,
90  bool LetImageContainerManageMemory);
91 
96  void SetRegion(const RegionType & region)
97  { if ( m_Region != region ) { m_Region = region; this->Modified(); } }
98 
103  const RegionType & GetRegion() const
104  { return m_Region; }
105 
108  itkSetMacro(Spacing, SpacingType);
109  itkGetConstReferenceMacro(Spacing, SpacingType);
110  itkSetVectorMacro(Spacing, const float, VImageDimension);
112 
115  itkSetMacro(Origin, OriginType);
116  itkGetConstReferenceMacro(Origin, OriginType);
117  itkSetVectorMacro(Origin, const float, VImageDimension);
119 
121 
124  virtual void SetDirection(const DirectionType & direction);
125 
128  itkGetConstReferenceMacro(Direction, DirectionType);
129 
130 protected:
132  ~ImportImageFilter() ITK_OVERRIDE;
133  void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
134 
137  virtual void GenerateData() ITK_OVERRIDE;
138 
142  virtual void GenerateOutputInformation() ITK_OVERRIDE;
143 
151  virtual void EnlargeOutputRequestedRegion(DataObject *output) ITK_OVERRIDE;
152 
153 private:
154  ITK_DISALLOW_COPY_AND_ASSIGN(ImportImageFilter);
155 
156  RegionType m_Region;
157  SpacingType m_Spacing;
158  OriginType m_Origin;
159  DirectionType m_Direction;
160 
161  typename ImportImageContainerType::Pointer m_ImportImageContainer;
163 };
164 } // end namespace itk
165 
166 #ifndef ITK_MANUAL_INSTANTIATION
167 #include "itkImportImageFilter.hxx"
168 #endif
169 
170 #endif
ImportImageContainer< SizeValueType, TPixel > ImportImageContainerType
SmartPointer< const Self > ConstPointer
Represent the size (bounds) of a n-dimensional image.
Definition: itkSize.h:52
Matrix< SpacePrecisionType, VImageDimension, VImageDimension > DirectionType
Index< VImageDimension > IndexType
An image region represents a structured region of data.
Base class for all process objects that output image data.
unsigned long SizeValueType
Definition: itkIntTypes.h:143
Import data from a standard C array into an itk::Image.
ImageRegion< VImageDimension > RegionType
ImageSource< OutputImageType > Superclass
OutputImageType::Pointer OutputImagePointer
const RegionType & GetRegion() const
Control indentation during Print() invocation.
Definition: itkIndent.h:49
void SetRegion(const RegionType &region)
OutputImageType::SpacingType SpacingType
OutputImageType::PointType OriginType
Image< TPixel, VImageDimension > OutputImageType
Size< VImageDimension > SizeType
Base class for all data objects in ITK.
Templated n-dimensional image class.
Definition: itkImage.h:75
SmartPointer< Self > Pointer
Defines an itk::Image front-end to a standard C-array.