ITK  4.13.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:
138 
140  itkNewMacro( Self );
141 
145 
147  itkStaticConstMacro( ImageDimension, unsigned int,
148  TOutputImage::ImageDimension );
149 
150  typedef TOutputImage ImageType;
151  typedef TInputPointSet PointSetType;
152 
154  typedef typename ImageType::PixelType PixelType;
155  typedef typename ImageType::RegionType RegionType;
156  typedef typename ImageType::SizeType SizeType;
157  typedef typename ImageType::IndexType IndexType;
158 
161  typedef typename PointSetType::Pointer PointSetPointer;
162  typedef typename PointSetType::PixelType PointDataType;
163  typedef typename PointSetType::PointDataContainer PointDataContainerType;
164 
166  typedef float RealType;
168 
170  typedef Image<PointDataType,
171  itkGetStaticConstMacro( ImageDimension )> PointDataImageType;
172  typedef Image<RealType,
173  itkGetStaticConstMacro( ImageDimension )> RealImageType;
176  typedef FixedArray<unsigned,
177  itkGetStaticConstMacro( ImageDimension )> ArrayType;
178  typedef FixedArray<RealType,
179  itkGetStaticConstMacro( ImageDimension )> RealArrayType;
181 
188 
189 
194  void SetSplineOrder( unsigned int );
195 
199  void SetSplineOrder( const ArrayType & );
200 
204  itkGetConstReferenceMacro( SplineOrder, ArrayType );
205 
210  itkSetMacro( NumberOfControlPoints, ArrayType );
211  itkGetConstReferenceMacro( NumberOfControlPoints, ArrayType );
213 
218  itkGetConstReferenceMacro( CurrentNumberOfControlPoints, ArrayType );
219 
224  void SetNumberOfLevels( unsigned int );
225 
230  void SetNumberOfLevels( const ArrayType & );
231 
236  itkGetConstReferenceMacro( NumberOfLevels, ArrayType );
237 
244  itkSetMacro( BSplineEpsilon, RealType );
245  itkGetConstMacro( BSplineEpsilon, RealType );
247 
264  itkSetMacro( CloseDimension, ArrayType );
265  itkGetConstReferenceMacro( CloseDimension, ArrayType );
267 
270  void SetPointWeights( WeightsContainerType *weights );
271 
275  itkSetMacro( GenerateOutputImage, bool );
276  itkGetConstReferenceMacro( GenerateOutputImage, bool );
277  itkBooleanMacro( GenerateOutputImage );
279 
282  {
283  return static_cast<PointDataImageType *>( this->ProcessObject::GetOutput( 1 ) );
284  }
285 
286 protected:
288  virtual ~BSplineScatteredDataPointSetToImageFilter() ITK_OVERRIDE;
289 
290  void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
291 
292  void ThreadedGenerateData( const RegionType &, ThreadIdType ) ITK_OVERRIDE;
293 
294  void BeforeThreadedGenerateData() ITK_OVERRIDE;
295 
296  void AfterThreadedGenerateData() ITK_OVERRIDE;
297 
298  unsigned int SplitRequestedRegion( unsigned int, unsigned int, RegionType & ) ITK_OVERRIDE;
299 
300  void GenerateData() ITK_OVERRIDE;
301 
302 private:
303 
304  ITK_DISALLOW_COPY_AND_ASSIGN(BSplineScatteredDataPointSetToImageFilter);
305 
308  void RefineControlPointLattice();
309 
311  void UpdatePointSet();
312 
315  void GenerateOutputImage();
316 
318  void ThreadedGenerateDataForFitting( const RegionType &, ThreadIdType );
319 
321  void ThreadedGenerateDataForReconstruction( const RegionType &, ThreadIdType );
322 
325  void CollapsePhiLattice( PointDataImageType *, PointDataImageType *,
326  const RealType, const unsigned int );
327 
330  void SetPhiLatticeParametricDomainParameters();
331 
334  IndexType NumberToIndex( const unsigned int, const SizeType );
335 
336  bool m_DoMultilevel;
337  bool m_GenerateOutputImage;
338  bool m_UsePointWeights;
339  unsigned int m_MaximumNumberOfLevels;
340  unsigned int m_CurrentLevel;
341  ArrayType m_NumberOfControlPoints;
342  ArrayType m_CurrentNumberOfControlPoints;
343  ArrayType m_CloseDimension;
344  ArrayType m_SplineOrder;
345  ArrayType m_NumberOfLevels;
346 
347  typename WeightsContainerType::Pointer m_PointWeights;
348 
349  typename PointDataImageType::Pointer m_PhiLattice;
350  typename PointDataImageType::Pointer m_PsiLattice;
351 
352  vnl_matrix<RealType> m_RefinedLatticeCoefficients[ImageDimension];
353 
354  typename PointDataContainerType::Pointer m_InputPointData;
355  typename PointDataContainerType::Pointer m_OutputPointData;
356 
357  typename KernelType::Pointer m_Kernel[ImageDimension];
358 
359  typename KernelOrder0Type::Pointer m_KernelOrder0;
360  typename KernelOrder1Type::Pointer m_KernelOrder1;
361  typename KernelOrder2Type::Pointer m_KernelOrder2;
362  typename KernelOrder3Type::Pointer m_KernelOrder3;
363 
364  std::vector<RealImagePointer> m_OmegaLatticePerThread;
365  std::vector<PointDataImagePointer> m_DeltaLatticePerThread;
366 
367  RealType m_BSplineEpsilon;
368  bool m_IsFittingComplete;
369 };
370 } // end namespace itk
371 
372 #ifndef ITK_MANUAL_INSTANTIATION
373 #include "itkBSplineScatteredDataPointSetToImageFilter.hxx"
374 #endif
375 
376 #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)