1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213
|
/*=========================================================================
Program: Insight Segmentation & Registration Toolkit
Module: $RCSfile: itkFastMarchingTest.cxx,v $
Language: C++
Date: $Date: 2008-02-03 04:05:34 $
Version: $Revision: 1.29 $
Copyright (c) Insight Software Consortium. All rights reserved.
See ITKCopyright.txt or http://www.itk.org/HTML/Copyright.htm for details.
This software is distributed WITHOUT ANY WARRANTY; without even
the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
PURPOSE. See the above copyright notices for more information.
=========================================================================*/
#if defined(_MSC_VER)
#pragma warning ( disable : 4786 )
#endif
#include "itkFastMarchingImageFilter.h"
#include "itkImage.h"
#include "itkImageRegionIterator.h"
#include "itkTextOutput.h"
#include "itkCommand.h"
#include "vnl/vnl_math.h"
namespace{
// The following class is used to support callbacks
// on the filter in the pipeline that follows later
class ShowProgressObject
{
public:
ShowProgressObject(itk::ProcessObject* o)
{m_Process = o;}
void ShowProgress()
{std::cout << "Progress " << m_Process->GetProgress() << std::endl;}
itk::ProcessObject::Pointer m_Process;
};
}
int itkFastMarchingTest(int, char* [] )
{
itk::OutputWindow::SetInstance(itk::TextOutput::New().GetPointer());
// create a fastmarching object
typedef float PixelType;
typedef itk::Image<PixelType,2> FloatImage;
typedef itk::FastMarchingImageFilter<FloatImage,FloatImage> FloatFMType;
FloatFMType::Pointer marcher = FloatFMType::New();
ShowProgressObject progressWatch(marcher);
itk::SimpleMemberCommand<ShowProgressObject>::Pointer command;
command = itk::SimpleMemberCommand<ShowProgressObject>::New();
command->SetCallbackFunction(&progressWatch,
&ShowProgressObject::ShowProgress);
marcher->AddObserver( itk::ProgressEvent(), command);
typedef FloatFMType::NodeType NodeType;
typedef FloatFMType::NodeContainer NodeContainer;
// setup alive points
NodeContainer::Pointer alivePoints = NodeContainer::New();
NodeType node;
FloatImage::OffsetType offset0 = {{28,35}};
itk::Index<2> index;
index.Fill(0);
node.SetValue( 0.0 );
node.SetIndex( index + offset0 );
alivePoints->InsertElement(0, node);
node.SetValue( 42.0 );
index.Fill( 200 );
node.SetIndex( index ); // this node is out of range
alivePoints->InsertElement(1, node);
marcher->SetAlivePoints( alivePoints );
// setup trial points
NodeContainer::Pointer trialPoints = NodeContainer::New();
node.SetValue( 1.0 );
index.Fill(0);
index += offset0;
index[0] += 1;
node.SetIndex( index );
trialPoints->InsertElement(0, node);
index[0] -= 1;
index[1] += 1;
node.SetIndex( index );
trialPoints->InsertElement(1, node);
index[0] -= 1;
index[1] -= 1;
node.SetIndex( index );
trialPoints->InsertElement(2, node);
index[0] += 1;
index[1] -= 1;
node.SetIndex( index );
trialPoints->InsertElement(3, node);
node.SetValue( 42.0 );
index.Fill( 300 ); // this node is out of ranage
node.SetIndex( index );
trialPoints->InsertElement(4, node);
marcher->SetTrialPoints( trialPoints );
// specify the size of the output image
FloatImage::SizeType size = {{64,64}};
marcher->SetOutputSize( size );
// setup a speed image of ones
FloatImage::Pointer speedImage = FloatImage::New();
FloatImage::RegionType region;
region.SetSize( size );
speedImage->SetLargestPossibleRegion( region );
speedImage->SetBufferedRegion( region );
speedImage->Allocate();
itk::ImageRegionIterator<FloatImage>
speedIter( speedImage, speedImage->GetBufferedRegion() );
for ( ; !speedIter.IsAtEnd(); ++speedIter )
{
speedIter.Set( 1.0 );
}
speedImage->Print( std::cout );
marcher->SetInput( speedImage );
marcher->SetStoppingValue( 100.0 );
// turn on debugging
marcher->DebugOn();
// update the marcher
marcher->Update();
// check the results
FloatImage::Pointer output = marcher->GetOutput();
itk::ImageRegionIterator<FloatImage>
iterator( output, output->GetBufferedRegion() );
bool passed = true;
for ( ; !iterator.IsAtEnd(); ++iterator )
{
FloatImage::IndexType tempIndex;
double distance;
float outputValue;
tempIndex = iterator.GetIndex();
tempIndex -= offset0;
distance = 0.0;
for ( int j = 0; j < 2; j++ )
{
distance += tempIndex[j] * tempIndex[j];
}
distance = vcl_sqrt( distance );
outputValue = (float) iterator.Get();
if (distance == 0)
{
continue;
}
if ( vnl_math_abs( outputValue ) / distance > 1.42 )
{
std::cout << iterator.GetIndex() << " ";
std::cout << vnl_math_abs( outputValue ) / distance << " ";
std::cout << vnl_math_abs( outputValue ) << " " << distance << std::endl;
passed = false;
}
}
// Exercise other member functions
std::cout << "SpeedConstant: " << marcher->GetSpeedConstant() << std::endl;
std::cout << "StoppingValue: " << marcher->GetStoppingValue() << std::endl;
std::cout << "CollectPoints: " << marcher->GetCollectPoints() << std::endl;
marcher->SetNormalizationFactor( 2.0 );
std::cout << "NormalizationFactor: " << marcher->GetNormalizationFactor();
std::cout << std::endl;
std::cout << "SpeedImage: " << marcher->GetInput();
std::cout << std::endl;
marcher->Print( std::cout );
if ( passed )
{
std::cout << "Fast Marching test passed" << std::endl;
return EXIT_SUCCESS;
}
else
{
std::cout << "Fast Marching test failed" << std::endl;
return EXIT_FAILURE;
}
}
|