<html xmlns:v="urn:schemas-microsoft-com:vml" xmlns:o="urn:schemas-microsoft-com:office:office" xmlns:w="urn:schemas-microsoft-com:office:word" xmlns:m="http://schemas.microsoft.com/office/2004/12/omml" xmlns="http://www.w3.org/TR/REC-html40">
<head>
<meta http-equiv="Content-Type" content="text/html; charset=utf-8">
<meta name="Generator" content="Microsoft Word 15 (filtered medium)">
<style><!--
/* Font Definitions */
@font-face
        {font-family:"Cambria Math";
        panose-1:2 4 5 3 5 4 6 3 2 4;}
@font-face
        {font-family:Calibri;
        panose-1:2 15 5 2 2 2 4 3 2 4;}
/* Style Definitions */
p.MsoNormal, li.MsoNormal, div.MsoNormal
        {margin:0in;
        margin-bottom:.0001pt;
        font-size:11.0pt;
        font-family:"Calibri",sans-serif;}
a:link, span.MsoHyperlink
        {mso-style-priority:99;
        color:#0563C1;
        text-decoration:underline;}
a:visited, span.MsoHyperlinkFollowed
        {mso-style-priority:99;
        color:#954F72;
        text-decoration:underline;}
p.msonormal0, li.msonormal0, div.msonormal0
        {mso-style-name:msonormal;
        mso-style-priority:99;
        mso-margin-top-alt:auto;
        margin-right:0in;
        mso-margin-bottom-alt:auto;
        margin-left:0in;
        font-size:11.0pt;
        font-family:"Calibri",sans-serif;}
span.EmailStyle19
        {mso-style-type:personal;
        font-family:"Calibri",sans-serif;
        color:windowtext;}
span.EmailStyle20
        {mso-style-type:personal;
        font-family:"Calibri",sans-serif;
        color:windowtext;}
span.EmailStyle22
        {mso-style-type:personal-reply;
        font-family:"Calibri",sans-serif;
        color:windowtext;}
.MsoChpDefault
        {mso-style-type:export-only;
        font-size:10.0pt;}
@page WordSection1
        {size:8.5in 11.0in;
        margin:1.0in 1.0in 1.0in 1.0in;}
div.WordSection1
        {page:WordSection1;}
--></style><!--[if gte mso 9]><xml>
<o:shapedefaults v:ext="edit" spidmax="1026" />
</xml><![endif]--><!--[if gte mso 9]><xml>
<o:shapelayout v:ext="edit">
<o:idmap v:ext="edit" data="1" />
</o:shapelayout></xml><![endif]-->
</head>
<body bgcolor="white" lang="EN-US" link="#0563C1" vlink="#954F72">
<div class="WordSection1">
<p class="MsoNormal">Hi Brad,<o:p></o:p></p>
<p class="MsoNormal"><o:p> </o:p></p>
<p class="MsoNormal">Thank you for your suggestions. The performance profiling I did was rather rudimentary, just inserting chrono timings in ITK code. I created a small example program to illustrate the situation, but I’ve seen this in many cases of multiresolution
 registration, especially if the coarse levels of the registration require a lot of smoothing.<o:p></o:p></p>
<p class="MsoNormal"><o:p> </o:p></p>
<p class="MsoNormal">In my example program (below), I’m performing 3D-to-3D registration using itk::CenteredEuler3DTransform, itk::CorrelationImageToImageMetricv4, and itk::RegularStepGradientDescentOptimizerv4. My fixed image is a 3D mouse brain image from
 the Allen Institute, with size [456, 320, 528] and 25 um isotropic spacing; my moving image is just a synthetically rotated and translated version of the same image. I run registration with a single resolution level, smoothing the image with a standard deviation
 of 100 um and shrinking by a factor of 8.<o:p></o:p></p>
<p class="MsoNormal"><o:p> </o:p></p>
<p class="MsoNormal">Sample timings from this program show that the smoothing that occurs in itk::ImageRegistrationMethodv4::InitializeRegistrationAtEachLevel() dominates the level initialization (15 s), the moving image’s gradient image computation dominates
 the metric initialization (4 s), and the registration itself runs pretty speedily (1.6 s). (The chrono timing part is not in my example code, so I’ve indicated which method the timings have been added to.)<o:p></o:p></p>
<p class="MsoNormal"><o:p> </o:p></p>
<p class="MsoNormal">Fixed image:<o:p></o:p></p>
<p class="MsoNormal">    size [456, 320, 528]<o:p></o:p></p>
<p class="MsoNormal">    spacing [25, 25, 25]<o:p></o:p></p>
<p class="MsoNormal">Initial transform:<o:p></o:p></p>
<p class="MsoNormal">    parameters [0, 0, 0, 5695.034988478505, 3726.9738602938187, 7309.779123793645, 279.5341527061055, -247.51563557593, -1.3564556700512185]<o:p></o:p></p>
<p class="MsoNormal">    center     [5695.03, 3726.97, 7309.78]<o:p></o:p></p>
<p class="MsoNormal">ImageRegistrationMethodv4::GenerateData() level: 0<o:p></o:p></p>
<p class="MsoNormal">InitializeRegistrationAtEachLevel() smoothing: 15.4068<o:p></o:p></p>
<p class="MsoNormal">this->InitializeRegistrationAtEachLevel(): 15.4224<o:p></o:p></p>
<p class="MsoNormal">ImageToImageMetricv4::ComputeMovingImageGradientFilterImage(): 3.9064<o:p></o:p></p>
<p class="MsoNormal">this->m_Metric->Initialize(): 3.9064<o:p></o:p></p>
<p class="MsoNormal">this->m_Optimizer->StartOptimization(): 1.62503<o:p></o:p></p>
<p class="MsoNormal">Final:  -0.987401 @ [0.0007390244255403159, 0.03171806411143246, 0.25914684071268346, 5693.978269101054, 3725.58761384107, 7309.850962050369, 278.47743332865645, -248.90188202867876, -1.2846174133290214]<o:p></o:p></p>
<p class="MsoNormal"><o:p> </o:p></p>
<p class="MsoNormal">I fully understand that smoothing takes time and that it’s a critically important part of multiresolution registration. My example program will not converge with lower smoothing values. It’s just that there are situations where the smoothed
 images could be faster to compute another way or may be readily available from elsewhere, so I was hoping there would be a way to exert more control over the multiresolution component of the framework. In fact, switching to itk::SmoothingRecursiveGaussianImageFilter
 reduces the initial smoothing time significantly (1.4 s), but it has different features that may have minor knock-on effects to ImageRegistrtaionMethodv4. I may also take your suggestion to separately configure and run each level, but as I understand it there
 will always be some smoothing applied to the fixed and moving images in itk::ImageRegistrationMethodv4, no matter how those images have been preprocessed upstream of registration.<o:p></o:p></p>
<p class="MsoNormal"><o:p> </o:p></p>
<p class="MsoNormal">Thank you,<o:p></o:p></p>
<p class="MsoNormal">Brian Eastwood<o:p></o:p></p>
<div style="mso-element:para-border-div;border:none;border-bottom:solid windowtext 1.0pt;padding:0in 0in 1.0pt 0in">
<p class="MsoNormal" style="border:none;padding:0in"><o:p> </o:p></p>
</div>
<p class="MsoNormal"><o:p> </o:p></p>
<p class="MsoNormal">/**<o:p></o:p></p>
<p class="MsoNormal">* \file TestRegistrationv4.cpp<o:p></o:p></p>
<p class="MsoNormal">* \brief ITKv4 registration test.<o:p></o:p></p>
<p class="MsoNormal">* \author Brian Eastwood<o:p></o:p></p>
<p class="MsoNormal">* \date 9/26/2017<o:p></o:p></p>
<p class="MsoNormal">*<o:p></o:p></p>
<p class="MsoNormal">* Example using ITK registration framework. This is a place to try out several<o:p></o:p></p>
<p class="MsoNormal">* different registration procedures on different data sets.<o:p></o:p></p>
<p class="MsoNormal">*/<o:p></o:p></p>
<p class="MsoNormal">#include <iostream><o:p></o:p></p>
<p class="MsoNormal"><o:p> </o:p></p>
<p class="MsoNormal">#include "itkCenteredEuler3DTransform.h"<o:p></o:p></p>
<p class="MsoNormal">#include "itkCenteredTransformInitializer.h"<o:p></o:p></p>
<p class="MsoNormal">#include "itkCommand.h"<o:p></o:p></p>
<p class="MsoNormal">#include "itkCorrelationImageToImageMetricv4.h"<o:p></o:p></p>
<p class="MsoNormal">#include "itkImage.h"<o:p></o:p></p>
<p class="MsoNormal">#include "itkImageFileReader.h"<o:p></o:p></p>
<p class="MsoNormal">#include "itkImageRegistrationMethodv4.h"<o:p></o:p></p>
<p class="MsoNormal">#include "itkLinearInterpolateImageFunction.h"<o:p></o:p></p>
<p class="MsoNormal">#include "itkObjectToObjectOptimizerBase.h"<o:p></o:p></p>
<p class="MsoNormal">#include "itkRegularStepGradientDescentOptimizerv4.h"<o:p></o:p></p>
<p class="MsoNormal">#include "itkVersorRigid3DTransform.h"<o:p></o:p></p>
<p class="MsoNormal"><o:p> </o:p></p>
<p class="MsoNormal"><o:p> </o:p></p>
<p class="MsoNormal">/**<o:p></o:p></p>
<p class="MsoNormal">* 3D-to-3D volume registration.<o:p></o:p></p>
<p class="MsoNormal">*/<o:p></o:p></p>
<p class="MsoNormal">void testVolumeRegistration(int argc, char* argv[])<o:p></o:p></p>
<p class="MsoNormal">{<o:p></o:p></p>
<p class="MsoNormal">                // define image types<o:p></o:p></p>
<p class="MsoNormal">                const unsigned int Dimension = 3;<o:p></o:p></p>
<p class="MsoNormal">                typedef float PixelType;<o:p></o:p></p>
<p class="MsoNormal">                typedef itk::Image< PixelType, Dimension > ImageType;<o:p></o:p></p>
<p class="MsoNormal">                typedef itk::ImageFileReader< ImageType > ReaderType;<o:p></o:p></p>
<p class="MsoNormal"><o:p> </o:p></p>
<p class="MsoNormal">                // define registration types<o:p></o:p></p>
<p class="MsoNormal">                typedef itk::CenteredEuler3DTransform<double> TransformType;<o:p></o:p></p>
<p class="MsoNormal">                typedef itk::CenteredTransformInitializer< TransformType, ImageType, ImageType > InitializerType;<o:p></o:p></p>
<p class="MsoNormal">                typedef itk::CorrelationImageToImageMetricv4< ImageType, ImageType > MetricType;<o:p></o:p></p>
<p class="MsoNormal">                typedef itk::RegularStepGradientDescentOptimizerv4< double > OptimizerType;<o:p></o:p></p>
<p class="MsoNormal">                typedef itk::ImageRegistrationMethodv4< ImageType, ImageType > RegistrationType;<o:p></o:p></p>
<p class="MsoNormal"><o:p> </o:p></p>
<p class="MsoNormal">                if (argc < 3)<o:p></o:p></p>
<p class="MsoNormal">                {<o:p></o:p></p>
<p class="MsoNormal">                                std::cerr << "Usage " << argv[0] << " fixedImage movingImage" << std::endl;<o:p></o:p></p>
<p class="MsoNormal">                                exit(1);<o:p></o:p></p>
<p class="MsoNormal">                }<o:p></o:p></p>
<p class="MsoNormal"><o:p> </o:p></p>
<p class="MsoNormal">                // create image pre-processing pipeline<o:p></o:p></p>
<p class="MsoNormal">                std::string fixedIn(argv[1]);<o:p></o:p></p>
<p class="MsoNormal">                std::string movingIn(argv[2]);<o:p></o:p></p>
<p class="MsoNormal">                ReaderType::Pointer fixedReader = ReaderType::New();<o:p></o:p></p>
<p class="MsoNormal">                fixedReader->SetFileName(fixedIn);<o:p></o:p></p>
<p class="MsoNormal">                fixedReader->Update();<o:p></o:p></p>
<p class="MsoNormal">                ReaderType::Pointer movingReader = ReaderType::New();<o:p></o:p></p>
<p class="MsoNormal">                movingReader->SetFileName(movingIn);<o:p></o:p></p>
<p class="MsoNormal">                movingReader->Update();<o:p></o:p></p>
<p class="MsoNormal">                std::cout << "Fixed image:" << std::endl;<o:p></o:p></p>
<p class="MsoNormal">                std::cout << "    size " << fixedReader->GetOutput()->GetLargestPossibleRegion().GetSize() << std::endl;<o:p></o:p></p>
<p class="MsoNormal">                std::cout << "    spacing " << fixedReader->GetOutput()->GetSpacing() << std::endl;<o:p></o:p></p>
<p class="MsoNormal"><o:p> </o:p></p>
<p class="MsoNormal">                // initialize transform parameters<o:p></o:p></p>
<p class="MsoNormal">                TransformType::Pointer transform = TransformType::New();<o:p></o:p></p>
<p class="MsoNormal">                InitializerType::Pointer initializer = InitializerType::New();<o:p></o:p></p>
<p class="MsoNormal">                initializer->SetFixedImage(fixedReader->GetOutput());<o:p></o:p></p>
<p class="MsoNormal">                initializer->SetMovingImage(movingReader->GetOutput());<o:p></o:p></p>
<p class="MsoNormal">                initializer->SetTransform(transform);<o:p></o:p></p>
<p class="MsoNormal">                initializer->MomentsOn();<o:p></o:p></p>
<p class="MsoNormal">                initializer->InitializeTransform();<o:p></o:p></p>
<p class="MsoNormal">                std::cout << "Initial transform:" << std::endl;<o:p></o:p></p>
<p class="MsoNormal">                std::cout << "    parameters " << transform->GetParameters() << std::endl;<o:p></o:p></p>
<p class="MsoNormal">                std::cout << "    center     " << transform->GetCenter() << std::endl;<o:p></o:p></p>
<p class="MsoNormal"><o:p> </o:p></p>
<p class="MsoNormal">                // create registration objects<o:p></o:p></p>
<p class="MsoNormal">                OptimizerType::Pointer optimizer = OptimizerType::New();<o:p></o:p></p>
<p class="MsoNormal">                MetricType::Pointer metric = MetricType::New();<o:p></o:p></p>
<p class="MsoNormal">                RegistrationType::Pointer registration = RegistrationType::New();<o:p></o:p></p>
<p class="MsoNormal">                registration->SetFixedImage(fixedReader->GetOutput());<o:p></o:p></p>
<p class="MsoNormal">                registration->SetMovingImage(movingReader->GetOutput());<o:p></o:p></p>
<p class="MsoNormal">                registration->SetOptimizer(optimizer);<o:p></o:p></p>
<p class="MsoNormal">                registration->SetMetric(metric);<o:p></o:p></p>
<p class="MsoNormal">                registration->SetInitialTransform(transform);<o:p></o:p></p>
<p class="MsoNormal">                registration->InPlaceOn();<o:p></o:p></p>
<p class="MsoNormal"><o:p> </o:p></p>
<p class="MsoNormal">                // configure resolution levels<o:p></o:p></p>
<p class="MsoNormal">                itk::SizeValueType levels = 1;<o:p></o:p></p>
<p class="MsoNormal">                ImageType::SpacingType spacing = fixedReader->GetOutput()->GetSpacing();<o:p></o:p></p>
<p class="MsoNormal">                RegistrationType::SmoothingSigmasArrayType sigmas(levels);<o:p></o:p></p>
<p class="MsoNormal">                RegistrationType::ShrinkFactorsArrayType shrink(levels);<o:p></o:p></p>
<p class="MsoNormal">                sigmas[0] = spacing[0] * 4;<o:p></o:p></p>
<p class="MsoNormal">                shrink[0] = 8;<o:p></o:p></p>
<p class="MsoNormal">                registration->SetNumberOfLevels(levels);<o:p></o:p></p>
<p class="MsoNormal">                registration->SetSmoothingSigmasPerLevel(sigmas);<o:p></o:p></p>
<p class="MsoNormal">                registration->SetSmoothingSigmasAreSpecifiedInPhysicalUnits(true);<o:p></o:p></p>
<p class="MsoNormal">                registration->SetShrinkFactorsPerLevel(shrink);<o:p></o:p></p>
<p class="MsoNormal"><o:p> </o:p></p>
<p class="MsoNormal">                // configure optimizer<o:p></o:p></p>
<p class="MsoNormal">                OptimizerType::ScalesType scales(transform->GetNumberOfParameters());<o:p></o:p></p>
<p class="MsoNormal">                scales.fill(1e-4);<o:p></o:p></p>
<p class="MsoNormal">                scales[0] = scales[1] = scales[2] = 1;<o:p></o:p></p>
<p class="MsoNormal">                optimizer->SetScales(scales);<o:p></o:p></p>
<p class="MsoNormal">                optimizer->SetLearningRate(1);<o:p></o:p></p>
<p class="MsoNormal">                optimizer->SetMinimumStepLength(1e-3);<o:p></o:p></p>
<p class="MsoNormal">                optimizer->SetMinimumConvergenceValue(1e-4);<o:p></o:p></p>
<p class="MsoNormal">                optimizer->SetNumberOfIterations(60);<o:p></o:p></p>
<p class="MsoNormal"><o:p> </o:p></p>
<p class="MsoNormal">                // run registration<o:p></o:p></p>
<p class="MsoNormal">                try<o:p></o:p></p>
<p class="MsoNormal">                {<o:p></o:p></p>
<p class="MsoNormal">                                registration->Update();<o:p></o:p></p>
<p class="MsoNormal">                }<o:p></o:p></p>
<p class="MsoNormal">                catch (itk::ExceptionObject& err)<o:p></o:p></p>
<p class="MsoNormal">                {<o:p></o:p></p>
<p class="MsoNormal">                                std::cerr << "Exception occurred during registration:\n" << err << std::endl;<o:p></o:p></p>
<p class="MsoNormal">                }<o:p></o:p></p>
<p class="MsoNormal"><o:p> </o:p></p>
<p class="MsoNormal">                // print result of registration<o:p></o:p></p>
<p class="MsoNormal">                std::cout << "Final:  " << optimizer->GetCurrentMetricValue() << " @ " << transform->GetParameters() << std::endl;            
<o:p></o:p></p>
<p class="MsoNormal">}<o:p></o:p></p>
<p class="MsoNormal"><o:p> </o:p></p>
<p class="MsoNormal">int main(int argc, char* argv[])<o:p></o:p></p>
<p class="MsoNormal">{<o:p></o:p></p>
<p class="MsoNormal">                testVolumeRegistration(argc, argv);<o:p></o:p></p>
<p class="MsoNormal">                return 0;<o:p></o:p></p>
<p class="MsoNormal">}<o:p></o:p></p>
<div style="mso-element:para-border-div;border:none;border-bottom:solid windowtext 1.0pt;padding:0in 0in 1.0pt 0in">
<p class="MsoNormal" style="border:none;padding:0in"><o:p> </o:p></p>
</div>
<p class="MsoNormal"><o:p> </o:p></p>
<p class="MsoNormal"><o:p> </o:p></p>
<div>
<div style="border:none;border-top:solid #E1E1E1 1.0pt;padding:3.0pt 0in 0in 0in">
<p class="MsoNormal"><b>From:</b> Lowekamp, Bradley (NIH/NLM/LHC) [C] [mailto:blowekamp@mail.nih.gov]
<br>
<b>Sent:</b> Monday, September 25, 2017 10:02 AM<br>
<b>To:</b> Brian Eastwood <brian@mbfbioscience.com>; community@itk.org<br>
<b>Cc:</b> Nick Tustison <ntustison@gmail.com><br>
<b>Subject:</b> Re: [ITK] Image Registration Initialization<o:p></o:p></p>
</div>
</div>
<p class="MsoNormal"><o:p> </o:p></p>
<p class="MsoNormal">Hello,<o:p></o:p></p>
<p class="MsoNormal"><o:p> </o:p></p>
<p class="MsoNormal">Did you do any specific profiling to determine these bottlenecks? If so could you please share the profiling code along with the timing you got?<o:p></o:p></p>
<p class="MsoNormal"><o:p> </o:p></p>
<p class="MsoNormal">Also you did not provide any details about the situation where you found these performance issues. What transform you are using? What is the size of the image? What optimizer are you using?<o:p></o:p></p>
<p class="MsoNormal"><o:p> </o:p></p>
<p class="MsoNormal">With the v4 framework, in some situations the initialization, and iterations are slower, however for many registration setups with a large number of parameters good results can be achieved with fewer iterations. Additionally, when the framework
 is combined with the parameters estimators the results are more robust. And with the built in multi-resolution and metric sampling efficient results can be achieved by using all these features of the framework.<o:p></o:p></p>
<p class="MsoNormal"><o:p> </o:p></p>
<p class="MsoNormal">Based on your comments, I looked through the itk::ImageRegistrationMethodv4. I see that it is using the itk::DiscreteGaussianImageFilter. The itk::SmoothingRecursiveGuassianImageFilter is generally known to be more efficient. I would suggest
 hacking the ImageRegistrationMethodv4 to use the itk::SmoothingRecursiveGaussianImageFilter, and determine if there is a performance improvement.<o:p></o:p></p>
<p class="MsoNormal"><o:p> </o:p></p>
<p class="MsoNormal">You cannot swap the fixed and moving inputs images at each resolution scale iteration. However, you could take control of scale iterations yourself, by just running the ImageRegistrationMethod at one level, where you choose the appropriate
 scale for setting the fixed and moving input images.<o:p></o:p></p>
<p class="MsoNormal"><o:p> </o:p></p>
<p class="MsoNormal">HTH,<o:p></o:p></p>
<p class="MsoNormal">Brad<o:p></o:p></p>
<p class="MsoNormal"><o:p> </o:p></p>
<p class="MsoNormal"><o:p> </o:p></p>
<div style="border:none;border-top:solid #B5C4DF 1.0pt;padding:3.0pt 0in 0in 0in">
<p class="MsoNormal" style="margin-left:.5in"><b><span style="font-size:12.0pt;color:black">From:
</span></b><span style="font-size:12.0pt;color:black">Brian Eastwood <<a href="mailto:brian@mbfbioscience.com">brian@mbfbioscience.com</a>><br>
<b>Date: </b>Friday, September 22, 2017 at 4:55 PM<br>
<b>To: </b>"<a href="mailto:community@itk.org">community@itk.org</a>" <<a href="mailto:community@itk.org">community@itk.org</a>><br>
<b>Subject: </b>[ITK] Image Registration Initialization<o:p></o:p></span></p>
</div>
<div>
<p class="MsoNormal" style="margin-left:.5in"><o:p> </o:p></p>
</div>
<p style="mso-margin-top-alt:0in;margin-right:0in;margin-bottom:0in;margin-left:.5in;margin-bottom:.0001pt">
<span style="font-size:14.0pt">Hi Folks,</span><o:p></o:p></p>
<p style="mso-margin-top-alt:0in;margin-right:0in;margin-bottom:0in;margin-left:.5in;margin-bottom:.0001pt">
<span style="font-size:14.0pt"> </span><o:p></o:p></p>
<p style="mso-margin-top-alt:0in;margin-right:0in;margin-bottom:0in;margin-left:.5in;margin-bottom:.0001pt">
<span style="font-size:14.0pt">I really like the changes implemented for ITK v4 registration. The composite transform makes multistage registration much easier to implement. Nice work!</span><o:p></o:p></p>
<p style="mso-margin-top-alt:0in;margin-right:0in;margin-bottom:0in;margin-left:.5in;margin-bottom:.0001pt">
<span style="font-size:14.0pt"> </span><o:p></o:p></p>
<p style="mso-margin-top-alt:0in;margin-right:0in;margin-bottom:0in;margin-left:.5in;margin-bottom:.0001pt">
<span style="font-size:14.0pt">Is there a way to explicitly set the fixed and moving images to be used at each level of a multiresolution registration? How about the gradient images used by the image metric? I see you can set gradient filters on the metric,
 but what if you already have the gradients precomputed?</span><o:p></o:p></p>
<p style="mso-margin-top-alt:0in;margin-right:0in;margin-bottom:0in;margin-left:.5in;margin-bottom:.0001pt">
<span style="font-size:14.0pt"> </span><o:p></o:p></p>
<p style="mso-margin-top-alt:0in;margin-right:0in;margin-bottom:0in;margin-left:.5in;margin-bottom:.0001pt">
<span style="font-size:14.0pt">I get the feeling that v4 registration runs slower than similar v3 registration, but I think it’s mainly down to the initialization where the fixed and moving images are smoothed and the metric computes the image gradients.</span><o:p></o:p></p>
<p style="mso-margin-top-alt:0in;margin-right:0in;margin-bottom:0in;margin-left:.5in;margin-bottom:.0001pt">
<span style="font-size:14.0pt"> </span><o:p></o:p></p>
<p style="mso-margin-top-alt:0in;margin-right:0in;margin-bottom:0in;margin-left:.5in;margin-bottom:.0001pt">
<span style="font-size:14.0pt">It seems like there could be some real advantages to using precomputed smoothed and gradient images. For example, the smoothing could be implemented with a cascade of Gaussians where the standard deviation used at a level k is
 given by s^2 = s_k^2 - s_{k-1}^2. Or the smoothing could be implemented on the GPU. Or we could be reusing a resolution level in a different registration stage, perhaps optimizing a different set of transform parameters.</span><o:p></o:p></p>
<p style="mso-margin-top-alt:0in;margin-right:0in;margin-bottom:0in;margin-left:.5in;margin-bottom:.0001pt">
<span style="font-size:14.0pt"> </span><o:p></o:p></p>
<p style="mso-margin-top-alt:0in;margin-right:0in;margin-bottom:0in;margin-left:.5in;margin-bottom:.0001pt">
<span style="font-size:14.0pt">Have I missed something in the API where it is possible to explicitly set the images and gradients to use at each level of the registration? I could probably create an itk::Command that swaps the fixed and moving images before
 each registration iteration, but is there a way to disable the smoothing step?</span><o:p></o:p></p>
<p style="mso-margin-top-alt:0in;margin-right:0in;margin-bottom:0in;margin-left:.5in;margin-bottom:.0001pt">
<span style="font-size:14.0pt"> </span><o:p></o:p></p>
<p style="mso-margin-top-alt:0in;margin-right:0in;margin-bottom:0in;margin-left:.5in;margin-bottom:.0001pt">
<span style="font-size:14.0pt">Regards,</span><o:p></o:p></p>
<p style="mso-margin-top-alt:0in;margin-right:0in;margin-bottom:0in;margin-left:.5in;margin-bottom:.0001pt">
<span style="font-size:14.0pt">Brian</span><o:p></o:p></p>
</div>
</body>
</html>