[Insight-users] no. of iterations in 3D registration...

Luis Ibanez luis.ibanez at kitware.com
Mon Aug 14 08:47:44 EDT 2006


Hi Kajetan,


1) The RegularStepGradientescentOptimizer strictly speaking does not
    find the minimum of the cost function. It finds a place where the
    gradient is null.  That's why it can be trapped in local minima,
    as many other gradient-based optimizers do.

    It was suggested a while ago to add to the optimizers an extra
    member variable that will keep the best value that the optimizer
    found during its path. This option would make sure that you keep
    the value at iteration 5 instead of the value that is found by
    the optimizer.


2) This may or may not qualify as a bug. Here is what happens.

    When reporting Values and CurrentPosition during normal iterations
    the Value reported by the optimizer is computed by the metric
    and grabbed in line 112 of the file

      Insight/Code/Numerics/
         itkRegularStepGradientDescentBaseOptimizer.cxx

     This value is computed by passing to the Metrics the
     "CurrentPosition" array. After the m_Value is computed,
     the AdvanceOneStep() method is invoked, where the *new*
     current position is computed, and *after* that, the
     IterationEvent is invoked.

     The result is that by the time the iteration event is invoked,
     the m_CurrentPosition and the m_Value do not match anymore.
     The m_Value is that corresponding to the "previous" current
     position.

     One easy solution is to revert the order of the calls:


        this->StepAlongGradient()              in line 263
        this->InvokeEvent( IterationEvent() )  in line  267

     of the file

        Insight/Code/Numerics/
            itkRegularStepGradientDescentBaseOptimizer.cxx


     In this way, the iteration will be reported at a moment
     when the m_Value and the CurrentPosition still match.


     All this, means that the reports that you get so far
     as shifted between the Value and the Position, while you
     are looking at the iterations during the progress of the
     optimizer. That is, the report of Iteration 11 is showing
     you the next "current position" along with the current value.

     Since the next "current position" happens to be the one
     where the optimizer converges, then the X,Y translations
     in Iteration 11 report and the final solution report are
     the same. However the m_value of the final solution is the
     Metric that *actually* corresponds to that position.


     We just added a bug entry for this issue.

     http://public.kitware.com/Bug/bug.php?op=show&bugid=3606&pos=0

     We should take care of it before the next release.




   Please let us know if you have further questions,



      Thanks



         Luis




--------------------------
Kajetan Berlinger wrote:
> Hi Luis,
> 
> I have a question with respect to registration/optimization as well.
> In the example below, I use the regualar step gradient descent method and 
> correlation coefficient for a simple rigid 2D-2D registration with 
> translation only. Visual inspection of the result is quite satisfying but I 
> don't understand the way of the optimizer through its parameter space.
> 
> (1) Why is the final result not the translation of iteration 5? This iteration 
> obviously yielded the minimum value.
> 
> (2) Why is the metric value of iteration 11 (last output of the observer) 
> different to the final metric, although they are both based on the same 
> translation?
> 
> ************************************************
> Current Iterations: 0
> Translation X: -0.826109
> Translation Y: 3.91376
> Similarity Value: -0.92855
> Pixels considered: 46767
> ************************************************
> 
> 
> ************************************************
> Current Iterations: 1
> Translation X: -1.454
> Translation Y: 7.86417
> Similarity Value: -0.948779
> Pixels considered: 46767
> ************************************************
> 
> 
> ************************************************
> Current Iterations: 2
> Translation X: -1.75847
> Translation Y: 11.8526
> Similarity Value: -0.956318
> Pixels considered: 46767
> ************************************************
> 
> 
> ************************************************
> Current Iterations: 3
> Translation X: 0.0892366
> Translation Y: 11.0871
> Similarity Value: -0.953901
> Pixels considered: 46767
> ************************************************
> 
> 
> ************************************************
> Current Iterations: 4
> Translation X: 0.302107
> Translation Y: 9.09845
> Similarity Value: -0.955402
> Pixels considered: 46767
> ************************************************
> 
> 
> ************************************************
> Current Iterations: 5
> Translation X: 0.856397
> Translation Y: 9.93077
> Similarity Value: -0.956675      <---------------------------------------(1)
> Pixels considered: 46767
> ************************************************
> 
> 
> ************************************************
> Current Iterations: 6
> Translation X: 0.410382
> Translation Y: 9.70478
> Similarity Value: -0.956346
> Pixels considered: 46767
> ************************************************
> 
> 
> ************************************************
> Current Iterations: 7
> Translation X: 0.474056
> Translation Y: 9.20885
> Similarity Value: -0.956433
> Pixels considered: 46767
> ************************************************
> 
> 
> ************************************************
> Current Iterations: 8
> Translation X: 0.201551
> Translation Y: 8.78964
> Similarity Value: -0.956553
> Pixels considered: 46767
> ************************************************
> 
> 
> ************************************************
> Current Iterations: 9
> Translation X: 0.44661
> Translation Y: 8.74017
> Similarity Value: -0.956645
> Pixels considered: 46767
> ************************************************
> 
> 
> ************************************************
> Current Iterations: 10
> Translation X: 0.391267
> Translation Y: 8.85226
> Similarity Value: -0.956657
> Pixels considered: 46767
> ************************************************
> 
> 
> ************************************************
> Current Iterations: 11
> Translation X: 0.397056
> Translation Y: 8.79002
> Similarity Value: -0.956625      <----------------------------- (2)
> Pixels considered: 46767
> ************************************************
> 
> ************************************************
> Image B:
> Current DRR Pair: 5
> ************************************************
> Final Result:
> Translation X: 0.397056
> Translation Y: 8.79002
> Iterations: 13
> Similarity Value: -0.956572   <-------------------------------- (2)
> Pixels considered: 46767
> Optimizer: Regular Step Gradient Descent
> Metric: Correlation Coefficient
> ************************************************
> ************************************************
> 
> 
> Thanks for helping me!
> 
> Kaj
> 
> 
> On Thursday 10 August 2006 09:39, Luis Ibanez wrote:
> 
>>Hi Shahab,
>>
>>As long as you set up the optimizer to minimize, the metric should
>>never increase for more than one consecutive iteration. That is,
>>the metric can increase, but then it should decrease in the next
>>iteration.
>>
>>There is a very very small chance that a metric can increase in
>>two consecutive iterations, but it requires a convoluted configuration
>>of the cost function.
>>
>>Can you post a plot of the metric values versus iteration number,
>>just as it is done in the ITK Software Guide examples ?
>>
>>     http://www.itk.org/ItkSoftwareGuide.pdf
>>
>>How many iterations to use depends on what kind of optimizer you
>>are using and how large you define the step length (or its equivalent).
>>
>>You may want to look at the tutorial on image registration:
>>
>>at
>>
>>       http://www.itk.org/HTML/Tutorials.htm
>>
>>more specifically at
>>
>>
>>http://www.itk.org/CourseWare/Training/RegistrationMethodsOverview.pdf
>>
>>
>>The effect and step length versus number of iterations is also
>>discussed in the Image Registration chapter of the Software Guide.
>>
>>
>>
>>Regards
>>
>>
>>
>>    Luis
>>
>>
>>
>>--------------------------
>>
>>Shahabuddin Ansari wrote:
>>
>>>Hello,
>>>
>>>In 3D rigid registration, I am using versor transform, mean square
>>>metric, and linear interpolator. The observation of the final parameters
>>>shows that the metric value first decreases to a minimum and then start
>>>back to increase. In my understanding, it should fluctuate around the
>>>optimal value even the number of iterations are much higher (say 200)
>>>then the required. Please commnets (say 50).
>>>
>>>Thanks
>>>Shahab
>>
>>_______________________________________________
>>Insight-users mailing list
>>Insight-users at itk.org
>>http://www.itk.org/mailman/listinfo/insight-users
> 
> _______________________________________________
> Insight-users mailing list
> Insight-users at itk.org
> http://www.itk.org/mailman/listinfo/insight-users
> 
> 




More information about the Insight-users mailing list