ITK  4.13.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:
63  typedef ImageFunction<TInputImage,
64  typename TInputImage::PixelType, TCoordRep> Superclass;
67 
69  itkNewMacro( Self );
70 
73 
75  itkStaticConstMacro( ImageDimension, unsigned int, TInputImage::ImageDimension );
76 
78  typedef TInputImage ControlPointLatticeType;
79  typedef TInputImage InputImageType;
80  typedef TCoordRep CoordRepType;
81  typedef typename InputImageType::PixelType PixelType;
82  typedef typename InputImageType::RegionType RegionType;
84  typedef typename Superclass::PointType PointType;
85  typedef typename InputImageType::RegionType InputImageRegionType;
86 
87  typedef typename InputImageType::SpacingType SpacingType;
90 
95 
100  typedef typename Superclass::ContinuousIndexType ContinuousIndexType;
101  typedef float RealType;
102 
109 
114  virtual void SetInputImage( const InputImageType * ) ITK_OVERRIDE;
115 
120  void SetSplineOrder( const unsigned int );
121 
126  void SetSplineOrder( const ArrayType & );
127 
131  itkGetConstReferenceMacro( SplineOrder, ArrayType );
132 
149  itkSetMacro( CloseDimension, ArrayType );
150 
154  itkGetConstReferenceMacro( CloseDimension, ArrayType );
155 
159  itkSetMacro( Spacing, SpacingType );
160  itkGetConstMacro( Spacing, SpacingType );
162 
166  itkSetMacro( Origin, OriginType );
167  itkGetConstMacro( Origin, OriginType );
169 
173  itkSetMacro( Size, SizeType );
174  itkGetConstMacro( Size, SizeType );
176 
185  itkSetMacro( BSplineEpsilon, RealType );
186  itkGetConstMacro( BSplineEpsilon, RealType );
188 
193  OutputType EvaluateAtParametricPoint( const PointType & ) const;
194 
199  virtual OutputType EvaluateAtIndex( const IndexType & ) const ITK_OVERRIDE;
200 
205  virtual OutputType EvaluateAtContinuousIndex(
206  const ContinuousIndexType & ) const ITK_OVERRIDE;
207 
213  virtual OutputType Evaluate( const PointType & ) const ITK_OVERRIDE;
214 
219  GradientType EvaluateGradientAtParametricPoint( const PointType & ) const;
220 
225  GradientType EvaluateGradientAtIndex( const IndexType & ) const;
226 
231  GradientType EvaluateGradientAtContinuousIndex(
232  const ContinuousIndexType & ) const;
233 
239  GradientType EvaluateGradient( const PointType & ) const;
240 
246  HessianComponentType EvaluateHessianAtParametricPoint(
247  const PointType &, const unsigned int ) const;
248 
254  HessianComponentType EvaluateHessianAtIndex(
255  const IndexType &, const unsigned int ) const;
256 
262  HessianComponentType EvaluateHessianAtContinuousIndex(
263  const ContinuousIndexType &, const unsigned int ) const;
264 
270  HessianComponentType EvaluateHessian(
271  const PointType &, const unsigned int ) const;
272 
273 protected:
275  virtual ~BSplineControlPointImageFunction() ITK_OVERRIDE;
276  void PrintSelf( std::ostream& os, Indent indent ) const ITK_OVERRIDE;
277 
278 private:
279  ITK_DISALLOW_COPY_AND_ASSIGN(BSplineControlPointImageFunction);
280 
282  SizeType m_Size;
283  SpacingType m_Spacing;
284  OriginType m_Origin;
285 
286  ArrayType m_NumberOfControlPoints;
287  ArrayType m_CloseDimension;
288  ArrayType m_SplineOrder;
289 
290  RealImagePointer m_NeighborhoodWeightImage;
291 
292  typename KernelType::Pointer m_Kernel[ImageDimension];
293  typename KernelOrder0Type::Pointer m_KernelOrder0;
294  typename KernelOrder1Type::Pointer m_KernelOrder1;
295  typename KernelOrder2Type::Pointer m_KernelOrder2;
296  typename KernelOrder3Type::Pointer m_KernelOrder3;
297 
298  CoordRepType m_BSplineEpsilon;
299 };
300 
301 } // end namespace itk
302 
303 #ifndef ITK_MANUAL_INSTANTIATION
304 #include "itkBSplineControlPointImageFunction.hxx"
305 #endif
306 
307 #endif
Light weight base class for most itk classes.
Represent the size (bounds) of a n-dimensional image.
Definition: itkSize.h:52
A templated class holding a M x N size Matrix.
BSpline kernel used for density estimation and nonparameteric regression.
Image< CoordRepType, ImageDimension > RealImageType
BSpline kernel used for density estimation and nonparameteric regression.
Evaluate a B-spline object given a grid of control points.
ImageFunction< TInputImage, typename TInputImage::PixelType, TCoordRep > Superclass
FixedArray< unsigned, ImageDimension > ArrayType
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
VariableSizeMatrix< CoordRepType > HessianComponentType