ITK  5.0.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:
64  ITK_DISALLOW_COPY_AND_ASSIGN(BSplineControlPointImageFilter);
65 
70 
72  itkNewMacro(Self);
73 
76 
78  static constexpr unsigned int ImageDimension = TInputImage::ImageDimension;
79 
80  using ControlPointLatticeType = TInputImage;
81  using OutputImageType = TOutputImage;
82 
84  using PixelType = typename OutputImageType::PixelType;
89 
90  using SpacingType = typename OutputImageType::SpacingType;
94 
96  using RealType = float;
98  Self::ImageDimension >;
100 
101  using ArrayType = FixedArray<unsigned,
102  Self::ImageDimension >;
104  Self::ImageDimension >;
105 
108  Self::ImageDimension >;
112  Self::ImageDimension >;
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  ~BSplineControlPointImageFilter() override = default;
214  void PrintSelf( std::ostream& os, Indent indent ) const override;
215 
217  void DynamicThreadedGenerateData( const OutputImageRegionType & ) override;
218 
219 
220 private:
221 
226  void BeforeThreadedGenerateData() override;
227 
232  unsigned int SplitRequestedRegion( unsigned int, unsigned int, OutputImageRegionType & ) override;
233 
238  void CollapsePhiLattice( PointDataImageType *, PointDataImageType *,
239  const RealType, const unsigned int );
240 
244  void SetNumberOfLevels( ArrayType );
245 
251 
252  bool m_DoMultilevel{ false };
253  unsigned int m_MaximumNumberOfLevels{ 1 };
258 
259  vnl_matrix<RealType> m_RefinedLatticeCoefficients[ImageDimension];
260 
261  typename KernelType::Pointer m_Kernel[ImageDimension];
266 
267  RealType m_BSplineEpsilon{ static_cast< RealType >( 1e-3 ) };
268 
269  inline typename RealImageType::IndexType
270  NumberToIndex( unsigned int number, typename RealImageType::SizeType size )
271  {
272  typename RealImageType::IndexType k;
273  k[0] = 1;
274 
275  for ( unsigned int i = 1; i < ImageDimension; i++ )
276  {
277  k[i] = size[ImageDimension-i-1]*k[i-1];
278  }
279  typename RealImageType::IndexType index;
280  for ( unsigned int i = 0; i < ImageDimension; i++ )
281  {
282  index[ImageDimension-i-1]
283  = static_cast<unsigned int>( number/k[ImageDimension-i-1] );
284  number %= k[ImageDimension-i-1];
285  }
286  return index;
287  }
288 
289 };
290 
291 } // end namespace itk
292 
293 #ifndef ITK_MANUAL_INSTANTIATION
294 #include "itkBSplineControlPointImageFilter.hxx"
295 #endif
296 
297 #endif
typename OutputImageType::PointType PointType
Process a given a B-spline grid of control points.
typename OutputImageType::IndexType IndexType
typename MeshTraits::PixelType PixelType
Definition: itkPointSet.h:103
typename PointSetType::PointDataContainer PointDataContainerType
typename PointDataImageType::Pointer PointDataImagePointer
typename OutputImageType::RegionType RegionType
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Base class for all process objects that output image data.
typename OutputImageType::PixelType PixelType
BSpline kernel used for density estimation and nonparameteric regression.
BSpline kernel used for density estimation and nonparameteric regression.
typename OutputImageType::RegionType OutputImageRegionType
RealImageType::IndexType NumberToIndex(unsigned int number, typename RealImageType::SizeType size)
TOutputImage OutputImageType
Represent a n-dimensional size (bounds) of a n-dimensional image.
Definition: itkSize.h:68
A superclass of the N-dimensional mesh structure; supports point (geometric coordinate and attribute)...
Definition: itkPointSet.h:84
typename Superclass::SizeType SizeType
Definition: itkImage.h:128
typename Superclass::IndexType IndexType
Definition: itkImage.h:121
typename MeshTraits::PointDataContainer PointDataContainer
Definition: itkPointSet.h:110
typename OutputImageType::SpacingType SpacingType
static constexpr double e
The base of the natural logarithm or Euler&#39;s number
Definition: itkMath.h:53
Base class for filters that take an image as input and produce an image as output.
Control indentation during Print() invocation.
Definition: itkIndent.h:49
typename OutputImageType::PointType OriginType
Templated n-dimensional image class.
Definition: itkImage.h:75
typename OutputImageType::DirectionType DirectionType