ITK  5.0.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>
44 class ITK_TEMPLATE_EXPORT DisplacementFieldToBSplineImageFilter
45  : public ImageToImageFilter<TInputImage, TOutputImage>
46 {
47 public:
48  ITK_DISALLOW_COPY_AND_ASSIGN(DisplacementFieldToBSplineImageFilter);
49 
54 
56  itkNewMacro( Self );
57 
59  static constexpr unsigned int ImageDimension = TInputImage::ImageDimension;
60 
61  using InputFieldType = TInputImage;
62  using InputPointSetType = TInputPointSet;
63  using OutputFieldType = TOutputImage;
64 
68 
70  using PixelType = typename OutputFieldType::PixelType;
71  using VectorType = typename OutputFieldType::PixelType;
74 
75  using SpacingType = typename OutputFieldType::SpacingType;
79 
80  using RealType = typename VectorType::RealValueType;
82 
85  using PointDataType = typename InputPointSetType::PixelType;
86  using PointsContainerType = typename InputPointSetType::PointsContainer;
87  using PointDataContainerType = typename InputPointSetType::PointDataContainer;
88 
95 
97  void SetDisplacementField( const InputFieldType * field )
98  {
99  this->SetInput( 0, field );
100  }
101 
104  {
105  return this->GetInput( 0 );
106  }
107 
113  void SetConfidenceImage( const RealImageType *image )
114  {
115  this->SetNthInput( 1, const_cast<RealImageType *>( image ) );
116  }
117  void SetInput1( const RealImageType *image ) { this->SetConfidenceImage( image ); }
119 
122  {
123  return static_cast<const RealImageType*>( this->ProcessObject::GetInput( 1 ) );
124  }
125 
127  void SetPointSet( const InputPointSetType * points )
128  {
129  this->SetNthInput( 2, const_cast<InputPointSetType *>( points ) );
130  }
131  void SetInput2( const InputPointSetType * points ) { this->SetPointSet( points ); }
133 
136  {
137  return static_cast<const InputPointSetType *>( this->ProcessObject::GetInput( 2 ) );
138  }
139 
141  void SetPointSetConfidenceWeights( WeightsContainerType *weights );
142 
145  {
146  return static_cast<const DisplacementFieldControlPointLatticeType*>( this->GetOutput( 1 ) );
147  }
148 
150  void SetBSplineDomainFromImage( RealImageType * );
151 
154  { this->SetBSplineDomainFromImage( const_cast<RealImageType *>( image ) ); }
155 
157  void SetBSplineDomainFromImage( InputFieldType * );
158 
161  { this->SetBSplineDomainFromImage( const_cast<InputFieldType *>( field ) ); }
162 
164  void SetBSplineDomain( OriginType, SpacingType, SizeType, DirectionType );
165 
166  /* Set/Get b-spline domain origin. */
167  itkGetConstMacro( BSplineDomainOrigin, OriginType );
168 
169  /* Set/Get b-spline domain spacing. */
170  itkGetConstMacro( BSplineDomainSpacing, SpacingType );
171 
172  /* Set/Get b-spline domain size. */
173  itkGetConstMacro( BSplineDomainSize, SizeType );
174 
175  /* Set/Get b-spline domain direction. */
176  itkGetConstMacro( BSplineDomainDirection, DirectionType );
177 
178  /* Use input field to define the B-spline doain. */
179  itkSetMacro( UseInputFieldToDefineTheBSplineDomain, bool );
180  itkGetConstMacro( UseInputFieldToDefineTheBSplineDomain, bool )
181  itkBooleanMacro( UseInputFieldToDefineTheBSplineDomain );
182 
186  itkSetMacro( SplineOrder, unsigned int );
187 
191  itkGetConstMacro( SplineOrder, unsigned int );
192 
200  itkSetMacro( NumberOfControlPoints, ArrayType );
201 
209  itkGetConstMacro( NumberOfControlPoints, ArrayType );
210 
217  itkSetMacro( NumberOfFittingLevels, ArrayType );
218 
225  void SetNumberOfFittingLevels( unsigned int n )
226  {
227  ArrayType nlevels;
228 
229  nlevels.Fill( n );
230  this->SetNumberOfFittingLevels( nlevels );
231  }
232 
239  itkGetConstMacro( NumberOfFittingLevels, ArrayType );
240 
244  itkBooleanMacro( EstimateInverse );
245  itkSetMacro( EstimateInverse, bool );
246  itkGetConstMacro( EstimateInverse, bool );
248 
252  itkBooleanMacro( EnforceStationaryBoundary );
253  itkSetMacro( EnforceStationaryBoundary, bool );
254  itkGetConstMacro( EnforceStationaryBoundary, bool );
256 
257 protected:
258 
261 
263  ~DisplacementFieldToBSplineImageFilter() override = default;
264 
266  void PrintSelf( std::ostream& os, Indent indent ) const override;
267 
269  void GenerateData() override;
270 
271 private:
272  bool m_EstimateInverse{ false };
273  bool m_EnforceStationaryBoundary{ true };
274  unsigned int m_SplineOrder{ 3 };
277 
278  typename WeightsContainerType::Pointer m_PointWeights;
279  bool m_UsePointWeights{ false };
280 
285 
286  bool m_BSplineDomainIsDefined{ true };
287  bool m_UseInputFieldToDefineTheBSplineDomain{ false };
288 };
289 
290 } // end namespace itk
291 
292 #ifndef ITK_MANUAL_INSTANTIATION
293 #include "itkDisplacementFieldToBSplineImageFilter.hxx"
294 #endif
295 
296 #endif
typename InputPointSetType::PointDataContainer PointDataContainerType
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
const DisplacementFieldControlPointLatticeType * GetDisplacementFieldControlPointLattice() const
Base class for all process objects that output image data.
typename BSplineFilterType::WeightsContainerType WeightsContainerType
Image filter which provides a B-spline output approximation.
DataObject * GetInput(const DataObjectIdentifierType &key)
Return an input.
typename BSplineFilterType::PointDataImageType DisplacementFieldControlPointLatticeType
typename InputPointSetType::PointsContainer PointsContainerType
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.
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...
Templated n-dimensional image class.
Definition: itkImage.h:75