[Insight-users] RE: increasing bspline control point density

Luis Ibanez luis.ibanez at kitware.com
Sat Jul 24 11:50:37 EDT 2004


Hi Stefan,

Setting the low resolution spline as a bulk transform
of the high resolution one will be an elegant way of
solving this problem with a single line of code.  :-)
As you suggested, it will be enough to initialize the
high resolution metric coefficients to zero.

However,

The computational load of the BSpline transform is quite
heavy. The transform and the interpolators are the two
bottlenecks of the image registration metric. At every
evaluation of the Image Metric, pixel coordinates are
mapped using the Transform, and image values are computed
with the interpolator.

For a typical 3D scan, with 512 x 512 x 200 pixels you
will be transforming 5.2e7 points per metric evaluation.

If you enchain two BSpline transforms, that will certainly
double the computational time of your registration.

At this point the question is:


    What is more critical for you ?

    1) The programming time required to get
       this application up and running

                   or

    2) The time it will take for running the
       application every time you register two
       images.


Note also that if you go with the addition of the Bulk Transform
you will have to use the same couple of transforms for resampling
the final registered image.



  Regards,


    Luis



-----------------------
Stefan Klein wrote:
> Luis,
> 
> Thanks for the example; actually i didn't realise it would be so simple. 
> Soon i will adapt my code as well,  instead of the BSplineUpsampleFilter 
> that i mentioned earlier, which is a little crude i have to admit. :-)
> 
> One question now 'raises in my head' though: would setting the 
> low-resolution BSplineTransform as a BulkTransform in the 
> high-resolution transform not give the exact result as well (and then 
> setting the initial Coefficients of the high-resolution BSpline to zero)?
> 
> It would be slightly slower maybe (because everytime the TransformPoint 
> method is called, the TransformPoint method of the BulkTransform will 
> also be called). But it might be an easy way of getting the exact solution.
> 
> Or am i making a big mistake now?
> 
> Stefan.
> 
> 
> 
> 
> 
> At 10:15 24/07/04, you wrote:
> 
>> Hi Karl,
>>
>> 1) In theory you are right.
>>
>>    The formal way to do this would be to compute the
>>    control points of the high resolution grid directly
>>    from the control points of the low resolution grid.
>>    That's certainly an interesting exercise for an exam
>>    in a numerical methods course, or for a friendly
>>    homework.
>>
>>
>> 2) A semi-pragmatic way to get the control points for
>>    the high resolution spline could be to sample the
>>    deformation field produced by the first stage of
>>    the registration process, then estimate the control
>>    points in the high resolution BSpline that would
>>    generate a good approximation to those sampled
>>    deformations. A classical option is to do a
>>    least-squares minimization in order to find the
>>    control point values.
>>
>>
>> 3) In practice, However...
>>    You just need to initialize the high resolution
>>    BSpline with reasonable values. Once you initialize
>>    the BSpline transform, the *next* thing you do is
>>    to run an iterative optimization method that is
>>    going to refine the BSpline.  In that context it
>>    doesn't seem to be too critical to adhere to the
>>    formality of computing the perfect control points
>>    at initialization.
>>
>>
>> The full formal estimation of the high resolution
>> control points will be critical if you were going
>> to use it directly without further refinement...
>> but then,... given that this high resolution BSpline
>> will be fully equivalent to the low resolution one,
>> .... there will be no point in the whole computation
>> because you could have used the low resolution grid
>> from the beginning.     :-)
>>
>>
>> If your goal is to publish a paper, you probably want
>> to do the nice formal theoretical approach in (1).
>> That will populate your paper with beautiful equations
>> that will certainly help to impress the reviewers and
>> the readers with your intelligence and math skills.
>>
>> If want you want is to get these two images registered,
>> then, you better invest your energy and time in fine
>> tunning the parameters of the optimizer and in evaluating
>> the quality of the registration, since after all, the whole
>> purpose of the second stage of registration is to change
>> the deformation field represented by the BSpline in order
>> to better register the two images.
>>
>>
>>  "Flexibility increases your chances for survival"
>>
>>
>>
>> Regards,
>>
>>
>>
>>    Luis
>>
>>
>>
>> ------------------
>> Karl Baum wrote:
>>
>>> Luis,
>>> Thanks for your help.
>>> Could you please clarify something for me though.  In the sample when
>>> initializing the parameters for the second registration, you are 
>>> finding the
>>> displacement(deformation) from the first iteration at the coordinates 
>>> of the
>>> new control points.  You are then using this as the start for the new
>>> control points.  However doesn't this change the deformation.  
>>> Because the
>>> control points are not actually on the spline, don't you need to find a
>>> control point that causes this same displacement, not use this as the
>>> displacement of our new control point.
>>> Sorry this was hard to explain and sounds confusing.
>>> Thanks,
>>> Karl
>>> -----Original Message-----
>>> From: Luis Ibanez [mailto:luis.ibanez at kitware.com] Sent: Friday, July 
>>> 23, 2004 11:26 AM
>>> To: Karl Baum
>>> Cc: insight-users at itk.org
>>> Subject: Re: [Insight-users] RE: increasing bspline control point density
>>>
>>> Hi Karl,
>>> For your convenience we just added an example to ITK
>>> illustrating how you can increase the resolution of
>>> the BSpline transform grid in order to refine a
>>> deformable registration.
>>>
>>> You will find this new example under:
>>>    Insight/Examples/Registration/
>>>              DeformableRegistration6.cxx
>>> http://www.itk.org/cgi-bin/viewcvs.cgi/Examples/Registration/DeformableRegis
>>> tration6.cxx?root=Insight
>>> You will have to update your CVS checkout in order
>>> to get this new file.
>>>
>>> The example is mostly based on DeformableRegistration4.cxx,
>>> and a second stage was added in which another BSpline
>>> transform with higher resolution grid is configured.
>>> Please let us know if you find any problems or
>>> have further questions.
>>>
>>>     Thanks
>>>
>>>       Luis
>>>
>>> ----------------
>>> Karl Baum wrote:
>>>
>>>> Jim,
>>>> You are correct in stating that I would like to "upsample" the B-Spline
>>>
>>> and
>>>
>>>> not the image.
>>>> Thanks,
>>>> Karl
>>>>
>>>>
>>>> Date: Fri, 23 Jul 2004 10:25:43 -0400
>>>> From: "Miller, James V (Research)" <millerjv at crd.ge.com>
>>>> Subject: RE: [Insight-users] increasing bspline control point density
>>>> To: "'Stefan Klein'" <stefan at isi.uu.nl>, insight-users at itk.org
>>>>
>>>> I think the BSplineUpsampleImageFilter upsamples an image via B-spline
>>>> interpolation. I think Karl wants to "upsample" his B-spline not an 
>>>> image.
>>>> Jim
>>>>
>>>> -----Original Message-----
>>>> From: Stefan Klein [mailto:stefan at isi.uu.nl]
>>>> Sent: Friday, July 23, 2004 9:21 AM
>>>> To: insight-users at itk.org
>>>> Subject: RE: [Insight-users] increasing bspline control point density
>>>>
>>>>
>>>> Hi Karl,
>>>>
>>>> You may use the itk::BSplineUpsampleImageFilter.
>>>> Look at the help of itk::BSplineUpsampleImageFilterBase. There it's
>>>> explained how it works.
>>>>
>>>> Regards,
>>>> Stefan.
>>>>
>>>>
>>>>
>>>>
>>>> At 08:46 23/07/04, Miller, James V (Research) wrote:
>>>>
>>>>
>>>> I do not think ITK provides a mechanism for increasing the number of 
>>>> control points in a B-spline grid.  However, the mathematics of 
>>>> B-splines
>>>> do support this operation.  There are formulas for adding a control 
>>>> point
>>>> to a B-spline that maintain the trace of the B-spline (shape of the
>>>> B-spline).
>>>>
>>>> Any good book on B-splines will provide the formulas.
>>>> Jim
>>>>
>>>>
>>>>
>>>> -----Original Message-----
>>>> From: Karl Baum [ mailto:kgbaum at syr.edu <mailto:kgbaum at syr.edu> ]
>>>> Sent: Thursday, July 22, 2004 5:22 PM
>>>> To: insight-users at itk.org
>>>> Subject: [Insight-users] increasing bspline control point density
>>>>
>>>>
>>>> Hi,
>>>> I am using BSplines to provide the transformations necessary for my
>>>> registration.  I would like to take an iterative approach where I 
>>>> first do
>>>> registration with a low control point density, and then use the 
>>>> results as
>>>> input to a registration with a higher control point density.  The 
>>>> issue I
>>>
>>> am
>>>
>>>> facing at this point in time is how to initialize the parameters of the
>>>> second iteration with the results of the first iteration.  In other 
>>>> words
>>>
>>> I
>>>
>>>> need to take the grid defined by the deformation coefficients of the 
>>>> first
>>>> iteration, and use it as input into the second.  To do this I need to
>>>> calculate the starting deformation coefficients on the higher density
>>>
>>> grid.
>>>
>>>> Does itk provide anything to do this?  I was unable to find any
>>>> documentation addressing this.
>>>>
>>>> Thanks,
>>>> Karl
>>>>
>>>>
>>>> _______________________________________________
>>>> 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
>>
>>
>>
>>
>> _______________________________________________
>> 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