ITK  4.10.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 
142  itkStaticConstMacro( ImageDimension, unsigned int,
143  TOutputImage::ImageDimension );
144 
145  typedef TOutputImage ImageType;
146  typedef TInputPointSet PointSetType;
147 
149  typedef typename ImageType::PixelType PixelType;
150  typedef typename ImageType::RegionType RegionType;
151  typedef typename ImageType::SizeType SizeType;
152  typedef typename ImageType::IndexType IndexType;
153 
155  typedef typename PointSetType::PointType PointType;
156  typedef typename PointSetType::Pointer PointSetPointer;
157  typedef typename PointSetType::PixelType PointDataType;
158  typedef typename PointSetType::PointDataContainer PointDataContainerType;
159 
161  typedef float RealType;
163 
165  typedef Image<PointDataType,
166  itkGetStaticConstMacro( ImageDimension )> PointDataImageType;
167  typedef Image<RealType,
168  itkGetStaticConstMacro( ImageDimension )> RealImageType;
171  typedef FixedArray<unsigned,
172  itkGetStaticConstMacro( ImageDimension )> ArrayType;
173  typedef FixedArray<RealType,
174  itkGetStaticConstMacro( ImageDimension )> RealArrayType;
176 
185 
186  // Helper functions
187 
194  void SetSplineOrder( unsigned int );
195 
201  void SetSplineOrder( const ArrayType & );
202 
208  itkGetConstReferenceMacro( SplineOrder, ArrayType );
209 
215  itkSetMacro( NumberOfControlPoints, ArrayType );
216 
222  itkGetConstReferenceMacro( NumberOfControlPoints, ArrayType );
223 
229  itkGetConstReferenceMacro( CurrentNumberOfControlPoints, ArrayType );
230 
237  void SetNumberOfLevels( unsigned int );
238 
245  void SetNumberOfLevels( const ArrayType & );
246 
255  itkSetMacro( BSplineEpsilon, RealType );
256  itkGetConstMacro( BSplineEpsilon, RealType );
258 
265  itkGetConstReferenceMacro( NumberOfLevels, ArrayType );
266 
283  itkSetMacro( CloseDimension, ArrayType );
284 
301  itkGetConstReferenceMacro( CloseDimension, ArrayType );
302 
307  void SetPointWeights( WeightsContainerType *weights );
308 
314  itkSetMacro( GenerateOutputImage, bool );
315 
321  itkGetConstReferenceMacro( GenerateOutputImage, bool );
322 
328  itkBooleanMacro( GenerateOutputImage );
329 
333  PointDataImagePointer GetPhiLattice()
334  {
335  return static_cast<PointDataImageType *>( this->ProcessObject::GetOutput( 1 ) );
336  }
337 
338 protected:
341 
342  void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
343 
344  void ThreadedGenerateData( const RegionType &, ThreadIdType ) ITK_OVERRIDE;
345 
346  void BeforeThreadedGenerateData() ITK_OVERRIDE;
347 
348  void AfterThreadedGenerateData() ITK_OVERRIDE;
349 
350  unsigned int SplitRequestedRegion( unsigned int, unsigned int, RegionType & ) ITK_OVERRIDE;
351 
352  void GenerateData() ITK_OVERRIDE;
353 
354 private:
355 
356  BSplineScatteredDataPointSetToImageFilter( const Self & ) ITK_DELETE_FUNCTION;
357  void operator=( const Self & ) ITK_DELETE_FUNCTION;
358 
364 
368  void UpdatePointSet();
369 
374  void GenerateOutputImage();
375 
379  void ThreadedGenerateDataForFitting( const RegionType &, ThreadIdType );
380 
384  void ThreadedGenerateDataForReconstruction( const RegionType &, ThreadIdType );
385 
391  const RealType, const unsigned int );
392 
398 
403  IndexType NumberToIndex( const unsigned int, const SizeType );
404 
409  unsigned int m_CurrentLevel;
415 
417 
420 
422 
423  typename PointDataContainerType::Pointer m_InputPointData;
424  typename PointDataContainerType::Pointer m_OutputPointData;
425 
427 
432 
433  std::vector<RealImagePointer> m_OmegaLatticePerThread;
434  std::vector<PointDataImagePointer> m_DeltaLatticePerThread;
435 
438 };
439 } // end namespace itk
440 
441 #ifndef ITK_MANUAL_INSTANTIATION
442 #include "itkBSplineScatteredDataPointSetToImageFilter.hxx"
443 #endif
444 
445 #endif
void CollapsePhiLattice(PointDataImageType *, PointDataImageType *, const RealType, const unsigned int)
Image< RealType, itkGetStaticConstMacro(ImageDimension)> RealImageType
void ThreadedGenerateDataForFitting(const RegionType &, ThreadIdType)
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes...
Definition: itkArray.h:30
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
void ThreadedGenerateData(const RegionType &, ThreadIdType) override
IndexType NumberToIndex(const unsigned int, const SizeType)
Image< PointDataType, itkGetStaticConstMacro(ImageDimension)> PointDataImageType
BSpline kernel used for density estimation and nonparameteric regression.
FixedArray< unsigned, itkGetStaticConstMacro(ImageDimension)> ArrayType
unsigned int SplitRequestedRegion(unsigned int, unsigned int, RegionType &) override
Image filter which provides a B-spline output approximation.
unsigned int ThreadIdType
Definition: itkIntTypes.h:159
void PrintSelf(std::ostream &os, Indent indent) const override
void ThreadedGenerateDataForReconstruction(const RegionType &, ThreadIdType)
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 "vector" container that conforms to the IndexedContainerInterface.
Control indentation during Print() invocation.
Definition: itkIndent.h:49
void SetPointWeights(WeightsContainerType *weights)
PointSetToImageFilter< TInputPointSet, TOutputImage > Superclass
Templated n-dimensional image class.
Definition: itkImage.h:75
void operator=(const Self &) ITK_DELETE_FUNCTION
DataObject * GetOutput(const DataObjectIdentifierType &key)