ITK  5.4.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 TCoordRep = double>
58 class ITK_TEMPLATE_EXPORT BSplineControlPointImageFunction
59  : public ImageFunction<TInputImage, typename TInputImage::PixelType, TCoordRep>
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 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
130  SetSplineOrder(const ArrayType &);
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);
197  OutputType
198  EvaluateAtParametricPoint(const PointType &) const;
199 
204  OutputType
205  EvaluateAtIndex(const IndexType &) const override;
206 
211  OutputType
212  EvaluateAtContinuousIndex(const ContinuousIndexType &) const override;
213 
219  OutputType
220  Evaluate(const PointType &) const override;
221 
227  EvaluateGradientAtParametricPoint(const PointType &) const;
228 
234  EvaluateGradientAtIndex(const IndexType &) const;
235 
241  EvaluateGradientAtContinuousIndex(const ContinuousIndexType &) const;
242 
249  EvaluateGradient(const PointType &) const;
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 
283 protected:
285  ~BSplineControlPointImageFunction() override = default;
286  void
287  PrintSelf(std::ostream & os, Indent indent) const override;
288 
289 private:
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
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::RealImagePointer
typename RealImageType::Pointer RealImagePointer
Definition: itkBSplineControlPointImageFunction.h:100
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:71
itk::BSplineControlPointImageFunction::RegionType
typename InputImageType::RegionType RegionType
Definition: itkBSplineControlPointImageFunction.h:83
itk::GTest::TypedefsAndConstructors::Dimension2::PointType
ImageBaseType::PointType PointType
Definition: itkGTestTypedefsAndConstructors.h:51
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itk::BSplineControlPointImageFunction::InputImageType
TInputImage InputImageType
Definition: itkBSplineControlPointImageFunction.h:80
itkImage.h
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::VariableSizeMatrix
A templated class holding a M x N size Matrix.
Definition: itkVariableSizeMatrix.h:45
itk::BSplineControlPointImageFunction::ControlPointLatticeType
TInputImage ControlPointLatticeType
Definition: itkBSplineControlPointImageFunction.h:79
itk::BSplineControlPointImageFunction::SizeType
typename InputImageType::SizeType SizeType
Definition: itkBSplineControlPointImageFunction.h:90
itk::BSplineControlPointImageFunction::IndexType
typename InputImageType::IndexType IndexType
Definition: itkBSplineControlPointImageFunction.h:84
itk::BSplineControlPointImageFunction::InputImageRegionType
typename InputImageType::RegionType InputImageRegionType
Definition: itkBSplineControlPointImageFunction.h:86
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
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::BSplineControlPointImageFunction::SpacingType
typename InputImageType::SpacingType SpacingType
Definition: itkBSplineControlPointImageFunction.h:88
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itkFixedArray.h
itk::BSplineControlPointImageFunction::OutputType
PixelType OutputType
Definition: itkBSplineControlPointImageFunction.h:93
itk::FixedArray< unsigned int, ImageDimension >
itkVectorContainer.h
itkCoxDeBoorBSplineKernelFunction.h
itk::BSplineControlPointImageFunction::RealType
float RealType
Definition: itkBSplineControlPointImageFunction.h:102
itkBSplineKernelFunction.h
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itk::ContinuousIndex< TCoordRep, 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::BSplineControlPointImageFunction::PixelType
typename InputImageType::PixelType PixelType
Definition: itkBSplineControlPointImageFunction.h:82
itk::CoxDeBoorBSplineKernelFunction
BSpline kernel used for density estimation and nonparametric regression.
Definition: itkCoxDeBoorBSplineKernelFunction.h:57
itk::BSplineControlPointImageFunction::OriginType
typename InputImageType::PointType OriginType
Definition: itkBSplineControlPointImageFunction.h:89
itk::BSplineControlPointImageFunction::CoordRepType
TCoordRep CoordRepType
Definition: itkBSplineControlPointImageFunction.h:81