ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
itkPointSetToImageFilter.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 itkPointSetToImageFilter_h
19#define itkPointSetToImageFilter_h
20
21#include "itkImageSource.h"
22#include "itkConceptChecking.h"
23
24namespace itk
25{
34template <typename TInputPointSet, typename TOutputImage>
35class ITK_TEMPLATE_EXPORT PointSetToImageFilter : public ImageSource<TOutputImage>
36{
37public:
38 ITK_DISALLOW_COPY_AND_MOVE(PointSetToImageFilter);
39
45 using SizeType = typename TOutputImage::SizeType;
46 using OutputImageType = TOutputImage;
47 using OutputImagePointer = typename OutputImageType::Pointer;
48 using ValueType = typename OutputImageType::ValueType;
49
51 itkNewMacro(Self);
52
54 itkOverrideGetNameOfClassMacro(PointSetToImageFilter);
55
58
60 using InputPointSetType = TInputPointSet;
61 using InputPointSetPointer = typename InputPointSetType::Pointer;
62 using InputPointSetConstPointer = typename InputPointSetType::ConstPointer;
63
65 static constexpr unsigned int InputPointSetDimension = InputPointSetType::PointDimension;
66 static constexpr unsigned int OutputImageDimension = TOutputImage::ImageDimension;
67
69 using SpacingType = typename TOutputImage::SpacingType;
70 using DirectionType = typename TOutputImage::DirectionType;
71 using PointType = typename TOutputImage::PointType;
72
75 virtual void
77
78 virtual void
79 SetInput(unsigned int, const InputPointSetType * pointset);
80
81 const InputPointSetType *
83
84 const InputPointSetType *
85 GetInput(unsigned int idx);
86
91 itkSetMacro(Spacing, SpacingType);
92 virtual void
93 SetSpacing(const double * v);
95
96 virtual void
97 SetSpacing(const float * v);
98
103 itkGetConstReferenceMacro(Spacing, SpacingType);
104
108 itkSetMacro(Direction, DirectionType);
109 itkGetConstReferenceMacro(Direction, DirectionType);
111
116 itkSetMacro(Origin, PointType);
117 virtual void
118 SetOrigin(const double * v);
120
121 virtual void
122 SetOrigin(const float * v);
123
128 itkGetConstReferenceMacro(Origin, PointType);
129
136 itkSetMacro(InsideValue, ValueType);
137 itkGetConstMacro(InsideValue, ValueType);
139
146 itkSetMacro(OutsideValue, ValueType);
147 itkGetConstMacro(OutsideValue, ValueType);
149
151 itkSetMacro(Size, SizeType);
152 itkGetConstMacro(Size, SizeType);
154
155protected:
157 ~PointSetToImageFilter() override = default;
158
159 void
161 {} // do nothing
162 void
163 GenerateData() override;
164
166
168
170
172
175
176 void
177 PrintSelf(std::ostream & os, Indent indent) const override;
178};
179} // end namespace itk
180
181#ifndef ITK_MANUAL_INSTANTIATION
182# include "itkPointSetToImageFilter.hxx"
183#endif
184
185#endif
typename OutputImageType::RegionType OutputImageRegionType
Control indentation during Print() invocation.
Definition itkIndent.h:50
virtual void SetSpacing(const double *v)
typename TOutputImage::SizeType SizeType
void PrintSelf(std::ostream &os, Indent indent) const override
virtual void SetInput(unsigned int, const InputPointSetType *pointset)
typename InputPointSetType::ConstPointer InputPointSetConstPointer
const InputPointSetType * GetInput()
virtual void SetSpacing(const float *v)
static constexpr unsigned int InputPointSetDimension
void GenerateData() override
SmartPointer< const Self > ConstPointer
static constexpr unsigned int OutputImageDimension
typename OutputImageType::ValueType ValueType
typename OutputImageType::Pointer OutputImagePointer
typename InputPointSetType::Pointer InputPointSetPointer
typename TOutputImage::PointType PointType
ImageSource< TOutputImage > Superclass
const InputPointSetType * GetInput(unsigned int idx)
virtual void SetInput(const InputPointSetType *input)
virtual void SetOrigin(const double *v)
typename TOutputImage::DirectionType DirectionType
virtual void SetOrigin(const float *v)
~PointSetToImageFilter() override=default
typename TOutputImage::SpacingType SpacingType
virtual void SetInput(const DataObjectIdentifierType &key, DataObject *input)
Protected method for setting indexed and named inputs.
Implements transparent reference counting.
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Represent a n-dimensional size (bounds) of a n-dimensional image.
Definition itkSize.h:70