File: itkSimpleImageRegistrationTest4.cxx

package info (click to toggle)
insighttoolkit5 5.4.5-2
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 704,592 kB
  • sloc: cpp: 784,579; ansic: 628,724; xml: 44,704; fortran: 34,250; python: 22,934; sh: 4,078; pascal: 2,636; lisp: 2,158; makefile: 460; yacc: 328; asm: 205; perl: 203; lex: 146; tcl: 132; javascript: 98; csh: 81
file content (177 lines) | stat: -rw-r--r-- 5,598 bytes parent folder | download | duplicates (2)
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
/*=========================================================================
 *
 *  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.
 *
 *=========================================================================*/

#include "itkImageFileReader.h"
#include "itkImageFileWriter.h"


#include "itkEuler2DTransform.h"
#include "itkMeanSquaresImageToImageMetricv4.h"
#include "itkTestingMacros.h"
#include "itkImageRegistrationMethodv4.h"
#include "itkConjugateGradientLineSearchOptimizerv4.h"

#include <iomanip>

namespace
{
template <typename TOptimizer>
class CommandIterationUpdate : public itk::Command
{
public:
  using Self = CommandIterationUpdate;
  using Superclass = itk::Command;
  using Pointer = itk::SmartPointer<Self>;
  itkNewMacro(Self);

protected:
  CommandIterationUpdate()
  {
    // mark used to avoid warnings
    (void)&Self::Clone;
  };

public:
  void
  Execute(itk::Object * caller, const itk::EventObject & event) override
  {
    Execute((const itk::Object *)caller, event);
  }

  void
  Execute(const itk::Object * object, const itk::EventObject & event) override
  {

    const auto * optimizer = dynamic_cast<const TOptimizer *>(object);

    if (typeid(event) != typeid(itk::IterationEvent) || !optimizer)
    {
      return;
    }

    // stash the stream state
    std::ios state(nullptr);
    state.copyfmt(std::cout);
    std::cout << std::fixed << std::setfill(' ') << std::setprecision(5);
    std::cout << std::setw(3) << optimizer->GetCurrentIteration();
    std::cout << " = " << std::setw(10) << optimizer->GetCurrentMetricValue();
    std::cout << " : " << optimizer->GetCurrentPosition() << std::endl;
    std::cout << "\nScales: " << optimizer->GetScales() << std::endl;
  }
};

template <unsigned int TDimension>
int
ImageRegistration(int itkNotUsed(argc), char * argv[])
{
  const unsigned int ImageDimension = TDimension;


  using PixelType = float;
  using FixedImageType = itk::Image<PixelType, ImageDimension>;
  using MovingImageType = itk::Image<PixelType, ImageDimension>;

  using ImageReaderType = itk::ImageFileReader<FixedImageType>;

  auto fixedImageReader = ImageReaderType::New();
  fixedImageReader->SetFileName(argv[2]);
  fixedImageReader->Update();
  typename FixedImageType::Pointer fixedImage = fixedImageReader->GetOutput();
  fixedImage->Update();
  fixedImage->DisconnectPipeline();

  auto movingImageReader = ImageReaderType::New();
  movingImageReader->SetFileName(argv[3]);
  movingImageReader->Update();
  typename MovingImageType::Pointer movingImage = movingImageReader->GetOutput();
  movingImage->Update();
  movingImage->DisconnectPipeline();

  // Set up the centered transform initializer
  using TransformType = itk::Euler2DTransform<double>;
  auto initialTransform = TransformType::New();


  using MetricType = itk::MeanSquaresImageToImageMetricv4<FixedImageType, MovingImageType>;
  auto metric = MetricType::New();

  using RegistrationType = itk::ImageRegistrationMethodv4<FixedImageType, MovingImageType, TransformType>;
  auto registration = RegistrationType::New();

  registration->SetFixedImage(fixedImage);
  registration->SetMovingImage(movingImage);
  registration->SetMetric(metric);
  registration->SetMovingInitialTransform(initialTransform);
  registration->SetNumberOfLevels(1);

  using Optimizerv4Type = itk::ConjugateGradientLineSearchOptimizerv4;
  auto optimizer = Optimizerv4Type::New();

  optimizer->SetLearningRate(1.0);
  optimizer->SetNumberOfIterations(100);
  optimizer->SetMinimumConvergenceValue(1e-5);
  optimizer->SetConvergenceWindowSize(2);

  double                                           scaleData[] = { 200000, 1.0, 1.0 };
  typename Optimizerv4Type::ScalesType::Superclass scales(scaleData, 3);
  optimizer->SetScales(scales);

  registration->SetOptimizer(optimizer);

  using CommandType = CommandIterationUpdate<Optimizerv4Type>;
  auto observer = CommandType::New();
  optimizer->AddObserver(itk::IterationEvent(), observer);

  ITK_TRY_EXPECT_NO_EXCEPTION(registration->Update());


  registration->GetTransform()->Print(std::cout);
  std::cout << optimizer->GetStopConditionDescription() << std::endl;
  typename TransformType::ParametersType results = registration->GetTransform()->GetParameters();

  std::cout << "Expecting close (+/- 0.3) to: ( 0.0, -2.8, 9.5 )" << std::endl;
  std::cout << "Parameters: " << results << std::endl;

  std::cout << "Number Of Iterations: " << optimizer->GetCurrentIteration();
  ITK_TEST_EXPECT_TRUE(optimizer->GetCurrentIteration() > 5);

  return EXIT_SUCCESS;
}
} // namespace

int
itkSimpleImageRegistrationTest4(int argc, char * argv[])
{
  if (argc < 4)
  {
    std::cerr << "Missing parameters." << std::endl;
    std::cerr << "Usage: " << itkNameOfTestExecutableMacro(argv);
    std::cerr << " imageDimension fixedImage movingImage" << std::endl;
    return EXIT_FAILURE;
  }

  switch (std::stoi(argv[1]))
  {
    case 2:
      return ImageRegistration<2>(argc, argv);

    default:
      std::cerr << "Unsupported dimension" << std::endl;
      return EXIT_FAILURE;
  }
}