File: itkPointSetToImageRegistrationTest.cxx

package info (click to toggle)
insighttoolkit4 4.13.3withdata-dfsg2-4
  • links: PTS, VCS
  • area: main
  • in suites: bookworm
  • size: 491,256 kB
  • sloc: cpp: 557,600; ansic: 180,546; fortran: 34,788; python: 16,572; sh: 2,187; lisp: 2,070; tcl: 993; java: 362; perl: 200; makefile: 133; csh: 81; pascal: 69; xml: 19; ruby: 10
file content (214 lines) | stat: -rw-r--r-- 7,298 bytes parent folder | download | duplicates (4)
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
/*=========================================================================
 *
 *  Copyright Insight Software Consortium
 *
 *  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
 *
 *         http://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 "itkTranslationTransform.h"
#include "itkLinearInterpolateImageFunction.h"
#include "itkNormalizedCorrelationPointSetToImageMetric.h"
#include "itkRegularStepGradientDescentOptimizer.h"
#include "itkPointSet.h"
#include "itkPointSetToImageRegistrationMethod.h"
#include "itkImageRegistrationMethodImageSource.h"
#include "itkTestingMacros.h"

#include <iostream>

/**
 *
 *  This program tests the registration of a PointSet against an image.
 *
 *  This test uses two 2D-Gaussians (standard deviation RegionSize/2)
 *  One is shifted by 5 pixels from the other.
 *
 */

int itkPointSetToImageRegistrationTest( int, char* [] )
{
  const unsigned int ImageDimension = 2;

  typedef double                   PixelType;

  typedef double                   CoordinateRepresentationType;

  typedef itk::Image< PixelType, ImageDimension >   MovingImageType;
  typedef itk::Image< PixelType, ImageDimension >   FixedImageType;

  typedef itk::testhelper::ImageRegistrationMethodImageSource<
    PixelType,
    PixelType,
    ImageDimension > ImageSourceType;

  ImageSourceType::Pointer imageSource = ImageSourceType::New();

  itk::Size< ImageDimension > size;
  size[0] = 100;
  size[1] = 100;

  imageSource->GenerateImages( size );

  // Create the two images
  MovingImageType::ConstPointer movingImage = imageSource->GetMovingImage();
  FixedImageType::ConstPointer  fixedImage  = imageSource->GetFixedImage();

  // Create the point set and load it with data by sampling
  // the fixed image.
  //
  typedef itk::PointSet< float, ImageDimension > FixedPointSetType;
  FixedPointSetType::Pointer fixedPointSet = FixedPointSetType::New();

  const unsigned int numberOfPoints = 10000;

  fixedPointSet->SetPointData( FixedPointSetType::PointDataContainer::New() );

  fixedPointSet->GetPoints()->Reserve( numberOfPoints );
  fixedPointSet->GetPointData()->Reserve( numberOfPoints );

  typedef itk::ImageRegionConstIterator< FixedImageType > ImageIteratorType;

  ImageIteratorType it( fixedImage, fixedImage->GetBufferedRegion() );

  const unsigned int skip =
    fixedImage->GetBufferedRegion().GetNumberOfPixels() / numberOfPoints;

  unsigned int counter = 0;

  FixedPointSetType::PointIdentifier pointId = 0;
  FixedPointSetType::PointType  point;

  it.GoToBegin();
  while( !it.IsAtEnd() )
    {
    if( counter == 0 )
      {
      fixedImage->TransformIndexToPhysicalPoint( it.GetIndex(), point );
      fixedPointSet->SetPoint( pointId, point );
      fixedPointSet->SetPointData( pointId, it.Get() );
      ++pointId;
      if( pointId == numberOfPoints )
        {
        break;
        }
      counter = skip;
      }
    --counter;
    ++it;
    }

  // Set up the Metric
  typedef itk::NormalizedCorrelationPointSetToImageMetric< FixedPointSetType,
    MovingImageType > MetricType;

  typedef MetricType::TransformType                 TransformBaseType;
  typedef TransformBaseType::ParametersType         ParametersType;

  MetricType::Pointer metric = MetricType::New();

  // Set up the Transform
  typedef itk::TranslationTransform< CoordinateRepresentationType,
    ImageDimension > TransformType;
  TransformType::Pointer transform = TransformType::New();

  // Set up the Interpolator
  typedef itk::LinearInterpolateImageFunction< MovingImageType, double >
    InterpolatorType;
  InterpolatorType::Pointer interpolator = InterpolatorType::New();

  // Set up the Optimizer
  typedef itk::RegularStepGradientDescentOptimizer OptimizerType;
  OptimizerType::Pointer optimizer = OptimizerType::New();

  // Set up the Registration method
  typedef itk::PointSetToImageRegistrationMethod< FixedPointSetType,
    MovingImageType > RegistrationType;

  RegistrationType::Pointer registration = RegistrationType::New();

  EXERCISE_BASIC_OBJECT_METHODS( registration, PointSetToImageRegistrationMethod,
    ProcessObject );

  typedef itk::CommandIterationUpdate< OptimizerType > CommandIterationType;

  // Instantiate an Observer to report the progress of the Optimization
  CommandIterationType::Pointer iterationCommand = CommandIterationType::New();
  iterationCommand->SetOptimizer( optimizer.GetPointer() );

  // Scale the translation components of the Transform in the Optimizer
  OptimizerType::ScalesType scales( transform->GetNumberOfParameters() );
  scales.Fill( 1.0 );

  unsigned long   numberOfIterations =   50;
  double          maximumStepLenght  =  1.0;  // no step will be larger than this
  double          minimumStepLenght  =  0.01;
  double          gradientTolerance  =  1e-6; // convergence criterion

  optimizer->SetScales( scales );
  optimizer->SetNumberOfIterations( numberOfIterations );
  optimizer->SetMinimumStepLength( minimumStepLenght );
  optimizer->SetMaximumStepLength( maximumStepLenght );
  optimizer->SetGradientMagnitudeTolerance( gradientTolerance );
  optimizer->MinimizeOn();

  // Connect all the components required for the registration
  registration->SetMetric( metric );
  TEST_SET_GET_VALUE( metric, registration->GetMetric() );

  registration->SetOptimizer( optimizer );
  TEST_SET_GET_VALUE( optimizer, registration->GetOptimizer() );

  registration->SetTransform( transform );
  TEST_SET_GET_VALUE( transform, registration->GetTransform() );

  registration->SetFixedPointSet( fixedPointSet );
  TEST_SET_GET_VALUE( fixedPointSet, registration->GetFixedPointSet() );

  registration->SetMovingImage( movingImage );
  TEST_SET_GET_VALUE( movingImage, registration->GetMovingImage() );

  registration->SetInterpolator( interpolator );
  TEST_SET_GET_VALUE( interpolator, registration->GetInterpolator() );

  // Set up transform parameters
  ParametersType parameters( transform->GetNumberOfParameters() );

  // Initialize the offset/vector part
  for( unsigned int k = 0; k < parameters.size(); k++ )
    {
    parameters[k] = 0.0f;
    }
  transform->SetParameters( parameters );
  registration->SetInitialTransformParameters( transform->GetParameters() );
  TEST_SET_GET_VALUE( transform->GetParameters(),
    registration->GetInitialTransformParameters() );


  TRY_EXPECT_NO_EXCEPTION( registration->Update() );


  // Print the last transform parameters to improve coverage
  //
  ParametersType finalParameters = registration->GetLastTransformParameters();

  const unsigned int numberOfParameters = parameters.Size();

  std::cout << "Last Transform Parameters: " << std::endl;
  for( unsigned int i = 0; i < numberOfParameters; ++i )
    {
    std::cout << finalParameters[i] << std::endl;
    }

  return EXIT_SUCCESS;
}