File: SmoothDisplacementField.cxx

package info (click to toggle)
ants 2.5.4%2Bdfsg-1
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid, trixie
  • size: 11,672 kB
  • sloc: cpp: 85,685; sh: 15,850; perl: 863; xml: 115; python: 111; makefile: 68
file content (305 lines) | stat: -rw-r--r-- 10,408 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
/*=========================================================================

  Program:   Advanced Normalization Tools

  Copyright (c) ConsortiumOfANTS. All rights reserved.
  See accompanying COPYING.txt or
 https://github.com/stnava/ANTs/blob/master/ANTSCopyright.txt for details.

     This software is distributed WITHOUT ANY WARRANTY; without even
     the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
     PURPOSE.  See the above copyright notices for more information.

=========================================================================*/

#include "antsUtilities.h"
#include <algorithm>
#include "ReadWriteData.h"

#include "itkDisplacementFieldToBSplineImageFilter.h"
#include "itkGaussianOperator.h"
#include "itkImage.h"
#include "itkImageDuplicator.h"
#include "itkTimeProbe.h"
#include "itkImageRegionIteratorWithIndex.h"
#include "itkImageRegionConstIteratorWithIndex.h"
#include "itkVectorNeighborhoodOperatorImageFilter.h"

namespace ants
{

template <unsigned int ImageDimension>
int
SmoothDisplacementField(int argc, char * argv[])
{

  using RealType = float;
  using RealImageType = itk::Image<RealType, ImageDimension>;
  using VectorType = itk::Vector<RealType, ImageDimension>;
  using DisplacementFieldType = itk::Image<VectorType, ImageDimension>;

  typename DisplacementFieldType::Pointer field = nullptr;
  ReadImage<DisplacementFieldType>(field, argv[2]);

  typename DisplacementFieldType::Pointer smoothField = DisplacementFieldType::New();

  float elapsedTime = 0.0;

  std::vector<float> var = ConvertVector<float>(std::string(argv[4]));
  if (var.size() == 1)
  {
    float variance = var[0];

    using GaussianSmoothingOperatorType = itk::GaussianOperator<float, ImageDimension>;
    using GaussianSmoothingSmootherType =
      itk::VectorNeighborhoodOperatorImageFilter<DisplacementFieldType, DisplacementFieldType>;

    GaussianSmoothingOperatorType gaussianSmoothingOperator;

    using DuplicatorType = itk::ImageDuplicator<DisplacementFieldType>;
    typename DuplicatorType::Pointer duplicator = DuplicatorType::New();
    duplicator->SetInputImage(field);
    duplicator->Update();

    smoothField = duplicator->GetOutput();

    itk::TimeProbe timer;
    timer.Start();
    typename GaussianSmoothingSmootherType::Pointer smoother = GaussianSmoothingSmootherType::New();

    for (unsigned int dimension = 0; dimension < ImageDimension; ++dimension)
    {
      // smooth along this dimension
      gaussianSmoothingOperator.SetDirection(dimension);
      gaussianSmoothingOperator.SetVariance(variance);
      gaussianSmoothingOperator.SetMaximumError(0.001);
      gaussianSmoothingOperator.SetMaximumKernelWidth(smoothField->GetRequestedRegion().GetSize()[dimension]);
      gaussianSmoothingOperator.CreateDirectional();

      // todo: make sure we only smooth within the buffered region
      smoother->SetOperator(gaussianSmoothingOperator);
      smoother->SetInput(smoothField);
      smoother->Update();

      smoothField = smoother->GetOutput();
      smoothField->Update();
      smoothField->DisconnectPipeline();
    }

    const VectorType zeroVector(0.0);

    // make sure boundary does not move
    float weight1 = itk::NumericTraits<float>::OneValue();
    if (variance < static_cast<float>(0.5))
    {
      weight1 = itk::NumericTraits<float>::OneValue() -
                itk::NumericTraits<float>::OneValue() * (variance / static_cast<float>(0.5));
    }
    float weight2 = itk::NumericTraits<float>::OneValue() - weight1;

    const typename DisplacementFieldType::RegionType region = field->GetLargestPossibleRegion();
    const typename DisplacementFieldType::SizeType   size = region.GetSize();
    const typename DisplacementFieldType::IndexType  startIndex = region.GetIndex();

    itk::ImageRegionConstIteratorWithIndex<DisplacementFieldType> fieldIt(field, field->GetLargestPossibleRegion());
    itk::ImageRegionIteratorWithIndex<DisplacementFieldType>      smoothedFieldIt(smoothField,
                                                                             smoothField->GetLargestPossibleRegion());
    for (fieldIt.GoToBegin(), smoothedFieldIt.GoToBegin(); !fieldIt.IsAtEnd(); ++fieldIt, ++smoothedFieldIt)
    {
      typename DisplacementFieldType::IndexType index = fieldIt.GetIndex();
      bool                                      isOnBoundary = false;
      for (unsigned int dimension = 0; dimension < ImageDimension; ++dimension)
      {
        if (index[dimension] == startIndex[dimension] ||
            index[dimension] == static_cast<int>(size[dimension]) - startIndex[dimension] - 1)
        {
          isOnBoundary = true;
          break;
        }
      }
      if (isOnBoundary)
      {
        smoothedFieldIt.Set(zeroVector);
      }
      else
      {
        smoothedFieldIt.Set(smoothedFieldIt.Get() * weight1 + fieldIt.Get() * weight2);
      }
    }
    timer.Stop();
    elapsedTime = timer.GetMean();
  }
  else if (var.size() == ImageDimension)
  {
    using BSplineFilterType = itk::DisplacementFieldToBSplineImageFilter<DisplacementFieldType>;

    unsigned int numberOfLevels = 1;
    if (argc > 5)
    {
      numberOfLevels = std::stoi(argv[5]);
    }

    unsigned int splineOrder = 3;
    if (argc > 6)
    {
      splineOrder = std::stoi(argv[6]);
    }

    typename BSplineFilterType::ArrayType ncps;
    for (unsigned int d = 0; d < ImageDimension; d++)
    {
      ncps[d] = var[d] + splineOrder;
    }

    typename BSplineFilterType::Pointer bspliner = BSplineFilterType::New();
    bspliner->SetDisplacementField(field);
    bspliner->SetBSplineDomainFromImage(field);
    bspliner->SetNumberOfControlPoints(ncps);
    bspliner->SetSplineOrder(splineOrder);
    bspliner->SetNumberOfFittingLevels(numberOfLevels);
    bspliner->SetEnforceStationaryBoundary(true);
    bspliner->SetEstimateInverse(false);

    if (argc > 7)
    {
      bspliner->SetEstimateInverse(static_cast<bool>(std::stoi(argv[7])));
    }

    if (argc > 8)
    {
      typename BSplineFilterType::RealImageType::Pointer confidenceImage = nullptr;
      ReadImage<typename BSplineFilterType::RealImageType>(confidenceImage, argv[8]);
      bspliner->SetConfidenceImage(confidenceImage);
    }

    itk::TimeProbe timer;
    timer.Start();
    bspliner->Update();
    timer.Stop();
    elapsedTime = timer.GetMean();

    smoothField = bspliner->GetOutput();
    smoothField->DisconnectPipeline();
  }
  else
  {
    std::cerr << "Error: unexpected variance format." << std::endl;
    return EXIT_FAILURE;
  }

  typename RealImageType::Pointer rmseImage = RealImageType::New();
  rmseImage->SetOrigin(field->GetOrigin());
  rmseImage->SetDirection(field->GetDirection());
  rmseImage->SetSpacing(field->GetSpacing());
  rmseImage->SetRegions(field->GetLargestPossibleRegion());
  rmseImage->AllocateInitialized();

  itk::ImageRegionConstIterator<DisplacementFieldType> fieldIt(field, field->GetLargestPossibleRegion());
  itk::ImageRegionConstIterator<DisplacementFieldType> smoothedFieldIt(smoothField,
                                                                       smoothField->GetLargestPossibleRegion());
  itk::ImageRegionIterator<RealImageType>              ItR(rmseImage, rmseImage->GetLargestPossibleRegion());

  float             rmse = 0.0;
  vnl_vector<float> rmse_comp(2);
  rmse_comp.fill(0.0);
  float N = 0.0;
  for (fieldIt.GoToBegin(), smoothedFieldIt.GoToBegin(), ItR.GoToBegin(); !fieldIt.IsAtEnd();
       ++fieldIt, ++smoothedFieldIt, ++ItR)
  {
    ItR.Set((fieldIt.Get() - smoothedFieldIt.Get()).GetSquaredNorm());
    for (unsigned int d = 0; d < ImageDimension; d++)
    {
      rmse_comp[d] += itk::Math::sqr(fieldIt.Get()[d] - smoothedFieldIt.Get()[d]);
    }
    rmse += static_cast<float>((fieldIt.Get() - smoothedFieldIt.Get()).GetSquaredNorm());
    N += itk::NumericTraits<float>::OneValue();
  }
  rmse = std::sqrt(rmse / N);

  std::cout << "Elapsed time: " << elapsedTime << std::endl;
  std::cout << "RMSE = " << rmse << std::endl;
  for (unsigned int d = 0; d < ImageDimension; d++)
  {
    std::cout << "  rmse[" << d << "] = " << std::sqrt(rmse_comp[d] / N) << std::endl;
  }

  ANTs::WriteImage<DisplacementFieldType>(smoothField, argv[3]);

  return EXIT_SUCCESS;
}

// entry point for the library; parameter 'args' is equivalent to 'argv' in (argc,argv) of commandline parameters to
// 'main()'
int
SmoothDisplacementField(std::vector<std::string> args, std::ostream * /*out_stream = nullptr */)
{
  // put the arguments coming in as 'args' into standard (argc,argv) format;
  // 'args' doesn't have the command name as first, argument, so add it manually;
  // 'args' may have adjacent arguments concatenated into one argument,
  // which the parser should handle
  args.insert(args.begin(), "SmoothDisplacementField");

  int     argc = args.size();
  char ** argv = new char *[args.size() + 1];
  for (unsigned int i = 0; i < args.size(); ++i)
  {
    // allocate space for the string plus a null character
    argv[i] = new char[args[i].length() + 1];
    std::strncpy(argv[i], args[i].c_str(), args[i].length());
    // place the null character in the end
    argv[i][args[i].length()] = '\0';
  }
  argv[argc] = nullptr;
  // class to automatically cleanup argv upon destruction
  class Cleanup_argv
  {
  public:
    Cleanup_argv(char ** argv_, int argc_plus_one_)
      : argv(argv_)
      , argc_plus_one(argc_plus_one_)
    {}

    ~Cleanup_argv()
    {
      for (unsigned int i = 0; i < argc_plus_one; ++i)
      {
        delete[] argv[i];
      }
      delete[] argv;
    }

  private:
    char **      argv;
    unsigned int argc_plus_one;
  };
  Cleanup_argv cleanup_argv(argv, argc + 1);

  // antscout->set_stream( out_stream );

  if (argc < 5)
  {
    std::cout << argv[0] << " imageDimension inputField outputField variance_or_mesh_size_base_level "
              << "[numberOfevels=1] [splineOrder=3] [estimateInverse=0] [confidenceImage]" << std::endl;

    return EXIT_FAILURE;
  }

  switch (std::stoi(argv[1]))
  {
    case 2:
      return SmoothDisplacementField<2>(argc, argv);
      break;
    case 3:
      return SmoothDisplacementField<3>(argc, argv);
      break;
    case 4:
      return SmoothDisplacementField<4>(argc, argv);
      break;
    default:
      std::cerr << "Unsupported dimension" << std::endl;
      return EXIT_FAILURE;
  }
  return EXIT_SUCCESS;
}

} // namespace ants