ITK  6.0.0
Insight Toolkit
itkBSplineControlPointImageFunction.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 itkBSplineControlPointImageFunction_h
19#define itkBSplineControlPointImageFunction_h
20
21#include "itkImageFunction.h"
22
25#include "itkFixedArray.h"
26#include "itkImage.h"
27#include "itkPointSet.h"
29#include "itkVector.h"
30#include "itkVectorContainer.h"
31
32namespace itk
33{
57template <typename TInputImage, typename TCoordinate = double>
58class ITK_TEMPLATE_EXPORT BSplineControlPointImageFunction
59 : public ImageFunction<TInputImage, typename TInputImage::PixelType, TCoordinate>
60{
61public:
62 ITK_DISALLOW_COPY_AND_MOVE(BSplineControlPointImageFunction);
63
68
70 itkNewMacro(Self);
71
73 itkOverrideGetNameOfClassMacro(BSplineControlPointImageFunction);
74
76 static constexpr unsigned int ImageDimension = TInputImage::ImageDimension;
77
79 using ControlPointLatticeType = TInputImage;
80 using InputImageType = TInputImage;
81 using CoordinateType = TCoordinate;
82#ifndef ITK_FUTURE_LEGACY_REMOVE
83 using CoordRepType ITK_FUTURE_DEPRECATED(
84 "ITK 6 discourages using `CoordRepType`. Please use `CoordinateType` instead!") = CoordinateType;
85#endif
86 using PixelType = typename InputImageType::PixelType;
89 using typename Superclass::PointType;
91
92 using SpacingType = typename InputImageType::SpacingType;
95
100
105 using typename Superclass::ContinuousIndexType;
106 using RealType = float;
107
114
119 void
120 SetInputImage(const InputImageType *) override;
121
126 void
127 SetSplineOrder(const unsigned int);
128
133 void
135
139 itkGetConstReferenceMacro(SplineOrder, ArrayType);
140
157 itkSetMacro(CloseDimension, ArrayType);
158
162 itkGetConstReferenceMacro(CloseDimension, ArrayType);
163
167 itkSetMacro(Spacing, SpacingType);
168 itkGetConstMacro(Spacing, SpacingType);
174 itkSetMacro(Origin, OriginType);
175 itkGetConstMacro(Origin, OriginType);
181 itkSetMacro(Size, SizeType);
182 itkGetConstMacro(Size, SizeType);
193 itkSetMacro(BSplineEpsilon, RealType);
194 itkGetConstMacro(BSplineEpsilon, RealType);
203
209 EvaluateAtIndex(const IndexType &) const override;
210
217
224 Evaluate(const PointType &) const override;
225
232
239
246
254
261 EvaluateHessianAtParametricPoint(const PointType &, const unsigned int) const;
262
269 EvaluateHessianAtIndex(const IndexType &, const unsigned int) const;
270
277 EvaluateHessianAtContinuousIndex(const ContinuousIndexType &, const unsigned int) const;
278
285 EvaluateHessian(const PointType &, const unsigned int) const;
286
287protected:
290 void
291 PrintSelf(std::ostream & os, Indent indent) const override;
292
293private:
295 SizeType m_Size{};
296 SpacingType m_Spacing{};
297 OriginType m_Origin{};
298
299 ArrayType m_NumberOfControlPoints{};
300 ArrayType m_CloseDimension{};
301 ArrayType m_SplineOrder{};
302
303 RealImagePointer m_NeighborhoodWeightImage{};
304
305 typename KernelType::Pointer m_Kernel[ImageDimension]{};
306 typename KernelOrder0Type::Pointer m_KernelOrder0{};
307 typename KernelOrder1Type::Pointer m_KernelOrder1{};
308 typename KernelOrder2Type::Pointer m_KernelOrder2{};
309 typename KernelOrder3Type::Pointer m_KernelOrder3{};
310
311 CoordinateType m_BSplineEpsilon{};
312};
313
314} // end namespace itk
315
316#ifndef ITK_MANUAL_INSTANTIATION
317# include "itkBSplineControlPointImageFunction.hxx"
318#endif
319
320#endif
Evaluate a B-spline object given a grid of control points.
void PrintSelf(std::ostream &os, Indent indent) const override
~BSplineControlPointImageFunction() override=default
GradientType EvaluateGradientAtContinuousIndex(const ContinuousIndexType &) const
HessianComponentType EvaluateHessianAtParametricPoint(const PointType &, const unsigned int) const
HessianComponentType EvaluateHessian(const PointType &, const unsigned int) const
void SetInputImage(const InputImageType *) override
GradientType EvaluateGradientAtParametricPoint(const PointType &) const
HessianComponentType EvaluateHessianAtIndex(const IndexType &, const unsigned int) const
GradientType EvaluateGradient(const PointType &) const
OutputType EvaluateAtParametricPoint(const PointType &) const
OutputType EvaluateAtContinuousIndex(const ContinuousIndexType &) const override
HessianComponentType EvaluateHessianAtContinuousIndex(const ContinuousIndexType &, const unsigned int) const
typename InputImageType::RegionType InputImageRegionType
OutputType Evaluate(const PointType &) const override
typename InputImageType::SpacingType SpacingType
OutputType EvaluateAtIndex(const IndexType &) const override
void SetSplineOrder(const ArrayType &)
GradientType EvaluateGradientAtIndex(const IndexType &) const
void SetSplineOrder(const unsigned int)
BSpline kernel used for density estimation and nonparametric regression.
BSpline kernel used for density estimation and nonparametric regression.
Evaluates a function of an image at specified position.
Templated n-dimensional image class.
Definition: itkImage.h:89
Control indentation during Print() invocation.
Definition: itkIndent.h:50
Light weight base class for most itk classes.
A templated class holding a M x N size Matrix.
SmartPointer< Self > Pointer
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