ITK  5.4.0
Insight Toolkit
itkDisplacementFieldToBSplineImageFilter.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright NumFOCUS
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  * https://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 : public ImageToImageFilter<TInputImage, TOutputImage>
45 {
46 public:
47  ITK_DISALLOW_COPY_AND_MOVE(DisplacementFieldToBSplineImageFilter);
48 
53 
55  itkNewMacro(Self);
56 
58  itkOverrideGetNameOfClassMacro(DisplacementFieldToBSplineImageFilter);
59 
61  static constexpr unsigned int ImageDimension = TInputImage::ImageDimension;
62 
63  using InputFieldType = TInputImage;
64  using InputPointSetType = TInputPointSet;
65  using OutputFieldType = TOutputImage;
66 
70 
72  using PixelType = typename OutputFieldType::PixelType;
73  using VectorType = typename OutputFieldType::PixelType;
76 
77  using SpacingType = typename OutputFieldType::SpacingType;
81 
82  using RealType = typename VectorType::RealValueType;
84 
87  using PointDataType = typename InputPointSetType::PixelType;
88  using PointsContainerType = typename InputPointSetType::PointsContainer;
89  using PointDataContainerType = typename InputPointSetType::PointDataContainer;
90 
96 
98  void
100  {
101  this->SetInput(0, field);
102  }
103 
105  const InputFieldType *
107  {
108  return this->GetInput(0);
109  }
110 
116  void
118  {
119  this->SetNthInput(1, const_cast<RealImageType *>(image));
120  }
121  void
122  SetInput1(const RealImageType * image)
123  {
124  this->SetConfidenceImage(image);
125  }
129  const RealImageType *
131  {
132  return static_cast<const RealImageType *>(this->ProcessObject::GetInput(1));
133  }
134 
136  void
138  {
139  this->SetNthInput(2, const_cast<InputPointSetType *>(points));
140  }
141  void
142  SetInput2(const InputPointSetType * points)
143  {
144  this->SetPointSet(points);
145  }
149  const InputPointSetType *
150  GetPointSet() const
151  {
152  return static_cast<const InputPointSetType *>(this->ProcessObject::GetInput(2));
153  }
154 
156  void
157  SetPointSetConfidenceWeights(WeightsContainerType * weights);
158 
160  const DisplacementFieldControlPointLatticeType *
162  {
163  return static_cast<const DisplacementFieldControlPointLatticeType *>(this->GetOutput(1));
164  }
165 
167  void
168  SetBSplineDomainFromImage(RealImageType *);
169 
171  void
173  {
174  this->SetBSplineDomainFromImage(const_cast<RealImageType *>(image));
175  }
176 
178  void
179  SetBSplineDomainFromImage(InputFieldType *);
180 
182  void
184  {
185  this->SetBSplineDomainFromImage(const_cast<InputFieldType *>(field));
186  }
187 
189  void SetBSplineDomain(OriginType, SpacingType, SizeType, DirectionType);
190 
191  /* Set/Get b-spline domain origin. */
192  itkGetConstMacro(BSplineDomainOrigin, OriginType);
193 
194  /* Set/Get b-spline domain spacing. */
195  itkGetConstMacro(BSplineDomainSpacing, SpacingType);
196 
197  /* Set/Get b-spline domain size. */
198  itkGetConstMacro(BSplineDomainSize, SizeType);
199 
200  /* Set/Get b-spline domain direction. */
201  itkGetConstMacro(BSplineDomainDirection, DirectionType);
202 
203  /* Use input field to define the B-spline domain. */
204  itkSetMacro(UseInputFieldToDefineTheBSplineDomain, bool);
205  itkGetConstMacro(UseInputFieldToDefineTheBSplineDomain, bool);
206  itkBooleanMacro(UseInputFieldToDefineTheBSplineDomain);
207 
211  itkSetMacro(SplineOrder, unsigned int);
212 
216  itkGetConstMacro(SplineOrder, unsigned int);
217 
225  itkSetMacro(NumberOfControlPoints, ArrayType);
226 
234  itkGetConstMacro(NumberOfControlPoints, ArrayType);
235 
242  itkSetMacro(NumberOfFittingLevels, ArrayType);
243 
250  void
251  SetNumberOfFittingLevels(unsigned int n)
252  {
253  ArrayType nlevels;
254 
255  nlevels.Fill(n);
256  this->SetNumberOfFittingLevels(nlevels);
257  }
258 
265  itkGetConstMacro(NumberOfFittingLevels, ArrayType);
266 
270  itkBooleanMacro(EstimateInverse);
271  itkSetMacro(EstimateInverse, bool);
272  itkGetConstMacro(EstimateInverse, bool);
278  itkBooleanMacro(EnforceStationaryBoundary);
279  itkSetMacro(EnforceStationaryBoundary, bool);
280  itkGetConstMacro(EnforceStationaryBoundary, bool);
283 protected:
286 
288  ~DisplacementFieldToBSplineImageFilter() override = default;
289 
291  void
292  PrintSelf(std::ostream & os, Indent indent) const override;
293 
295  void
296  GenerateData() override;
297 
298 private:
299  bool m_EstimateInverse{ false };
300  bool m_EnforceStationaryBoundary{ true };
301  unsigned int m_SplineOrder{ 3 };
302  ArrayType m_NumberOfControlPoints{};
303  ArrayType m_NumberOfFittingLevels{};
304 
305  typename WeightsContainerType::Pointer m_PointWeights{};
306  bool m_UsePointWeights{ false };
307 
308  OriginType m_BSplineDomainOrigin{};
309  SpacingType m_BSplineDomainSpacing{};
310  SizeType m_BSplineDomainSize{};
311  DirectionType m_BSplineDomainDirection{};
312 
313  bool m_BSplineDomainIsDefined{ true };
314  bool m_UseInputFieldToDefineTheBSplineDomain{ false };
315 };
316 
317 } // end namespace itk
318 
319 #ifndef ITK_MANUAL_INSTANTIATION
320 # include "itkDisplacementFieldToBSplineImageFilter.hxx"
321 #endif
322 
323 #endif
itk::DisplacementFieldToBSplineImageFilter::PointType
typename InputPointSetType::PointType PointType
Definition: itkDisplacementFieldToBSplineImageFilter.h:86
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
itk::DisplacementFieldToBSplineImageFilter::IndexType
typename OutputFieldType::IndexType IndexType
Definition: itkDisplacementFieldToBSplineImageFilter.h:75
itk::DisplacementFieldToBSplineImageFilter::RegionType
typename OutputFieldType::RegionType RegionType
Definition: itkDisplacementFieldToBSplineImageFilter.h:74
itk::DisplacementFieldToBSplineImageFilter::SpacingType
typename OutputFieldType::SpacingType SpacingType
Definition: itkDisplacementFieldToBSplineImageFilter.h:77
itk::GTest::TypedefsAndConstructors::Dimension2::DirectionType
ImageBaseType::DirectionType DirectionType
Definition: itkGTestTypedefsAndConstructors.h:52
itk::DisplacementFieldToBSplineImageFilter::SetBSplineDomainFromImage
void SetBSplineDomainFromImage(const RealImageType *image)
Definition: itkDisplacementFieldToBSplineImageFilter.h:172
itk::DisplacementFieldToBSplineImageFilter::SetInput1
void SetInput1(const RealImageType *image)
Definition: itkDisplacementFieldToBSplineImageFilter.h:122
itk::DisplacementFieldToBSplineImageFilter::InputFieldPointType
typename InputFieldType::PointType InputFieldPointType
Definition: itkDisplacementFieldToBSplineImageFilter.h:69
itk::DisplacementFieldToBSplineImageFilter::SetPointSet
void SetPointSet(const InputPointSetType *points)
Definition: itkDisplacementFieldToBSplineImageFilter.h:137
itk::DisplacementFieldToBSplineImageFilter::SetInput2
void SetInput2(const InputPointSetType *points)
Definition: itkDisplacementFieldToBSplineImageFilter.h:142
itk::DisplacementFieldToBSplineImageFilter::OriginType
typename OutputFieldType::PointType OriginType
Definition: itkDisplacementFieldToBSplineImageFilter.h:78
itk::DisplacementFieldToBSplineImageFilter::PointsContainerType
typename InputPointSetType::PointsContainer PointsContainerType
Definition: itkDisplacementFieldToBSplineImageFilter.h:88
itk::DisplacementFieldToBSplineImageFilter::SizeType
typename OutputFieldType::SizeType SizeType
Definition: itkDisplacementFieldToBSplineImageFilter.h:79
itk::GTest::TypedefsAndConstructors::Dimension2::PointType
ImageBaseType::PointType PointType
Definition: itkGTestTypedefsAndConstructors.h:51
itk::BSplineScatteredDataPointSetToImageFilter
Image filter which provides a B-spline output approximation.
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:132
itk::DisplacementFieldToBSplineImageFilter::GetDisplacementField
const InputFieldType * GetDisplacementField() const
Definition: itkDisplacementFieldToBSplineImageFilter.h:106
itk::DisplacementFieldToBSplineImageFilter
Class which takes a dense displacement field image and/or a set of points with associated displacemen...
Definition: itkDisplacementFieldToBSplineImageFilter.h:44
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itk::DisplacementFieldToBSplineImageFilter::GetPointSet
const InputPointSetType * GetPointSet() const
Definition: itkDisplacementFieldToBSplineImageFilter.h:150
itk::DisplacementFieldToBSplineImageFilter::DisplacementFieldControlPointLatticeType
typename BSplineFilterType::PointDataImageType DisplacementFieldControlPointLatticeType
Definition: itkDisplacementFieldToBSplineImageFilter.h:94
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::DisplacementFieldToBSplineImageFilter::PointDataType
typename InputPointSetType::PixelType PointDataType
Definition: itkDisplacementFieldToBSplineImageFilter.h:87
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itk::ImageToImageFilter
Base class for filters that take an image as input and produce an image as output.
Definition: itkImageToImageFilter.h:108
itk::ImageSource
Base class for all process objects that output image data.
Definition: itkImageSource.h:67
itk::DisplacementFieldToBSplineImageFilter::GetDisplacementFieldControlPointLattice
const DisplacementFieldControlPointLatticeType * GetDisplacementFieldControlPointLattice() const
Definition: itkDisplacementFieldToBSplineImageFilter.h:161
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itk::DisplacementFieldToBSplineImageFilter::InputFieldType
TInputImage InputFieldType
Definition: itkDisplacementFieldToBSplineImageFilter.h:63
itk::DisplacementFieldToBSplineImageFilter::PixelType
typename OutputFieldType::PixelType PixelType
Definition: itkDisplacementFieldToBSplineImageFilter.h:72
itkBSplineScatteredDataPointSetToImageFilter.h
itk::DisplacementFieldToBSplineImageFilter::DirectionType
typename OutputFieldType::DirectionType DirectionType
Definition: itkDisplacementFieldToBSplineImageFilter.h:80
itk::DisplacementFieldToBSplineImageFilter::InputPointSetType
TInputPointSet InputPointSetType
Definition: itkDisplacementFieldToBSplineImageFilter.h:64
itkImageToImageFilter.h
itk::FixedArray< unsigned int, Self::ImageDimension >
itk::DisplacementFieldToBSplineImageFilter::SetBSplineDomainFromImage
void SetBSplineDomainFromImage(const InputFieldType *field)
Definition: itkDisplacementFieldToBSplineImageFilter.h:183
itk::DisplacementFieldToBSplineImageFilter::RealType
typename VectorType::RealValueType RealType
Definition: itkDisplacementFieldToBSplineImageFilter.h:82
itk::DisplacementFieldToBSplineImageFilter::WeightsContainerType
typename BSplineFilterType::WeightsContainerType WeightsContainerType
Definition: itkDisplacementFieldToBSplineImageFilter.h:93
itk::DisplacementFieldToBSplineImageFilter::SetConfidenceImage
void SetConfidenceImage(const RealImageType *image)
Definition: itkDisplacementFieldToBSplineImageFilter.h:117
itk::DisplacementFieldToBSplineImageFilter::GetConfidenceImage
const RealImageType * GetConfidenceImage() const
Definition: itkDisplacementFieldToBSplineImageFilter.h:130
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itk::DisplacementFieldToBSplineImageFilter::VectorType
typename OutputFieldType::PixelType VectorType
Definition: itkDisplacementFieldToBSplineImageFilter.h:73
itk::DisplacementFieldToBSplineImageFilter::PointDataContainerType
typename InputPointSetType::PointDataContainer PointDataContainerType
Definition: itkDisplacementFieldToBSplineImageFilter.h:89
itk::DisplacementFieldToBSplineImageFilter::InverseDisplacementFieldType
OutputFieldType InverseDisplacementFieldType
Definition: itkDisplacementFieldToBSplineImageFilter.h:68
itk::ProcessObject
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Definition: itkProcessObject.h:139
itkVector.h
itk::DisplacementFieldToBSplineImageFilter::SetDisplacementField
void SetDisplacementField(const InputFieldType *field)
Definition: itkDisplacementFieldToBSplineImageFilter.h:99
itk::DisplacementFieldToBSplineImageFilter::ArrayType
typename BSplineFilterType::ArrayType ArrayType
Definition: itkDisplacementFieldToBSplineImageFilter.h:95
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:88
itk::DisplacementFieldToBSplineImageFilter::OutputFieldType
TOutputImage OutputFieldType
Definition: itkDisplacementFieldToBSplineImageFilter.h:65
itkPointSet.h
itk::DisplacementFieldToBSplineImageFilter::SetNumberOfFittingLevels
void SetNumberOfFittingLevels(unsigned int n)
Definition: itkDisplacementFieldToBSplineImageFilter.h:251
itk::DisplacementFieldToBSplineImageFilter::DisplacementFieldType
InputFieldType DisplacementFieldType
Definition: itkDisplacementFieldToBSplineImageFilter.h:67
itk::VectorContainer
Define a front-end to the STL "vector" container that conforms to the IndexedContainerInterface.
Definition: itkVectorContainer.h:48
itk::ProcessObject::GetInput
DataObject * GetInput(const DataObjectIdentifierType &key)
Return an input.