[Insight-users] setting the spacing when writing an image to disk
Lucas Lorenzo
lucas at cvrti.utah.edu
Thu, 22 Apr 2004 10:49:18 -0600
--Apple-Mail-3--881340016
Content-Type: multipart/alternative;
boundary=Apple-Mail-4--881340015
--Apple-Mail-4--881340015
Content-Transfer-Encoding: 7bit
Content-Type: text/plain;
charset=US-ASCII;
format=flowed
Hi Luis,
I'm attaching the file GeoACPSLS.cxx (which is a modified version of
itkGeodesicActiveContourShapePriorLevelSetImageFilterTest.cxx).
When compiling against released version 1.6 the output after running
one single iteration is as follows (as you see all image files are in
vtk format):
lucas:GeodesicAcContPrSh> GeoACPSLS PrincipalModes_fr06 5 0.75 -0.4 12
0.2 2.0 0.15 0.8 1 slc5_fr_filt09.vtk 112 120
Potential map created successfully !
../DATAP091703/Anatomical/SLICE_05/PrincipalModes_fr06/mean.vtk
shape components set successfully !
Filter parameters set successfully !
Thresholder parameters set successfully !
1: [0.00025, 0, 0, 0, 0, 0, 0] 0.185567
The input file header is:
# vtk DataFile Version 3.0
VTK File Generated by Insight Segmentation and Registration Toolkit
(ITK)
BINARY
DATASET STRUCTURED_POINTS
DIMENSIONS 256 256 1
SPACING 1.33 1.33 1.0
ORIGIN 0 0 0.0
POINT_DATA 65536
SCALARS scalars unsigned_short 1
LOOKUP_TABLE default
And the header file for the vtk image resulting from the fast marching
algorithm is like this:
# vtk DataFile Version 3.0
VTK File Generated by Insight Segmentation and Registration Toolkit
(ITK)
BINARY
DATASET STRUCTURED_POINTS
DIMENSIONS 256 256 1
SPACING 1 1 1.0
ORIGIN 0 0 0.0
POINT_DATA 65536
SCALARS scalars float 1
LOOKUP_TABLE default
When compiling against the CVS copy I checked out about one week and a
half ago the output is like this:
lucas:GeodesicAcContPrSh> GeoACPSLS PrincipalModes_fr06 5 0.75 -0.4 12
0.2 2.0 0.15 0.8 1 slc5_fr_filt09.vtk 112 120
Potential map created successfully !
../DATAP091703/Anatomical/SLICE_05/PrincipalModes_fr06/mean.vtk
shape components set successfully !
Filter parameters set successfully !
Thresholder parameters set successfully !
inputRequestedRegion: ImageRegion (0xbfffeb80)
Dimension: 2
Index: [-1, -1]
Size: [258, 258]
largestPossibleRegion: ImageRegion (0x1101f1c)
Dimension: 2
Index: [0, 0]
Size: [256, 256]
inputRequestedRegion: ImageRegion (0xbfffeb80)
Dimension: 2
Index: [-1, -1]
Size: [258, 258]
largestPossibleRegion: ImageRegion (0x1101f1c)
Dimension: 2
Index: [0, 0]
Size: [256, 256]
1: [0.00025, 0, 0, 0, 0, 0, 0] 0.185567
And the header files are still as before (that is, the spacing problem
still exists).
I hope that this information is clear enough. If not, please let me
know.
Thanks for your help,
Lucas
--Apple-Mail-4--881340015
Content-Transfer-Encoding: 7bit
Content-Type: text/enriched;
charset=US-ASCII
Hi Luis,
I'm attaching the file GeoACPSLS.cxx (which is a modified version of
itkGeodesicActiveContourShapePriorLevelSetImageFilterTest.cxx).
When compiling against released version 1.6 the output after running
one single iteration is as follows (as you see all image files are in
vtk format):
<fontfamily><param>Hei</param>lucas:GeodesicAcContPrSh> GeoACPSLS
PrincipalModes_fr06 5 0.75 -0.4 12 0.2 2.0 0.15 0.8 1
slc5_fr_filt09.vtk 112 120
Potential map created successfully !
../DATAP091703/Anatomical/SLICE_05/PrincipalModes_fr06/mean.vtk
shape components set successfully !
Filter parameters set successfully !
Thresholder parameters set successfully !
1: [0.00025, 0, 0, 0, 0, 0, 0] 0.185567</fontfamily>
The input file header is:
<fontfamily><param>Hei</param># vtk DataFile Version 3.0
VTK File Generated by Insight Segmentation and Registration Toolkit
(ITK)
BINARY
DATASET STRUCTURED_POINTS
DIMENSIONS 256 256 1
SPACING 1.33 1.33 1.0
ORIGIN 0 0 0.0
POINT_DATA 65536
SCALARS scalars unsigned_short 1
LOOKUP_TABLE default
</fontfamily>
And the header file for the vtk image resulting from the fast marching
algorithm is like this:
<fontfamily><param>Hei</param># vtk DataFile Version 3.0
VTK File Generated by Insight Segmentation and Registration Toolkit
(ITK)
BINARY
DATASET STRUCTURED_POINTS
DIMENSIONS 256 256 1
SPACING 1 1 1.0
ORIGIN 0 0 0.0
POINT_DATA 65536
SCALARS scalars float 1
LOOKUP_TABLE default</fontfamily>
When compiling against the CVS copy I checked out about one week and a
half ago the output is like this:
<fontfamily><param>Hei</param>lucas:GeodesicAcContPrSh> GeoACPSLS
PrincipalModes_fr06 5 0.75 -0.4 12 0.2 2.0 0.15 0.8 1
slc5_fr_filt09.vtk 112 120
Potential map created successfully !
../DATAP091703/Anatomical/SLICE_05/PrincipalModes_fr06/mean.vtk
shape components set successfully !
Filter parameters set successfully !
Thresholder parameters set successfully !
inputRequestedRegion: ImageRegion (0xbfffeb80)
Dimension: 2
Index: [-1, -1]
Size: [258, 258]
largestPossibleRegion: ImageRegion (0x1101f1c)
Dimension: 2
Index: [0, 0]
Size: [256, 256]
inputRequestedRegion: ImageRegion (0xbfffeb80)
Dimension: 2
Index: [-1, -1]
Size: [258, 258]
largestPossibleRegion: ImageRegion (0x1101f1c)
Dimension: 2
Index: [0, 0]
Size: [256, 256]
1: [0.00025, 0, 0, 0, 0, 0, 0] 0.185567</fontfamily>
And the header files are still as before (that is, the spacing problem
still exists).
I hope that this information is clear enough. If not, please let me
know.
Thanks for your help,
Lucas
--Apple-Mail-4--881340015--
--Apple-Mail-3--881340016
Content-Transfer-Encoding: 7bit
Content-Type: application/octet-stream;
x-unix-mode=0644;
name="GeoACPSLS.cxx"
Content-Disposition: attachment;
filename=GeoACPSLS.cxx
#include "itkImage.h"
#include "itkImageFileReader.h"
#include "itkImageFileWriter.h"
#include "vnl/vnl_matrix_fixed.h"
#include "vnl/vnl_math.h"
#include <cstdio>
#include <cstring>
#include <cstdlib>
#include "itkGeodesicActiveContourShapePriorLevelSetImageFilter.h"
#include "itkPCAShapeSignedDistanceFunction.h"
#include "itkShapePriorMAPCostFunction.h"
#include "itkAmoebaOptimizer.h"
#include "itkCastImageFilter.h"
#include "itkGradientMagnitudeRecursiveGaussianImageFilter.h"
#include "itkSigmoidImageFilter.h"
#include "itkFastMarchingImageFilter.h"
#include "itkNumericSeriesFileNames.h"
#include "itkBinaryThresholdImageFilter.h"
#include "itkCommand.h"
#include "itkGradientMagnitudeImageFilter.h"
/** This class is used to support callbacks
* on the segmentation filter in this test. */
namespace {
template<class TFilter>
class ShowIterationObject
{
public:
ShowIterationObject( TFilter * filter )
{ m_Filter = filter; }
void ShowIteration()
{
std::cout << m_Filter->GetElapsedIterations() << ": ";
std::cout << m_Filter->GetCurrentParameters() << " ";
std::cout << m_Filter->GetRMSChange() << std::endl;
}
typename TFilter::Pointer m_Filter;
};
}
int main( int argc, char * argv[] )
{
if( argc < 14 )
{
std::cerr << "Usage: " << argv[0];
std::cerr << " Basefilename NumberOfPrincipalComponents gradMagnitudesigma sigmoid_alpha sigmoid_beta propagation_scaling advec_esc curve_scaling shape_scaling filter_num_iter input_image ini_x ini_y";
//GeoACPSLS ../DATAP091703/Anatomical/SLICE_05/PrincipalModes_fr06 5 0.75 -0.4 12 0.2 2.0 0.15 0.8 1 ../DATAP091703/Anatomical/SLICE_05/slc5_fr_filt09.vtk 112 120
std::cerr << std::endl;
return 1;
}
// ********************************************************************* //
/**
* Read the numPComp eigenvalues of the principal components
*/
long numPComp;
char * pEnd;
numPComp = std::strtol(argv[2],&pEnd,0);
vnl_vector<double> eigenValues(numPComp);
std::FILE *infile;
char *eigen_base = argv[1];
std::string eigenvalues = eigen_base;
eigenvalues += "/eigenvalues.dat";
infile = fopen(eigenvalues.c_str(),"r");
if(infile==0) {
std::printf("Error: can't access eigenvalues.dat.\n");
return 1;
}
for(unsigned int i= 0; i<numPComp ; i++ )
{
std::fseek (infile,(i-1) * sizeof(double),SEEK_SET);
std::fread( static_cast<void*>(&(eigenValues[i])) , sizeof(double) , 1 , infile);
}
std::fclose(infile);
// ********************************************************************* //
// ********************************************************************* //
/* Typedefs of components. */
const unsigned int ImageDimension = 2;
typedef unsigned short InputPixelType;
typedef unsigned char BinaryPixelType;
typedef double PCAPixelType;
typedef float InternalPixelType;
typedef itk::Image<InputPixelType,ImageDimension> OriginalImageType;
typedef itk::Image<PCAPixelType,ImageDimension> PCAImageType;
typedef itk::Image<InternalPixelType,ImageDimension> InternalImageType;
typedef itk::Image<BinaryPixelType,ImageDimension> BinaryImageType;
typedef itk::GeodesicActiveContourShapePriorLevelSetImageFilter<InternalImageType,InternalImageType> FilterType;
typedef itk::PCAShapeSignedDistanceFunction<double,ImageDimension> ShapeFunctionType;
typedef itk::ShapePriorMAPCostFunction<InternalImageType,InternalPixelType> CostFunctionType;
typedef itk::AmoebaOptimizer OptimizerType;
typedef FilterType::ParametersType ParametersType;
FilterType::Pointer filter = FilterType::New();
ShapeFunctionType::Pointer shape = ShapeFunctionType::New();
CostFunctionType::Pointer costFunction = CostFunctionType::New();
OptimizerType::Pointer optimizer = OptimizerType::New();
//create a reader for the image to be segmented:
typedef itk::ImageFileReader< OriginalImageType > OriginalReaderType;
OriginalReaderType::Pointer original_reader = OriginalReaderType::New();
// ********************************************************************* //
// ********************************************************************* //
char *input_image = argv[11];
original_reader->SetFileName(input_image);
/**
* Create an edge potential map.
* First compute the image gradient magnitude using a derivative of gaussian filter.
* Then apply a sigmoid function to the gradient magnitude.
*/
typedef itk::CastImageFilter< OriginalImageType, InternalImageType > CastFilterType;
CastFilterType::Pointer caster = CastFilterType::New();
caster->SetInput(original_reader->GetOutput() );
typedef itk::GradientMagnitudeImageFilter< InternalImageType,InternalImageType > GradientImageType;
GradientImageType::Pointer gradMagnitude = GradientImageType::New();
gradMagnitude->SetInput( caster->GetOutput() );
// write the gradient magnitude image to a file
typedef itk::ImageFileWriter< InternalImageType > magWriterType;
magWriterType::Pointer mag_writer = magWriterType::New();
mag_writer->SetFileName("grad_mag.vtk");
mag_writer->SetInput(gradMagnitude->GetOutput());
mag_writer->Update();
typedef itk::SigmoidImageFilter< InternalImageType, InternalImageType >
SigmoidFilterType;
SigmoidFilterType::Pointer sigmoid = SigmoidFilterType::New();
sigmoid->SetOutputMinimum( 0.0 );
sigmoid->SetOutputMaximum( 1.0 );
float alpha = std::atof(argv[4]);
sigmoid->SetAlpha( alpha );
float beta = std::atof(argv[5]);
sigmoid->SetBeta( beta );
sigmoid->SetInput( gradMagnitude->GetOutput() );
sigmoid->Update();// added by me
//the following lines are for writing the feature image to a file:
typedef itk::ImageFileWriter< InternalImageType > SigWriterType;
SigWriterType::Pointer sig_writer = SigWriterType::New();
sig_writer->SetFileName("feature.vtk");
sig_writer->SetInput(sigmoid->GetOutput());
sig_writer->Update();
/**
* Create an initial level.
* Use fast marching to create an signed distance from a seed point.
*/
typedef itk::FastMarchingImageFilter<InternalImageType> FastMarchingFilterType;
FastMarchingFilterType::Pointer fastMarching = FastMarchingFilterType::New();
typedef FastMarchingFilterType::NodeContainer NodeContainer;
typedef FastMarchingFilterType::NodeType NodeType;
NodeContainer::Pointer seeds = NodeContainer::New();
// Choose an initial contour that is within the left ventricle
// The initial contour is a circle centered at {x_ini,y_ini} with radius
// 2.0
int x_ini = std::atoi(argv[12]);
int y_ini = std::atoi(argv[13]);
InternalImageType::IndexType seedPosition;
seedPosition[0] = x_ini;
seedPosition[1] = y_ini;
NodeType node;
node.SetValue( -2.0 );
node.SetIndex( seedPosition );
seeds->Initialize();
seeds->InsertElement( 0, node );
fastMarching->SetTrialPoints( seeds );
fastMarching->SetSpeedConstant( 1.0 );
OriginalImageType::RegionType inputRegion = original_reader->GetOutput()->GetLargestPossibleRegion();
OriginalImageType::SizeType imageSize = inputRegion.GetSize();
fastMarching->SetOutputSize( imageSize );
try
{
fastMarching->Update();
}
catch( itk::ExceptionObject & exp )
{
std::cerr << "Exception caught ! fastMarching" << std::endl;
std::cerr << exp << std::endl;
return -1;
}
typedef itk::ImageFileWriter< InternalImageType > FMWriterType;
FMWriterType::Pointer FM_writer = FMWriterType::New();
FM_writer->SetFileName("fast_marching.vtk");
FM_writer->SetInput(fastMarching->GetOutput());
FM_writer->Update();
// ********************************************************************* //
std::cout << "Potential map created successfully !" << std::endl;
// ********************************************************************* //
/**
* Set up the components of the shape prior segmentation filter
*/
// Set up the shape function
//
typedef ShapeFunctionType::ImageType ComponentImageType;
//create a reader for the PCA images
typedef itk::ImageFileReader<ComponentImageType > PCAReaderType;
PCAReaderType::Pointer PCAreader = PCAReaderType::New();
// Load the PCA components:
//
typedef ShapeFunctionType::ImagePointerVector ImageVectorType;
ImageVectorType pca;
pca.resize(numPComp );
// generate the file names:
typedef itk::NumericSeriesFileNames NumSerFNM;
NumSerFNM::Pointer filenames_gen = NumSerFNM::New();
filenames_gen->SetStartIndex(1);
filenames_gen->SetEndIndex(numPComp);
char *base1 = argv[1];
std::string base = base1;
std::strcat(base1,"/mode%02d.vtk");
filenames_gen->SetSeriesFormat(base1);
std::vector<std::string > flname = filenames_gen->GetFileNames();
for(unsigned int i= 0; i<numPComp ; i++ )
{
PCAreader->SetFileName(flname[i].c_str());
PCAreader->Update();
pca[i] = PCAreader->GetOutput();
pca[i]->DisconnectPipeline();
}
// Set up a translation transform
typedef itk::TranslationTransform<double,ImageDimension> TransformType;
TransformType::Pointer transform = TransformType::New();
// Set up the standard deviations
TransformType::ParametersType pcaStdDev(numPComp );
for(unsigned int i= 0; i<numPComp ; i++ )
{
pcaStdDev[i] = sqrt(eigenValues[i]);
}
shape->SetNumberOfPrincipalComponents( numPComp );
base += "/mean.vtk";
std::cout << base << std::endl;
PCAreader->SetFileName(base.c_str());
PCAreader->Update();
shape->SetMeanImage(PCAreader->GetOutput() );
shape->SetPrincipalComponentImages( pca );
shape->SetTransform( transform );
shape->SetPrincipalComponentStandardDeviations( pcaStdDev );
try {
shape->Initialize();
}
catch( itk::ExceptionObject & err )
{
std::cout << "shape" << err << std::endl;
return EXIT_FAILURE;
}
// shape->Print(std::cout);
std::cout << "shape components set successfully !" << std::endl;
// Set up the cost function
CostFunctionType::ArrayType mean( shape->GetNumberOfShapeParameters() );
CostFunctionType::ArrayType stddev( shape->GetNumberOfShapeParameters() );
mean.Fill(0.0);
stddev.Fill(1.0);
costFunction->SetShapeParameterMeans( mean );
costFunction->SetShapeParameterStandardDeviations( stddev );
CostFunctionType::WeightsType weights;
weights.Fill( 1.0 );
weights[1] = 10.0;
costFunction->SetWeights( weights );
// Set up the optimizer
optimizer->SetFunctionConvergenceTolerance( 0.1 );
optimizer->SetParametersConvergenceTolerance( 0.5 );
optimizer->SetMaximumNumberOfIterations( 50 );
// optimizer->SetCostFunction(costFunction);
// Set up the initial parameters
ParametersType parameters( shape->GetNumberOfParameters() );
//let start with the mean shape and then let us set the pose to
//(0;0):
parameters.Fill(0);
typedef ParametersType::iterator Iterator;
Iterator it = parameters.end();
it--;
*it = 0;
it--;
*it = 0;
// Set up the scaling between the level set terms
float prop_esc = std::atof(argv[6]);
filter->SetPropagationScaling( prop_esc );
float advec_esc = std::atof(argv[7]);
filter->SetAdvectionScaling(advec_esc );
float curv_esc = std::atof(argv[8]);
filter->SetCurvatureScaling( curv_esc );
float shape_esc = std::atof(argv[9]);
filter->SetShapePriorScaling(shape_esc);
// Hook up components to the filter
filter->SetInput( fastMarching->GetOutput() ); // initial level set
filter->SetFeatureImage( sigmoid->GetOutput() ); // edge potential map
filter->SetShapeFunction( shape );
filter->SetCostFunction( costFunction );
filter->SetOptimizer( optimizer );
filter->SetInitialParameters( parameters );
filter->SetNumberOfLayers( 4 );
filter->SetMaximumRMSError( 0.01 );
long filter_num_iter = std::atol(argv[10]);
filter->SetNumberOfIterations(filter_num_iter );
// filter->Update();
std::cout << "Filter parameters set successfully !" << std::endl;
/**
* Connect an observer to the filter
*/
typedef ShowIterationObject<FilterType> WatcherType;
WatcherType iterationWatcher(filter);
itk::SimpleMemberCommand<WatcherType>::Pointer command =
itk::SimpleMemberCommand<WatcherType>::New();
command->SetCallbackFunction( &iterationWatcher,
&WatcherType::ShowIteration );
filter->AddObserver( itk::IterationEvent(), command );
/**
* Threshold the output level set to display the final contour.
*/
typedef itk::BinaryThresholdImageFilter< InternalImageType,BinaryImageType >
ThresholdFilterType;
ThresholdFilterType::Pointer thresholder = ThresholdFilterType::New()
;
typedef itk::ImageFileWriter<BinaryImageType > BinaryWriterType;
BinaryWriterType::Pointer binary_writer = BinaryWriterType::New();
thresholder->SetInput( filter->GetOutput() );
thresholder->SetLowerThreshold( -1e+10 );
thresholder->SetUpperThreshold( 0.0 );
thresholder->SetOutsideValue( 0 );
thresholder->SetInsideValue( 255 );
std::cout << "Thresholder parameters set successfully !" << std::endl;
//write the binary image to a file:
binary_writer->SetFileName("out.vtk");
binary_writer->SetInput(thresholder->GetOutput());
try
{
binary_writer->Update();
}
catch( itk::ExceptionObject & exp )
{
std::cerr << "Exception caught ! binary image writer" << std::endl;
std::cerr << exp << std::endl;
return -1;
}
return 0;
}
--Apple-Mail-3--881340016
Content-Type: multipart/alternative;
boundary=Apple-Mail-5--881340013
--Apple-Mail-5--881340013
Content-Transfer-Encoding: 7bit
Content-Type: text/plain;
charset=US-ASCII;
format=flowed
On Apr 21, 2004, at 7:51 PM, Luis Ibanez wrote:
>
> Hi Lucas,
>
> Can you please post a minimal (working) example
> of code that illustrates this behavior. We need
> some initial code for reproducing the problem
> that you report.
>
>
> Thanks
>
>
> Luis
>
>
> -----------------------
> Lucas Lorenzo wrote:
>
>> HI Luis,
>> I've checked out a CVS copy last week and I still have the same
>> problem.
>> Also now, when running my code the following messages appear:
>> inputRequestedRegion: ImageRegion (0xbfffeb70)
>> Dimension: 2
>> Index: [-1, -1]
>> Size: [258, 258]
>> largestPossibleRegion: ImageRegion (0x1101c5c)
>> Dimension: 2
>> Index: [0, 0]
>> Size: [256, 256]
>> inputRequestedRegion: ImageRegion (0xbfffeb70)
>> Dimension: 2
>> Index: [-1, -1]
>> Size: [258, 258]
>> largestPossibleRegion: ImageRegion (0x1101c5c)
>> Dimension: 2
>> Index: [0, 0]
>> Size: [256, 256]
>> Do you have any suggestions ?
>> Thanks,
>> Lucas
>> On Feb 19, 2004, at 12:05 AM, Luis Ibanez wrote:
>>>
>>> Hi Lucas,
>>>
>>> The changes to the FastMarchingImageFilter have
>>> been commited. The output image now uses the
>>> spacing and origin of the input image.
>>>
>>> Please let us know if you find any problem.
>>>
>>>
>>> Thanks
>>>
>>>
>>> Luis
>>>
>>>
>>> ------------------------
>>> Lucas Lorenzo wrote:
>>>
>>>> Hi Luis,
>>>> sorry for answering so late.
>>>> I've tried what you suggested but I have a run time error "Abort
>>>> trap" when trying to apply the GetOutput() method to my
>>>> FastMarchingImageFilter object.
>>>> Let me know if there are any other changes we could try.
>>>> Thanks,
>>>> Lucas
>>>> On Feb 11, 2004, at 5:58 AM, Luis Ibanez wrote:
>>>>
>>>>>
>>>>> Hi Lucas,
>>>>>
>>>>> Thanks for pointing this out.
>>>>>
>>>>>
>>>>> Please try the following:
>>>>>
>>>>> Edit the file:
>>>>>
>>>>> Insight/Code/Algorithms/
>>>>> itkFastMarchingImageFilter.txx
>>>>>
>>>>> Go to to line : 150
>>>>> and after the statement
>>>>>
>>>>>> output->Allocate();
>>>>>
>>>>>
>>>>>
>>>>> add
>>>>>
>>>>> output->CopyInformation( this->GetInput() );
>>>>>
>>>>> Then, go to line : 157
>>>>> and after the statement
>>>>>
>>>>>
>>>>>> m_LabelImage->Allocate();
>>>>>
>>>>>
>>>>>
>>>>>
>>>>> add
>>>>>
>>>>> m_LabelImage->CopyInformation( this->GetInput() );
>>>>>
>>>>>
>>>>> This should copy the origin and spacing
>>>>> of the input image into the output and
>>>>> label images.
>>>>>
>>>>>
>>>>> Then recompile your application and try
>>>>> running it to see if the correct spacing
>>>>> appears in the file.
>>>>>
>>>>>
>>>>> Please let us know what you find, so we proceed
>>>>> to do the same changes in the repository.
>>>>>
>>>>>
>>>>>
>>>>> Thanks
>>>>>
>>>>>
>>>>>
>>>>> Luis
>>>>>
>>>>>
>>>>>
>>>>> ----------------------
>>>>> Lucas Lorenzo wrote:
>>>>>
>>>>>> Hi Luis,
>>>>>> sorry for having such a mess in my code.
>>>>>> I'm using an application based on ITK/Testing/Code/Algorithms/
>>>>>> itkGeodesicActiveContourShapePriorLevelSetImageFilterTest_2.cxx
>>>>>> The spacing is being carried through the pipeline with no
>>>>>> problem, except when arriving to the point when I have to
>>>>>> generate my initial contour (signed distance map) from a seed
>>>>>> point using FastMarchingImageFilter. It is the output from this
>>>>>> filter the one that has the "default" spacing (and I think that
>>>>>> this new spacing is carried to the end to the output image) and
>>>>>> I can't find any method to set the correct spacing.
>>>>>> Thanks,
>>>>>> Lucas
>>>>>> On Tuesday, February 10, 2004, at 08:38 PM, Luis Ibanez wrote:
>>>>>> Hi Lucas,
>>>>>> Why are you setting the spacing on the ImageIO
>>>>>> object instead of the image itself ?
>>>>>> You should just carry the spacing through the
>>>>>> pipeline. Does your input image has an invalid
>>>>>> spacing ?
>>>>>> An option in that case is to use the
>>>>>> ChangeInformationImageFilter
>>>>>> http://www.itk.org/Insight/Doxygen/html/
>>>>>> classitk_1_1ChangeInformationImageFilter.html
>>>>>> This filter carries the data buffer of the
>>>>>> input image to the output image, and allows
>>>>>> you to alter the meta-data such as image
>>>>>> origin and spacing.
>>>>>> Please don't use this filter for processing
>>>>>> images of human beings or any other living
>>>>>> organisms, since chances are that you will
>>>>>> make somebody operate in a liver instead of
>>>>>> a lung.
>>>>>> In the long term the right thing to do is
>>>>>> to fix the source of your images which is
>>>>>> where the real spacing information should
>>>>>> be comming from.
>>>>>> Regards,
>>>>>> Luis
>>>>>> ======================================
>>>>>> -------------------
>>>>>> Lucas Lorenzo wrote:
>>>>>> Hi all,
>>>>>> I'm trying to write an image to disk in vtk format.
>>>>>> By default the spacing is set to 1 1 1. I'd like to
>>>>>> change it so
>>>>>> I'm doing the following:
>>>>>> /#include "itkVTKImageIO.h"
>>>>>> int main( int argc, char * argv[] )
>>>>>> {
>>>>>> /* /* Typedefs of components. */*/
>>>>>> const unsigned int ImageDimension = 2;
>>>>>> typedef unsigned char BinaryPixelType;
>>>>>> typedef itk::Image<BinaryPixelType,ImageDimension>
>>>>>> BinaryImageType;
>>>>>> /*// read the input image and get the spacing from it:*/
>>>>>> typedef itk::VTKImageIO ImageIOType;
>>>>>> ImageIOType::Pointer IO1 = ImageIOType::New();
>>>>>> original_reader->SetImageIO(IO1);
>>>>>> double dx,dy,dz;
>>>>>> original_reader->Update();
>>>>>> dx = IO1->GetSpacing(0);
>>>>>> dy = IO1->GetSpacing(1);
>>>>>> dz = IO1->GetSpacing(2);
>>>>>> /*// here I'm omitting when I process the input image
>>>>>> // write the image to a file but perviously set the
>>>>>> spacing
>>>>>> */
>>>>>> binary_writer->SetFileName("out.vtk");
>>>>>> binary_writer->SetInput(thresholder->GetOutput());
>>>>>> ImageIOType::Pointer IO2 = ImageIOType::New();
>>>>>> IO2->SetSpacing(0,dx);
>>>>>> IO2->SetSpacing(1,dy);
>>>>>> IO2->SetSpacing(2,dz);
>>>>>> binary_writer->SetImageIO(IO2);
>>>>>> try
>>>>>> {
>>>>>> binary_writer->Update();
>>>>>> }
>>>>>> catch( itk::ExceptionObject & exp )
>>>>>> {
>>>>>> std::cerr << "Exception caught ! binary image writer" <<
>>>>>> std::endl;
>>>>>> std::cerr << exp << std::endl;
>>>>>> return -1;
>>>>>> }
>>>>>> return 0;
>>>>>> }
>>>>>> /
>>>>>> When I execute this program I have a "segmentation fault"
>>>>>> run
>>>>>> time error.
>>>>>> If I ommit the line setting the spacing in the z axis ("/
>>>>>> IO2->SetSpacing(2,dz); /") the programs executes without
>>>>>> crashing but the spacing is not really set, that is, in
>>>>>> the
>>>>>> header of the output vtk file (out.vtk) you can still see
>>>>>> "SPACING 1 1 1.0" instead of "SPACING dx dy 1.0".
>>>>>> Any clue of what am I doing wrong ?
>>>>>> Thanks,
>>>>>> Lucas Lorenzo
>>>>>> University of Utah
>>>>>> Nora Eccles Harrison CardioVascular Research and Training
>>>>>> Institute
>>>>>> Fellows Room
>>>>>> 95 South 2000 East
>>>>>> Salt Lake City, UT 84112-5000
>>>>>> e-mail: lucas at cvrti.utah.edu
>>>>>> telephone: 801-587-9536
>>>>>> _______________________________________________
>>>>>> Insight-users mailing list
>>>>>> Insight-users at itk.org
>>>>>> http://www.itk.org/mailman/listinfo/insight-users
>>>>>> Lucas Lorenzo
>>>>>> University of Utah
>>>>>> Nora Eccles Harrison CardioVascular Research and Training
>>>>>> Institute
>>>>>> Fellows Room
>>>>>> 95 South 2000 East
>>>>>> Salt Lake City, UT 84112-5000
>>>>>> e-mail: lucas at cvrti.utah.edu
>>>>>> telephone: 801-587-9536
>>>>>
>>>>>
>>>>>
>>>>>
>>>>>
>>>>>
>>>> Lucas Lorenzo
>>>> University of Utah
>>>> Nora Eccles Harrison CardioVascular Research and Training Institute
>>>> Fellows Room
>>>> 95 South 2000 East
>>>> Salt Lake City, UT 84112-5000
>>>> e-mail: lucas at cvrti.utah.edu
>>>> telephone: 801-587-9536
>>>> _______________________________________________
>>>> Insight-users mailing list
>>>> Insight-users at itk.org
>>>> http://www.itk.org/mailman/listinfo/insight-users
>>>
>>>
>>>
>>>
>>>
>> Lucas Lorenzo
>> University of Utah
>> Nora Eccles Harrison CardioVascular Research and Training Institute
>> Fellows Room
>> 95 South 2000 East
>> Salt Lake City, UT 84112-5000
>> e-mail: lucas at cvrti.utah.edu
>> telephone: 801-587-9536
>
>
>
> _______________________________________________
> Insight-users mailing list
> Insight-users at itk.org
> http://www.itk.org/mailman/listinfo/insight-users
>
Lucas Lorenzo
University of Utah
Nora Eccles Harrison CardioVascular Research and Training Institute
Fellows Room
95 South 2000 East
Salt Lake City, UT 84112-5000
e-mail: lucas at cvrti.utah.edu
telephone: 801-587-9536
--Apple-Mail-5--881340013
Content-Transfer-Encoding: 7bit
Content-Type: text/enriched;
charset=US-ASCII
On Apr 21, 2004, at 7:51 PM, Luis Ibanez wrote:
<excerpt>
Hi Lucas,
Can you please post a minimal (working) example
of code that illustrates this behavior. We need
some initial code for reproducing the problem
that you report.
Thanks
Luis
-----------------------
Lucas Lorenzo wrote:
<excerpt>HI Luis,
I've checked out a CVS copy last week and I still have the same
problem.
Also now, when running my code the following messages appear:
inputRequestedRegion: ImageRegion (0xbfffeb70)
Dimension: 2
Index: [-1, -1]
Size: [258, 258]
largestPossibleRegion: ImageRegion (0x1101c5c)
Dimension: 2
Index: [0, 0]
Size: [256, 256]
inputRequestedRegion: ImageRegion (0xbfffeb70)
Dimension: 2
Index: [-1, -1]
Size: [258, 258]
largestPossibleRegion: ImageRegion (0x1101c5c)
Dimension: 2
Index: [0, 0]
Size: [256, 256]
Do you have any suggestions ?
Thanks,
Lucas
On Feb 19, 2004, at 12:05 AM, Luis Ibanez wrote:
<excerpt>
Hi Lucas,
The changes to the FastMarchingImageFilter have
been commited. The output image now uses the
spacing and origin of the input image.
Please let us know if you find any problem.
Thanks
Luis
------------------------
Lucas Lorenzo wrote:
<excerpt>Hi Luis,
sorry for answering so late.
I've tried what you suggested but I have a run time error "Abort trap"
when trying to apply the GetOutput() method to my
FastMarchingImageFilter object.
Let me know if there are any other changes we could try.
Thanks,
Lucas
On Feb 11, 2004, at 5:58 AM, Luis Ibanez wrote:
<excerpt>
Hi Lucas,
Thanks for pointing this out.
Please try the following:
Edit the file:
Insight/Code/Algorithms/
itkFastMarchingImageFilter.txx
Go to to line : 150
and after the statement
<excerpt> output->Allocate();
</excerpt>
add
output->CopyInformation( this->GetInput() );
Then, go to line : 157
and after the statement
<excerpt> m_LabelImage->Allocate();
</excerpt>
add
m_LabelImage->CopyInformation( this->GetInput() );
This should copy the origin and spacing
of the input image into the output and
label images.
Then recompile your application and try
running it to see if the correct spacing
appears in the file.
Please let us know what you find, so we proceed
to do the same changes in the repository.
Thanks
Luis
----------------------
Lucas Lorenzo wrote:
<excerpt>Hi Luis,
sorry for having such a mess in my code.
I'm using an application based on ITK/Testing/Code/Algorithms/
itkGeodesicActiveContourShapePriorLevelSetImageFilterTest_2.cxx The
spacing is being carried through the pipeline with no problem, except
when arriving to the point when I have to generate my initial contour
(signed distance map) from a seed point using
FastMarchingImageFilter. It is the output from this filter the one
that has the "default" spacing (and I think that this new spacing is
carried to the end to the output image) and I can't find any method
to set the correct spacing.
Thanks,
Lucas
On Tuesday, February 10, 2004, at 08:38 PM, Luis Ibanez wrote:
Hi Lucas,
Why are you setting the spacing on the ImageIO
object instead of the image itself ?
You should just carry the spacing through the
pipeline. Does your input image has an invalid
spacing ?
An option in that case is to use the
ChangeInformationImageFilter
http://www.itk.org/Insight/Doxygen/html/
classitk_1_1ChangeInformationImageFilter.html
This filter carries the data buffer of the
input image to the output image, and allows
you to alter the meta-data such as image
origin and spacing.
Please don't use this filter for processing
images of human beings or any other living
organisms, since chances are that you will
make somebody operate in a liver instead of
a lung.
In the long term the right thing to do is
to fix the source of your images which is
where the real spacing information should
be comming from.
Regards,
Luis
======================================
-------------------
Lucas Lorenzo wrote:
Hi all,
I'm trying to write an image to disk in vtk format.
By default the spacing is set to 1 1 1. I'd like to change it
so
I'm doing the following:
/#include "itkVTKImageIO.h"
int main( int argc, char * argv[] )
{
/* /* Typedefs of components. */*/
const unsigned int ImageDimension = 2;
typedef unsigned char BinaryPixelType;
typedef itk::Image<<BinaryPixelType,ImageDimension>
BinaryImageType;
/*// read the input image and get the spacing from it:*/
typedef itk::VTKImageIO ImageIOType;
ImageIOType::Pointer IO1 = ImageIOType::New();
original_reader->SetImageIO(IO1);
double dx,dy,dz;
original_reader->Update();
dx = IO1->GetSpacing(0);
dy = IO1->GetSpacing(1);
dz = IO1->GetSpacing(2);
/*// here I'm omitting when I process the input image
// write the image to a file but perviously set the spacing
*/
binary_writer->SetFileName("out.vtk");
binary_writer->SetInput(thresholder->GetOutput());
ImageIOType::Pointer IO2 = ImageIOType::New();
IO2->SetSpacing(0,dx);
IO2->SetSpacing(1,dy);
IO2->SetSpacing(2,dz);
binary_writer->SetImageIO(IO2);
try
{
binary_writer->Update();
}
catch( itk::ExceptionObject & exp )
{
std::cerr <<<< "Exception caught ! binary image writer" <<<<
std::endl;
std::cerr <<<< exp <<<< std::endl;
return -1;
}
return 0;
}
/
When I execute this program I have a "segmentation fault" run
time error.
If I ommit the line setting the spacing in the z axis ("/
IO2->SetSpacing(2,dz); /") the programs executes without
crashing but the spacing is not really set, that is, in the
header of the output vtk file (out.vtk) you can still see
"SPACING 1 1 1.0" instead of "SPACING dx dy 1.0".
Any clue of what am I doing wrong ?
Thanks,
Lucas Lorenzo
University of Utah
Nora Eccles Harrison CardioVascular Research and Training
Institute
Fellows Room
95 South 2000 East
Salt Lake City, UT 84112-5000
e-mail: lucas at cvrti.utah.edu
telephone: 801-587-9536
_______________________________________________
Insight-users mailing list
Insight-users at itk.org
http://www.itk.org/mailman/listinfo/insight-users
Lucas Lorenzo
University of Utah
Nora Eccles Harrison CardioVascular Research and Training Institute
Fellows Room
95 South 2000 East
Salt Lake City, UT 84112-5000
e-mail: lucas at cvrti.utah.edu
telephone: 801-587-9536
</excerpt>
</excerpt>Lucas Lorenzo
University of Utah
Nora Eccles Harrison CardioVascular Research and Training Institute
Fellows Room
95 South 2000 East
Salt Lake City, UT 84112-5000
e-mail: lucas at cvrti.utah.edu
telephone: 801-587-9536
_______________________________________________
Insight-users mailing list
Insight-users at itk.org
http://www.itk.org/mailman/listinfo/insight-users
</excerpt>
</excerpt>Lucas Lorenzo
University of Utah
Nora Eccles Harrison CardioVascular Research and Training Institute
Fellows Room
95 South 2000 East
Salt Lake City, UT 84112-5000
e-mail: lucas at cvrti.utah.edu
telephone: 801-587-9536
</excerpt>
_______________________________________________
Insight-users mailing list
Insight-users at itk.org
http://www.itk.org/mailman/listinfo/insight-users
</excerpt>Lucas Lorenzo
University of Utah
Nora Eccles Harrison CardioVascular Research and Training Institute
Fellows Room
95 South 2000 East
Salt Lake City, UT 84112-5000
e-mail: lucas at cvrti.utah.edu
telephone: 801-587-9536
--Apple-Mail-5--881340013--
--Apple-Mail-3--881340016--