ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
itkDisplacementFieldToBSplineImageFilter.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 itkDisplacementFieldToBSplineImageFilter_h
19#define itkDisplacementFieldToBSplineImageFilter_h
20
22
24#include "itkPointSet.h"
25#include "itkVector.h"
26
27namespace itk
28{
29
40
41template <typename TInputImage,
43 typename TOutputImage = TInputImage>
44class ITK_TEMPLATE_EXPORT DisplacementFieldToBSplineImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
45{
46public:
47 ITK_DISALLOW_COPY_AND_MOVE(DisplacementFieldToBSplineImageFilter);
48
53
55 itkNewMacro(Self);
56
58 itkOverrideGetNameOfClassMacro(DisplacementFieldToBSplineImageFilter);
59
61 static constexpr unsigned int ImageDimension = TInputImage::ImageDimension;
62
63 using InputFieldType = TInputImage;
64 using InputPointSetType = TInputPointSet;
65 using OutputFieldType = TOutputImage;
66
69 using InputFieldPointType = typename InputFieldType::PointType;
70
72 using PixelType = typename OutputFieldType::PixelType;
73 using VectorType = typename OutputFieldType::PixelType;
74 using RegionType = typename OutputFieldType::RegionType;
75 using IndexType = typename OutputFieldType::IndexType;
76
77 using SpacingType = typename OutputFieldType::SpacingType;
78 using OriginType = typename OutputFieldType::PointType;
79 using SizeType = typename OutputFieldType::SizeType;
80 using DirectionType = typename OutputFieldType::DirectionType;
81
82 using RealType = typename VectorType::RealValueType;
84
90
96
98 void
100 {
101 this->SetInput(0, field);
102 }
103
105 const InputFieldType *
107 {
108 return this->GetInput(0);
109 }
110
117 void
119 {
120 this->SetNthInput(1, const_cast<RealImageType *>(image));
121 }
122 void
124 {
125 this->SetConfidenceImage(image);
126 }
127
128
130 const RealImageType *
132 {
133 return static_cast<const RealImageType *>(this->ProcessObject::GetInput(1));
134 }
135
138 void
140 {
141 this->SetNthInput(2, const_cast<InputPointSetType *>(points));
142 }
143 void
145 {
146 this->SetPointSet(points);
147 }
148
149
151 const InputPointSetType *
153 {
154 return static_cast<const InputPointSetType *>(this->ProcessObject::GetInput(2));
155 }
156
158 void
160
164 {
165 return static_cast<const DisplacementFieldControlPointLatticeType *>(this->GetOutput(1));
166 }
167
169 void
171
173 void
175 {
176 this->SetBSplineDomainFromImage(const_cast<RealImageType *>(image));
177 }
178
180 void
182
184 void
186 {
187 this->SetBSplineDomainFromImage(const_cast<InputFieldType *>(field));
188 }
189
192
193 /* Set/Get b-spline domain origin. */
194 itkGetConstMacro(BSplineDomainOrigin, OriginType);
195
196 /* Set/Get b-spline domain spacing. */
197 itkGetConstMacro(BSplineDomainSpacing, SpacingType);
198
199 /* Set/Get b-spline domain size. */
200 itkGetConstMacro(BSplineDomainSize, SizeType);
201
202 /* Set/Get b-spline domain direction. */
203 itkGetConstMacro(BSplineDomainDirection, DirectionType);
204
205 /* Use input field to define the B-spline domain. */
206 itkSetMacro(UseInputFieldToDefineTheBSplineDomain, bool);
207 itkGetConstMacro(UseInputFieldToDefineTheBSplineDomain, bool);
208 itkBooleanMacro(UseInputFieldToDefineTheBSplineDomain);
209
213 itkSetMacro(SplineOrder, unsigned int);
214
218 itkGetConstMacro(SplineOrder, unsigned int);
219
227 itkSetMacro(NumberOfControlPoints, ArrayType);
228
236 itkGetConstMacro(NumberOfControlPoints, ArrayType);
237
244 itkSetMacro(NumberOfFittingLevels, ArrayType);
245
252 void
254 {
255 auto nlevels = MakeFilled<ArrayType>(n);
256 this->SetNumberOfFittingLevels(nlevels);
257 }
258
265 itkGetConstMacro(NumberOfFittingLevels, ArrayType);
266
271 itkBooleanMacro(EstimateInverse);
272 itkSetMacro(EstimateInverse, bool);
273 itkGetConstMacro(EstimateInverse, bool);
275
280 itkBooleanMacro(EnforceStationaryBoundary);
281 itkSetMacro(EnforceStationaryBoundary, bool);
282 itkGetConstMacro(EnforceStationaryBoundary, bool);
284
285protected:
288
291
293 void
294 PrintSelf(std::ostream & os, Indent indent) const override;
295
297 void
298 GenerateData() override;
299
300private:
301 bool m_EstimateInverse{ false };
303 unsigned int m_SplineOrder{ 3 };
306
307 typename WeightsContainerType::Pointer m_PointWeights{};
308 bool m_UsePointWeights{ false };
309
314
317};
318
319} // end namespace itk
320
321#ifndef ITK_MANUAL_INSTANTIATION
322# include "itkDisplacementFieldToBSplineImageFilter.hxx"
323#endif
324
325#endif
Image filter which provides a B-spline output approximation.
void SetPointSetConfidenceWeights(WeightsContainerType *weights)
const DisplacementFieldControlPointLatticeType * GetDisplacementFieldControlPointLattice() const
PointSet< typename ConstantVelocityFieldType::PixelType, ConstantVelocityFieldType::ImageDimension > InputPointSetType
ImageToImageFilter< ConstantVelocityFieldType, ConstantVelocityFieldType > Superclass
void PrintSelf(std::ostream &os, Indent indent) const override
void SetBSplineDomain(OriginType, SpacingType, SizeType, DirectionType)
void SetBSplineDomainFromImage(InputFieldType *)
void SetBSplineDomainFromImage(RealImageType *)
BSplineScatteredDataPointSetToImageFilter< InputPointSetType, OutputFieldType > BSplineFilterType
~DisplacementFieldToBSplineImageFilter() override=default
OutputImageType * GetOutput()
virtual void SetInput(const InputImageType *input)
const InputImageType * GetInput() const
Templated n-dimensional image class.
Definition itkImage.h:89
Control indentation during Print() invocation.
Definition itkIndent.h:50
TPointsContainer PointsContainer
A superclass of the N-dimensional mesh structure; supports point (geometric coordinate and attribute)...
Definition itkPointSet.h:82
typename MeshTraits::PixelType PixelType
typename MeshTraits::PointType PointType
typename MeshTraits::PointDataContainer PointDataContainer
virtual void SetNthInput(DataObjectPointerArraySizeType idx, DataObject *input)
DataObject * GetInput(const DataObjectIdentifierType &key)
Return an input.
Implements transparent reference counting.
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
constexpr TContainer MakeFilled(typename TContainer::const_reference value)