ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
itkSpatialObjectToImageFilter.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 itkSpatialObjectToImageFilter_h
19#define itkSpatialObjectToImageFilter_h
20
21#include "itkImageSource.h"
22#include "itkConceptChecking.h"
23
24namespace itk
25{
40template <typename TInputSpatialObject, typename TOutputImage>
41class ITK_TEMPLATE_EXPORT SpatialObjectToImageFilter : public ImageSource<TOutputImage>
42{
43public:
44 ITK_DISALLOW_COPY_AND_MOVE(SpatialObjectToImageFilter);
45
51
52 using OutputImageType = TOutputImage;
53 using IndexType = typename OutputImageType::IndexType;
54 using SizeType = typename OutputImageType::SizeType;
55 using PointType = typename OutputImageType::PointType;
56 using OutputImagePointer = typename OutputImageType::Pointer;
57 using ValueType = typename OutputImageType::ValueType;
58 using SpacingType = typename OutputImageType::SpacingType;
59 using DirectionType = typename OutputImageType::DirectionType;
60
62 itkNewMacro(Self);
63
65 itkOverrideGetNameOfClassMacro(SpatialObjectToImageFilter);
66
69
71 using InputSpatialObjectType = TInputSpatialObject;
72 using InputSpatialObjectPointer = typename InputSpatialObjectType::Pointer;
73 using InputSpatialObjectConstPointer = typename InputSpatialObjectType::ConstPointer;
74 using ChildrenListType = typename TInputSpatialObject::ChildrenListType;
75
77 static constexpr unsigned int ObjectDimension = InputSpatialObjectType::ObjectDimension;
78
79 static constexpr unsigned int OutputImageDimension = TOutputImage::ImageDimension;
80
83 virtual void
85
86 virtual void
87 SetInput(unsigned int, const InputSpatialObjectType * object);
88
91
93 GetInput(unsigned int idx);
94
97 template <class TReferenceImage>
98 void
99 SetReferenceImage(TReferenceImage * refImage)
100 {
101 this->SetOrigin(refImage->GetOrigin());
102 this->SetSpacing(refImage->GetSpacing());
103 this->SetDirection(refImage->GetDirection());
104 this->SetIndex(refImage->GetLargestPossibleRegion().GetIndex());
105 this->SetSize(refImage->GetLargestPossibleRegion().GetSize());
106 }
107
108
113 virtual void
114 SetSpacing(const SpacingType & spacing);
115
116 virtual void
117 SetSpacing(const double * spacing);
118
119 virtual void
120 SetSpacing(const float * spacing);
121
122 virtual const double *
123 GetSpacing() const;
124
125 virtual const SpacingType &
127
130 virtual void
132
133 virtual const DirectionType &
135
142 itkSetMacro(InsideValue, ValueType);
143 itkGetConstMacro(InsideValue, ValueType);
145
152 itkSetMacro(OutsideValue, ValueType);
153 itkGetConstMacro(OutsideValue, ValueType);
155
160 virtual void
161 SetOrigin(const PointType & origin);
162
163 virtual void
164 SetOrigin(const double * origin);
165
166 virtual void
167 SetOrigin(const float * origin);
168
169 virtual const double *
170 GetOrigin() const;
171
172 virtual const PointType &
174
178 itkSetMacro(Index, IndexType);
179 itkGetConstMacro(Index, IndexType);
181
186 itkSetMacro(ChildrenDepth, unsigned int);
187 itkGetConstMacro(ChildrenDepth, unsigned int);
189
191 itkSetMacro(Size, SizeType);
192 itkGetConstMacro(Size, SizeType);
194
197 itkSetMacro(UseObjectValue, bool);
198 itkGetConstMacro(UseObjectValue, bool);
199 itkBooleanMacro(UseObjectValue);
201
202protected:
204 ~SpatialObjectToImageFilter() override = default;
205
206 void
208 {} // do nothing
209 void
210 GenerateData() override;
211
217
220
221 unsigned int m_ChildrenDepth{};
222
225
227
228 void
229 PrintSelf(std::ostream & os, Indent indent) const override;
230
231private:
232};
233} // end namespace itk
234
235#ifndef ITK_MANUAL_INSTANTIATION
236# include "itkSpatialObjectToImageFilter.hxx"
237#endif
238
239#endif
typename OutputImageType::RegionType OutputImageRegionType
Control indentation during Print() invocation.
Definition itkIndent.h:50
virtual void SetInput(const DataObjectIdentifierType &key, DataObject *input)
Protected method for setting indexed and named inputs.
Implements transparent reference counting.
virtual void SetSize(SizeType _arg)
virtual void SetIndex(IndexType _arg)
virtual void SetOrigin(const PointType &origin)
typename OutputImageType::DirectionType DirectionType
virtual const PointType & GetOriginPoint() const
void SetReferenceImage(TReferenceImage *refImage)
virtual const double * GetSpacing() const
virtual void SetOrigin(const float *origin)
virtual const SpacingType & GetSpacingVector() const
static constexpr unsigned int ObjectDimension
typename OutputImageType::IndexType IndexType
virtual const double * GetOrigin() const
typename TInputSpatialObject::ChildrenListType ChildrenListType
virtual void SetSpacing(const SpacingType &spacing)
virtual void SetOrigin(const double *origin)
virtual void SetInput(unsigned int, const InputSpatialObjectType *object)
~SpatialObjectToImageFilter() override=default
const InputSpatialObjectType * GetInput()
virtual void SetInput(const InputSpatialObjectType *input)
static constexpr unsigned int OutputImageDimension
void PrintSelf(std::ostream &os, Indent indent) const override
virtual const DirectionType & GetDirection() const
typename OutputImageType::PointType PointType
virtual void SetSpacing(const double *spacing)
typename OutputImageType::SpacingType SpacingType
typename OutputImageType::ValueType ValueType
virtual void SetSpacing(const float *spacing)
typename InputSpatialObjectType::ConstPointer InputSpatialObjectConstPointer
typename OutputImageType::Pointer OutputImagePointer
typename OutputImageType::SizeType SizeType
typename InputSpatialObjectType::Pointer InputSpatialObjectPointer
const InputSpatialObjectType * GetInput(unsigned int idx)
virtual void SetDirection(const DirectionType &dir)
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Represent a n-dimensional index in a n-dimensional image.
Definition itkIndex.h:69
Represent a n-dimensional size (bounds) of a n-dimensional image.
Definition itkSize.h:70