ITK  5.1.0
Insight Toolkit
itkBSplineControlPointImageFilter.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright NumFOCUS
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 : public ImageToImageFilter<TInputImage, TOutputImage>
61 {
62 public:
63  ITK_DISALLOW_COPY_AND_ASSIGN(BSplineControlPointImageFilter);
64 
69 
71  itkNewMacro(Self);
72 
75 
77  static constexpr unsigned int ImageDimension = TInputImage::ImageDimension;
78 
79  using ControlPointLatticeType = TInputImage;
80  using OutputImageType = TOutputImage;
81 
83  using PixelType = typename OutputImageType::PixelType;
88 
89  using SpacingType = typename OutputImageType::SpacingType;
93 
95  using RealType = float;
98 
101 
108 
115 
120  void
121  SetSplineOrder(unsigned int);
122 
127  void SetSplineOrder(ArrayType);
128 
132  itkGetConstReferenceMacro(SplineOrder, ArrayType);
133 
150  itkSetMacro(CloseDimension, ArrayType);
151  itkGetConstReferenceMacro(CloseDimension, ArrayType);
153 
157  itkSetMacro(Spacing, SpacingType);
158  itkGetConstMacro(Spacing, SpacingType);
160 
164  itkSetMacro(Origin, OriginType);
165  itkGetConstMacro(Origin, OriginType);
167 
171  itkSetMacro(Size, SizeType);
172  itkGetConstMacro(Size, SizeType);
174 
189  itkSetMacro(Direction, DirectionType);
190 
194  itkGetConstMacro(Direction, DirectionType);
195 
203  typename ControlPointLatticeType::Pointer RefineControlPointLattice(ArrayType);
204 
205 protected:
207  ~BSplineControlPointImageFilter() override = default;
208  void
209  PrintSelf(std::ostream & os, Indent indent) const override;
210 
212  void
213  DynamicThreadedGenerateData(const OutputImageRegionType &) override;
214 
215 
216 private:
221  void
222  BeforeThreadedGenerateData() override;
223 
228  unsigned int
229  SplitRequestedRegion(unsigned int, unsigned int, OutputImageRegionType &) override;
230 
235  void
236  CollapsePhiLattice(PointDataImageType *, PointDataImageType *, const RealType, const unsigned int);
237 
241  void SetNumberOfLevels(ArrayType);
242 
248 
249  bool m_DoMultilevel{ false };
250  unsigned int m_MaximumNumberOfLevels{ 1 };
255 
256  vnl_matrix<RealType> m_RefinedLatticeCoefficients[ImageDimension];
257 
258  typename KernelType::Pointer m_Kernel[ImageDimension];
263 
264  RealType m_BSplineEpsilon{ static_cast<RealType>(1e-3) };
265 
266  inline typename RealImageType::IndexType
267  NumberToIndex(unsigned int number, typename RealImageType::SizeType size)
268  {
269  typename RealImageType::IndexType k;
270  k[0] = 1;
271 
272  for (unsigned int i = 1; i < ImageDimension; i++)
273  {
274  k[i] = size[ImageDimension - i - 1] * k[i - 1];
275  }
276  typename RealImageType::IndexType index;
277  for (unsigned int i = 0; i < ImageDimension; i++)
278  {
279  index[ImageDimension - i - 1] = static_cast<unsigned int>(number / k[ImageDimension - i - 1]);
280  number %= k[ImageDimension - i - 1];
281  }
282  return index;
283  }
284 };
285 
286 } // end namespace itk
287 
288 #ifndef ITK_MANUAL_INSTANTIATION
289 # include "itkBSplineControlPointImageFilter.hxx"
290 #endif
291 
292 #endif
itk::BSplineControlPointImageFilter::RealImagePointer
typename RealImageType::Pointer RealImagePointer
Definition: itkBSplineControlPointImageFilter.h:97
itk::BSplineControlPointImageFilter::m_KernelOrder1
KernelOrder1Type::Pointer m_KernelOrder1
Definition: itkBSplineControlPointImageFilter.h:260
itk::BSplineControlPointImageFilter::m_KernelOrder0
KernelOrder0Type::Pointer m_KernelOrder0
Definition: itkBSplineControlPointImageFilter.h:259
itk::BSplineKernelFunction
BSpline kernel used for density estimation and nonparametric regression.
Definition: itkBSplineKernelFunction.h:43
itk::GTest::TypedefsAndConstructors::Dimension2::DirectionType
ImageBaseType::DirectionType DirectionType
Definition: itkGTestTypedefsAndConstructors.h:52
itk::PointSet
A superclass of the N-dimensional mesh structure; supports point (geometric coordinate and attribute)...
Definition: itkPointSet.h:82
itkVariableSizeMatrix.h
itk::Size
Represent a n-dimensional size (bounds) of a n-dimensional image.
Definition: itkSize.h:69
itk::BSplineControlPointImageFilter::ControlPointLatticeType
TInputImage ControlPointLatticeType
Definition: itkBSplineControlPointImageFilter.h:79
itk::BSplineControlPointImageFilter::m_NumberOfLevels
ArrayType m_NumberOfLevels
Definition: itkBSplineControlPointImageFilter.h:254
itk::GTest::TypedefsAndConstructors::Dimension2::PointType
ImageBaseType::PointType PointType
Definition: itkGTestTypedefsAndConstructors.h:51
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itk::BSplineControlPointImageFilter::m_Origin
OriginType m_Origin
Definition: itkBSplineControlPointImageFilter.h:246
itk::BSplineControlPointImageFilter::RegionType
typename OutputImageType::RegionType RegionType
Definition: itkBSplineControlPointImageFilter.h:84
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::BSplineControlPointImageFilter::m_CloseDimension
ArrayType m_CloseDimension
Definition: itkBSplineControlPointImageFilter.h:252
itk::BSplineControlPointImageFilter::PointType
typename OutputImageType::PointType PointType
Definition: itkBSplineControlPointImageFilter.h:86
itk::BSplineControlPointImageFilter::PixelType
typename OutputImageType::PixelType PixelType
Definition: itkBSplineControlPointImageFilter.h:83
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itk::ImageToImageFilter
Base class for filters that take an image as input and produce an image as output.
Definition: itkImageToImageFilter.h:108
itk::ImageSource
Base class for all process objects that output image data.
Definition: itkImageSource.h:67
itk::BSplineControlPointImageFilter::OriginType
typename OutputImageType::PointType OriginType
Definition: itkBSplineControlPointImageFilter.h:90
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itk::BSplineControlPointImageFilter::NumberToIndex
RealImageType::IndexType NumberToIndex(unsigned int number, typename RealImageType::SizeType size)
Definition: itkBSplineControlPointImageFilter.h:267
itk::PointSet::PixelType
typename MeshTraits::PixelType PixelType
Definition: itkPointSet.h:101
itkFixedArray.h
itk::BSplineControlPointImageFilter::PointDataContainerType
typename PointSetType::PointDataContainer PointDataContainerType
Definition: itkBSplineControlPointImageFilter.h:105
itkImageToImageFilter.h
itk::BSplineControlPointImageFilter::m_Size
SizeType m_Size
Definition: itkBSplineControlPointImageFilter.h:244
itk::BSplineControlPointImageFilter::m_Spacing
SpacingType m_Spacing
Definition: itkBSplineControlPointImageFilter.h:245
itk::FixedArray< unsigned, Self::ImageDimension >
itk::BSplineControlPointImageFilter::m_KernelOrder2
KernelOrder2Type::Pointer m_KernelOrder2
Definition: itkBSplineControlPointImageFilter.h:261
itk::BSplineControlPointImageFilter::SpacingType
typename OutputImageType::SpacingType SpacingType
Definition: itkBSplineControlPointImageFilter.h:89
itk::BSplineControlPointImageFilter::PointDataImagePointer
typename PointDataImageType::Pointer PointDataImagePointer
Definition: itkBSplineControlPointImageFilter.h:107
itk::BSplineControlPointImageFilter::IndexType
typename OutputImageType::IndexType IndexType
Definition: itkBSplineControlPointImageFilter.h:85
itk::ImageSource::OutputImageRegionType
typename OutputImageType::RegionType OutputImageRegionType
Definition: itkImageSource.h:92
itkVectorContainer.h
itkCoxDeBoorBSplineKernelFunction.h
itk::BSplineControlPointImageFilter::m_Direction
DirectionType m_Direction
Definition: itkBSplineControlPointImageFilter.h:247
itkBSplineKernelFunction.h
itk::BSplineControlPointImageFilter::PointDataType
typename PointSetType::PixelType PointDataType
Definition: itkBSplineControlPointImageFilter.h:104
itk::Image::IndexType
typename Superclass::IndexType IndexType
Definition: itkImage.h:132
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkArray.h:26
itk::ProcessObject
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Definition: itkProcessObject.h:138
itkVector.h
itk::BSplineControlPointImageFilter::m_NumberOfControlPoints
ArrayType m_NumberOfControlPoints
Definition: itkBSplineControlPointImageFilter.h:251
itk::PointSet::PointDataContainer
typename MeshTraits::PointDataContainer PointDataContainer
Definition: itkPointSet.h:108
itk::Math::e
static constexpr double e
Definition: itkMath.h:54
itk::BSplineControlPointImageFilter::SizeType
typename OutputImageType::SizeType SizeType
Definition: itkBSplineControlPointImageFilter.h:91
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:86
itk::BSplineControlPointImageFilter::m_KernelOrder3
KernelOrder3Type::Pointer m_KernelOrder3
Definition: itkBSplineControlPointImageFilter.h:262
itkPointSet.h
itk::BSplineControlPointImageFilter::m_SplineOrder
ArrayType m_SplineOrder
Definition: itkBSplineControlPointImageFilter.h:253
itk::CoxDeBoorBSplineKernelFunction
BSpline kernel used for density estimation and nonparametric regression.
Definition: itkCoxDeBoorBSplineKernelFunction.h:58
itk::BSplineControlPointImageFilter::RealType
float RealType
Definition: itkBSplineControlPointImageFilter.h:95
itk::BSplineControlPointImageFilter::DirectionType
typename OutputImageType::DirectionType DirectionType
Definition: itkBSplineControlPointImageFilter.h:92
itk::BSplineControlPointImageFilter
Process a given a B-spline grid of control points.
Definition: itkBSplineControlPointImageFilter.h:60
itk::Image::SizeType
typename Superclass::SizeType SizeType
Definition: itkImage.h:139
itk::ImageSource::OutputImageType
TOutputImage OutputImageType
Definition: itkImageSource.h:90