template <class TOptimizer>
{
public:
using Self = IterationCallback;
itkOverrideGetNameOfClassMacro(IterationCallback);
using OptimizerType = TOptimizer;
void
SetOptimizer(OptimizerType * optimizer)
{
m_Optimizer = optimizer;
m_Optimizer->AddObserver(itk::IterationEvent(), this);
}
void
{
}
void
{
if (typeid(event) == typeid(itk::StartEvent))
{
std::cout << std::endl << "Position Value";
std::cout << std::endl << std::endl;
}
else if (typeid(event) == typeid(itk::IterationEvent))
{
std::cout << m_Optimizer->GetCurrentIteration() << " ";
std::cout << m_Optimizer->GetValue() << " ";
std::cout << m_Optimizer->GetCurrentPosition() << std::endl;
}
else if (typeid(event) == typeid(itk::EndEvent))
{
std::cout << std::endl << std::endl;
std::cout << "After " << m_Optimizer->GetCurrentIteration();
std::cout << " iterations " << std::endl;
std::cout << "Solution is = " << m_Optimizer->GetCurrentPosition();
std::cout << std::endl;
}
}
protected:
IterationCallback() = default;
};
template <typename TFixedImage, typename TMovingSpatialObject>
class SimpleImageToSpatialObjectMetric
{
public:
using Self = SimpleImageToSpatialObjectMetric;
using PointListType = std::list<PointType>;
using MovingSpatialObjectType = TMovingSpatialObject;
using typename Superclass::ParametersType;
using typename Superclass::DerivativeType;
using typename Superclass::MeasureType;
itkOverrideGetNameOfClassMacro(SimpleImageToSpatialObjectMetric);
static constexpr unsigned int ParametricSpaceDimension = 3;
void
{
if (!this->m_FixedImage)
{
std::cout << "Please set the image before the moving spatial object"
<< std::endl;
return;
}
this->m_MovingSpatialObject = object;
m_PointList.clear();
using myIteratorType =
itk::ImageRegionConstIteratorWithIndex<TFixedImage>;
myIteratorType it(this->m_FixedImage,
this->m_FixedImage->GetBufferedRegion());
itk::Point<double, 2>
point;
while (!it.IsAtEnd())
{
this->m_FixedImage->TransformIndexToPhysicalPoint(it.GetIndex(), point);
if (this->m_MovingSpatialObject->IsInsideInWorldSpace(point, 99999))
{
m_PointList.push_back(point);
}
++it;
}
std::cout << "Number of points in the metric = "
<< static_cast<unsigned long>(m_PointList.size()) << std::endl;
}
void
GetDerivative(
const ParametersType &, DerivativeType &)
const override
{
return;
}
MeasureType
GetValue(
const ParametersType & parameters)
const override
{
double value;
this->m_Transform->SetParameters(parameters);
value = 0;
for (auto point : m_PointList)
{
this->m_Transform->TransformPoint(point);
if (this->m_Interpolator->IsInsideBuffer(transformedPoint))
{
value += this->m_Interpolator->Evaluate(transformedPoint);
}
}
return value;
}
void
MeasureType & Value,
DerivativeType & Derivative) const override
{
}
private:
PointListType m_PointList;
};
int
main(int argc, char * argv[])
{
if (argc > 1)
{
std::cerr << "Too many parameters " << std::endl;
std::cerr << "Usage: " << argv[0] << std::endl;
}
auto ellipse1 = EllipseType::New();
auto ellipse2 = EllipseType::New();
auto ellipse3 = EllipseType::New();
ellipse1->SetRadiusInObjectSpace(10.0);
ellipse2->SetRadiusInObjectSpace(10.0);
ellipse3->SetRadiusInObjectSpace(10.0);
EllipseType::TransformType::OffsetType offset;
offset[0] = 100.0;
offset[1] = 40.0;
ellipse1->GetModifiableObjectToParentTransform()->SetOffset(offset);
ellipse1->Update();
offset[0] = 40.0;
offset[1] = 150.0;
ellipse2->GetModifiableObjectToParentTransform()->SetOffset(offset);
ellipse2->Update();
offset[0] = 150.0;
offset[1] = 150.0;
ellipse3->GetModifiableObjectToParentTransform()->SetOffset(offset);
ellipse3->Update();
auto group = GroupType::New();
group->AddChild(ellipse1);
group->AddChild(ellipse2);
group->AddChild(ellipse3);
using SpatialObjectToImageFilterType =
auto imageFilter = SpatialObjectToImageFilterType::New();
imageFilter->SetInput(group);
ImageType::SizeType size;
size[0] = 200;
size[1] = 200;
imageFilter->SetSize(size);
imageFilter->Update();
using GaussianFilterType =
auto gaussianFilter = GaussianFilterType::New();
gaussianFilter->SetInput(imageFilter->GetOutput());
constexpr double variance = 20;
gaussianFilter->SetVariance(variance);
gaussianFilter->Update();
using RegistrationType =
auto registration = RegistrationType::New();
using MetricType = SimpleImageToSpatialObjectMetric<ImageType, GroupType>;
auto metric = MetricType::New();
using InterpolatorType =
auto interpolator = InterpolatorType::New();
auto optimizer = OptimizerType::New();
auto transform = TransformType::New();
generator->Initialize(12345);
optimizer->SetNormalVariateGenerator(generator);
optimizer->Initialize(10);
optimizer->SetMaximumIteration(400);
TransformType::ParametersType parametersScale;
parametersScale.set_size(3);
parametersScale[0] = 1000;
for (unsigned int i = 1; i < 3; ++i)
{
parametersScale[i] = 2;
}
optimizer->SetScales(parametersScale);
using IterationCallbackType = IterationCallback<OptimizerType>;
auto callback = IterationCallbackType::New();
callback->SetOptimizer(optimizer);
registration->SetFixedImage(gaussianFilter->GetOutput());
registration->SetMovingSpatialObject(group);
registration->SetTransform(transform);
registration->SetInterpolator(interpolator);
registration->SetOptimizer(optimizer);
registration->SetMetric(metric);
TransformType::ParametersType initialParameters(
transform->GetNumberOfParameters());
initialParameters[0] = 0.2;
initialParameters[1] = 7.0;
initialParameters[2] = 6.0;
registration->SetInitialTransformParameters(initialParameters);
std::cout << "Initial Parameters : " << initialParameters << std::endl;
optimizer->MaximizeOn();
try
{
registration->Update();
std::cout << "Optimizer stop condition: "
<< registration->GetOptimizer()->GetStopConditionDescription()
<< std::endl;
}
{
std::cerr << "Exception caught ! " << std::endl;
std::cerr << exp << std::endl;
}
const RegistrationType::ParametersType finalParameters =
registration->GetLastTransformParameters();
std::cout << "Final Solution is : " << finalParameters << std::endl;
return EXIT_SUCCESS;
}
Superclass for callback/observer methods.
virtual void Execute(Object *caller, const EventObject &event)=0
Blurs an image by separable convolution with discrete gaussian kernels. This filter performs Gaussian...
Abstraction of the Events used to communicating among filters and with GUIs.
Standard exception handling object.
Representation of a group based on the spatial object classes.
Computes similarity between a moving spatial object and an Image to be registered.
virtual void SetMovingSpatialObject(const MovingSpatialObjectType *_arg)
void GetValueAndDerivative(const ParametersType ¶meters, MeasureType &Value, DerivativeType &Derivative) const override=0
Base class for Image Registration Methods.
Templated n-dimensional image class.
Linearly interpolate an image at specified positions.
Base class for most ITK classes.
1+1 evolutionary strategy optimizer
A templated class holding a geometric point in n-Dimensional space.
virtual MeasureType GetValue(const ParametersType ¶meters) const =0
virtual void GetDerivative(const ParametersType ¶meters, DerivativeType &derivative) const =0
Implements transparent reference counting.
Base class for filters that take a SpatialObject as input and produce an image as output....
Implements a weak reference to an object.
SmartPointer< const Self > ConstPointer
BinaryGeneratorImageFilter< TInputImage1, TInputImage2, TOutputImage > Superclass
SmartPointer< Self > Pointer
ImageBaseType::PointType PointType
*par Constraints *The filter image with at least two dimensions and a vector *length of at least The theory supports extension to scalar but *the implementation of the itk vector classes do not **The template parameter TRealType must be floating point(float or double) or *a user-defined "real" numerical type with arithmetic operations defined *sufficient to compute derivatives. **\par Performance *This filter will automatically multithread if run with *SetUsePrincipleComponents