[Insight-users] if there are somethings wrong with deform image using deform fields?
Luis Ibanez
luis.ibanez at kitware.com
Mon Feb 2 14:25:50 EST 2009
Hi Han,
Your code looks good, except for the fact that you
are computing the displacement backwards
Instead of:
displacement = fixedPoint -movingPoint;
You should do:
displacement = movingPoint -fixedPoint;
This is clearer if you look at the source code
of the WarpImageFilter, to see how the displacement
is going to be used during resampling:
Insight/Code/BasicFilters/itkWarpImageFilter.txx
lines 202-231:
while( !outputIt.IsAtEnd() )
{
// get the output image index
index = outputIt.GetIndex();
outputPtr->TransformIndexToPhysicalPoint( index, point );
// get the required displacement
displacement = fieldIt.Get();
// compute the required input image point
for(unsigned int j = 0; j < ImageDimension; j++ )
{
point[j] += displacement[j];
}
// get the interpolated value
if( m_Interpolator->IsInsideBuffer( point ) )
{
PixelType value = static_cast<PixelType>(
m_Interpolator->Evaluate( point ) );
outputIt.Set( value );
}
else
{
outputIt.Set( m_EdgePaddingValue );
}
++outputIt;
++fieldIt;
progress.CompletedPixel();
}
Regards,
Luis
----------------
Han D. F wrote:
> I use the bspline deform field (which is computed from the bapline
> registration algorithm) to wrap the image to the new image.
> When I use the the bspline parameters, the result is ok.
> when using the defor field, the result is wrong.
> Below is the code, id there is some thing wrong?
> thnanks
>
> assuming the bsplineTransform is computed and I will compute the deform
> field and warp to a new image.
> help me check the code and many thanks to your help
>
> // define the types for computing deform field
> typedef itk::Point< float, ImageDimension > PointType;
> typedef itk::Vector< float, ImageDimension > VectorType;
> typedef itk::Image< VectorType, ImageDimension > DeformationFieldType;
> DeformationFieldType::Pointer field = DeformationFieldType::New();
> field->SetRegions( fixedRegion );
> field->SetOrigin( fixedOrigin );
> field->SetSpacing( fixedSpacing );
> field->SetDirection( fixedDirection );
> field->Allocate();
> typedef itk::ImageRegionIterator< DeformationFieldType > FieldIterator;
> FieldIterator fi( field, fixedRegion );
> fi.GoToBegin();
> TransformType::InputPointType fixedPoint;
> TransformType::OutputPointType movingPoint;
> DeformationFieldType::IndexType index;
> VectorType displacement;
> // getting the deform field, I have see the result the computed
> displacement is right
> while( ! fi.IsAtEnd() )
> {
> index = fi.GetIndex();
> field->TransformIndexToPhysicalPoint( index, fixedPoint );
> movingPoint = bsplineTransform->TransformPoint( fixedPoint );
> displacement = fixedPoint -movingPoint ;
> fi.Set( displacement );
> //std::cout<<displacement<<std::endl;
> ++fi;
>
> }
>
> typedef itk::WarpImageFilter< FixedImageType,
> FixedImageType,
> DeformationFieldType > DeformFilterType;
> DeformFilterType::Pointer filter = DeformFilterType::New();
> // fiil the deform filter,such as SetInterpolator,SetOutputSpacing,....
>
> typedef itk::LinearInterpolateImageFunction<
> FixedImageType, double > InterpolatorType1;
> InterpolatorType1::Pointer interpolator1 = InterpolatorType1::New();
> filter->SetInterpolator( interpolator1 );
> filter->SetOutputSpacing( field->GetSpacing() );
> filter->SetOutputOrigin( field->GetOrigin() );
> filter->SetDeformationField( field );
> filter->SetInput( fixedReader->GetOutput() );
>
> //define the writer
>
> typedef itk::ImageFileWriter< FixedImageType > WriterType;
> WriterType::Pointer writer = WriterType::New();
> writer->SetFileName( argv[5] );
> writer->SetInput( filter->GetOutput() );
> // try the wrap filter and get the result
> try
> {
> writer->Update();
> }
> catch( itk::ExceptionObject & excp )
> {
> std::cerr << "Exception thrown " << std::endl;
> std::cerr << excp << std::endl;
> }
>
>
>
>
>
> ------------------------------------------------------------------------
> Ãâ·ÑËÍÄãºÍ°®ÈËȥмÓƹýÇéÈ˽ڣ¡
> <http://love.mail.163.com/valentine/main.do>
>
>
> ------------------------------------------------------------------------
>
> _______________________________________________
> 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