ITK  4.13.0
Insight Segmentation and Registration Toolkit
itkBSplineControlPointImageFilter.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 itkBSplineControlPointImageFilter_h
19 #define itkBSplineControlPointImageFilter_h
20 
21 #include "itkImageToImageFilter.h"
24 #include "itkFixedArray.h"
25 #include "itkPointSet.h"
26 #include "itkVariableSizeMatrix.h"
27 #include "itkVector.h"
28 #include "itkVectorContainer.h"
29 
30 #include "vnl/vnl_matrix.h"
31 #include "vnl/vnl_vector.h"
32 
33 namespace itk
34 {
35 
59 template <typename TInputImage, typename TOutputImage = TInputImage>
60 class ITK_TEMPLATE_EXPORT BSplineControlPointImageFilter
61  : public ImageToImageFilter<TInputImage, TOutputImage>
62 {
63 public:
68 
70  itkNewMacro(Self);
71 
74 
76  itkStaticConstMacro( ImageDimension, unsigned int,
77  TInputImage::ImageDimension );
78 
79  typedef TInputImage ControlPointLatticeType;
80  typedef TOutputImage OutputImageType;
81 
83  typedef typename OutputImageType::PixelType PixelType;
84  typedef typename OutputImageType::RegionType RegionType;
87  typedef typename OutputImageType::RegionType OutputImageRegionType;
88 
89  typedef typename OutputImageType::SpacingType SpacingType;
93 
95  typedef float RealType;
96  typedef Image<RealType,
97  itkGetStaticConstMacro( ImageDimension )> RealImageType;
99 
100  typedef FixedArray<unsigned,
101  itkGetStaticConstMacro( ImageDimension )> ArrayType;
102  typedef FixedArray<RealType,
103  itkGetStaticConstMacro( ImageDimension )> RealArrayType;
104 
106  typedef PointSet<PixelType,
107  itkGetStaticConstMacro( ImageDimension )> PointSetType;
110  typedef Image<PointDataType,
111  itkGetStaticConstMacro( ImageDimension )> PointDataImageType;
114 
121 
126  void SetSplineOrder( unsigned int );
127 
132  void SetSplineOrder( ArrayType );
133 
137  itkGetConstReferenceMacro( SplineOrder, ArrayType );
138 
155  itkSetMacro( CloseDimension, ArrayType );
156  itkGetConstReferenceMacro( CloseDimension, ArrayType );
158 
162  itkSetMacro( Spacing, SpacingType );
163  itkGetConstMacro( Spacing, SpacingType );
165 
169  itkSetMacro( Origin, OriginType );
170  itkGetConstMacro( Origin, OriginType );
172 
176  itkSetMacro( Size, SizeType );
177  itkGetConstMacro( Size, SizeType );
179 
194  itkSetMacro( Direction, DirectionType );
195 
199  itkGetConstMacro( Direction, DirectionType );
200 
208  typename ControlPointLatticeType::Pointer
209  RefineControlPointLattice( ArrayType );
210 
211 protected:
213  virtual ~BSplineControlPointImageFilter() ITK_OVERRIDE;
214  void PrintSelf( std::ostream& os, Indent indent ) const ITK_OVERRIDE;
215 
219  void ThreadedGenerateData( const OutputImageRegionType &, ThreadIdType ) ITK_OVERRIDE;
220 
221 private:
222  ITK_DISALLOW_COPY_AND_ASSIGN(BSplineControlPointImageFilter);
223 
224 
229  void BeforeThreadedGenerateData() ITK_OVERRIDE;
230 
235  virtual unsigned int SplitRequestedRegion( unsigned int, unsigned int, OutputImageRegionType & ) ITK_OVERRIDE;
236 
241  void CollapsePhiLattice( PointDataImageType *, PointDataImageType *,
242  const RealType, const unsigned int );
243 
247  void SetNumberOfLevels( ArrayType );
248 
250  SizeType m_Size;
251  SpacingType m_Spacing;
252  OriginType m_Origin;
253  DirectionType m_Direction;
254 
255  bool m_DoMultilevel;
256  unsigned int m_MaximumNumberOfLevels;
257  ArrayType m_NumberOfControlPoints;
258  ArrayType m_CloseDimension;
259  ArrayType m_SplineOrder;
260  ArrayType m_NumberOfLevels;
261 
262  vnl_matrix<RealType> m_RefinedLatticeCoefficients[ImageDimension];
263 
264  typename KernelType::Pointer m_Kernel[ImageDimension];
265  typename KernelOrder0Type::Pointer m_KernelOrder0;
266  typename KernelOrder1Type::Pointer m_KernelOrder1;
267  typename KernelOrder2Type::Pointer m_KernelOrder2;
268  typename KernelOrder3Type::Pointer m_KernelOrder3;
269 
270  RealType m_BSplineEpsilon;
271 
272  inline typename RealImageType::IndexType
273  NumberToIndex( unsigned int number, typename RealImageType::SizeType size )
274  {
275  typename RealImageType::IndexType k;
276  k[0] = 1;
277 
278  for ( unsigned int i = 1; i < ImageDimension; i++ )
279  {
280  k[i] = size[ImageDimension-i-1]*k[i-1];
281  }
282  typename RealImageType::IndexType index;
283  for ( unsigned int i = 0; i < ImageDimension; i++ )
284  {
285  index[ImageDimension-i-1]
286  = static_cast<unsigned int>( number/k[ImageDimension-i-1] );
287  number %= k[ImageDimension-i-1];
288  }
289  return index;
290  }
291 
292 };
293 
294 } // end namespace itk
295 
296 #ifndef ITK_MANUAL_INSTANTIATION
297 #include "itkBSplineControlPointImageFilter.hxx"
298 #endif
299 
300 #endif
Image< PointDataType, itkGetStaticConstMacro(ImageDimension)> PointDataImageType
CoxDeBoorBSplineKernelFunction< 3 > KernelType
Process a given a B-spline grid of control points.
Represent the size (bounds) of a n-dimensional image.
Definition: itkSize.h:52
PointSetType::PointDataContainer PointDataContainerType
Image< RealType, itkGetStaticConstMacro(ImageDimension)> RealImageType
FixedArray< RealType, itkGetStaticConstMacro(ImageDimension)> RealArrayType
Base class for all process objects that output image data.
MeshTraits::PixelType PixelType
Definition: itkPointSet.h:101
Simulate a standard C array with copy semnatics.
Definition: itkFixedArray.h:50
BSpline kernel used for density estimation and nonparameteric regression.
ImageToImageFilter< TInputImage, TOutputImage > Superclass
BSpline kernel used for density estimation and nonparameteric regression.
FixedArray< unsigned, itkGetStaticConstMacro(ImageDimension)> ArrayType
A superclass of the N-dimensional mesh structure; supports point (geometric coordinate and attribute)...
Definition: itkPointSet.h:84
unsigned int ThreadIdType
Definition: itkIntTypes.h:159
Base class for filters that take an image as input and produce an image as output.
OutputImageType::RegionType OutputImageRegionType
Control indentation during Print() invocation.
Definition: itkIndent.h:49
MeshTraits::PointDataContainer PointDataContainer
Definition: itkPointSet.h:108
PointSet< PixelType, itkGetStaticConstMacro(ImageDimension)> PointSetType
Templated n-dimensional image class.
Definition: itkImage.h:75