ITK  4.11.0
Insight Segmentation and Registration Toolkit
itkBSplineScatteredDataPointSetToImageFilter.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 itkBSplineScatteredDataPointSetToImageFilter_h
19 #define itkBSplineScatteredDataPointSetToImageFilter_h
20 
22 
25 #include "itkVectorContainer.h"
26 
27 #include "vnl/vnl_matrix.h"
28 
29 namespace itk
30 {
128 template< typename TInputPointSet, typename TOutputImage >
130  public PointSetToImageFilter< TInputPointSet, TOutputImage >
131 {
132 public:
137 
139  itkNewMacro( Self );
140 
144 
146  itkStaticConstMacro( ImageDimension, unsigned int,
147  TOutputImage::ImageDimension );
148 
149  typedef TOutputImage ImageType;
150  typedef TInputPointSet PointSetType;
151 
153  typedef typename ImageType::PixelType PixelType;
154  typedef typename ImageType::RegionType RegionType;
155  typedef typename ImageType::SizeType SizeType;
156  typedef typename ImageType::IndexType IndexType;
157 
159  typedef typename PointSetType::PointType PointType;
160  typedef typename PointSetType::Pointer PointSetPointer;
161  typedef typename PointSetType::PixelType PointDataType;
162  typedef typename PointSetType::PointDataContainer PointDataContainerType;
163 
165  typedef float RealType;
167 
169  typedef Image<PointDataType,
170  itkGetStaticConstMacro( ImageDimension )> PointDataImageType;
171  typedef Image<RealType,
172  itkGetStaticConstMacro( ImageDimension )> RealImageType;
175  typedef FixedArray<unsigned,
176  itkGetStaticConstMacro( ImageDimension )> ArrayType;
177  typedef FixedArray<RealType,
178  itkGetStaticConstMacro( ImageDimension )> RealArrayType;
180 
189 
190  // Helper functions
191 
198  void SetSplineOrder( unsigned int );
199 
205  void SetSplineOrder( const ArrayType & );
206 
212  itkGetConstReferenceMacro( SplineOrder, ArrayType );
213 
219  itkSetMacro( NumberOfControlPoints, ArrayType );
220 
226  itkGetConstReferenceMacro( NumberOfControlPoints, ArrayType );
227 
233  itkGetConstReferenceMacro( CurrentNumberOfControlPoints, ArrayType );
234 
241  void SetNumberOfLevels( unsigned int );
242 
249  void SetNumberOfLevels( const ArrayType & );
250 
259  itkSetMacro( BSplineEpsilon, RealType );
260  itkGetConstMacro( BSplineEpsilon, RealType );
262 
269  itkGetConstReferenceMacro( NumberOfLevels, ArrayType );
270 
287  itkSetMacro( CloseDimension, ArrayType );
288 
305  itkGetConstReferenceMacro( CloseDimension, ArrayType );
306 
311  void SetPointWeights( WeightsContainerType *weights );
312 
318  itkSetMacro( GenerateOutputImage, bool );
319 
325  itkGetConstReferenceMacro( GenerateOutputImage, bool );
326 
332  itkBooleanMacro( GenerateOutputImage );
333 
338  {
339  return static_cast<PointDataImageType *>( this->ProcessObject::GetOutput( 1 ) );
340  }
341 
342 protected:
345 
346  void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
347 
348  void ThreadedGenerateData( const RegionType &, ThreadIdType ) ITK_OVERRIDE;
349 
350  void BeforeThreadedGenerateData() ITK_OVERRIDE;
351 
352  void AfterThreadedGenerateData() ITK_OVERRIDE;
353 
354  unsigned int SplitRequestedRegion( unsigned int, unsigned int, RegionType & ) ITK_OVERRIDE;
355 
356  void GenerateData() ITK_OVERRIDE;
357 
358 private:
359 
360  ITK_DISALLOW_COPY_AND_ASSIGN(BSplineScatteredDataPointSetToImageFilter);
361 
366  void RefineControlPointLattice();
367 
371  void UpdatePointSet();
372 
377  void GenerateOutputImage();
378 
382  void ThreadedGenerateDataForFitting( const RegionType &, ThreadIdType );
383 
387  void ThreadedGenerateDataForReconstruction( const RegionType &, ThreadIdType );
388 
393  void CollapsePhiLattice( PointDataImageType *, PointDataImageType *,
394  const RealType, const unsigned int );
395 
400  void SetPhiLatticeParametricDomainParameters();
401 
406  IndexType NumberToIndex( const unsigned int, const SizeType );
407 
408  bool m_DoMultilevel;
409  bool m_GenerateOutputImage;
410  bool m_UsePointWeights;
411  unsigned int m_MaximumNumberOfLevels;
412  unsigned int m_CurrentLevel;
413  ArrayType m_NumberOfControlPoints;
414  ArrayType m_CurrentNumberOfControlPoints;
415  ArrayType m_CloseDimension;
416  ArrayType m_SplineOrder;
417  ArrayType m_NumberOfLevels;
418 
419  typename WeightsContainerType::Pointer m_PointWeights;
420 
421  typename PointDataImageType::Pointer m_PhiLattice;
422  typename PointDataImageType::Pointer m_PsiLattice;
423 
424  vnl_matrix<RealType> m_RefinedLatticeCoefficients[ImageDimension];
425 
426  typename PointDataContainerType::Pointer m_InputPointData;
427  typename PointDataContainerType::Pointer m_OutputPointData;
428 
429  typename KernelType::Pointer m_Kernel[ImageDimension];
430 
431  typename KernelOrder0Type::Pointer m_KernelOrder0;
432  typename KernelOrder1Type::Pointer m_KernelOrder1;
433  typename KernelOrder2Type::Pointer m_KernelOrder2;
434  typename KernelOrder3Type::Pointer m_KernelOrder3;
435 
436  std::vector<RealImagePointer> m_OmegaLatticePerThread;
437  std::vector<PointDataImagePointer> m_DeltaLatticePerThread;
438 
439  RealType m_BSplineEpsilon;
440  bool m_IsFittingComplete;
441 };
442 } // end namespace itk
443 
444 #ifndef ITK_MANUAL_INSTANTIATION
445 #include "itkBSplineScatteredDataPointSetToImageFilter.hxx"
446 #endif
447 
448 #endif
Image< RealType, itkGetStaticConstMacro(ImageDimension)> RealImageType
Base class for all process objects that output image data.
Simulate a standard C array with copy semnatics.
Definition: itkFixedArray.h:50
BSpline kernel used for density estimation and nonparameteric regression.
FixedArray< RealType, itkGetStaticConstMacro(ImageDimension)> RealArrayType
Image< PointDataType, itkGetStaticConstMacro(ImageDimension)> PointDataImageType
BSpline kernel used for density estimation and nonparameteric regression.
FixedArray< unsigned, itkGetStaticConstMacro(ImageDimension)> ArrayType
Image filter which provides a B-spline output approximation.
unsigned int ThreadIdType
Definition: itkIntTypes.h:159
Base class for filters that take a PointSet as input and produce an image as output. By default, if the user does not specify the size of the output image, the maximum size of the point-set&#39;s bounding box is used.
Define a front-end to the STL &quot;vector&quot; container that conforms to the IndexedContainerInterface.
Control indentation during Print() invocation.
Definition: itkIndent.h:49
PointSetToImageFilter< TInputPointSet, TOutputImage > Superclass
Templated n-dimensional image class.
Definition: itkImage.h:75
DataObject * GetOutput(const DataObjectIdentifierType &key)