ITK  4.8.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 #include "vnl/vnl_matrix.h"
33 #include "vnl/vnl_vector.h"
34 
35 namespace itk
36 {
65 template <typename TInputImage, typename TCoordRep = double>
67 : public ImageFunction<TInputImage, typename TInputImage::PixelType, TCoordRep>
68 {
69 public:
71  typedef ImageFunction<TInputImage,
72  typename TInputImage::PixelType, TCoordRep> Superclass;
75 
77  itkNewMacro(Self);
78 
80  itkStaticConstMacro( ImageDimension, unsigned int, TInputImage::ImageDimension );
81 
83  typedef TInputImage ControlPointLatticeType;
84  typedef TInputImage InputImageType;
85  typedef TCoordRep CoordRepType;
86  typedef typename InputImageType::PixelType PixelType;
87  typedef typename InputImageType::RegionType RegionType;
88  typedef typename InputImageType::IndexType IndexType;
89  typedef typename Superclass::PointType PointType;
90  typedef typename InputImageType::RegionType InputImageRegionType;
91 
92  typedef typename InputImageType::SpacingType SpacingType;
93  typedef typename InputImageType::PointType OriginType;
94  typedef typename InputImageType::SizeType SizeType;
95 
100 
106 
113 
118  virtual void SetInputImage( const InputImageType * ) ITK_OVERRIDE;
119 
124  void SetSplineOrder( const unsigned int );
125 
130  void SetSplineOrder( const ArrayType & );
131 
135  itkGetConstReferenceMacro( SplineOrder, ArrayType );
136 
153  itkSetMacro( CloseDimension, ArrayType );
154 
158  itkGetConstReferenceMacro( CloseDimension, ArrayType );
159 
163  itkSetMacro( Spacing, SpacingType );
164 
168  itkGetConstMacro( Spacing, SpacingType );
169 
173  itkSetMacro( Origin, OriginType );
174 
178  itkGetConstMacro( Origin, OriginType );
179 
183  itkSetMacro( Size, SizeType );
184 
188  itkGetConstMacro( Size, SizeType );
189 
195 
200  virtual OutputType EvaluateAtIndex( const IndexType & ) const ITK_OVERRIDE;
201 
207  const ContinuousIndexType & ) const ITK_OVERRIDE;
208 
214  virtual OutputType Evaluate( const PointType & ) const ITK_OVERRIDE;
215 
221 
227 
233  const ContinuousIndexType & ) const;
234 
240  GradientType EvaluateGradient( const PointType & ) const;
241 
248  const PointType &, const unsigned int ) const;
249 
256  const IndexType &, const unsigned int ) const;
257 
264  const ContinuousIndexType &, const unsigned int ) const;
265 
272  const PointType &, const unsigned int ) const;
273 
274 protected:
276  virtual ~BSplineControlPointImageFunction();
277  void PrintSelf( std::ostream& os, Indent indent ) const ITK_OVERRIDE;
278 
279 private:
280  BSplineControlPointImageFunction( const Self& ); //purposely not implemented
281  void operator=( const Self& ); //purposely not implemented
282 
284  SizeType m_Size;
285  SpacingType m_Spacing;
286  OriginType m_Origin;
287 
289  ArrayType m_CloseDimension;
290  ArrayType m_SplineOrder;
291 
293 
299 
301 
302 };
303 
304 } // end namespace itk
305 
306 #ifndef ITK_MANUAL_INSTANTIATION
307 #include "itkBSplineControlPointImageFunction.hxx"
308 #endif
309 
310 #endif
void SetSplineOrder(const unsigned int)
HessianComponentType EvaluateHessianAtContinuousIndex(const ContinuousIndexType &, const unsigned int) const
Light weight base class for most itk classes.
Point< TCoordRep, itkGetStaticConstMacro(ImageDimension) > PointType
virtual OutputType EvaluateAtContinuousIndex(const ContinuousIndexType &) const override
virtual void SetInputImage(const InputImageType *) override
GradientType EvaluateGradientAtParametricPoint(const PointType &) const
Represent the size (bounds) of a n-dimensional image.
Definition: itkSize.h:52
void PrintSelf(std::ostream &os, Indent indent) const override
A templated class holding a M x N size Matrix.
GradientType EvaluateGradient(const PointType &) const
HessianComponentType EvaluateHessianAtParametricPoint(const PointType &, const unsigned int) const
BSpline kernel used for density estimation and nonparameteric regression.
HessianComponentType EvaluateHessianAtIndex(const IndexType &, const unsigned int) const
HessianComponentType EvaluateHessian(const PointType &, const unsigned int) const
Image< CoordRepType, ImageDimension > RealImageType
BSpline kernel used for density estimation and nonparameteric regression.
Evaluate a B-spline object given a grid of control points.
ContinuousIndex< TCoordRep, itkGetStaticConstMacro(ImageDimension) > ContinuousIndexType
ImageFunction< TInputImage, typename TInputImage::PixelType, TCoordRep > Superclass
FixedArray< unsigned, ImageDimension > ArrayType
virtual OutputType Evaluate(const PointType &) const override
Control indentation during Print() invocation.
Definition: itkIndent.h:49
OutputType EvaluateAtParametricPoint(const PointType &) const
GradientType EvaluateGradientAtContinuousIndex(const ContinuousIndexType &) const
virtual OutputType EvaluateAtIndex(const IndexType &) const override
Evaluates a function of an image at specified position.
Templated n-dimensional image class.
Definition: itkImage.h:75
GradientType EvaluateGradientAtIndex(const IndexType &) const
VariableSizeMatrix< CoordRepType > HessianComponentType