File: itkRegularStepGradientDescentOptimizerv4Test.cxx

package info (click to toggle)
insighttoolkit4 4.10.1-dfsg1-1.1
  • links: PTS, VCS
  • area: main
  • in suites: stretch
  • size: 416,780 kB
  • ctags: 104,347
  • sloc: cpp: 553,142; ansic: 142,389; fortran: 34,788; python: 16,392; lisp: 2,070; sh: 1,862; tcl: 993; java: 362; perl: 200; makefile: 111; csh: 81; pascal: 69; xml: 19; ruby: 10
file content (364 lines) | stat: -rw-r--r-- 11,357 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
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
/*=========================================================================
 *
 *  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 "itkRegularStepGradientDescentOptimizerv4.h"
#include "itkMath.h"

/**
 *  The objectif function is the quadratic form:
 *
 *  1/2 x^T A x - b^T x
 *
 *  Where A is a matrix and b is a vector
 *  The system in this example is:
 *
 *     | 3  2 ||x|   | 2|   |0|
 *     | 2  6 ||y| + |-8| = |0|
 *
 *
 *   the solution is the vector | 2 -2 |
 *
 * \class RSGv4TestMetric
 */
class RSGv4TestMetric : public itk::ObjectToObjectMetricBase
{
public:

  typedef RSGv4TestMetric               Self;
  typedef itk::ObjectToObjectMetricBase Superclass;
  typedef itk::SmartPointer<Self>       Pointer;
  typedef itk::SmartPointer<const Self> ConstPointer;
  itkNewMacro( Self );

  enum { SpaceDimension=2 };

  typedef Superclass::ParametersType      ParametersType;
  typedef Superclass::DerivativeType      DerivativeType;
  typedef Superclass::ParametersValueType ParametersValueType;
  typedef Superclass::MeasureType         MeasureType;


  RSGv4TestMetric()
  {
    m_Parameters.SetSize( SpaceDimension );
    m_Parameters.Fill( 0 );
  }

  virtual void Initialize(void) throw ( itk::ExceptionObject ) ITK_OVERRIDE {}

  virtual void GetDerivative( DerivativeType & derivative ) const ITK_OVERRIDE
  {
    MeasureType value;
    GetValueAndDerivative( value, derivative );
  }

  void GetValueAndDerivative( MeasureType & value,
                              DerivativeType & derivative ) const ITK_OVERRIDE
  {
    if( derivative.Size() != 2 )
      {
      derivative.SetSize(2);
      }

    double x = m_Parameters[0];
    double y = m_Parameters[1];

    std::cout << "GetValueAndDerivative( ";
    std::cout << x << " ";
    std::cout << y << ") = ";

    value = 0.5*(3*x*x+4*x*y+6*y*y) - 2*x + 8*y;

    std::cout << "value: " << value << std::endl;

    /* The optimizer simply takes the derivative from the metric
     * and adds it to the transform after scaling. So instead of
     * setting a 'minimize' option in the gradient, we return
     * a minimizing derivative. */
    derivative[0] = -( 3 * x + 2 * y -2 );
    derivative[1] = -( 2 * x + 6 * y +8 );

    std::cout << "derivative: " << derivative << std::endl;
  }

  virtual MeasureType  GetValue() const ITK_OVERRIDE
  {
    return 0.0;
  }

  virtual void UpdateTransformParameters( const DerivativeType & update, ParametersValueType factor ) ITK_OVERRIDE
  {
    m_Parameters += update * factor;
  }

  virtual unsigned int GetNumberOfParameters(void) const ITK_OVERRIDE
  {
    return SpaceDimension;
  }

  virtual bool HasLocalSupport() const ITK_OVERRIDE
  {
    return false;
  }

  virtual unsigned int GetNumberOfLocalParameters() const ITK_OVERRIDE
  {
    return SpaceDimension;
  }

  /* These Set/Get methods are only needed for this test derivation that
   * isn't using a transform */
  virtual void SetParameters( ParametersType & parameters ) ITK_OVERRIDE
  {
    m_Parameters = parameters;
  }

  virtual const ParametersType & GetParameters() const ITK_OVERRIDE
  {
    return m_Parameters;
  }

private:
  ParametersType m_Parameters;
};

///////////////////////////////////////////////////////////
int itkRegularStepGradientDescentOptimizerv4Test(int, char* [] )
{
  std::cout << "RegularStepGradientDescentOptimizerv4 Test ";
  std::cout << std::endl << std::endl;

  typedef  itk::RegularStepGradientDescentOptimizerv4<double>  OptimizerType;

  typedef  OptimizerType::ScalesType  ScalesType;


  // Declaration of a itkOptimizer
  OptimizerType::Pointer  itkOptimizer = OptimizerType::New();


  // Declaration of the metric
  RSGv4TestMetric::Pointer metric = RSGv4TestMetric::New();


  itkOptimizer->SetMetric( metric );


  typedef RSGv4TestMetric::ParametersType    ParametersType;


  const unsigned int spaceDimension =
                      metric->GetNumberOfParameters();

  // We start not so far from  | 2 -2 |
  ParametersType  initialPosition( spaceDimension );
  initialPosition[0] =  100;
  initialPosition[1] = -100;
  metric->SetParameters( initialPosition );

  ScalesType    parametersScale( spaceDimension );
  parametersScale[0] = 1.0;
  parametersScale[1] = 1.0;

  itkOptimizer->SetScales( parametersScale );
  itkOptimizer->SetGradientMagnitudeTolerance( 1e-6 );
  itkOptimizer->SetLearningRate( 100 );
  itkOptimizer->SetMinimumStepLength( 1e-6 );
  itkOptimizer->SetNumberOfIterations( 900 );

  try
    {
    std::cout << "currentPosition before optimization: " << itkOptimizer->GetMetric()->GetParameters() << std::endl;
    itkOptimizer->StartOptimization();
    std::cout << "currentPosition after optimization: " << itkOptimizer->GetMetric()->GetParameters() << std::endl;
    std::cout << " Stop Condition  = " << itkOptimizer->GetStopConditionDescription() << std::endl;
    }
  catch( itk::ExceptionObject & e )
    {
    std::cout << "Exception thrown ! " << std::endl;
    std::cout << "An error occurred during Optimization" << std::endl;
    std::cout << "Location    = " << e.GetLocation()    << std::endl;
    std::cout << "Description = " << e.GetDescription() << std::endl;
    return EXIT_FAILURE;
    }


  ParametersType finalPosition = itkOptimizer->GetMetric()->GetParameters();
  std::cout << "Solution        = (";
  std::cout << finalPosition[0] << ",";
  std::cout << finalPosition[1] << ")" << std::endl;

  //
  // check results to see if it is within range
  //
  bool pass = true;
  double trueParameters[2] = { 2, -2 };
  for( unsigned int j = 0; j < 2; j++ )
    {
    if( itk::Math::FloatAlmostEqual( finalPosition[j], trueParameters[j] ) )
      {
      pass = false;
      }
    }

  if( !pass )
    {
    std::cout << "Test failed." << std::endl;
    return EXIT_FAILURE;
    }


  // Run now with a different relaxation factor
  std::cout << "\nRun test with a different relaxation factor: 0.8, instead of default value: 0.5." << std::endl;
  {
  metric->SetParameters( initialPosition );

  itkOptimizer->SetRelaxationFactor( 0.8 );
  try
    {
    std::cout << "currentPosition before optimization: " << itkOptimizer->GetMetric()->GetParameters() << std::endl;
    itkOptimizer->StartOptimization();
    std::cout << "currentPosition after optimization: " << itkOptimizer->GetMetric()->GetParameters() << std::endl;
    std::cout << " Stop Condition  = " << itkOptimizer->GetStopConditionDescription() << std::endl;
    }
  catch( itk::ExceptionObject & e )
    {
    std::cout << "Exception thrown ! " << std::endl;
    std::cout << "An error occurred during Optimization" << std::endl;
    std::cout << "Location    = " << e.GetLocation()    << std::endl;
    std::cout << "Description = " << e.GetDescription() << std::endl;
    return EXIT_FAILURE;
    }

  finalPosition = itkOptimizer->GetMetric()->GetParameters();
  std::cout << "Solution        = (";
  std::cout << finalPosition[0] << ",";
  std::cout << finalPosition[1] << ")" << std::endl;

  //
  // check results to see if it is within range
  //
  pass = true;
  for( unsigned int j = 0; j < 2; j++ )
    {
    if( itk::Math::FloatAlmostEqual( finalPosition[j], trueParameters[j] ) )
      {
      pass = false;
      }
    }

  if( !pass )
    {
    std::cout << "Test failed." << std::endl;
    return EXIT_FAILURE;
    }

  }

  //
  // Verify that the optimizer doesn't run if the
  // number of iterations is set to zero.
  //
  std::cout << "\nCheck the optimizer when number of iterations is set to zero:" << std::endl;
  {
  itkOptimizer->SetNumberOfIterations( 0 );
  metric->SetParameters( initialPosition );

  try
    {
    std::cout << "currentPosition before optimization: " << itkOptimizer->GetMetric()->GetParameters() << std::endl;
    itkOptimizer->StartOptimization();
    std::cout << "currentPosition after optimization: " << itkOptimizer->GetMetric()->GetParameters() << std::endl;
    std::cout << " Stop Condition  = " << itkOptimizer->GetStopConditionDescription() << std::endl;
    }
  catch( itk::ExceptionObject & excp )
    {
    std::cout << excp << std::endl;
    return EXIT_FAILURE;
    }

  if( itkOptimizer->GetCurrentIteration() > 0 )
    {
    std::cerr << "The optimizer is running iterations despite of ";
    std::cerr << "having a maximum number of iterations set to zero" << std::endl;
    return EXIT_FAILURE;
    }
  }

  //
  // Test the Exception if the GradientMagnitudeTolerance is set to a negative value
  //
  std::cout << "\nTest the Exception if the GradientMagnitudeTolerance is set to a negative value:" << std::endl;
  itkOptimizer->SetGradientMagnitudeTolerance( -1.0 );
  bool expectedExceptionReceived = false;
  try
    {
    std::cout << "currentPosition before optimization: " << itkOptimizer->GetMetric()->GetParameters() << std::endl;
    itkOptimizer->StartOptimization();
    std::cout << "currentPosition after optimization: " << itkOptimizer->GetMetric()->GetParameters() << std::endl;
    std::cout << " Stop Condition  = " << itkOptimizer->GetStopConditionDescription() << std::endl;
    }
  catch( itk::ExceptionObject & excp )
    {
    expectedExceptionReceived = true;
    std::cout << "Expected Exception " << std::endl;
    std::cout << excp << std::endl;
    }

  if( !expectedExceptionReceived )
    {
    std::cerr << "Failure to produce an exception when";
    std::cerr << "the GradientMagnitudeTolerance is negative " << std::endl;
    std::cerr << "TEST FAILED !" << std::endl;
    return EXIT_FAILURE;
    }

  //
  // Test the Exception if the RelaxationFactor is set to a value more than one.
  //
  std::cout << "\nTest the Exception if the RelaxationFactor is set to a value more than one:" << std::endl;
  itkOptimizer->SetNumberOfIterations( 100 );
  itkOptimizer->SetGradientMagnitudeTolerance( 0.01 );
  itkOptimizer->SetRelaxationFactor( 1.1 );
  expectedExceptionReceived = false;
  try
  {
  std::cout << "currentPosition before optimization: " << itkOptimizer->GetMetric()->GetParameters() << std::endl;
  itkOptimizer->StartOptimization();
  std::cout << "currentPosition after optimization: " << itkOptimizer->GetMetric()->GetParameters() << std::endl;
  std::cout << " Stop Condition  = " << itkOptimizer->GetStopConditionDescription() << std::endl;
  }
  catch( itk::ExceptionObject & excp )
  {
  expectedExceptionReceived = true;
  std::cout << "Expected Exception " << std::endl;
  std::cout << excp << std::endl;
  }

  if( !expectedExceptionReceived )
    {
    std::cerr << "Failure to produce an exception when";
    std::cerr << "the RelaxationFactor is more than one " << std::endl;
    std::cerr << "TEST FAILED !" << std::endl;
    return EXIT_FAILURE;
    }

  std::cout << "TEST PASSED !" << std::endl;
  return EXIT_SUCCESS;

}