ITK  4.8.0
Insight Segmentation and Registration Toolkit
itkDisplacementFieldToBSplineImageFilter.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 itkDisplacementFieldToBSplineImageFilter_h
19 #define itkDisplacementFieldToBSplineImageFilter_h
20 
21 #include "itkImageToImageFilter.h"
22 
24 #include "itkPointSet.h"
25 #include "itkVector.h"
26 
27 namespace itk
28 {
29 
41 template <typename TInputImage,
42  typename TInputPointSet = PointSet<typename TInputImage::PixelType, TInputImage::ImageDimension>,
43  typename TOutputImage = TInputImage>
45  : public ImageToImageFilter<TInputImage, TOutputImage>
46 {
47 public:
52 
54  itkNewMacro( Self );
55 
57  itkStaticConstMacro( ImageDimension, unsigned int, TInputImage::ImageDimension );
58 
59  typedef TInputImage InputFieldType;
60  typedef TInputPointSet InputPointSetType;
61  typedef TOutputImage OutputFieldType;
62 
65  typedef typename InputFieldType::PointType InputFieldPointType;
66 
68  typedef typename OutputFieldType::PixelType PixelType;
69  typedef typename OutputFieldType::PixelType VectorType;
70  typedef typename OutputFieldType::RegionType RegionType;
71  typedef typename OutputFieldType::IndexType IndexType;
72 
73  typedef typename OutputFieldType::SpacingType SpacingType;
74  typedef typename OutputFieldType::PointType OriginType;
75  typedef typename OutputFieldType::SizeType SizeType;
76  typedef typename OutputFieldType::DirectionType DirectionType;
77 
78  typedef typename VectorType::RealValueType RealType;
80 
82  typedef typename InputPointSetType::PointType PointType;
83  typedef typename InputPointSetType::PixelType PointDataType;
84  typedef typename InputPointSetType::PointsContainer PointsContainerType;
85  typedef typename InputPointSetType::PointDataContainer PointDataContainerType;
86 
93 
95  void SetDisplacementField( const InputFieldType * field )
96  {
97  this->SetInput( 0, field );
98  }
99 
102  {
103  return this->GetInput( 0 );
104  }
105 
111  void SetConfidenceImage( const RealImageType *image )
112  {
113  this->SetNthInput( 1, const_cast<RealImageType *>( image ) );
114  }
115  void SetInput1( const RealImageType *image ) { this->SetConfidenceImage( image ); }
117 
120  {
121  return static_cast<const RealImageType*>( this->ProcessObject::GetInput( 1 ) );
122  }
123 
125  void SetPointSet( const InputPointSetType * points )
126  {
127  this->SetNthInput( 2, const_cast<InputPointSetType *>( points ) );
128  }
129  void SetInput2( const InputPointSetType * points ) { this->SetPointSet( points ); }
131 
134  {
135  return static_cast<const InputPointSetType *>( this->ProcessObject::GetInput( 2 ) );
136  }
137 
140 
143  {
144  return static_cast<const DisplacementFieldControlPointLatticeType*>( this->GetOutput( 1 ) );
145  }
146 
149 
152  { this->SetBSplineDomainFromImage( const_cast<RealImageType *>( image ) ); }
153 
156 
159  { this->SetBSplineDomainFromImage( const_cast<InputFieldType *>( field ) ); }
160 
163 
164  /* Set/Get b-spline domain origin. */
165  itkGetConstMacro( BSplineDomainOrigin, OriginType );
166 
167  /* Set/Get b-spline domain spacing. */
168  itkGetConstMacro( BSplineDomainSpacing, SpacingType );
169 
170  /* Set/Get b-spline domain size. */
171  itkGetConstMacro( BSplineDomainSize, SizeType );
172 
173  /* Set/Get b-spline domain direction. */
174  itkGetConstMacro( BSplineDomainDirection, DirectionType );
175 
176  /* Use input field to define the B-spline doain. */
177  itkSetMacro( UseInputFieldToDefineTheBSplineDomain, bool );
178  itkGetConstMacro( UseInputFieldToDefineTheBSplineDomain, bool )
179  itkBooleanMacro( UseInputFieldToDefineTheBSplineDomain );
180 
184  itkSetMacro( SplineOrder, unsigned int );
185 
189  itkGetConstMacro( SplineOrder, unsigned int );
190 
198  itkSetMacro( NumberOfControlPoints, ArrayType );
199 
207  itkGetConstMacro( NumberOfControlPoints, ArrayType );
208 
215  itkSetMacro( NumberOfFittingLevels, ArrayType );
216 
223  void SetNumberOfFittingLevels( unsigned int n )
224  {
225  ArrayType nlevels;
226 
227  nlevels.Fill( n );
228  this->SetNumberOfFittingLevels( nlevels );
229  }
230 
237  itkGetConstMacro( NumberOfFittingLevels, ArrayType );
238 
242  itkBooleanMacro( EstimateInverse );
243  itkSetMacro( EstimateInverse, bool );
244  itkGetConstMacro( EstimateInverse, bool );
246 
250  itkBooleanMacro( EnforceStationaryBoundary );
251  itkSetMacro( EnforceStationaryBoundary, bool );
252  itkGetConstMacro( EnforceStationaryBoundary, bool );
254 
255 protected:
256 
259 
262 
264  void PrintSelf( std::ostream& os, Indent indent ) const ITK_OVERRIDE;
265 
267  void GenerateData() ITK_OVERRIDE;
268 
269 private:
270  DisplacementFieldToBSplineImageFilter( const Self& ); //purposely not implemented
271  void operator=( const Self& ); //purposely not implemented
272 
275  unsigned int m_SplineOrder;
278 
281 
286 
289 };
290 
291 } // end namespace itk
292 
293 #ifndef ITK_MANUAL_INSTANTIATION
294 #include "itkDisplacementFieldToBSplineImageFilter.hxx"
295 #endif
296 
297 #endif
void PrintSelf(std::ostream &os, Indent indent) const override
const DisplacementFieldControlPointLatticeType * GetDisplacementFieldControlPointLattice() const
virtual void SetNumberOfFittingLevels(ArrayType _arg)
Base class for all process objects that output image data.
void SetBSplineDomain(OriginType, SpacingType, SizeType, DirectionType)
void Fill(const ValueType &)
void SetPointSetConfidenceWeights(WeightsContainerType *weights)
BSplineScatteredDataPointSetToImageFilter< InputPointSetType, OutputFieldType > BSplineFilterType
BSplineFilterType::PointDataImageType DisplacementFieldControlPointLatticeType
virtual void SetInput(const InputImageType *image)
const InputImageType * GetInput() const
Image filter which provides a B-spline output approximation.
DataObject * GetInput(const DataObjectIdentifierType &key)
Return an input.
Base class for filters that take an image as input and produce an image as output.
Define a front-end to the STL &quot;vector&quot; container that conforms to the IndexedContainerInterface.
OutputImageType * GetOutput()
Control indentation during Print() invocation.
Definition: itkIndent.h:49
Class which takes a dense displacement field image and/or a set of points with associated displacemen...
virtual void SetNthInput(DataObjectPointerArraySizeType num, DataObject *input)
Templated n-dimensional image class.
Definition: itkImage.h:75
ImageToImageFilter< TInputImage, TOutputImage > Superclass
void SetBSplineDomainFromImage(RealImageType *)