ITK  4.9.0
Insight Segmentation and Registration Toolkit
itkPointSetToImageFilter.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 itkPointSetToImageFilter_h
19 #define itkPointSetToImageFilter_h
20 
21 #include "itkImageSource.h"
22 #include "itkConceptChecking.h"
23 
24 namespace itk
25 {
33 template< typename TInputPointSet, typename TOutputImage >
34 class PointSetToImageFilter:public ImageSource< TOutputImage >
35 {
36 public:
42  typedef typename TOutputImage::SizeType SizeType;
43  typedef TOutputImage OutputImageType;
44  typedef typename OutputImageType::Pointer OutputImagePointer;
45  typedef typename OutputImageType::ValueType ValueType;
46 
48  itkNewMacro(Self);
49 
52 
55 
57  typedef TInputPointSet InputPointSetType;
58  typedef typename InputPointSetType::Pointer InputPointSetPointer;
59  typedef typename InputPointSetType::ConstPointer InputPointSetConstPointer;
60 
62  itkStaticConstMacro(InputPointSetDimension, unsigned int,
63  InputPointSetType::PointDimension);
64  itkStaticConstMacro(OutputImageDimension, unsigned int,
65  TOutputImage::ImageDimension);
67 
69  typedef typename TOutputImage::SpacingType SpacingType;
70  typedef typename TOutputImage::DirectionType DirectionType;
71  typedef typename TOutputImage::PointType PointType;
72 
75  virtual void SetInput(const InputPointSetType *pointset);
76 
77  virtual void SetInput(unsigned int, const InputPointSetType *pointset);
78 
79  const InputPointSetType * GetInput();
80 
81  const InputPointSetType * GetInput(unsigned int idx);
82 
87  itkSetMacro(Spacing, SpacingType);
88  virtual void SetSpacing(const double *spacing);
90 
91  virtual void SetSpacing(const float *spacing);
92 
97  itkGetConstReferenceMacro(Spacing, SpacingType);
98 
102  itkSetMacro(Direction, DirectionType);
103  itkGetConstReferenceMacro(Direction, DirectionType);
105 
110  itkSetMacro(Origin, PointType);
111  virtual void SetOrigin(const double *origin);
113 
114  virtual void SetOrigin(const float *origin);
115 
120  itkGetConstReferenceMacro(Origin, PointType);
121 
128  itkSetMacro(InsideValue, ValueType);
129  itkGetConstMacro(InsideValue, ValueType);
131 
138  itkSetMacro(OutsideValue, ValueType);
139  itkGetConstMacro(OutsideValue, ValueType);
141 
143  itkSetMacro(Size, SizeType);
144  itkGetConstMacro(Size, SizeType);
146 
147 protected:
150 
151  virtual void GenerateOutputInformation() ITK_OVERRIDE {} // do nothing
152  virtual void GenerateData() ITK_OVERRIDE;
153 
155 
157 
159 
161 
164 
165  virtual void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
166 
167 private:
168  PointSetToImageFilter(const Self &) ITK_DELETE_FUNCTION;
169  void operator=(const Self &) ITK_DELETE_FUNCTION;
170 };
171 } // end namespace itk
172 
173 #ifndef ITK_MANUAL_INSTANTIATION
174 #include "itkPointSetToImageFilter.hxx"
175 #endif
176 
177 #endif
const InputPointSetType * GetInput()
OutputImageType::ValueType ValueType
virtual void SetOrigin(PointType _arg)
Represent the size (bounds) of a n-dimensional image.
Definition: itkSize.h:52
SmartPointer< const Self > ConstPointer
OutputImageType::Pointer OutputImagePointer
InputPointSetType::Pointer InputPointSetPointer
Base class for all process objects that output image data.
static const unsigned int InputPointSetDimension
virtual void SetInput(const InputPointSetType *pointset)
virtual void PrintSelf(std::ostream &os, Indent indent) const override
TOutputImage::DirectionType DirectionType
ImageSource< TOutputImage > Superclass
TOutputImage::SpacingType SpacingType
Superclass::OutputImageRegionType OutputImageRegionType
virtual void SetInput(const DataObjectIdentifierType &key, DataObject *input)
Protected method for setting indexed and named inputs.
static const unsigned int OutputImageDimension
Base class for filters that take a PointSet as input and produce an image as output. By default, if the user does not specify the size of the output image, the maximum size of the point-set&#39;s bounding box is used.
virtual void SetSpacing(SpacingType _arg)
virtual void GenerateOutputInformation() override
OutputImageType::RegionType OutputImageRegionType
Control indentation during Print() invocation.
Definition: itkIndent.h:49
TOutputImage::PointType PointType
virtual void GenerateData() override
InputPointSetType::ConstPointer InputPointSetConstPointer