[Insight-users] Seg. fault with ReleaseDataFlag and JoinImageFilter
Benjamin King
king.benjamin at mh-hannover.de
Fri, 30 Jan 2004 15:28:23 +0100
------------tZSJ85ntYvaiKtUqSG8IVH
Content-Type: text/plain; format=flowed; charset=iso-8859-15
Content-Transfer-Encoding: 8bit
Hello,
the attached code implements the following pipeline:
reader-->gradient magnitude-->JoinImageFilter
| ^
| |
|-----------------------------------|
When I use reader->ReleaseDataFlagOn(), the program produces a
segmentation fault on joiner->Update() _unless_ preceeded by a
gradient->Update(). Please see the comments in the code. I use a fairly
recent ITK checkout and Linux.
My questions are:
1) Why does this happen? Maybe during Update(), the gradient magnitude
filter makes the erroneous assumption that the reader's output is still
valid...
2) Is there a better way to reduce the memory footprint than using
ReleaseDataFlagOn()?
Thanks for the help,
Benjamin
#include <iostream>
#include <stdexcept>
#include "itkImage.h"
#include "itkImageFileReader.h"
#include "itkGradientMagnitudeImageFilter.h"
#include "itkJoinImageFilter.h"
using std::exception;
int main(int argc, char *argv[])
{
if (argc != 2)
{
std::cout <<
"Makes a two component image of scalar value vs. gradient
magnitude\n"
"\n"
"USAGE: " << argv[0] << " <volume>\n"
" <volume>: ITK volume file, e.g. Analyze. Example: foo.hdr\n"
<< std::endl;
return 1;
}
typedef itk::Image<unsigned short, 3> InputImageType;
typedef itk::Image<float, 3> DerivativeImageType;
typedef itk::ImageFileReader<InputImageType> ImageReaderType;
ImageReaderType::Pointer reader = ImageReaderType::New();
reader->SetFileName(argv[1]);
reader->ReleaseDataFlagOn(); // Either insert comment to prevent seg.
fault...
typedef itk::GradientMagnitudeImageFilter<InputImageType,
DerivativeImageType> GradientFilterType;
GradientFilterType::Pointer gradient = GradientFilterType::New();
gradient->ReleaseDataFlagOn();
gradient->SetInput(reader->GetOutput());
typedef itk::JoinImageFilter<InputImageType, DerivativeImageType>
JoinFilterType;
JoinFilterType::Pointer combined = JoinFilterType::New();
combined->SetInput1(reader->GetOutput());
combined->SetInput2(gradient->GetOutput());
try
{
// gradient->Update(); // ...or remove this comment to prevent seg.
fault
combined->Update();
}
catch (itk::ExceptionObject & e)
{
std::cerr << "ITK exception: " << e << std::endl;
}
catch (std::exception & e)
{
std::cerr << "std exception: " << e.what() << std::endl;
}
return 0;
}
--
Benjamin King
Institut für Medizinische Informatik
Medizinische Hochschule Hannover
Tel.: +49 511 532-2663
------------tZSJ85ntYvaiKtUqSG8IVH
Content-Disposition: attachment; filename=DebugSegFault.cpp
Content-Type: application/octet-stream; name=DebugSegFault.cpp
Content-Transfer-Encoding: 8bit
#include <stdexcept>
#include "itkImage.h"
#include "itkImageFileReader.h"
#include "itkGradientMagnitudeImageFilter.h"
#include "itkJoinImageFilter.h"
using std::exception;
int main(int argc, char *argv[])
{
if (argc != 2)
{
std::cout <<
"Makes a two component image of scalar value vs. gradient magnitude\n"
"\n"
"USAGE: " << argv[0] << " <volume>\n"
" <volume>: ITK volume file, e.g. Analyze. Example: foo.hdr\n"
<< std::endl;
return 1;
}
// Define image types
typedef itk::Image<unsigned short, 3> InputImageType;
typedef itk::Image<float, 3> DerivativeImageType;
// Read Volume
typedef itk::ImageFileReader<InputImageType> ImageReaderType;
ImageReaderType::Pointer reader = ImageReaderType::New();
reader->SetFileName(argv[1]);
reader->ReleaseDataFlagOn(); // Either insert comment to prevent seg. fault...
// Compute gradient magnitude
typedef itk::GradientMagnitudeImageFilter<InputImageType, DerivativeImageType> GradientFilterType;
GradientFilterType::Pointer gradient = GradientFilterType::New();
gradient->ReleaseDataFlagOn();
gradient->SetInput(reader->GetOutput());
// Join images
typedef itk::JoinImageFilter<InputImageType, DerivativeImageType> JoinFilterType;
JoinFilterType::Pointer combined = JoinFilterType::New();
combined->SetInput1(reader->GetOutput());
combined->SetInput2(gradient->GetOutput());
try
{
// gradient->Update(); // ...or remove this comment to prevent seg. fault
combined->Update();
}
catch (itk::ExceptionObject & e)
{
std::cerr << "ITK exception: " << e << std::endl;
}
catch (std::exception & e)
{
std::cerr << "std exception: " << e.what() << std::endl;
}
return 0;
}
------------tZSJ85ntYvaiKtUqSG8IVH--