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"
28 #include "itkVariableSizeMatrix.h"
29 #include "itkVector.h"
30 #include "itkVectorContainer.h"
31 
32 namespace itk
33 {
57 template <typename TInputImage, typename TCoordinate = double>
58 class ITK_TEMPLATE_EXPORT BSplineControlPointImageFunction
59  : public ImageFunction<TInputImage, typename TInputImage::PixelType, TCoordinate>
60 {
61 public:
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
134  SetSplineOrder(const ArrayType &);
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);
201  OutputType
202  EvaluateAtParametricPoint(const PointType &) const;
203 
208  OutputType
209  EvaluateAtIndex(const IndexType &) const override;
210 
215  OutputType
216  EvaluateAtContinuousIndex(const ContinuousIndexType &) const override;
217 
223  OutputType
224  Evaluate(const PointType &) const override;
225 
231  EvaluateGradientAtParametricPoint(const PointType &) const;
232 
238  EvaluateGradientAtIndex(const IndexType &) const;
239 
245  EvaluateGradientAtContinuousIndex(const ContinuousIndexType &) const;
246 
253  EvaluateGradient(const PointType &) const;
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 
287 protected:
289  ~BSplineControlPointImageFunction() override = default;
290  void
291  PrintSelf(std::ostream & os, Indent indent) const override;
292 
293 private:
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
itk::BSplineControlPointImageFunction::RegionType
typename InputImageType::RegionType RegionType
Definition: itkBSplineControlPointImageFunction.h:87
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
itkImageFunction.h
itk::BSplineKernelFunction
BSpline kernel used for density estimation and nonparametric regression.
Definition: itkBSplineKernelFunction.h:43
itk::BSplineControlPointImageFunction::PixelType
typename InputImageType::PixelType PixelType
Definition: itkBSplineControlPointImageFunction.h:86
itk::BSplineControlPointImageFunction::OutputType
PixelType OutputType
Definition: itkBSplineControlPointImageFunction.h:97
itk::BSplineControlPointImageFunction
Evaluate a B-spline object given a grid of control points.
Definition: itkBSplineControlPointImageFunction.h:58
itkVariableSizeMatrix.h
itk::Size
Represent a n-dimensional size (bounds) of a n-dimensional image.
Definition: itkSize.h:69
itk::GTest::TypedefsAndConstructors::Dimension2::PointType
ImageBaseType::PointType PointType
Definition: itkGTestTypedefsAndConstructors.h:51
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itkImage.h
itk::BSplineControlPointImageFunction::RealImagePointer
typename RealImageType::Pointer RealImagePointer
Definition: itkBSplineControlPointImageFunction.h:104
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::BSplineControlPointImageFunction::OriginType
typename InputImageType::PointType OriginType
Definition: itkBSplineControlPointImageFunction.h:93
itk::VariableSizeMatrix
A templated class holding a M x N size Matrix.
Definition: itkVariableSizeMatrix.h:45
itk::BSplineControlPointImageFunction::InputImageType
TInputImage InputImageType
Definition: itkBSplineControlPointImageFunction.h:80
itk::BSplineControlPointImageFunction::SpacingType
typename InputImageType::SpacingType SpacingType
Definition: itkBSplineControlPointImageFunction.h:92
itk::BSplineControlPointImageFunction::IndexType
typename InputImageType::IndexType IndexType
Definition: itkBSplineControlPointImageFunction.h:88
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itk::BSplineControlPointImageFunction::CoordinateType
TCoordinate CoordinateType
Definition: itkBSplineControlPointImageFunction.h:81
itk::ImageFunction
Evaluates a function of an image at specified position.
Definition: itkImageFunction.h:55
itk::LightObject
Light weight base class for most itk classes.
Definition: itkLightObject.h:55
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itkFixedArray.h
itk::BSplineControlPointImageFunction::ControlPointLatticeType
TInputImage ControlPointLatticeType
Definition: itkBSplineControlPointImageFunction.h:79
itk::FixedArray< unsigned int, ImageDimension >
itkVectorContainer.h
itkCoxDeBoorBSplineKernelFunction.h
itkBSplineKernelFunction.h
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnatomicalOrientation.h:29
itk::ContinuousIndex< TCoordinate, Self::ImageDimension >
itkVector.h
itk::Point
A templated class holding a geometric point in n-Dimensional space.
Definition: itkPoint.h:53
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:88
itkPointSet.h
itk::CoxDeBoorBSplineKernelFunction
BSpline kernel used for density estimation and nonparametric regression.
Definition: itkCoxDeBoorBSplineKernelFunction.h:57
itk::BSplineControlPointImageFunction::SizeType
typename InputImageType::SizeType SizeType
Definition: itkBSplineControlPointImageFunction.h:94
itk::BSplineControlPointImageFunction::RealType
float RealType
Definition: itkBSplineControlPointImageFunction.h:106
itk::BSplineControlPointImageFunction::InputImageRegionType
typename InputImageType::RegionType InputImageRegionType
Definition: itkBSplineControlPointImageFunction.h:90