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 TCoordRep = double>
58class ITK_TEMPLATE_EXPORT BSplineControlPointImageFunction
59 : public ImageFunction<TInputImage, typename TInputImage::PixelType, TCoordRep>
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 CoordRepType = TCoordRep;
82 using PixelType = typename InputImageType::PixelType;
85 using typename Superclass::PointType;
87
88 using SpacingType = typename InputImageType::SpacingType;
91
96
101 using typename Superclass::ContinuousIndexType;
102 using RealType = float;
103
110
115 void
116 SetInputImage(const InputImageType *) override;
117
122 void
123 SetSplineOrder(const unsigned int);
124
129 void
131
135 itkGetConstReferenceMacro(SplineOrder, ArrayType);
136
153 itkSetMacro(CloseDimension, ArrayType);
154
158 itkGetConstReferenceMacro(CloseDimension, ArrayType);
159
163 itkSetMacro(Spacing, SpacingType);
164 itkGetConstMacro(Spacing, SpacingType);
170 itkSetMacro(Origin, OriginType);
171 itkGetConstMacro(Origin, OriginType);
177 itkSetMacro(Size, SizeType);
178 itkGetConstMacro(Size, SizeType);
189 itkSetMacro(BSplineEpsilon, RealType);
190 itkGetConstMacro(BSplineEpsilon, RealType);
199
205 EvaluateAtIndex(const IndexType &) const override;
206
213
220 Evaluate(const PointType &) const override;
221
228
235
242
250
257 EvaluateHessianAtParametricPoint(const PointType &, const unsigned int) const;
258
265 EvaluateHessianAtIndex(const IndexType &, const unsigned int) const;
266
273 EvaluateHessianAtContinuousIndex(const ContinuousIndexType &, const unsigned int) const;
274
281 EvaluateHessian(const PointType &, const unsigned int) const;
282
283protected:
286 void
287 PrintSelf(std::ostream & os, Indent indent) const override;
288
289private:
291 SizeType m_Size{};
292 SpacingType m_Spacing{};
293 OriginType m_Origin{};
294
295 ArrayType m_NumberOfControlPoints{};
296 ArrayType m_CloseDimension{};
297 ArrayType m_SplineOrder{};
298
299 RealImagePointer m_NeighborhoodWeightImage{};
300
301 typename KernelType::Pointer m_Kernel[ImageDimension]{};
302 typename KernelOrder0Type::Pointer m_KernelOrder0{};
303 typename KernelOrder1Type::Pointer m_KernelOrder1{};
304 typename KernelOrder2Type::Pointer m_KernelOrder2{};
305 typename KernelOrder3Type::Pointer m_KernelOrder3{};
306
307 CoordRepType m_BSplineEpsilon{};
308};
309
310} // end namespace itk
311
312#ifndef ITK_MANUAL_INSTANTIATION
313# include "itkBSplineControlPointImageFunction.hxx"
314#endif
315
316#endif
Evaluate a B-spline object given a grid of control points.
HessianComponentType EvaluateHessianAtContinuousIndex(const ContinuousIndexType &, const unsigned int) const
~BSplineControlPointImageFunction() override=default
void SetInputImage(const InputImageType *) override
void SetSplineOrder(const ArrayType &)
OutputType EvaluateAtContinuousIndex(const ContinuousIndexType &) const override
HessianComponentType EvaluateHessianAtParametricPoint(const PointType &, const unsigned int) const
GradientType EvaluateGradientAtContinuousIndex(const ContinuousIndexType &) const
void PrintSelf(std::ostream &os, Indent indent) const override
GradientType EvaluateGradientAtIndex(const IndexType &) const
GradientType EvaluateGradientAtParametricPoint(const PointType &) const
OutputType Evaluate(const PointType &) const override
OutputType EvaluateAtParametricPoint(const PointType &) const
typename InputImageType::SpacingType SpacingType
void SetSplineOrder(const unsigned int)
HessianComponentType EvaluateHessian(const PointType &, const unsigned int) const
HessianComponentType EvaluateHessianAtIndex(const IndexType &, const unsigned int) const
GradientType EvaluateGradient(const PointType &) const
OutputType EvaluateAtIndex(const IndexType &) const override
typename InputImageType::RegionType InputImageRegionType
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