ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkBSplineControlPointImageFunction.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright Insight Software Consortium
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  * http://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_ASSIGN(BSplineControlPointImageFunction);
63 
65  using Superclass = ImageFunction<TInputImage,
66  typename TInputImage::PixelType, TCoordRep>;
69 
71  itkNewMacro( Self );
72 
75 
77  static constexpr unsigned int ImageDimension = TInputImage::ImageDimension;
78 
80  using ControlPointLatticeType = TInputImage;
81  using InputImageType = TInputImage;
82  using CoordRepType = TCoordRep;
83  using PixelType = typename InputImageType::PixelType;
86  using PointType = typename Superclass::PointType;
88 
89  using SpacingType = typename InputImageType::SpacingType;
92 
97 
102  using ContinuousIndexType = typename Superclass::ContinuousIndexType;
103  using RealType = float;
104 
111 
116  void SetInputImage( const InputImageType * ) override;
117 
122  void SetSplineOrder( const unsigned int );
123 
128  void SetSplineOrder( const ArrayType & );
129 
133  itkGetConstReferenceMacro( SplineOrder, ArrayType );
134 
151  itkSetMacro( CloseDimension, ArrayType );
152 
156  itkGetConstReferenceMacro( CloseDimension, ArrayType );
157 
161  itkSetMacro( Spacing, SpacingType );
162  itkGetConstMacro( Spacing, SpacingType );
164 
168  itkSetMacro( Origin, OriginType );
169  itkGetConstMacro( Origin, OriginType );
171 
175  itkSetMacro( Size, SizeType );
176  itkGetConstMacro( Size, SizeType );
178 
187  itkSetMacro( BSplineEpsilon, RealType );
188  itkGetConstMacro( BSplineEpsilon, RealType );
190 
195  OutputType EvaluateAtParametricPoint( const PointType & ) const;
196 
201  OutputType EvaluateAtIndex( const IndexType & ) const override;
202 
207  OutputType EvaluateAtContinuousIndex(
208  const ContinuousIndexType & ) const override;
209 
215  OutputType Evaluate( const PointType & ) const override;
216 
221  GradientType EvaluateGradientAtParametricPoint( const PointType & ) const;
222 
227  GradientType EvaluateGradientAtIndex( const IndexType & ) const;
228 
233  GradientType EvaluateGradientAtContinuousIndex(
234  const ContinuousIndexType & ) const;
235 
241  GradientType EvaluateGradient( const PointType & ) const;
242 
248  HessianComponentType EvaluateHessianAtParametricPoint(
249  const PointType &, const unsigned int ) const;
250 
256  HessianComponentType EvaluateHessianAtIndex(
257  const IndexType &, const unsigned int ) const;
258 
264  HessianComponentType EvaluateHessianAtContinuousIndex(
265  const ContinuousIndexType &, const unsigned int ) const;
266 
272  HessianComponentType EvaluateHessian(
273  const PointType &, const unsigned int ) const;
274 
275 protected:
277  ~BSplineControlPointImageFunction() override = default;
278  void PrintSelf( std::ostream& os, Indent indent ) const override;
279 
280 private:
285 
289 
291 
292  typename KernelType::Pointer m_Kernel[ImageDimension];
297 
299 };
300 
301 } // end namespace itk
302 
303 #ifndef ITK_MANUAL_INSTANTIATION
304 #include "itkBSplineControlPointImageFunction.hxx"
305 #endif
306 
307 #endif
typename InputImageType::SpacingType SpacingType
Light weight base class for most itk classes.
A templated class holding a M x N size Matrix.
typename Superclass::ContinuousIndexType ContinuousIndexType
typename InputImageType::RegionType InputImageRegionType
BSpline kernel used for density estimation and nonparameteric regression.
BSpline kernel used for density estimation and nonparameteric regression.
Evaluate a B-spline object given a grid of control points.
Represent a n-dimensional size (bounds) of a n-dimensional image.
Definition: itkSize.h:68
Control indentation during Print() invocation.
Definition: itkIndent.h:49
Evaluates a function of an image at specified position.
Templated n-dimensional image class.
Definition: itkImage.h:75