ITK  5.1.0
Insight Toolkit
itkBSplineScatteredDataPointSetToImageFilter.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 itkBSplineScatteredDataPointSetToImageFilter_h
19 #define itkBSplineScatteredDataPointSetToImageFilter_h
20 
22 
25 #include "itkVectorContainer.h"
26 
27 #include "vnl/vnl_matrix.h"
28 
29 namespace itk
30 {
132 template <typename TInputPointSet, typename TOutputImage>
134  : public PointSetToImageFilter<TInputPointSet, TOutputImage>
135 {
136 public:
137  ITK_DISALLOW_COPY_AND_ASSIGN(BSplineScatteredDataPointSetToImageFilter);
138 
144 
146  itkNewMacro(Self);
147 
150 
152  static constexpr unsigned int ImageDimension = TOutputImage::ImageDimension;
153 
154  using ImageType = TOutputImage;
155  using PointSetType = TInputPointSet;
156 
158  using PixelType = typename ImageType::PixelType;
161  using SizeType = typename ImageType::SizeType;
162  using IndexType = typename ImageType::IndexType;
163 
166  using PointSetPointer = typename PointSetType::Pointer;
167  using PointDataType = typename PointSetType::PixelType;
168  using PointDataContainerType = typename PointSetType::PointDataContainer;
169 
171  using RealType = float;
173 
181 
188 
189 
194  void
195  SetSplineOrder(unsigned int);
196 
200  void
201  SetSplineOrder(const ArrayType &);
202 
206  itkGetConstReferenceMacro(SplineOrder, ArrayType);
207 
212  itkSetMacro(NumberOfControlPoints, ArrayType);
213  itkGetConstReferenceMacro(NumberOfControlPoints, ArrayType);
215 
220  itkGetConstReferenceMacro(CurrentNumberOfControlPoints, ArrayType);
221 
226  void
227  SetNumberOfLevels(unsigned int);
228 
233  void
234  SetNumberOfLevels(const ArrayType &);
235 
240  itkGetConstReferenceMacro(NumberOfLevels, ArrayType);
241 
248  itkSetMacro(BSplineEpsilon, RealType);
249  itkGetConstMacro(BSplineEpsilon, RealType);
251 
268  itkSetMacro(CloseDimension, ArrayType);
269  itkGetConstReferenceMacro(CloseDimension, ArrayType);
271 
274  void
275  SetPointWeights(WeightsContainerType * weights);
276 
280  itkSetMacro(GenerateOutputImage, bool);
281  itkGetConstReferenceMacro(GenerateOutputImage, bool);
282  itkBooleanMacro(GenerateOutputImage);
284 
288  {
289  return static_cast<PointDataImageType *>(this->ProcessObject::GetOutput(1));
290  }
291 
292 protected:
294  ~BSplineScatteredDataPointSetToImageFilter() override = default;
295 
296  void
297  PrintSelf(std::ostream & os, Indent indent) const override;
298 
299  void
300  ThreadedGenerateData(const RegionType &, ThreadIdType) override;
301 
302  void
304  {
305  itkExceptionMacro("This class requires threadId so it must use classic multi-threading model");
306  }
307 
308  void
309  BeforeThreadedGenerateData() override;
310 
311  void
312  AfterThreadedGenerateData() override;
313 
314  unsigned int
315  SplitRequestedRegion(unsigned int, unsigned int, RegionType &) override;
316 
317  void
318  GenerateData() override;
319 
320 private:
323  void
324  RefineControlPointLattice();
325 
327  void
328  UpdatePointSet();
329 
332  void
333  GenerateOutputImage();
334 
336  void
337  ThreadedGenerateDataForFitting(const RegionType &, ThreadIdType);
338 
340  void
341  ThreadedGenerateDataForReconstruction(const RegionType &, ThreadIdType);
342 
345  void
346  CollapsePhiLattice(PointDataImageType *, PointDataImageType *, const RealType, const unsigned int);
347 
350  void
351  SetPhiLatticeParametricDomainParameters();
352 
355  IndexType
356  NumberToIndex(const unsigned int, const SizeType);
357 
358  bool m_DoMultilevel{ false };
359  bool m_GenerateOutputImage{ true };
360  bool m_UsePointWeights{ false };
361  unsigned int m_MaximumNumberOfLevels{ 1 };
362  unsigned int m_CurrentLevel{ 0 };
368 
370 
373 
374  vnl_matrix<RealType> m_RefinedLatticeCoefficients[ImageDimension];
375 
376  typename PointDataContainerType::Pointer m_InputPointData;
377  typename PointDataContainerType::Pointer m_OutputPointData;
378 
379  typename KernelType::Pointer m_Kernel[ImageDimension];
380 
385 
386  std::vector<RealImagePointer> m_OmegaLatticePerThread;
387  std::vector<PointDataImagePointer> m_DeltaLatticePerThread;
388 
389  RealType m_BSplineEpsilon{ static_cast<RealType>(1e-3) };
390  bool m_IsFittingComplete{ false };
391 };
392 } // end namespace itk
393 
394 #ifndef ITK_MANUAL_INSTANTIATION
395 # include "itkBSplineScatteredDataPointSetToImageFilter.hxx"
396 #endif
397 
398 #endif
itk::PointSetToImageFilter::SizeType
typename TOutputImage::SizeType SizeType
Definition: itkPointSetToImageFilter.h:45
itk::BSplineKernelFunction
BSpline kernel used for density estimation and nonparametric regression.
Definition: itkBSplineKernelFunction.h:43
itk::BSplineScatteredDataPointSetToImageFilter::m_NumberOfControlPoints
ArrayType m_NumberOfControlPoints
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:363
itk::PointSetToImageFilter::PointType
typename TOutputImage::PointType PointType
Definition: itkPointSetToImageFilter.h:71
itk::BSplineScatteredDataPointSetToImageFilter::ImageType
TOutputImage ImageType
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:154
itk::BSplineScatteredDataPointSetToImageFilter::m_SplineOrder
ArrayType m_SplineOrder
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:366
itk::PointSetToImageFilter
Base class for filters that take a PointSet as input and produce an image as output....
Definition: itkPointSetToImageFilter.h:35
itk::BSplineScatteredDataPointSetToImageFilter::GetPhiLattice
PointDataImagePointer GetPhiLattice()
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:287
itk::BSplineScatteredDataPointSetToImageFilter::m_CurrentNumberOfControlPoints
ArrayType m_CurrentNumberOfControlPoints
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:364
itk::BSplineScatteredDataPointSetToImageFilter::PointDataType
typename PointSetType::PixelType PointDataType
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:167
itk::BSplineScatteredDataPointSetToImageFilter::PointSetType
TInputPointSet PointSetType
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:155
itk::GTest::TypedefsAndConstructors::Dimension2::PointType
ImageBaseType::PointType PointType
Definition: itkGTestTypedefsAndConstructors.h:51
itk::BSplineScatteredDataPointSetToImageFilter
Image filter which provides a B-spline output approximation.
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:133
itk::BSplineScatteredDataPointSetToImageFilter::m_OutputPointData
PointDataContainerType::Pointer m_OutputPointData
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:377
itk::BSplineScatteredDataPointSetToImageFilter::m_KernelOrder1
KernelOrder1Type::Pointer m_KernelOrder1
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:382
itk::ProcessObject::GetOutput
DataObject * GetOutput(const DataObjectIdentifierType &key)
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::BSplineScatteredDataPointSetToImageFilter::m_DeltaLatticePerThread
std::vector< PointDataImagePointer > m_DeltaLatticePerThread
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:387
itk::BSplineScatteredDataPointSetToImageFilter::m_PhiLattice
PointDataImageType::Pointer m_PhiLattice
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:371
itk::BSplineScatteredDataPointSetToImageFilter::RealImagePointer
typename RealImageType::Pointer RealImagePointer
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:177
itk::BSplineScatteredDataPointSetToImageFilter::m_PointWeights
WeightsContainerType::Pointer m_PointWeights
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:369
itk::ThreadIdType
unsigned int ThreadIdType
Definition: itkIntTypes.h:99
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itk::ImageSource
Base class for all process objects that output image data.
Definition: itkImageSource.h:67
itk::BSplineScatteredDataPointSetToImageFilter::PointDataContainerType
typename PointSetType::PointDataContainer PointDataContainerType
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:168
itk::BSplineScatteredDataPointSetToImageFilter::DynamicThreadedGenerateData
void DynamicThreadedGenerateData(const RegionType &) override
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:303
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itk::BSplineScatteredDataPointSetToImageFilter::m_InputPointData
PointDataContainerType::Pointer m_InputPointData
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:376
itk::BSplineScatteredDataPointSetToImageFilter::RealType
float RealType
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:171
itk::FixedArray< unsigned, Self::ImageDimension >
itk::BSplineScatteredDataPointSetToImageFilter::m_NumberOfLevels
ArrayType m_NumberOfLevels
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:367
itk::BSplineScatteredDataPointSetToImageFilter::PointDataImagePointer
typename PointDataImageType::Pointer PointDataImagePointer
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:178
itk::BSplineScatteredDataPointSetToImageFilter::PointSetPointer
typename PointSetType::Pointer PointSetPointer
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:166
itk::BSplineScatteredDataPointSetToImageFilter::IndexType
typename ImageType::IndexType IndexType
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:162
itk::BSplineScatteredDataPointSetToImageFilter::m_KernelOrder3
KernelOrder3Type::Pointer m_KernelOrder3
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:384
itk::ImageSource::OutputImageRegionType
typename OutputImageType::RegionType OutputImageRegionType
Definition: itkImageSource.h:92
itk::BSplineScatteredDataPointSetToImageFilter::m_PsiLattice
PointDataImageType::Pointer m_PsiLattice
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:372
itkVectorContainer.h
itkCoxDeBoorBSplineKernelFunction.h
itkBSplineKernelFunction.h
itk::BSplineScatteredDataPointSetToImageFilter::m_KernelOrder2
KernelOrder2Type::Pointer m_KernelOrder2
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:383
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
itk::Math::e
static constexpr double e
Definition: itkMath.h:54
itk::BSplineScatteredDataPointSetToImageFilter::m_CloseDimension
ArrayType m_CloseDimension
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:365
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:86
itk::BSplineScatteredDataPointSetToImageFilter::PixelType
typename ImageType::PixelType PixelType
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:158
itk::CoxDeBoorBSplineKernelFunction
BSpline kernel used for density estimation and nonparametric regression.
Definition: itkCoxDeBoorBSplineKernelFunction.h:58
itk::BSplineScatteredDataPointSetToImageFilter::m_KernelOrder0
KernelOrder0Type::Pointer m_KernelOrder0
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:381
itk::BSplineScatteredDataPointSetToImageFilter::m_OmegaLatticePerThread
std::vector< RealImagePointer > m_OmegaLatticePerThread
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:386
itkPointSetToImageFilter.h
itk::BSplineScatteredDataPointSetToImageFilter::RegionType
typename ImageType::RegionType RegionType
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:159
itk::VectorContainer
Define a front-end to the STL "vector" container that conforms to the IndexedContainerInterface.
Definition: itkVectorContainer.h:48