ITK  4.6.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;
174 
183 
184  // Helper functions
185 
192  void SetSplineOrder( unsigned int );
193 
199  void SetSplineOrder( const ArrayType & );
200 
206  itkGetConstReferenceMacro( SplineOrder, ArrayType );
207 
213  itkSetMacro( NumberOfControlPoints, ArrayType );
214 
220  itkGetConstReferenceMacro( NumberOfControlPoints, ArrayType );
221 
227  itkGetConstReferenceMacro( CurrentNumberOfControlPoints, ArrayType );
228 
235  void SetNumberOfLevels( unsigned int );
236 
243  void SetNumberOfLevels( const ArrayType & );
244 
251  itkGetConstReferenceMacro( NumberOfLevels, ArrayType );
252 
269  itkSetMacro( CloseDimension, ArrayType );
270 
287  itkGetConstReferenceMacro( CloseDimension, ArrayType );
288 
293  void SetPointWeights( WeightsContainerType *weights );
294 
300  itkSetMacro( GenerateOutputImage, bool );
301 
307  itkGetConstReferenceMacro( GenerateOutputImage, bool );
308 
314  itkBooleanMacro( GenerateOutputImage );
315 
319  itkGetConstMacro( PhiLattice, PointDataImagePointer );
320 
321 protected:
324 
325  void PrintSelf(std::ostream & os, Indent indent) const;
326 
328 
330 
332 
333  unsigned int SplitRequestedRegion( unsigned int, unsigned int, RegionType & );
334 
335  void GenerateData();
336 
337 private:
338 
339  //purposely not implemented
341  void operator=( const Self & );
342 
348 
352  void UpdatePointSet();
353 
358  void GenerateOutputImage();
359 
364 
369 
375  const RealType, const unsigned int );
376 
382 
387  IndexType NumberToIndex( const unsigned int, const SizeType );
388 
393  unsigned int m_CurrentLevel;
399 
401 
404 
406 
407  typename PointDataContainerType::Pointer m_InputPointData;
408  typename PointDataContainerType::Pointer m_OutputPointData;
409 
411 
416 
417  std::vector<RealImagePointer> m_OmegaLatticePerThread;
418  std::vector<PointDataImagePointer> m_DeltaLatticePerThread;
419 
422 };
423 } // end namespace itk
424 
425 #ifndef ITK_MANUAL_INSTANTIATION
426 #include "itkBSplineScatteredDataPointSetToImageFilter.hxx"
427 #endif
428 
429 #endif
void PrintSelf(std::ostream &os, Indent indent) const
void CollapsePhiLattice(PointDataImageType *, PointDataImageType *, const RealType, const unsigned int)
Image< RealType, itkGetStaticConstMacro(ImageDimension)> RealImageType
void ThreadedGenerateDataForFitting(const RegionType &, ThreadIdType)
void ThreadedGenerateData(const RegionType &, ThreadIdType)
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.
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
Image filter which provides a B-spline output approximation.
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.
unsigned int SplitRequestedRegion(unsigned int, unsigned int, RegionType &)
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
void SetPointWeights(WeightsContainerType *weights)
PointSetToImageFilter< TInputPointSet, TOutputImage > Superclass
Templated n-dimensional image class.
Definition: itkImage.h:75
unsigned int ThreadIdType
Definition: itkIntTypes.h:159