ITK  5.2.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 {
131 template <typename TInputPointSet, typename TOutputImage>
133  : public PointSetToImageFilter<TInputPointSet, TOutputImage>
134 {
135 public:
136  ITK_DISALLOW_COPY_AND_MOVE(BSplineScatteredDataPointSetToImageFilter);
137 
143 
145  itkNewMacro(Self);
146 
149 
151  static constexpr unsigned int ImageDimension = TOutputImage::ImageDimension;
152 
153  using ImageType = TOutputImage;
154  using PointSetType = TInputPointSet;
155 
157  using PixelType = typename ImageType::PixelType;
160  using SizeType = typename ImageType::SizeType;
161  using IndexType = typename ImageType::IndexType;
162 
165  using PointSetPointer = typename PointSetType::Pointer;
166  using PointDataType = typename PointSetType::PixelType;
167  using PointDataContainerType = typename PointSetType::PointDataContainer;
168 
170  using RealType = float;
172 
180 
187 
188 
193  void
194  SetSplineOrder(unsigned int);
195 
199  void
200  SetSplineOrder(const ArrayType &);
201 
205  itkGetConstReferenceMacro(SplineOrder, ArrayType);
206 
211  itkSetMacro(NumberOfControlPoints, ArrayType);
212  itkGetConstReferenceMacro(NumberOfControlPoints, ArrayType);
214 
219  itkGetConstReferenceMacro(CurrentNumberOfControlPoints, ArrayType);
220 
225  void
226  SetNumberOfLevels(unsigned int);
227 
232  void
233  SetNumberOfLevels(const ArrayType &);
234 
239  itkGetConstReferenceMacro(NumberOfLevels, ArrayType);
240 
247  itkSetMacro(BSplineEpsilon, RealType);
248  itkGetConstMacro(BSplineEpsilon, RealType);
250 
267  itkSetMacro(CloseDimension, ArrayType);
268  itkGetConstReferenceMacro(CloseDimension, ArrayType);
270 
273  void
274  SetPointWeights(WeightsContainerType * weights);
275 
279  itkSetMacro(GenerateOutputImage, bool);
280  itkGetConstReferenceMacro(GenerateOutputImage, bool);
281  itkBooleanMacro(GenerateOutputImage);
283 
287  {
288  return static_cast<PointDataImageType *>(this->ProcessObject::GetOutput(1));
289  }
290 
291 protected:
293  ~BSplineScatteredDataPointSetToImageFilter() override = default;
294 
295  void
296  PrintSelf(std::ostream & os, Indent indent) const override;
297 
298  void
299  ThreadedGenerateData(const RegionType &, ThreadIdType) override;
300 
301  void
303  {
304  itkExceptionMacro("This class requires threadId so it must use classic multi-threading model");
305  }
306 
307  void
308  BeforeThreadedGenerateData() override;
309 
310  void
311  AfterThreadedGenerateData() override;
312 
313  unsigned int
314  SplitRequestedRegion(unsigned int, unsigned int, RegionType &) override;
315 
316  void
317  GenerateData() override;
318 
319 private:
322  void
323  RefineControlPointLattice();
324 
326  void
327  UpdatePointSet();
328 
331  void
332  GenerateOutputImage();
333 
335  void
336  ThreadedGenerateDataForFitting(const RegionType &, ThreadIdType);
337 
339  void
340  ThreadedGenerateDataForReconstruction(const RegionType &, ThreadIdType);
341 
344  void
345  CollapsePhiLattice(PointDataImageType *, PointDataImageType *, const RealType, const unsigned int);
346 
349  void
350  SetPhiLatticeParametricDomainParameters();
351 
354  IndexType
355  NumberToIndex(const unsigned int, const SizeType);
356 
357  bool m_DoMultilevel{ false };
358  bool m_GenerateOutputImage{ true };
359  bool m_UsePointWeights{ false };
360  unsigned int m_MaximumNumberOfLevels{ 1 };
361  unsigned int m_CurrentLevel{ 0 };
367 
369 
372 
373  vnl_matrix<RealType> m_RefinedLatticeCoefficients[ImageDimension];
374 
375  typename PointDataContainerType::Pointer m_InputPointData;
376  typename PointDataContainerType::Pointer m_OutputPointData;
377 
378  typename KernelType::Pointer m_Kernel[ImageDimension];
379 
384 
385  std::vector<RealImagePointer> m_OmegaLatticePerThread;
386  std::vector<PointDataImagePointer> m_DeltaLatticePerThread;
387 
388  RealType m_BSplineEpsilon{ static_cast<RealType>(1e-3) };
389  bool m_IsFittingComplete{ false };
390 };
391 } // end namespace itk
392 
393 #ifndef ITK_MANUAL_INSTANTIATION
394 # include "itkBSplineScatteredDataPointSetToImageFilter.hxx"
395 #endif
396 
397 #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:362
itk::PointSetToImageFilter::PointType
typename TOutputImage::PointType PointType
Definition: itkPointSetToImageFilter.h:71
itk::BSplineScatteredDataPointSetToImageFilter::ImageType
TOutputImage ImageType
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:153
itk::BSplineScatteredDataPointSetToImageFilter::m_SplineOrder
ArrayType m_SplineOrder
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:365
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:286
itk::BSplineScatteredDataPointSetToImageFilter::m_CurrentNumberOfControlPoints
ArrayType m_CurrentNumberOfControlPoints
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:363
itk::BSplineScatteredDataPointSetToImageFilter::PointDataType
typename PointSetType::PixelType PointDataType
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:166
itk::BSplineScatteredDataPointSetToImageFilter::PointSetType
TInputPointSet PointSetType
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:154
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:132
itk::BSplineScatteredDataPointSetToImageFilter::m_OutputPointData
PointDataContainerType::Pointer m_OutputPointData
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:376
itk::BSplineScatteredDataPointSetToImageFilter::m_KernelOrder1
KernelOrder1Type::Pointer m_KernelOrder1
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:381
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:386
itk::BSplineScatteredDataPointSetToImageFilter::m_PhiLattice
PointDataImageType::Pointer m_PhiLattice
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:370
itk::BSplineScatteredDataPointSetToImageFilter::RealImagePointer
typename RealImageType::Pointer RealImagePointer
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:176
itk::BSplineScatteredDataPointSetToImageFilter::m_PointWeights
WeightsContainerType::Pointer m_PointWeights
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:368
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:167
itk::BSplineScatteredDataPointSetToImageFilter::DynamicThreadedGenerateData
void DynamicThreadedGenerateData(const RegionType &) override
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:302
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itk::BSplineScatteredDataPointSetToImageFilter::m_InputPointData
PointDataContainerType::Pointer m_InputPointData
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:375
itk::BSplineScatteredDataPointSetToImageFilter::RealType
float RealType
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:170
itk::FixedArray< unsigned, Self::ImageDimension >
itk::BSplineScatteredDataPointSetToImageFilter::m_NumberOfLevels
ArrayType m_NumberOfLevels
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:366
itk::BSplineScatteredDataPointSetToImageFilter::PointDataImagePointer
typename PointDataImageType::Pointer PointDataImagePointer
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:177
itk::BSplineScatteredDataPointSetToImageFilter::PointSetPointer
typename PointSetType::Pointer PointSetPointer
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:165
itk::BSplineScatteredDataPointSetToImageFilter::IndexType
typename ImageType::IndexType IndexType
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:161
itk::BSplineScatteredDataPointSetToImageFilter::m_KernelOrder3
KernelOrder3Type::Pointer m_KernelOrder3
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:383
itk::ImageSource::OutputImageRegionType
typename OutputImageType::RegionType OutputImageRegionType
Definition: itkImageSource.h:92
itk::BSplineScatteredDataPointSetToImageFilter::m_PsiLattice
PointDataImageType::Pointer m_PsiLattice
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:371
itkVectorContainer.h
itkCoxDeBoorBSplineKernelFunction.h
itkBSplineKernelFunction.h
itk::BSplineScatteredDataPointSetToImageFilter::m_KernelOrder2
KernelOrder2Type::Pointer m_KernelOrder2
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:382
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
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:364
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:86
itk::BSplineScatteredDataPointSetToImageFilter::PixelType
typename ImageType::PixelType PixelType
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:157
itk::CoxDeBoorBSplineKernelFunction
BSpline kernel used for density estimation and nonparametric regression.
Definition: itkCoxDeBoorBSplineKernelFunction.h:57
itk::BSplineScatteredDataPointSetToImageFilter::m_KernelOrder0
KernelOrder0Type::Pointer m_KernelOrder0
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:380
itk::BSplineScatteredDataPointSetToImageFilter::m_OmegaLatticePerThread
std::vector< RealImagePointer > m_OmegaLatticePerThread
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:385
itkPointSetToImageFilter.h
itk::BSplineScatteredDataPointSetToImageFilter::RegionType
typename ImageType::RegionType RegionType
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:158
itk::VectorContainer
Define a front-end to the STL "vector" container that conforms to the IndexedContainerInterface.
Definition: itkVectorContainer.h:48