File: itkDisplacementFieldToBSplineImageFilter.h

package info (click to toggle)
insighttoolkit5 5.4.3-5
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid, trixie
  • size: 704,384 kB
  • sloc: cpp: 783,592; ansic: 628,724; xml: 44,704; fortran: 34,250; python: 22,874; sh: 4,078; pascal: 2,636; lisp: 2,158; makefile: 464; yacc: 328; asm: 205; perl: 203; lex: 146; tcl: 132; javascript: 98; csh: 81
file content (319 lines) | stat: -rw-r--r-- 10,788 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
/*=========================================================================
 *
 *  Copyright NumFOCUS
 *
 *  Licensed under the Apache License, Version 2.0 (the "License");
 *  you may not use this file except in compliance with the License.
 *  You may obtain a copy of the License at
 *
 *         https://www.apache.org/licenses/LICENSE-2.0.txt
 *
 *  Unless required by applicable law or agreed to in writing, software
 *  distributed under the License is distributed on an "AS IS" BASIS,
 *  WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
 *  See the License for the specific language governing permissions and
 *  limitations under the License.
 *
 *=========================================================================*/
#ifndef itkDisplacementFieldToBSplineImageFilter_h
#define itkDisplacementFieldToBSplineImageFilter_h

#include "itkImageToImageFilter.h"

#include "itkBSplineScatteredDataPointSetToImageFilter.h"
#include "itkPointSet.h"
#include "itkVector.h"

namespace itk
{

/**
 * \class DisplacementFieldToBSplineImageFilter
 * \brief Class which takes a dense displacement field image and/or a set of points
 * with associated displacements and smooths them using B-splines.  The inverse
 * can also be estimated.
 *
 * \author Nick Tustison
 *
 * \ingroup ITKDisplacementField
 */

template <typename TInputImage,
          typename TInputPointSet = PointSet<typename TInputImage::PixelType, TInputImage::ImageDimension>,
          typename TOutputImage = TInputImage>
class ITK_TEMPLATE_EXPORT DisplacementFieldToBSplineImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
{
public:
  ITK_DISALLOW_COPY_AND_MOVE(DisplacementFieldToBSplineImageFilter);

  using Self = DisplacementFieldToBSplineImageFilter;
  using Superclass = ImageToImageFilter<TInputImage, TOutputImage>;
  using Pointer = SmartPointer<Self>;
  using ConstPointer = SmartPointer<const Self>;

  /** Method for creation through the object factory. */
  itkNewMacro(Self);

  /** \see LightObject::GetNameOfClass() */
  itkOverrideGetNameOfClassMacro(DisplacementFieldToBSplineImageFilter);

  /** Extract dimension from input image. */
  static constexpr unsigned int ImageDimension = TInputImage::ImageDimension;

  using InputFieldType = TInputImage;
  using InputPointSetType = TInputPointSet;
  using OutputFieldType = TOutputImage;

  using DisplacementFieldType = InputFieldType;
  using InverseDisplacementFieldType = OutputFieldType;
  using InputFieldPointType = typename InputFieldType::PointType;

  /** Image type alias support */
  using PixelType = typename OutputFieldType::PixelType;
  using VectorType = typename OutputFieldType::PixelType;
  using RegionType = typename OutputFieldType::RegionType;
  using IndexType = typename OutputFieldType::IndexType;

  using SpacingType = typename OutputFieldType::SpacingType;
  using OriginType = typename OutputFieldType::PointType;
  using SizeType = typename OutputFieldType::SizeType;
  using DirectionType = typename OutputFieldType::DirectionType;

  using RealType = typename VectorType::RealValueType;
  using RealImageType = Image<RealType, ImageDimension>;

  /** Point set type alias support */
  using PointType = typename InputPointSetType::PointType;
  using PointDataType = typename InputPointSetType::PixelType;
  using PointsContainerType = typename InputPointSetType::PointsContainer;
  using PointDataContainerType = typename InputPointSetType::PointDataContainer;

  /** B-sline filter type alias */
  using BSplineFilterType = BSplineScatteredDataPointSetToImageFilter<InputPointSetType, OutputFieldType>;
  using WeightsContainerType = typename BSplineFilterType::WeightsContainerType;
  using DisplacementFieldControlPointLatticeType = typename BSplineFilterType::PointDataImageType;
  using ArrayType = typename BSplineFilterType::ArrayType;

  /** Set the displacement field */
  void
  SetDisplacementField(const InputFieldType * field)
  {
    this->SetInput(0, field);
  }

  /** Get the input displacement field. */
  const InputFieldType *
  GetDisplacementField() const
  {
    return this->GetInput(0);
  }

  /**
   * Set confidence image function.  If a confidence image is specified,
   * estimation of the displacement field weights the contribution of each voxel
   * according the value of the corresponding voxel in the confidence image.
   */
  void
  SetConfidenceImage(const RealImageType * image)
  {
    this->SetNthInput(1, const_cast<RealImageType *>(image));
  }
  void
  SetInput1(const RealImageType * image)
  {
    this->SetConfidenceImage(image);
  }

  /** Get confidence image function. */
  const RealImageType *
  GetConfidenceImage() const
  {
    return static_cast<const RealImageType *>(this->ProcessObject::GetInput(1));
  }

  /** Set the input point set */
  void
  SetPointSet(const InputPointSetType * points)
  {
    this->SetNthInput(2, const_cast<InputPointSetType *>(points));
  }
  void
  SetInput2(const InputPointSetType * points)
  {
    this->SetPointSet(points);
  }

  /** Get the input point set. */
  const InputPointSetType *
  GetPointSet() const
  {
    return static_cast<const InputPointSetType *>(this->ProcessObject::GetInput(2));
  }

  /** Set the confidence weights associated with the input point set*/
  void
  SetPointSetConfidenceWeights(WeightsContainerType * weights);

  /** Get the displacement field control point lattice. */
  const DisplacementFieldControlPointLatticeType *
  GetDisplacementFieldControlPointLattice() const
  {
    return static_cast<const DisplacementFieldControlPointLatticeType *>(this->GetOutput(1));
  }

  /** Define the b-spline domain from an image */
  void
  SetBSplineDomainFromImage(RealImageType *);

  /** Define the b-spline domain from an image */
  void
  SetBSplineDomainFromImage(const RealImageType * image)
  {
    this->SetBSplineDomainFromImage(const_cast<RealImageType *>(image));
  }

  /** Define the b-spline domain from a displacement field */
  void
  SetBSplineDomainFromImage(InputFieldType *);

  /** Define the b-spline domain from a displacement field */
  void
  SetBSplineDomainFromImage(const InputFieldType * field)
  {
    this->SetBSplineDomainFromImage(const_cast<InputFieldType *>(field));
  }

  /** Define the b-spline domain explicitly. */
  void SetBSplineDomain(OriginType, SpacingType, SizeType, DirectionType);

  /* Set/Get b-spline domain origin. */
  itkGetConstMacro(BSplineDomainOrigin, OriginType);

  /* Set/Get b-spline domain spacing. */
  itkGetConstMacro(BSplineDomainSpacing, SpacingType);

  /* Set/Get b-spline domain size. */
  itkGetConstMacro(BSplineDomainSize, SizeType);

  /* Set/Get b-spline domain direction. */
  itkGetConstMacro(BSplineDomainDirection, DirectionType);

  /* Use input field to define the B-spline domain. */
  itkSetMacro(UseInputFieldToDefineTheBSplineDomain, bool);
  itkGetConstMacro(UseInputFieldToDefineTheBSplineDomain, bool);
  itkBooleanMacro(UseInputFieldToDefineTheBSplineDomain);

  /**
   * Set the spline order defining the bias field estimate.  Default = 3.
   */
  itkSetMacro(SplineOrder, unsigned int);

  /**
   * Get the spline order defining the bias field estimate.  Default = 3.
   */
  itkGetConstMacro(SplineOrder, unsigned int);

  /**
   * Set the control point grid size defining the B-spline estimate of the
   * scalar bias field.  In each dimension, the B-spline mesh size is equal
   * to the number of control points in that dimension minus the spline order.
   * Default = 4 control points in each dimension for a mesh size of 1 in each
   * dimension.
   */
  itkSetMacro(NumberOfControlPoints, ArrayType);

  /**
   * Get the control point grid size defining the B-spline estimate of the
   * scalar bias field.  In each dimension, the B-spline mesh size is equal
   * to the number of control points in that dimension minus the spline order.
   * Default = 4 control points in each dimension for a mesh size of 1 in each
   * dimension.
   */
  itkGetConstMacro(NumberOfControlPoints, ArrayType);

  /**
   * Set the number of fitting levels.  One of the contributions of N4 is the
   * introduction of a multi-scale approach to fitting. This allows one to
   * specify a B-spline mesh size for initial fitting followed by a doubling of
   * the mesh resolution for each subsequent fitting level.  Default = 1 level.
   */
  itkSetMacro(NumberOfFittingLevels, ArrayType);

  /**
   * Set the number of fitting levels.  One of the contributions of N4 is the
   * introduction of a multi-scale approach to fitting. This allows one to
   * specify a B-spline mesh size for initial fitting followed by a doubling of
   * the mesh resolution for each subsequent fitting level.  Default = 1 level.
   */
  void
  SetNumberOfFittingLevels(unsigned int n)
  {
    ArrayType nlevels;

    nlevels.Fill(n);
    this->SetNumberOfFittingLevels(nlevels);
  }

  /**
   * Get the number of fitting levels.  One of the contributions of N4 is the
   * introduction of a multi-scale approach to fitting. This allows one to
   * specify a B-spline mesh size for initial fitting followed by a doubling of
   * the mesh resolution for each subsequent fitting level.  Default = 1 level.
   */
  itkGetConstMacro(NumberOfFittingLevels, ArrayType);

  /**
   * Estimate the inverse field instead of the forward field.  Default = false.
   */
  itkBooleanMacro(EstimateInverse);
  itkSetMacro(EstimateInverse, bool);
  itkGetConstMacro(EstimateInverse, bool);

  /**
   * Enforce stationary boundary conditions.  Default = false.
   */
  itkBooleanMacro(EnforceStationaryBoundary);
  itkSetMacro(EnforceStationaryBoundary, bool);
  itkGetConstMacro(EnforceStationaryBoundary, bool);

protected:
  /** Constructor */
  DisplacementFieldToBSplineImageFilter();

  /** Deconstructor */
  ~DisplacementFieldToBSplineImageFilter() override = default;

  /** Standard print self function **/
  void
  PrintSelf(std::ostream & os, Indent indent) const override;

  /** preprocessing function */
  void
  GenerateData() override;

private:
  bool         m_EstimateInverse{ false };
  bool         m_EnforceStationaryBoundary{ true };
  unsigned int m_SplineOrder{ 3 };
  ArrayType    m_NumberOfControlPoints{};
  ArrayType    m_NumberOfFittingLevels{};

  typename WeightsContainerType::Pointer m_PointWeights{};
  bool                                   m_UsePointWeights{ false };

  OriginType    m_BSplineDomainOrigin{};
  SpacingType   m_BSplineDomainSpacing{};
  SizeType      m_BSplineDomainSize{};
  DirectionType m_BSplineDomainDirection{};

  bool m_BSplineDomainIsDefined{ true };
  bool m_UseInputFieldToDefineTheBSplineDomain{ false };
};

} // end namespace itk

#ifndef ITK_MANUAL_INSTANTIATION
#  include "itkDisplacementFieldToBSplineImageFilter.hxx"
#endif

#endif