SimpleITK  2.0.0
ImageRegistrationMethodDisplacement1/ImageRegistrationMethodDisplacement1.cxx
/*=========================================================================
*
* Copyright NumFOCUS
*
* 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.
*
*=========================================================================*/
// This one header will include all SimpleITK filters and external
// objects.
#include <SimpleITK.h>
#include <iostream>
#include <stdlib.h>
#include <iomanip>
namespace sitk = itk::simple;
class IterationUpdate
: public sitk::Command
{
public:
IterationUpdate( const sitk::ImageRegistrationMethod &m)
: m_Method(m)
{}
void Execute( ) override
{
// use sitk's output operator for std::vector etc..
using sitk::operator<<;
// stash the stream state
std::ios state(NULL);
state.copyfmt(std::cout);
std::cout << std::fixed << std::setfill(' ') << std::setprecision( 5 );
if ( m_Method.GetOptimizerIteration() == 0 )
{
std::cout << "\tLevel: " << std::setw(3) << m_Method.GetCurrentLevel() << std::endl;
std::cout << "\tScales: " << m_Method.GetOptimizerScales() << std::endl;
}
std::cout << '#' << m_Method.GetOptimizerIteration() << std::endl;
std::cout << "\tMetric Value: " << m_Method.GetMetricValue() << std::endl;
std::cout << "\tLearning Rate: " << m_Method.GetOptimizerLearningRate() << std::endl;
if (m_Method.GetOptimizerConvergenceValue() != std::numeric_limits<double>::max())
{
std::cout << "\tConvergence Value: " << std::scientific << m_Method.GetOptimizerConvergenceValue() << std::endl;
}
std::cout.copyfmt(state);
}
private:
};
class MultiResolutionIterationUpdate
: public sitk::Command
{
public:
MultiResolutionIterationUpdate( const sitk::ImageRegistrationMethod &m)
: m_Method(m)
{}
void Execute( ) override
{
// use sitk's output operator for std::vector etc..
using sitk::operator<<;
// stash the stream state
std::ios state(NULL);
state.copyfmt(std::cout);
std::cout << std::fixed << std::setfill(' ') << std::setprecision( 5 );
std::cout << "\tStop Condition: " << m_Method.GetOptimizerStopConditionDescription() << std::endl;
std::cout << "============= Resolution Change =============" << std::endl;
std::cout.copyfmt(state);
}
private:
};
int main(int argc, char *argv[])
{
if ( argc < 4 )
{
std::cerr << "Usage: " << argv[0] << " <fixedImageFilter> <movingImageFile> <outputTransformFile>" << std::endl;
return 1;
}
{
std::vector<unsigned int> shrinkFactors = { 3, 2, 1 };
std::vector<double> smoothingSigmas = { 2.0, 1.0, 1.0 };
R.SetShrinkFactorsPerLevel(shrinkFactors);
R.SetSmoothingSigmasPerLevel(smoothingSigmas);
}
{
double learningRate = 1.0;
unsigned int numberOfIterations = 100;
double convergenceMinimumValue = 1e-6;
unsigned int convergenceWindowSize = 10;
numberOfIterations,
convergenceMinimumValue,
convergenceWindowSize,
estimateLearningRate
);
}
R.SetInitialTransform(initialTx);
IterationUpdate cmd(R);
MultiResolutionIterationUpdate cmd2(R);
sitk::Transform outTx1 = R.Execute( fixed, moving );
std::cout << "-------" << std::endl;
std::cout << outTx1.ToString() << std::endl;
std::cout << "Optimizer stop condition: " << R.GetOptimizerStopConditionDescription() << std::endl;
std::cout << " Iteration: " << R.GetOptimizerIteration() << std::endl;
std::cout << " Metric value: " << R.GetMetricValue() << std::endl;
sitk::Image displacementField = sitk::Image(fixed.GetSize(), sitk::sitkVectorFloat64);
displacementField.CopyInformation(fixed);
sitk::DisplacementFieldTransform displacementTx(displacementField);
const double varianceForUpdateField = 0.0;
const double varianceForTotalField = 1.5;
displacementTx.SetSmoothingGaussianOnUpdate(varianceForUpdateField,
varianceForTotalField);
R.SetInitialTransform(displacementTx, true);
{
std::vector<unsigned int> shrinkFactors = { 3, 2, 1 };
std::vector<double> smoothingSigmas = { 2.0, 1.0, 1.0 };
R.SetShrinkFactorsPerLevel(shrinkFactors);
R.SetSmoothingSigmasPerLevel(smoothingSigmas);
}
{
double learningRate = 1.0;
unsigned int numberOfIterations = 300;
double convergenceMinimumValue = 1e-6;
unsigned int convergenceWindowSize = 10;
numberOfIterations,
convergenceMinimumValue,
convergenceWindowSize,
estimateLearningRate
);
}
R.Execute(fixed, moving);
std::cout << "-------" << std::endl;
std::cout << displacementTx.ToString() << std::endl;
std::cout << "Optimizer stop condition: " << R.GetOptimizerStopConditionDescription() << std::endl;
std::cout << " Iteration: " << R.GetOptimizerIteration() << std::endl;
std::cout << " Metric value: " << R.GetMetricValue() << std::endl;
sitk::CompositeTransform outTx( {outTx1, displacementTx} );
sitk::WriteTransform(outTx, argv[3]);
return 0;
}
itk::simple::Image
The Image class for SimpleITK.
Definition: sitkImage.h:75
itk::simple::ImageRegistrationMethod::SetOptimizerScalesFromPhysicalShift
Self & SetOptimizerScalesFromPhysicalShift(unsigned int centralRegionRadius=5, double smallParameterVariation=0.01)
Estimating scales of transform parameters a step sizes, from the maximum voxel shift in physical spac...
itk::simple::ImageRegistrationMethod::SetOptimizerAsGradientDescent
Self & SetOptimizerAsGradientDescent(double learningRate, unsigned int numberOfIterations, double convergenceMinimumValue=1e-6, unsigned int convergenceWindowSize=10, EstimateLearningRateType estimateLearningRate=Once, double maximumStepSizeInPhysicalUnits=0.0)
Gradient descent optimizer.
itk::simple::ProcessObject::AddCommand
virtual int AddCommand(itk::simple::EventEnum event, itk::simple::Command &cmd)
Add a Command Object to observer the event.
itk::simple::sitkLinear
@ sitkLinear
N-D linear interpolation.
Definition: sitkInterpolator.h:39
itk::simple::ImageRegistrationMethod::GetOptimizerStopConditionDescription
std::string GetOptimizerStopConditionDescription() const
itk::simple::WriteTransform
SITKCommon_EXPORT void WriteTransform(const Transform &transform, const std::string &filename)
itk::simple::DisplacementFieldTransform
A dense deformable transform over a bounded spatial domain for 2D or 3D coordinates space.
Definition: sitkDisplacementFieldTransform.h:36
itk::simple::Transform
A simplified wrapper around a variety of ITK transforms.
Definition: sitkTransform.h:80
itk::simple::ImageRegistrationMethod::EachIteration
@ EachIteration
Estimate learning rate at each iteration.
Definition: sitkImageRegistrationMethod.h:257
SimpleITK.h
itk::simple::ImageRegistrationMethod::MetricUseFixedImageGradientFilterOff
Self & MetricUseFixedImageGradientFilterOff()
Enable image gradient computation by a filter.
Definition: sitkImageRegistrationMethod.h:521
itk::simple::Image::GetSize
std::vector< unsigned int > GetSize() const
itk::simple::Command
An implementation of the Command design pattern for callback.
Definition: sitkCommand.h:44
itk::simple::ImageRegistrationMethod::SetMovingInitialTransform
Self & SetMovingInitialTransform(const Transform &transform)
Set a fixed transform component towards moving domain.
Definition: sitkImageRegistrationMethod.h:180
itk::simple::sitkVectorFloat64
@ sitkVectorFloat64
Multi-component of 64 bit float.
Definition: sitkPixelIDValues.h:113
itk::simple::ImageRegistrationMethod::GetMetricValue
double GetMetricValue() const
itk::simple::ImageRegistrationMethod::SetSmoothingSigmasPerLevel
Self & SetSmoothingSigmasPerLevel(const std::vector< double > &smoothingSigmas)
Set the sigmas of Gaussian used for smoothing.
itk::simple::Image::GetDimension
unsigned int GetDimension() const
itk::simple::ImageRegistrationMethod::SetShrinkFactorsPerLevel
Self & SetShrinkFactorsPerLevel(const std::vector< unsigned int > &shrinkFactors)
Set the isotropic shrink factors for each level.
itk::simple::sitkFloat32
@ sitkFloat32
32 bit float
Definition: sitkPixelIDValues.h:100
itk::simple::CompositeTransform
This class contains a stack of transforms and concatenates them by composition.
Definition: sitkCompositeTransform.h:58
itk::simple::ReadImage
SITKIO_EXPORT Image ReadImage(const std::string &filename, PixelIDValueEnum outputPixelType=sitkUnknown, const std::string &imageIO="")
ReadImage is a procedural interface to the ImageFileReader class which is convenient for most image r...
itk::simple::ImageRegistrationMethod::SetMetricAsANTSNeighborhoodCorrelation
Self & SetMetricAsANTSNeighborhoodCorrelation(unsigned int radius)
Use normalized cross correlation using a small neighborhood for each voxel between two images,...
itk::simple::ImageRegistrationMethod
An interface method to the modular ITKv4 registration framework.
Definition: sitkImageRegistrationMethod.h:89
itk::simple::ImageRegistrationMethod::SetInitialTransform
Self & SetInitialTransform(const Transform &transform)
Set the initial transform and parameters to optimize.
itk::simple::ImageRegistrationMethod::Execute
Transform Execute(const Image &fixed, const Image &moving)
Optimize the configured registration problem.
itk::simple::sitkIterationEvent
@ sitkIterationEvent
Occurs with some algorithms that run for a fixed or undetermined number of iterations.
Definition: sitkEvent.h:46
itk::simple::Command::Execute
virtual void Execute()
itk::simple::AffineTransform
An affine transformation about a fixed center with translation for a 2D or 3D coordinate.
Definition: sitkAffineTransform.h:35
itk::simple::ImageRegistrationMethod::SetInterpolator
Self & SetInterpolator(InterpolatorEnum Interpolator)
Set and get the interpolator to use.
Definition: sitkImageRegistrationMethod.h:117
itk::simple::ImageRegistrationMethod::SetMetricAsJointHistogramMutualInformation
Self & SetMetricAsJointHistogramMutualInformation(unsigned int numberOfHistogramBins=20, double varianceForJointPDFSmoothing=1.5)
Use mutual information between two images.
itk::simple::sitkMultiResolutionIterationEvent
@ sitkMultiResolutionIterationEvent
Occurs when some filters change processing to a different scale.
Definition: sitkEvent.h:59
itk::simple::ImageRegistrationMethod::EstimateLearningRateType
EstimateLearningRateType
Definition: sitkImageRegistrationMethod.h:255
itk::simple
Definition: sitkAdditionalProcedures.h:29
itk::simple::ImageRegistrationMethod::GetOptimizerIteration
unsigned int GetOptimizerIteration() const
itk::simple::Image::CopyInformation
void CopyInformation(const Image &srcImage)
Copy common meta-data from an image to this one.
itk::simple::CenteredTransformInitializer
Transform CenteredTransformInitializer(const Image &fixedImage, const Image &movingImage, const Transform &transform, CenteredTransformInitializerFilter::OperationModeType operationMode=itk::simple::CenteredTransformInitializerFilter::MOMENTS)
CenteredTransformInitializer is a helper class intended to initialize the center of rotation and the ...
itk::simple::Transform::ToString
std::string ToString() const