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: {BrainProtonDensitySliceR10X13Y17.png}
21 // OUTPUTS: {ImageRegistration6Output.png}
22 // OUTPUTS: {ImageRegistration6DifferenceBefore.png}
23 // OUTPUTS: {ImageRegistration6DifferenceAfter.png}
24
25
26 //
27 // This example illustrates the use of the \doxygen{CenteredRigid2DTransform}
28 // for performing registration. The example code is for the most part
29 // identical to the one presented in Section~\ref{sec:RigidRegistrationIn2D}.
30 // Even though this current example is done in $2D$, the class
31 // \doxygen{CenteredTransformInitializer} is quite generic and could be used
32 // in other dimensions. The objective of the initializer class is to simplify
33 // the computation of the center of rotation and the translation required to
34 // initialize certain transforms such as the
35 // CenteredRigid2DTransform. The initializer accepts two images and
36 // a transform as inputs. The images are considered to be the fixed and
37 // moving images of the registration problem, while the transform is the one
38 // used to register the images.
39 //
40 // The CenteredRigid2DTransform supports two modes of operation. In the first
41 // mode, the centers of the images are computed as space coordinates using the
42 // image origin, size and spacing. The center of the fixed image is assigned as
43 // the rotational center of the transform while the vector going from the fixed
44 // image center to the moving image center is passed as the initial translation
45 // of the transform. In the second mode, the image centers are not computed
46 // geometrically but by using the moments of the intensity gray levels. The
47 // center of mass of each image is computed using the helper class
48 // \doxygen{ImageMomentsCalculator}. The center of mass of the fixed image is
49 // passed as the rotational center of the transform while the vector going from
50 // the fixed image center of mass to the moving image center of mass is passed
51 // as the initial translation of the transform. This second mode of operation
52 // is quite convenient when the anatomical structures of interest are not
53 // centered in the image. In such cases the alignment of the centers of mass
54 // provides a better rough initial registration than the simple use of the
55 // geometrical centers. The validity of the initial registration should be
56 // questioned when the two images are acquired in different imaging modalities.
57 // In those cases, the center of mass of intensities in one modality does not
58 // necessarily matches the center of mass of intensities in the other imaging
59 // modality.
60 //
61 // \index{itk::CenteredRigid2DTransform}
62 // \index{itk::ImageMomentsCalculator}
63 //
64 //
65
66 #include "itkImageRegistrationMethod.h"
67 #include "itkMeanSquaresImageToImageMetric.h"
68 #include "itkRegularStepGradientDescentOptimizer.h"
69
70
71 //
72 // The following are the most relevant headers in this example.
73 //
74 // \index{itk::CenteredRigid2DTransform!header}
75 //
76
77 #include "itkCenteredRigid2DTransform.h"
78 #include "itkCenteredTransformInitializer.h"
79
80
81 #include "itkImageFileReader.h"
82 #include "itkImageFileWriter.h"
83
84 #include "itkResampleImageFilter.h"
85 #include "itkCastImageFilter.h"
86 #include "itkRescaleIntensityImageFilter.h"
87 #include "itkSubtractImageFilter.h"
88
89
90 //
91 // The following section of code implements a command observer
92 // that will monitor the evolution of the registration process.
93 //
94 #include "itkCommand.h"
95 class CommandIterationUpdate : public itk::Command
96 {
97 public:
98 using Self = CommandIterationUpdate;
99 using Superclass = itk::Command;
100 using Pointer = itk::SmartPointer<Self>;
101 itkNewMacro( Self );
102
103 protected:
104 CommandIterationUpdate() = default;
105
106 public:
107 using OptimizerType = itk::RegularStepGradientDescentOptimizer;
108 using OptimizerPointer = const OptimizerType *;
109
Execute(itk::Object * caller,const itk::EventObject & event)110 void Execute(itk::Object *caller, const itk::EventObject & event) override
111 {
112 Execute( (const itk::Object *)caller, event);
113 }
114
Execute(const itk::Object * object,const itk::EventObject & event)115 void Execute(const itk::Object * object, const itk::EventObject & event) override
116 {
117 auto optimizer = static_cast< OptimizerPointer >( object );
118 if( ! itk::IterationEvent().CheckEvent( &event ) )
119 {
120 return;
121 }
122 std::cout << optimizer->GetCurrentIteration() << " ";
123 std::cout << optimizer->GetValue() << " ";
124 std::cout << optimizer->GetCurrentPosition() << std::endl;
125 }
126 };
127
128 #include "itkTestDriverIncludeRequiredIOFactories.h"
main(int argc,char * argv[])129 int main( int argc, char *argv[] )
130 {
131 RegisterRequiredFactories();
132 if( argc < 4 )
133 {
134 std::cerr << "Missing Parameters " << std::endl;
135 std::cerr << "Usage: " << argv[0];
136 std::cerr << " fixedImageFile movingImageFile ";
137 std::cerr << " outputImagefile [differenceBeforeRegistration] ";
138 std::cerr << " [differenceAfterRegistration] "<< std::endl;
139 return EXIT_FAILURE;
140 }
141
142 constexpr unsigned int Dimension = 2;
143 using PixelType = float;
144
145 using FixedImageType = itk::Image< PixelType, Dimension >;
146 using MovingImageType = itk::Image< PixelType, Dimension >;
147
148
149 //
150 // The transform type is instantiated using the code below. The only
151 // template parameter of this class is the representation type of the
152 // space coordinates.
153 //
154 // \index{itk::CenteredRigid2DTransform!Instantiation}
155 //
156
157 using TransformType = itk::CenteredRigid2DTransform< double >;
158
159
160 using OptimizerType = itk::RegularStepGradientDescentOptimizer;
161 using MetricType = itk::MeanSquaresImageToImageMetric<
162 FixedImageType,
163 MovingImageType >;
164 using InterpolatorType = itk:: LinearInterpolateImageFunction<
165 MovingImageType,
166 double >;
167 using RegistrationType = itk::ImageRegistrationMethod<
168 FixedImageType,
169 MovingImageType >;
170
171 MetricType::Pointer metric = MetricType::New();
172 OptimizerType::Pointer optimizer = OptimizerType::New();
173 InterpolatorType::Pointer interpolator = InterpolatorType::New();
174 RegistrationType::Pointer registration = RegistrationType::New();
175
176
177 registration->SetMetric( metric );
178 registration->SetOptimizer( optimizer );
179 registration->SetInterpolator( interpolator );
180
181
182 //
183 // Like the previous section, a direct initialization method is used here.
184 // The transform object is constructed below. This transform will
185 // be initialized, and its initial parameters will be considered as
186 // the parameters to be used when the registration process begins.
187 //
188 // \index{itk::CenteredRigid2DTransform!New()}
189 // \index{itk::CenteredRigid2DTransform!Pointer}
190 // \index{itk::RegistrationMethod!SetTransform()}
191 //
192
193 TransformType::Pointer transform = TransformType::New();
194 registration->SetTransform( transform );
195
196 using FixedImageReaderType = itk::ImageFileReader< FixedImageType >;
197 using MovingImageReaderType = itk::ImageFileReader< MovingImageType >;
198 FixedImageReaderType::Pointer fixedImageReader = FixedImageReaderType::New();
199 MovingImageReaderType::Pointer movingImageReader = MovingImageReaderType::New();
200
201 fixedImageReader->SetFileName( argv[1] );
202 movingImageReader->SetFileName( argv[2] );
203
204
205 registration->SetFixedImage( fixedImageReader->GetOutput() );
206 registration->SetMovingImage( movingImageReader->GetOutput() );
207 fixedImageReader->Update();
208
209 registration->SetFixedImageRegion(
210 fixedImageReader->GetOutput()->GetBufferedRegion() );
211
212
213 //
214 // The input images are taken from readers. It is not necessary to
215 // explicitly call \code{Update()} on the readers since the
216 // CenteredTransformInitializer class will do it as part of its
217 // initialization. The following code instantiates the initializer. This
218 // class is templated over the fixed and moving images type as well as the
219 // transform type. An initializer is then constructed by calling the
220 // \code{New()} method and assigning the result to a
221 // \doxygen{SmartPointer}.
222 //
223 // \index{itk::CenteredRigid2DTransform!Instantiation}
224 // \index{itk::CenteredRigid2DTransform!New()}
225 // \index{itk::CenteredRigid2DTransform!SmartPointer}
226 //
227
228 using TransformInitializerType = itk::CenteredTransformInitializer<
229 TransformType,
230 FixedImageType,
231 MovingImageType >;
232
233 TransformInitializerType::Pointer initializer =
234 TransformInitializerType::New();
235
236 //
237 // The initializer is now connected to the transform and to the fixed and
238 // moving images.
239 //
240
241
242 initializer->SetTransform( transform );
243 initializer->SetFixedImage( fixedImageReader->GetOutput() );
244 initializer->SetMovingImage( movingImageReader->GetOutput() );
245
246 //
247 // The use of the geometrical centers is selected by calling
248 // \code{GeometryOn()} while the use of center of mass is selected by
249 // calling \code{MomentsOn()}. Below we select the center of mass mode.
250 //
251 // \index{CenteredTransformInitializer!MomentsOn()}
252 // \index{CenteredTransformInitializer!GeometryOn()}
253 //
254
255 initializer->MomentsOn();
256
257
258 //
259 // Finally, the computation of the center and translation is triggered by
260 // the \code{InitializeTransform()} method. The resulting values will be
261 // passed directly to the transform.
262 //
263
264
265 initializer->InitializeTransform();
266
267
268 //
269 // The remaining parameters of the transform are initialized as before.
270 //
271
272 transform->SetAngle( 0.0 );
273
274
275 //
276 // Now the initialized transform object will be set to the registration method,
277 // and the starting point of the registration is defined by its initial parameters.
278 //
279
280 registration->SetInitialTransformParameters( transform->GetParameters() );
281
282
283 using OptimizerScalesType = OptimizerType::ScalesType;
284 OptimizerScalesType optimizerScales( transform->GetNumberOfParameters() );
285 const double translationScale = 1.0 / 1000.0;
286
287 optimizerScales[0] = 1.0;
288 optimizerScales[1] = translationScale;
289 optimizerScales[2] = translationScale;
290 optimizerScales[3] = translationScale;
291 optimizerScales[4] = translationScale;
292
293 optimizer->SetScales( optimizerScales );
294
295 optimizer->SetMaximumStepLength( 0.1 );
296 optimizer->SetMinimumStepLength( 0.001 );
297 optimizer->SetNumberOfIterations( 200 );
298
299
300 // Create the Command observer and register it with the optimizer.
301 //
302 CommandIterationUpdate::Pointer observer = CommandIterationUpdate::New();
303 optimizer->AddObserver( itk::IterationEvent(), observer );
304
305 try
306 {
307 registration->Update();
308 std::cout << "Optimizer stop condition: "
309 << registration->GetOptimizer()->GetStopConditionDescription()
310 << std::endl;
311 }
312 catch( itk::ExceptionObject & err )
313 {
314 std::cerr << "ExceptionObject caught !" << std::endl;
315 std::cerr << err << std::endl;
316 return EXIT_FAILURE;
317 }
318
319 OptimizerType::ParametersType finalParameters =
320 registration->GetLastTransformParameters();
321
322
323 const double finalAngle = finalParameters[0];
324 const double finalRotationCenterX = finalParameters[1];
325 const double finalRotationCenterY = finalParameters[2];
326 const double finalTranslationX = finalParameters[3];
327 const double finalTranslationY = finalParameters[4];
328
329 const unsigned int numberOfIterations = optimizer->GetCurrentIteration();
330 const double bestValue = optimizer->GetValue();
331
332 // Print out results
333 //
334 const double finalAngleInDegrees = finalAngle * 180.0 / itk::Math::pi;
335
336 std::cout << "Result = " << std::endl;
337 std::cout << " Angle (radians) " << finalAngle << std::endl;
338 std::cout << " Angle (degrees) " << finalAngleInDegrees << std::endl;
339 std::cout << " Center X = " << finalRotationCenterX << std::endl;
340 std::cout << " Center Y = " << finalRotationCenterY << std::endl;
341 std::cout << " Translation X = " << finalTranslationX << std::endl;
342 std::cout << " Translation Y = " << finalTranslationY << std::endl;
343 std::cout << " Iterations = " << numberOfIterations << std::endl;
344 std::cout << " Metric value = " << bestValue << std::endl;
345
346
347 //
348 // Let's execute this example over some of the images provided in
349 // \code{Examples/Data}, for example:
350 //
351 // \begin{itemize}
352 // \item \code{BrainProtonDensitySliceBorder20.png}
353 // \item \code{BrainProtonDensitySliceR10X13Y17.png}
354 // \end{itemize}
355 //
356 // The second image is the result of intentionally rotating the first
357 // image by $10$ degrees and shifting it $13mm$ in $X$ and $17mm$ in
358 // $Y$. Both images have unit-spacing and are shown in Figure
359 // \ref{fig:FixedMovingImageRegistration5}. The registration takes $22$
360 // iterations and produces:
361 //
362 // \begin{center}
363 // \begin{verbatim}
364 // [0.174475, 111.177, 131.572, 12.4566, 16.0729]
365 // \end{verbatim}
366 // \end{center}
367 //
368 // These parameters are interpreted as
369 //
370 // \begin{itemize}
371 // \item Angle = $0.174475$ radians
372 // \item Center = $( 111.177 , 131.572 )$ millimeters
373 // \item Translation = $( 12.4566 , 16.0729 )$ millimeters
374 // \end{itemize}
375 //
376 // Note that the reported translation is not the translation of $(13,17)$
377 // that might be expected. The reason is that the five parameters of the
378 // CenteredRigid2DTransform are redundant. The actual movement
379 // in space is described by only $3$ parameters. This means that there are
380 // infinite combinations of rotation center and translations that will
381 // represent the same actual movement in space. It is more illustrative in
382 // this case to take a look at the actual rotation matrix and offset
383 // resulting form the five parameters.
384 //
385
386 transform->SetParameters( finalParameters );
387
388 TransformType::MatrixType matrix = transform->GetMatrix();
389 TransformType::OffsetType offset = transform->GetOffset();
390
391 std::cout << "Matrix = " << std::endl << matrix << std::endl;
392 std::cout << "Offset = " << std::endl << offset << std::endl;
393
394 //
395 // Which produces the following output.
396 //
397 // \begin{verbatim}
398 // Matrix =
399 // 0.984818 -0.173591
400 // 0.173591 0.984818
401 //
402 // Offset =
403 // [36.9843, -1.22896]
404 // \end{verbatim}
405 //
406 // This output illustrates how counter-intuitive the mix of center of
407 // rotation and translations can be. Figure
408 // \ref{fig:TranslationAndRotationCenter} will clarify this situation. The
409 // figure shows the original image on the left. A rotation of $10^{\circ}$
410 // around the center of the image is shown in the middle. The same rotation
411 // performed around the origin of coordinates is shown on the right. It can
412 // be seen here that changing the center of rotation introduces additional
413 // translations.
414 //
415 // Let's analyze what happens to the center of the image that we just
416 // registered. Under the point of view of rotating $10^{\circ}$ around the
417 // center and then applying a translation of $(13mm,17mm)$. The image has
418 // a size of $(221 \times 257)$ pixels and unit spacing. Hence its center
419 // has coordinates $(110.5,128.5)$. Since the rotation is done around this
420 // point, the center behaves as the fixed point of the transformation and
421 // remains unchanged. Then with the $(13mm,17mm)$ translation it is mapped
422 // to $(123.5,145.5)$ which becomes its final position.
423 //
424 // The matrix and offset that we obtained at the end of the registration
425 // indicate that this should be equivalent to a rotation of $10^{\circ}$
426 // around the origin, followed by a translations of $(36.98,-1.22)$. Let's
427 // compute this in detail. First the rotation of the image center by
428 // $10^{\circ}$ around the origin will move the point to
429 // $(86.52,147.97)$. Now, applying a translation of $(36.98,-1.22)$ maps
430 // this point to $(123.5,146.75)$. Which is close to the result of our
431 // previous computation.
432 //
433 // It is unlikely that we could have chosen such translations as the
434 // initial guess, since we tend to think about image in a coordinate
435 // system whose origin is in the center of the image.
436 //
437 // \begin{figure}
438 // \center
439 // \includegraphics[width=\textwidth]{TranslationAndRotationCenter}
440 // \itkcaption[Effect of changing the center of rotation]{Effect of changing
441 // the center of rotation.}
442 // \label{fig:TranslationAndRotationCenter}
443 // \end{figure}
444 //
445
446
447 //
448 // You may be wondering why the actual movement is represented by three
449 // parameters when we take the trouble of using five. In particular, why
450 // use a $5$-dimensional optimizer space instead of a $3$-dimensional
451 // one. The answer is that by using five parameters we have a much simpler
452 // way of initializing the transform with the rotation matrix and
453 // offset. Using the minimum three parameters it is not obvious how to
454 // determine what the initial rotation and translations should be.
455 //
456
457
458 //
459 // \begin{figure}
460 // \center
461 // \includegraphics[width=0.44\textwidth]{BrainProtonDensitySliceBorder20}
462 // \includegraphics[width=0.44\textwidth]{BrainProtonDensitySliceR10X13Y17}
463 // \itkcaption[CenteredTransformInitializer input images]{Fixed and moving
464 // images provided as input to the registration method using
465 // CenteredTransformInitializer.}
466 // \label{fig:FixedMovingImageRegistration6}
467 // \end{figure}
468 //
469 //
470 // \begin{figure}
471 // \center
472 // \includegraphics[width=0.32\textwidth]{ImageRegistration6Output}
473 // \includegraphics[width=0.32\textwidth]{ImageRegistration6DifferenceBefore}
474 // \includegraphics[width=0.32\textwidth]{ImageRegistration6DifferenceAfter}
475 // \itkcaption[CenteredTransformInitializer output images]{Resampled moving
476 // image (left). Differences between fixed and moving images, before
477 // registration (center) and after registration (right) with the
478 // CenteredTransformInitializer.}
479 // \label{fig:ImageRegistration6Outputs}
480 // \end{figure}
481 //
482 // Figure \ref{fig:ImageRegistration6Outputs} shows the output of the
483 // registration. The image on the right of this figure shows the differences
484 // between the fixed image and the resampled moving image after registration.
485 //
486 // \begin{figure}
487 // \center
488 // \includegraphics[height=0.32\textwidth]{ImageRegistration6TraceMetric}
489 // \includegraphics[height=0.32\textwidth]{ImageRegistration6TraceAngle}
490 // \includegraphics[height=0.32\textwidth]{ImageRegistration6TraceTranslations}
491 // \itkcaption[CenteredTransformInitializer output plots]{Plots of the Metric,
492 // rotation angle, center of rotation and translations during the
493 // registration using CenteredTransformInitializer.}
494 // \label{fig:ImageRegistration6Plots}
495 // \end{figure}
496 //
497 // Figure \ref{fig:ImageRegistration6Plots} plots the output parameters of
498 // the registration process. It includes, the metric values at every
499 // iteration, the angle values at every iteration, and the values of the
500 // translation components as the registration progress. Note that this is
501 // the complementary translation as used in the transform, not the actual
502 // total translation that is used in the transform offset. We could modify
503 // the observer to print the total offset instead of printing the array of
504 // parameters. Let's call that an exercise for the reader!
505 //
506
507
508 using ResampleFilterType = itk::ResampleImageFilter<
509 MovingImageType,
510 FixedImageType >;
511 TransformType::Pointer finalTransform = TransformType::New();
512 finalTransform->SetParameters( finalParameters );
513 finalTransform->SetFixedParameters( transform->GetFixedParameters() );
514
515 ResampleFilterType::Pointer resample = ResampleFilterType::New();
516
517 resample->SetTransform( finalTransform );
518 resample->SetInput( movingImageReader->GetOutput() );
519
520 FixedImageType::Pointer fixedImage = fixedImageReader->GetOutput();
521
522 resample->SetSize( fixedImage->GetLargestPossibleRegion().GetSize() );
523 resample->SetOutputOrigin( fixedImage->GetOrigin() );
524 resample->SetOutputSpacing( fixedImage->GetSpacing() );
525 resample->SetOutputDirection( fixedImage->GetDirection() );
526 resample->SetDefaultPixelValue( 100 );
527
528 using OutputPixelType = unsigned char;
529
530 using OutputImageType = itk::Image< OutputPixelType, Dimension >;
531
532 using CastFilterType = itk::CastImageFilter<
533 FixedImageType,
534 OutputImageType >;
535
536 using WriterType = itk::ImageFileWriter< OutputImageType >;
537
538
539 WriterType::Pointer writer = WriterType::New();
540 CastFilterType::Pointer caster = CastFilterType::New();
541
542
543 writer->SetFileName( argv[3] );
544
545
546 caster->SetInput( resample->GetOutput() );
547 writer->SetInput( caster->GetOutput() );
548 writer->Update();
549
550 // Now compute the difference between the images
551 // before and after registration.
552 //
553 using DifferenceImageType = itk::Image< float, Dimension >;
554
555 using DifferenceFilterType = itk::SubtractImageFilter<
556 FixedImageType,
557 FixedImageType,
558 DifferenceImageType >;
559
560 DifferenceFilterType::Pointer difference = DifferenceFilterType::New();
561
562 using OutputPixelType = unsigned char;
563
564 using OutputImageType = itk::Image< OutputPixelType, Dimension >;
565
566 using RescalerType = itk::RescaleIntensityImageFilter<
567 DifferenceImageType,
568 OutputImageType >;
569
570 RescalerType::Pointer intensityRescaler = RescalerType::New();
571
572 intensityRescaler->SetOutputMinimum( 0 );
573 intensityRescaler->SetOutputMaximum( 255 );
574
575 difference->SetInput1( fixedImageReader->GetOutput() );
576 difference->SetInput2( resample->GetOutput() );
577
578 resample->SetDefaultPixelValue( 1 );
579
580 intensityRescaler->SetInput( difference->GetOutput() );
581
582 using WriterType = itk::ImageFileWriter< OutputImageType >;
583
584 WriterType::Pointer writer2 = WriterType::New();
585
586 writer2->SetInput( intensityRescaler->GetOutput() );
587
588
589 try
590 {
591 // Compute the difference image between the
592 // fixed and moving image after registration.
593 if( argc > 5 )
594 {
595 writer2->SetFileName( argv[5] );
596 writer2->Update();
597 }
598
599 // Compute the difference image between the
600 // fixed and resampled moving image after registration.
601 TransformType::Pointer identityTransform = TransformType::New();
602 identityTransform->SetIdentity();
603 resample->SetTransform( identityTransform );
604 if( argc > 4 )
605 {
606 writer2->SetFileName( argv[4] );
607 writer2->Update();
608 }
609 }
610 catch( itk::ExceptionObject & excp )
611 {
612 std::cerr << "Error while writing difference images" << std::endl;
613 std::cerr << excp << std::endl;
614 return EXIT_FAILURE;
615 }
616
617 return EXIT_SUCCESS;
618 }
619