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
19 // INPUTS: {BrainProtonDensitySliceBorder20.png}
20 // INPUTS: {BrainProtonDensitySliceRotated10.png}
21 // OUTPUTS: {ImageRegistration5Output.png}
22 // OUTPUTS: {ImageRegistration5DifferenceAfter.png}
23 // OUTPUTS: {ImageRegistration5DifferenceBefore.png}
24 // ARGUMENTS: 0.1
25
26 // INPUTS: {BrainProtonDensitySliceBorder20.png}
27 // INPUTS: {BrainProtonDensitySliceR10X13Y17.png}
28 // OUTPUTS: {ImageRegistration5Output2.png}
29 // OUTPUTS: {ImageRegistration5DifferenceAfter2.png}
30 // OUTPUTS: {ImageRegistration5DifferenceBefore2.png}
31 // ARGUMENTS: 1.0
32
33
34 //
35 // This example illustrates the use of the \doxygen{CenteredRigid2DTransform}
36 // for performing rigid registration in $2D$. The example code is for the
37 // most part identical to that presented in Section
38 // \ref{sec:IntroductionImageRegistration}. The main difference is the use
39 // of the CenteredRigid2DTransform here instead of the
40 // \doxygen{TranslationTransform}.
41 //
42 // \index{itk::CenteredRigid2DTransform}
43 //
44
45 #include "itkImageRegistrationMethod.h"
46 #include "itkMeanSquaresImageToImageMetric.h"
47 #include "itkRegularStepGradientDescentOptimizer.h"
48
49
50 //
51 // In addition to the headers included in previous examples, the
52 // following header must also be included.
53 //
54 // \index{itk::CenteredRigid2DTransform!header}
55 //
56
57 #include "itkCenteredRigid2DTransform.h"
58
59
60 #include "itkImageFileReader.h"
61 #include "itkImageFileWriter.h"
62
63 #include "itkResampleImageFilter.h"
64 #include "itkSubtractImageFilter.h"
65 #include "itkRescaleIntensityImageFilter.h"
66
67
68 // The following section of code implements a Command observer
69 // that will monitor the evolution of the registration process.
70 //
71 #include "itkCommand.h"
72 class CommandIterationUpdate : public itk::Command
73 {
74 public:
75 using Self = CommandIterationUpdate;
76 using Superclass = itk::Command;
77 using Pointer = itk::SmartPointer<Self>;
78 itkNewMacro( Self );
79
80 protected:
81 CommandIterationUpdate() = default;
82
83 public:
84 using OptimizerType = itk::RegularStepGradientDescentOptimizer;
85 using OptimizerPointer = const OptimizerType *;
86
Execute(itk::Object * caller,const itk::EventObject & event)87 void Execute(itk::Object *caller, const itk::EventObject & event) override
88 {
89 Execute( (const itk::Object *)caller, event);
90 }
91
Execute(const itk::Object * object,const itk::EventObject & event)92 void Execute(const itk::Object * object, const itk::EventObject & event) override
93 {
94 auto optimizer = static_cast< OptimizerPointer >( object );
95 if( ! itk::IterationEvent().CheckEvent( &event ) )
96 {
97 return;
98 }
99 std::cout << optimizer->GetCurrentIteration() << " ";
100 std::cout << optimizer->GetValue() << " ";
101 std::cout << optimizer->GetCurrentPosition() << std::endl;
102 }
103 };
104
105 #include "itkTestDriverIncludeRequiredIOFactories.h"
main(int argc,char * argv[])106 int main( int argc, char *argv[] )
107 {
108 RegisterRequiredFactories();
109 if( argc < 4 )
110 {
111 std::cerr << "Missing Parameters " << std::endl;
112 std::cerr << "Usage: " << argv[0];
113 std::cerr << " fixedImageFile movingImageFile ";
114 std::cerr << " outputImagefile [differenceAfterRegistration] ";
115 std::cerr << " [differenceBeforeRegistration] ";
116 std::cerr << " [initialStepLength] "<< std::endl;
117 return EXIT_FAILURE;
118 }
119
120 constexpr unsigned int Dimension = 2;
121 using PixelType = unsigned char;
122
123 using FixedImageType = itk::Image< PixelType, Dimension >;
124 using MovingImageType = itk::Image< PixelType, Dimension >;
125
126
127 //
128 // The transform type is instantiated using the code below. The only
129 // template parameter for this class is the representation type of the
130 // space coordinates.
131 //
132 // \index{itk::CenteredRigid2DTransform!Instantiation}
133 //
134
135 using TransformType = itk::CenteredRigid2DTransform< double >;
136
137
138 using OptimizerType = itk::RegularStepGradientDescentOptimizer;
139 using MetricType = itk::MeanSquaresImageToImageMetric<
140 FixedImageType,
141 MovingImageType >;
142 using InterpolatorType = itk:: LinearInterpolateImageFunction<
143 MovingImageType,
144 double >;
145 using RegistrationType = itk::ImageRegistrationMethod<
146 FixedImageType,
147 MovingImageType >;
148
149 MetricType::Pointer metric = MetricType::New();
150 OptimizerType::Pointer optimizer = OptimizerType::New();
151 InterpolatorType::Pointer interpolator = InterpolatorType::New();
152 RegistrationType::Pointer registration = RegistrationType::New();
153
154 registration->SetMetric( metric );
155 registration->SetOptimizer( optimizer );
156 registration->SetInterpolator( interpolator );
157
158
159 //
160 // The transform object is constructed below and passed to the registration
161 // method.
162 //
163 // \index{itk::CenteredRigid2DTransform!New()}
164 // \index{itk::CenteredRigid2DTransform!Pointer}
165 // \index{itk::RegistrationMethod!SetTransform()}
166 //
167
168 TransformType::Pointer transform = TransformType::New();
169 registration->SetTransform( transform );
170
171
172 using FixedImageReaderType = itk::ImageFileReader< FixedImageType >;
173 using MovingImageReaderType = itk::ImageFileReader< MovingImageType >;
174
175 FixedImageReaderType::Pointer fixedImageReader = FixedImageReaderType::New();
176 MovingImageReaderType::Pointer movingImageReader = MovingImageReaderType::New();
177
178 fixedImageReader->SetFileName( argv[1] );
179 movingImageReader->SetFileName( argv[2] );
180
181
182 registration->SetFixedImage( fixedImageReader->GetOutput() );
183 registration->SetMovingImage( movingImageReader->GetOutput() );
184 fixedImageReader->Update();
185
186 registration->SetFixedImageRegion(
187 fixedImageReader->GetOutput()->GetBufferedRegion() );
188
189
190 //
191 // In this example, the input images are taken from readers. The code
192 // below updates the readers in order to ensure that the image parameters
193 // (size, origin and spacing) are valid when used to initialize the
194 // transform. We intend to use the center of the fixed image as the
195 // rotation center and then use the vector between the fixed image center
196 // and the moving image center as the initial translation to be applied
197 // after the rotation.
198 //
199
200 fixedImageReader->Update();
201 movingImageReader->Update();
202
203 using SpacingType = FixedImageType::SpacingType;
204 using OriginType = FixedImageType::PointType;
205 using RegionType = FixedImageType::RegionType;
206 using SizeType = FixedImageType::SizeType;
207
208 //
209 // The center of rotation is computed using the origin, size and spacing of
210 // the fixed image.
211 //
212
213 FixedImageType::Pointer fixedImage = fixedImageReader->GetOutput();
214
215 const SpacingType fixedSpacing = fixedImage->GetSpacing();
216 const OriginType fixedOrigin = fixedImage->GetOrigin();
217 const RegionType fixedRegion = fixedImage->GetLargestPossibleRegion();
218 const SizeType fixedSize = fixedRegion.GetSize();
219
220 TransformType::InputPointType centerFixed;
221
222 centerFixed[0] = fixedOrigin[0] + fixedSpacing[0] * fixedSize[0] / 2.0;
223 centerFixed[1] = fixedOrigin[1] + fixedSpacing[1] * fixedSize[1] / 2.0;
224
225
226 //
227 // The center of the moving image is computed in a similar way.
228 //
229
230 MovingImageType::Pointer movingImage = movingImageReader->GetOutput();
231
232 const SpacingType movingSpacing = movingImage->GetSpacing();
233 const OriginType movingOrigin = movingImage->GetOrigin();
234 const RegionType movingRegion = movingImage->GetLargestPossibleRegion();
235 const SizeType movingSize = movingRegion.GetSize();
236
237 TransformType::InputPointType centerMoving;
238
239 centerMoving[0] = movingOrigin[0] + movingSpacing[0] * movingSize[0] / 2.0;
240 centerMoving[1] = movingOrigin[1] + movingSpacing[1] * movingSize[1] / 2.0;
241
242
243 //
244 // Then, we initialize the transform by
245 // passing the center of the fixed image as the rotation center with the
246 // \code{SetCenter()} method. Also, the translation is set as the vector
247 // relating the center of the moving image to the center of the fixed
248 // image. This last vector is passed with the method
249 // \code{SetTranslation()}.
250 //
251
252 transform->SetCenter( centerFixed );
253 transform->SetTranslation( centerMoving - centerFixed );
254
255
256 //
257 // Let's finally initialize the rotation with a zero angle.
258 //
259
260 transform->SetAngle( 0.0 );
261
262
263 //
264 // Now we pass the current transform's parameters as the initial
265 // parameters to be used when the registration process starts.
266 //
267
268 registration->SetInitialTransformParameters( transform->GetParameters() );
269
270
271 //
272 // Keep in mind that the scale of units in rotation and translation is
273 // quite different, we take advantage of the scaling functionality provided
274 // by the optimizers. We know that the first element of the parameters array
275 // corresponds to the angle that is measured in radians, while the other
276 // parameters correspond to translations that are measured in millimeters.
277 // For this reason we use small factors in the scales associated with
278 // translations and the coordinates of the rotation center .
279 //
280
281 using OptimizerScalesType = OptimizerType::ScalesType;
282 OptimizerScalesType optimizerScales( transform->GetNumberOfParameters() );
283 const double translationScale = 1.0 / 1000.0;
284
285 optimizerScales[0] = 1.0;
286 optimizerScales[1] = translationScale;
287 optimizerScales[2] = translationScale;
288 optimizerScales[3] = translationScale;
289 optimizerScales[4] = translationScale;
290
291 optimizer->SetScales( optimizerScales );
292
293
294 //
295 // Next we set the normal parameters of the optimization method. In this
296 // case we are using an \doxygen{RegularStepGradientDescentOptimizer}.
297 // Below, we define the optimization parameters like the relaxation factor,
298 // initial step length, minimal step length and number of iterations. These
299 // last two act as stopping criteria for the optimization.
300 //
301 // \index{Regular\-Step\-Gradient\-Descent\-Optimizer!SetRelaxationFactor()}
302 // \index{Regular\-Step\-Gradient\-Descent\-Optimizer!SetMaximumStepLength()}
303 // \index{Regular\-Step\-Gradient\-Descent\-Optimizer!SetMinimumStepLength()}
304 // \index{Regular\-Step\-Gradient\-Descent\-Optimizer!SetNumberOfIterations()}
305 //
306
307 double initialStepLength = 0.1;
308
309 if( argc > 6 )
310 {
311 initialStepLength = std::stod( argv[6] );
312 }
313
314 optimizer->SetRelaxationFactor( 0.6 );
315 optimizer->SetMaximumStepLength( initialStepLength );
316 optimizer->SetMinimumStepLength( 0.001 );
317 optimizer->SetNumberOfIterations( 200 );
318
319
320 // Create the Command observer and register it with the optimizer.
321 //
322 CommandIterationUpdate::Pointer observer = CommandIterationUpdate::New();
323 optimizer->AddObserver( itk::IterationEvent(), observer );
324
325 try
326 {
327 registration->Update();
328 std::cout << "Optimizer stop condition: "
329 << registration->GetOptimizer()->GetStopConditionDescription()
330 << std::endl;
331 }
332 catch( itk::ExceptionObject & err )
333 {
334 std::cerr << "ExceptionObject caught !" << std::endl;
335 std::cerr << err << std::endl;
336 return EXIT_FAILURE;
337 }
338
339 OptimizerType::ParametersType finalParameters =
340 registration->GetLastTransformParameters();
341
342 const double finalAngle = finalParameters[0];
343 const double finalRotationCenterX = finalParameters[1];
344 const double finalRotationCenterY = finalParameters[2];
345 const double finalTranslationX = finalParameters[3];
346 const double finalTranslationY = finalParameters[4];
347
348 const unsigned int numberOfIterations = optimizer->GetCurrentIteration();
349
350 const double bestValue = optimizer->GetValue();
351
352
353 // Print out results
354 //
355 const double finalAngleInDegrees = finalAngle * 180.0 / itk::Math::pi;
356
357 std::cout << "Result = " << std::endl;
358 std::cout << " Angle (radians) = " << finalAngle << std::endl;
359 std::cout << " Angle (degrees) = " << finalAngleInDegrees << std::endl;
360 std::cout << " Translation X = " << finalTranslationX << std::endl;
361 std::cout << " Translation Y = " << finalTranslationY << std::endl;
362 std::cout << " Center X = " << finalRotationCenterX << std::endl;
363 std::cout << " Center Y = " << finalRotationCenterY << std::endl;
364 std::cout << " Iterations = " << numberOfIterations << std::endl;
365 std::cout << " Metric value = " << bestValue << std::endl;
366
367
368 //
369 // Let's execute this example over two of the images provided in
370 // \code{Examples/Data}:
371 //
372 // \begin{itemize}
373 // \item \code{BrainProtonDensitySliceBorder20.png}
374 // \item \code{BrainProtonDensitySliceRotated10.png}
375 // \end{itemize}
376 //
377 // The second image is the result of intentionally rotating the first image
378 // by $10$ degrees around the geometrical center of the image. Both images
379 // have unit-spacing and are shown in Figure
380 // \ref{fig:FixedMovingImageRegistration5}. The registration takes $20$
381 // iterations and produces the results:
382 //
383 // \begin{center}
384 // \begin{verbatim}
385 // [0.177458, 110.489, 128.488, 0.0106296, 0.00194103]
386 // \end{verbatim}
387 // \end{center}
388 //
389 // These results are interpreted as
390 //
391 // \begin{itemize}
392 // \item Angle = $0.177458$ radians
393 // \item Center = $( 110.489 , 128.488 )$ millimeters
394 // \item Translation = $( 0.0106296, 0.00194103 )$ millimeters
395 // \end{itemize}
396 //
397 // As expected, these values match the misalignment intentionally introduced
398 // into the moving image quite well, since $10$ degrees is about $0.174532$
399 // radians.
400 //
401 // \begin{figure}
402 // \center
403 // \includegraphics[width=0.44\textwidth]{BrainProtonDensitySliceBorder20}
404 // \includegraphics[width=0.44\textwidth]{BrainProtonDensitySliceRotated10}
405 // \itkcaption[Rigid2D Registration input images]{Fixed and moving images
406 // are provided as input to the registration method using the CenteredRigid2D
407 // transform.}
408 // \label{fig:FixedMovingImageRegistration5}
409 // \end{figure}
410 //
411 //
412 // \begin{figure}
413 // \center
414 // \includegraphics[width=0.32\textwidth]{ImageRegistration5Output}
415 // \includegraphics[width=0.32\textwidth]{ImageRegistration5DifferenceBefore}
416 // \includegraphics[width=0.32\textwidth]{ImageRegistration5DifferenceAfter}
417 // \itkcaption[Rigid2D Registration output images]{Resampled moving image
418 // (left). Differences between the fixed and moving images, before (center)
419 // and after (right) registration using the CenteredRigid2D transform.}
420 // \label{fig:ImageRegistration5Outputs}
421 // \end{figure}
422 //
423 // Figure \ref{fig:ImageRegistration5Outputs} shows from left to right the
424 // resampled moving image after registration, the difference between fixed
425 // and moving images before registration, and the difference between fixed
426 // and resampled moving image after registration. It can be seen from the
427 // last difference image that the rotational component has been solved but
428 // that a small centering misalignment persists.
429 //
430 // \begin{figure}
431 // \center
432 // \includegraphics[height=0.32\textwidth]{ImageRegistration5TraceMetric1}
433 // \includegraphics[height=0.32\textwidth]{ImageRegistration5TraceAngle1}
434 // \includegraphics[height=0.32\textwidth]{ImageRegistration5TraceTranslations1}
435 // \itkcaption[Rigid2D Registration output plots]{Metric values, rotation
436 // angle and translations during registration with the CenteredRigid2D
437 // transform.}
438 // \label{fig:ImageRegistration5Plots}
439 // \end{figure}
440 //
441 // Figure \ref{fig:ImageRegistration5Plots} shows plots of the main output
442 // parameters produced from the registration process. This includes the
443 // metric values at every iteration, the angle values at every iteration,
444 // and the translation components of the transform as the registration
445 // progress.
446 //
447
448
449 using ResampleFilterType = itk::ResampleImageFilter<
450 MovingImageType,
451 FixedImageType >;
452
453 TransformType::Pointer finalTransform = TransformType::New();
454
455 finalTransform->SetParameters( finalParameters );
456 finalTransform->SetFixedParameters( transform->GetFixedParameters() );
457
458 ResampleFilterType::Pointer resample = ResampleFilterType::New();
459
460 resample->SetTransform( finalTransform );
461 resample->SetInput( movingImageReader->GetOutput() );
462
463 resample->SetSize( fixedImage->GetLargestPossibleRegion().GetSize() );
464 resample->SetOutputOrigin( fixedImage->GetOrigin() );
465 resample->SetOutputSpacing( fixedImage->GetSpacing() );
466 resample->SetOutputDirection( fixedImage->GetDirection() );
467 resample->SetDefaultPixelValue( 100 );
468
469 using WriterFixedType = itk::ImageFileWriter< FixedImageType >;
470
471 WriterFixedType::Pointer writer = WriterFixedType::New();
472
473 writer->SetFileName( argv[3] );
474
475 writer->SetInput( resample->GetOutput() );
476
477 try
478 {
479 writer->Update();
480 }
481 catch( itk::ExceptionObject & excp )
482 {
483 std::cerr << "ExceptionObject while writing the resampled image !" << std::endl;
484 std::cerr << excp << std::endl;
485 return EXIT_FAILURE;
486 }
487
488 using DifferenceImageType = itk::Image< float, Dimension >;
489
490 using DifferenceFilterType = itk::SubtractImageFilter<
491 FixedImageType,
492 FixedImageType,
493 DifferenceImageType >;
494
495 DifferenceFilterType::Pointer difference = DifferenceFilterType::New();
496
497 using OutputPixelType = unsigned char;
498
499 using OutputImageType = itk::Image< OutputPixelType, Dimension >;
500
501 using RescalerType = itk::RescaleIntensityImageFilter<
502 DifferenceImageType,
503 OutputImageType >;
504
505 RescalerType::Pointer intensityRescaler = RescalerType::New();
506
507 intensityRescaler->SetOutputMinimum( 0 );
508 intensityRescaler->SetOutputMaximum( 255 );
509
510 difference->SetInput1( fixedImageReader->GetOutput() );
511 difference->SetInput2( resample->GetOutput() );
512
513 resample->SetDefaultPixelValue( 1 );
514
515 intensityRescaler->SetInput( difference->GetOutput() );
516
517 using WriterType = itk::ImageFileWriter< OutputImageType >;
518
519 WriterType::Pointer writer2 = WriterType::New();
520
521 writer2->SetInput( intensityRescaler->GetOutput() );
522
523
524 try
525 {
526 // Compute the difference image between the
527 // fixed and moving image after registration.
528 if( argc > 4 )
529 {
530 writer2->SetFileName( argv[4] );
531 writer2->Update();
532 }
533
534 // Compute the difference image between the
535 // fixed and resampled moving image after registration.
536 TransformType::Pointer identityTransform = TransformType::New();
537 identityTransform->SetIdentity();
538 resample->SetTransform( identityTransform );
539 if( argc > 5 )
540 {
541 writer2->SetFileName( argv[5] );
542 writer2->Update();
543 }
544 }
545 catch( itk::ExceptionObject & excp )
546 {
547 std::cerr << "Error while writing difference images" << std::endl;
548 std::cerr << excp << std::endl;
549 return EXIT_FAILURE;
550 }
551
552 //
553 // Let's now consider the case in which rotations and translations are
554 // present in the initial registration, as in the following pair
555 // of images:
556 //
557 // \begin{itemize}
558 // \item \code{BrainProtonDensitySliceBorder20.png}
559 // \item \code{BrainProtonDensitySliceR10X13Y17.png}
560 // \end{itemize}
561 //
562 // The second image is the result of intentionally rotating the first image
563 // by $10$ degrees and then translating it $13mm$ in $X$ and $17mm$ in $Y$.
564 // Both images have unit-spacing and are shown in Figure
565 // \ref{fig:FixedMovingImageRegistration5b}. In order to accelerate
566 // convergence it is convenient to use a larger step length as shown here.
567 //
568 // \code{optimizer->SetMaximumStepLength( 1.0 );}
569 //
570 // The registration now takes $46$ iterations and produces the following
571 // results:
572 //
573 // \begin{center}
574 // \begin{verbatim}
575 // [0.174454, 110.361, 128.647, 12.977, 15.9761]
576 // \end{verbatim}
577 // \end{center}
578 //
579 // These parameters are interpreted as
580 //
581 // \begin{itemize}
582 // \item Angle = $0.174454$ radians
583 // \item Center = $( 110.361 , 128.647 )$ millimeters
584 // \item Translation = $( 12.977 , 15.9761 )$ millimeters
585 // \end{itemize}
586 //
587 // These values approximately match the initial misalignment intentionally
588 // introduced into the moving image, since $10$ degrees is about $0.174532$
589 // radians. The horizontal translation is well resolved while the vertical
590 // translation ends up being off by about one millimeter.
591 //
592 // \begin{figure}
593 // \center
594 // \includegraphics[width=0.44\textwidth]{BrainProtonDensitySliceBorder20}
595 // \includegraphics[width=0.44\textwidth]{BrainProtonDensitySliceR10X13Y17}
596 // \itkcaption[Rigid2D Registration input images]{Fixed and moving images
597 // provided as input to the registration method using the CenteredRigid2D
598 // transform.}
599 // \label{fig:FixedMovingImageRegistration5b}
600 // \end{figure}
601 //
602 //
603 // \begin{figure}
604 // \center
605 // \includegraphics[width=0.32\textwidth]{ImageRegistration5Output2}
606 // \includegraphics[width=0.32\textwidth]{ImageRegistration5DifferenceBefore2}
607 // \includegraphics[width=0.32\textwidth]{ImageRegistration5DifferenceAfter2}
608 // \itkcaption[Rigid2D Registration output images]{Resampled moving image
609 // (left). Differences between the fixed and moving images, before (center)
610 // and after (right) registration with the CenteredRigid2D transform.}
611 // \label{fig:ImageRegistration5Outputs2}
612 // \end{figure}
613 //
614 // Figure \ref{fig:ImageRegistration5Outputs2} shows the output of the
615 // registration. The rightmost image of this figure shows the difference
616 // between the fixed image and the resampled moving image after registration.
617 //
618 // \begin{figure}
619 // \center
620 // \includegraphics[height=0.32\textwidth]{ImageRegistration5TraceMetric2}
621 // \includegraphics[height=0.32\textwidth]{ImageRegistration5TraceAngle2}
622 // \includegraphics[height=0.32\textwidth]{ImageRegistration5TraceTranslations2}
623 // \itkcaption[Rigid2D Registration output plots]{Metric values, rotation
624 // angle and translations during the registration using the CenteredRigid2D
625 // transform on an image with rotation and translation mis-registration.}
626 // \label{fig:ImageRegistration5Plots2}
627 // \end{figure}
628 //
629 // Figure \ref{fig:ImageRegistration5Plots2} shows plots of the main output
630 // registration parameters when the rotation and translations are combined.
631 // These results include the metric values at every iteration, the angle
632 // values at every iteration, and the translation components of the
633 // registration as the registration converges. It can be seen from the
634 // smoothness of these plots that a larger step length could have been
635 // supported easily by the optimizer. You may want to modify this value in
636 // order to get a better idea of how to tune the parameters.
637 //
638
639
640 return EXIT_SUCCESS;
641 }
642