ITK  4.9.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>
61  : public ImageToImageFilter<TInputImage, TOutputImage>
62 {
63 public:
68 
70  itkNewMacro(Self);
71 
73  itkStaticConstMacro( ImageDimension, unsigned int,
74  TInputImage::ImageDimension );
75 
76  typedef TInputImage ControlPointLatticeType;
77  typedef TOutputImage OutputImageType;
78 
80  typedef typename OutputImageType::PixelType PixelType;
81  typedef typename OutputImageType::RegionType RegionType;
82  typedef typename OutputImageType::IndexType IndexType;
83  typedef typename OutputImageType::PointType PointType;
84  typedef typename OutputImageType::RegionType OutputImageRegionType;
85 
86  typedef typename OutputImageType::SpacingType SpacingType;
87  typedef typename OutputImageType::PointType OriginType;
88  typedef typename OutputImageType::SizeType SizeType;
89  typedef typename OutputImageType::DirectionType DirectionType;
90 
92  typedef float RealType;
93  typedef Image<RealType,
94  itkGetStaticConstMacro( ImageDimension )> RealImageType;
96 
97  typedef FixedArray<unsigned,
98  itkGetStaticConstMacro( ImageDimension )> ArrayType;
99  typedef FixedArray<RealType,
100  itkGetStaticConstMacro( ImageDimension )> RealArrayType;
101 
103  typedef PointSet<PixelType,
104  itkGetStaticConstMacro( ImageDimension )> PointSetType;
107  typedef Image<PointDataType,
108  itkGetStaticConstMacro( ImageDimension )> PointDataImageType;
111 
118 
123  void SetSplineOrder( unsigned int );
124 
129  void SetSplineOrder( ArrayType );
130 
134  itkGetConstReferenceMacro( SplineOrder, ArrayType );
135 
152  itkSetMacro( CloseDimension, ArrayType );
153  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 
191  itkSetMacro( Direction, DirectionType );
192 
196  itkGetConstMacro( Direction, DirectionType );
197 
205  typename ControlPointLatticeType::Pointer
207 
208 protected:
211  void PrintSelf( std::ostream& os, Indent indent ) const ITK_OVERRIDE;
212 
216  void ThreadedGenerateData( const OutputImageRegionType &, ThreadIdType ) ITK_OVERRIDE;
217 
218 private:
219  BSplineControlPointImageFilter( const Self& ) ITK_DELETE_FUNCTION;
220  void operator=( const Self& ) ITK_DELETE_FUNCTION;
221 
222 
227  void BeforeThreadedGenerateData() ITK_OVERRIDE;
228 
233  virtual unsigned int SplitRequestedRegion( unsigned int, unsigned int, OutputImageRegionType & ) ITK_OVERRIDE;
234 
239  void CollapsePhiLattice( PointDataImageType *, PointDataImageType *,
240  const RealType, const unsigned int );
241 
246 
252 
256  ArrayType m_CloseDimension;
257  ArrayType m_SplineOrder;
258  ArrayType m_NumberOfLevels;
259 
261 
262  typename KernelType::Pointer m_Kernel[ImageDimension];
267 
269 
270  inline typename RealImageType::IndexType
271  NumberToIndex( unsigned int number, typename RealImageType::SizeType size )
272  {
273  typename RealImageType::IndexType k;
274  k[0] = 1;
275 
276  for ( unsigned int i = 1; i < ImageDimension; i++ )
277  {
278  k[i] = size[ImageDimension-i-1]*k[i-1];
279  }
280  typename RealImageType::IndexType index;
281  for ( unsigned int i = 0; i < ImageDimension; i++ )
282  {
283  index[ImageDimension-i-1]
284  = static_cast<unsigned int>( number/k[ImageDimension-i-1] );
285  number %= k[ImageDimension-i-1];
286  }
287  return index;
288  }
289 
290 };
291 
292 } // end namespace itk
293 
294 #ifndef ITK_MANUAL_INSTANTIATION
295 #include "itkBSplineControlPointImageFilter.hxx"
296 #endif
297 
298 #endif
Image< PointDataType, itkGetStaticConstMacro(ImageDimension)> PointDataImageType
CoxDeBoorBSplineKernelFunction< 3 > KernelType
Process a given a B-spline grid of control points.
KernelType::Pointer m_Kernel[ImageDimension]
Represent the size (bounds) of a n-dimensional image.
Definition: itkSize.h:52
virtual unsigned int SplitRequestedRegion(unsigned int, unsigned int, OutputImageRegionType &) override
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.
void PrintSelf(std::ostream &os, Indent indent) const override
ImageToImageFilter< TInputImage, TOutputImage > Superclass
vnl_matrix< RealType > m_RefinedLatticeCoefficients[ImageDimension]
BSpline kernel used for density estimation and nonparameteric regression.
ControlPointLatticeType::Pointer RefineControlPointLattice(ArrayType)
FixedArray< unsigned, itkGetStaticConstMacro(ImageDimension)> ArrayType
RealImageType::IndexType NumberToIndex(unsigned int number, typename RealImageType::SizeType size)
Superclass::IndexType IndexType
Definition: itkImage.h:119
void CollapsePhiLattice(PointDataImageType *, PointDataImageType *, const RealType, const unsigned int)
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
void ThreadedGenerateData(const OutputImageRegionType &, ThreadIdType) override
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
void BeforeThreadedGenerateData() override
Templated n-dimensional image class.
Definition: itkImage.h:75