#include <iostream>
#include <stdlib.h>
#include <iomanip>
{
public:
: m_Method(m)
{}
void
{
using sitk::operator<<;
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:
};
{
public:
: m_Method(m)
{}
void
{
using sitk::operator<<;
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:
const sitk::ImageRegistrationMethod & m_Method;
};
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 };
}
{
double learningRate = 1.0;
unsigned int numberOfIterations = 100;
double convergenceMinimumValue = 1e-6;
unsigned int convergenceWindowSize = 10;
learningRate, numberOfIterations, convergenceMinimumValue, convergenceWindowSize, estimateLearningRate);
}
IterationUpdate cmd(R);
MultiResolutionIterationUpdate cmd2(R);
std::cout << "-------" << std::endl;
std::cout << outTx1.
ToString() << std::endl;
const double varianceForUpdateField = 0.0;
const double varianceForTotalField = 1.5;
displacementTx.SetSmoothingGaussianOnUpdate(varianceForUpdateField, varianceForTotalField);
{
std::vector<unsigned int> shrinkFactors = { 3, 2, 1 };
std::vector<double> smoothingSigmas = { 2.0, 1.0, 1.0 };
}
{
double learningRate = 1.0;
unsigned int numberOfIterations = 300;
double convergenceMinimumValue = 1e-6;
unsigned int convergenceWindowSize = 10;
learningRate, numberOfIterations, convergenceMinimumValue, convergenceWindowSize, estimateLearningRate);
}
std::cout << "-------" << std::endl;
std::cout << displacementTx.
ToString() << std::endl;
return 0;
}
An implementation of the Command design pattern for callback.
An interface method to the modular ITKv4 registration framework.
Self & SetInterpolator(InterpolatorEnum Interpolator)
Set and get the interpolator to use.
Transform Execute(const Image &fixed, const Image &moving)
Optimize the configured registration problem.
Self & SetShrinkFactorsPerLevel(const std::vector< unsigned int > &shrinkFactors)
Set the isotropic shrink factors for each level.
Self & SetSmoothingSigmasPerLevel(const std::vector< double > &smoothingSigmas)
Set the sigmas of Gaussian used for smoothing.
std::string GetOptimizerStopConditionDescription() const
Self & SetMetricAsANTSNeighborhoodCorrelation(unsigned int radius)
Use normalized cross correlation using a small neighborhood for each voxel between two images,...
Self & SetMovingInitialTransform(const Transform &transform)
Set a fixed transform component towards moving domain.
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...
Self & SetMetricAsJointHistogramMutualInformation(unsigned int numberOfHistogramBins=20, double varianceForJointPDFSmoothing=1.5)
Use mutual information between two images.
@ EachIteration
Estimate learning rate at each iteration.
unsigned int GetOptimizerIteration() const
double GetMetricValue() const
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.
Self & SetInitialTransform(const Transform &transform)
Set the initial transform and parameters to optimize.
Self & MetricUseFixedImageGradientFilterOff()
Enable image gradient computation by a filter.
The Image class for SimpleITK.
unsigned int GetDimension() const
std::vector< unsigned int > GetSize() const
void CopyInformation(const Image &srcImage)
Copy common meta-data from an image to this one.
virtual int AddCommand(itk::simple::EventEnum event, itk::simple::Command &cmd)
Add a Command Object to observer the event.
SITKIO_EXPORT Image ReadImage(const PathType &filename, PixelIDValueEnum outputPixelType=sitkUnknown, const std::string &imageIO="")
ReadImage is a procedural interface to the ImageFileReader class which is convenient for most image r...
@ sitkLinear
N-D linear interpolation.
@ sitkMultiResolutionIterationEvent
Occurs when some filters change processing to a different scale.
@ sitkIterationEvent
Occurs with some algorithms that run for a fixed or undetermined number of iterations.
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 ...
SITKCommon_EXPORT void WriteTransform(const Transform &transform, const PathType &filename)