[ITK] [ITK-users] How to make ITK-binaries NOT ignore kernel out of memory messages?
Matt McCormick
matt.mccormick at kitware.com
Mon Jul 14 20:55:01 EDT 2014
Hi Roman,
The CMake configuration variable:
ITK_COMPUTER_MEMORY_SIZE
is only used to choose which unit tests to run, and has no effect on
your program.
It may be worth a try to run
export ITK_GLOBAL_DEFAULT_NUMBER_OF_THREADS=1
before starting your program on the server to see the effect. It is
possible that too many threads or being spawned or that the amount of
memory consumed is scaling (poorly) with the number of threads used.
Hope this helps,
Matt
On Mon, Jun 23, 2014 at 6:36 AM, Dr. Roman Grothausmann
<grothausmann.roman at mh-hannover.de> wrote:
> Dear Luis,
>
>
> Many thanks for Your reply. It took me some time to reproduce the problem,
> because Your test program works as expected on our server:
>
> ~$ ~/itk/simple/build_itk-4.5.1/memtest_01 1000
> ITK Hello World !
>
> ~$ ~/itk/simple/build_itk-4.5.1/memtest_01 10000
> /opt/itk-4.5.1/include/ITK-4.5/itkImportImageContainer.hxx:192:
>
> Failed to allocate memory for image.
>
> One program of mine for which the problem exists is below. The problem
> occurs when BinaryMask3DMeshSource acquires more memory after running for
> some minutes. The very strange thing to me is, that the program is killed by
> the kernel as expected on my desktop PC:
>
> ~/itk/simple/build_itk4_pc3G48243/euler-characteristic_01 13-511_EAA_ot.mha
> 255
> component type is: unsigned_char
> numDimensions: 3
> component size: 1
> pixel type (string): scalar
> pixel type: 1
> _________________________
>
> Executing ImageFileReader done. Took 2.78 seconds.
> Executing BinaryMask3DMeshSource Killed
>
>
> The header of 13-511_EAA_ot.mha:
>
> ObjectType = Image
> NDims = 3
> BinaryData = True
> BinaryDataByteOrderMSB = False
> CompressedData = True
> CompressedDataSize = 20474538
> TransformMatrix = 1 0 0 0 1 0 0 0 1
> Offset = 0 0 0
> CenterOfRotation = 0 0 0
> AnatomicalOrientation = RAI
> ElementSpacing = 1 1 1
> DimSize = 998 998 1034
> ElementType = MET_UCHAR
> ElementDataFile = LOCAL
>
>
> However on the server it is not terminated when the RAM runs out, instead
> many processes go nuts that had been mostly sleeping and finally the system
> does not respond any more at all and a cold start is necessary. The
> itk-4.5.1 is configured on the server with (see attached
> CMakeCache_itk-4.5.1.txt):
>
> ITK_COMPUTER_MEMORY_SIZE 190
>
> The server has 256GB of RAM an no Swap enabled:
>
> free
> total used free shared buffers cached
> Mem: 264498496 1341020 263157476 0 285760 289120
> -/+ buffers/cache: 766140 263732356
> Swap: 0 0 0
>
>
> Has 12x2 Intel(R) Xeon(R) CPUs:
>
> cat /proc/cpuinfo
> processor : 0
> vendor_id : GenuineIntel
> cpu family : 6
> model : 45
> model name : Intel(R) Xeon(R) CPU E5-2630L 0 @ 2.00GHz
> stepping : 7
> microcode : 0x710
> cpu MHz : 1200.000
> cache size : 15360 KB
> physical id : 0
> siblings : 12
> core id : 0
> cpu cores : 6
> apicid : 0
> initial apicid : 0
> fpu : yes
> fpu_exception : yes
> cpuid level : 13
> wp : yes
> flags : fpu vme de pse tsc msr pae mce cx8 apic sep mtrr pge mca
> cmov pat pse36 clflush dts acpi mmx fxsr sse sse2 ss ht tm pbe syscall nx
> pdpe1gb rdtscp lm constant_tsc arch_perfmon pebs bts rep_good nopl xtopology
> nonstop_tsc aperfmperf pni pclmulqdq dtes64 monitor ds_cpl vmx smx est tm2
> ssse3 cx16 xtpr pdcm pcid dca sse4_1 sse4_2 x2apic popcnt tsc_deadline_timer
> aes xsave avx lahf_lm ida arat epb xsaveopt pln pts dtherm tpr_shadow vnmi
> flexpriority ept vpid
> bogomips : 4000.56
> clflush size : 64
> cache_alignment : 64
> address sizes : 46 bits physical, 48 bits virtual
> power management:
> .
> .
> .
>
>
> The server runs a Debian Wheezy (7.2):
> uname -a
> Linux img1-serv 3.2.0-4-amd64 #1 SMP Debian 3.2.57-3+deb7u2 x86_64 GNU/Linux
>
> My desktop PC has the same Debian (also same update state):
> uname -a
> Linux pc3G48243 3.2.0-4-amd64 #1 SMP Debian 3.2.57-3+deb7u2 x86_64 GNU/Linux
>
> free
> total used free shared buffers cached
> Mem: 8180356 1185364 6994992 0 23348 216084
> -/+ buffers/cache: 945932 7234424
> Swap: 0 0 0
>
>
> and 4 AMD FX(tm)-4100 Quad-Core Processors:
>
> cat /proc/cpuinfo
> processor : 0
> vendor_id : AuthenticAMD
> cpu family : 21
> model : 1
> model name : AMD FX(tm)-4100 Quad-Core Processor
> stepping : 2
> microcode : 0x6000626
> cpu MHz : 1400.000
> cache size : 2048 KB
> physical id : 0
> siblings : 4
> core id : 0
> cpu cores : 2
> apicid : 16
> initial apicid : 0
> fpu : yes
> fpu_exception : yes
> cpuid level : 13
> wp : yes
> flags : fpu vme de pse tsc msr pae mce cx8 apic sep mtrr pge mca
> cmov pat pse36 clflush mmx fxsr sse sse2 ht syscall nx mmxext fxsr_opt
> pdpe1gb rdtscp lm constant_tsc rep_good nopl nonstop_tsc extd_apicid
> aperfmperf pni pclmulqdq monitor ssse3 cx16 sse4_1 sse4_2 popcnt aes xsave
> avx lahf_lm cmp_legacy svm extapic cr8_legacy abm sse4a misalignsse
> 3dnowprefetch osvw ibs xop skinit wdt lwp fma4 nodeid_msr topoext
> perfctr_core arat cpb hw_pstate npt lbrv svm_lock nrip_save tsc_scale
> vmcb_clean flushbyasid decodeassists pausefilter pfthreshold
> bogomips : 7233.41
> TLB size : 1536 4K pages
> clflush size : 64
> cache_alignment : 64
> address sizes : 48 bits physical, 48 bits virtual
> power management: ts ttp tm 100mhzsteps hwpstate [9]
>
>
> Knowing now, that itk checks for available RAM, the only sources for this
> problem I can think of would be an odd configuration of ITK itself or its
> use in programs compiled on the server, or some configuration differences in
> the (Debian) system, or some hardware problem.
>
> Do You have any idea what it could be?
> Or do You have any suggestions how I should go on searching for the cause?
>
>
> Any help or hints are very much appriciated.
> Many thanks for Your help again!
>
> Roman
>
>
> _____________________________________________
>
>
>
> ////program to calculate the Euler characteristic of a binary image
> converted to itkQuadEdgeMesh
> ////THEREFORE does NOT take filled regions into account only covering
> surfaces!
> ////if fillings should be accounted use instead a vtk program taking 3D
> cells into account
>
>
>
> #include <itkImageFileReader.h>
> //#include <itkImage.h>
> #include <itkBinaryMask3DMeshSource.h> //itk marching cubes (mc) very slow
> and resource hungry :-(
> #include <itkQuadEdgeMeshTopologyChecker.h>
> #include <itkQuadEdgeMeshBoundaryEdgesMeshFunction.h>
>
>
>
> #include "itkFilterWatcher2.h"
>
>
>
>
> void dispatch_cT(itk::ImageIOBase::IOPixelType,
> itk::ImageIOBase::IOComponentType, size_t, int, char **);
>
> template<typename InputComponentType>
> void dispatch_pT(itk::ImageIOBase::IOPixelType pixelType, size_t, int, char
> **);
>
> template<typename InputComponentType, typename InputPixelType>
> void dispatch_D(size_t, int, char **);
>
> template<typename InputComponentType, typename InputPixelType, size_t
> Dimension>
> int DoIt(int, char *argv[]);
>
>
>
>
>
> template<typename InputComponentType, typename InputPixelType, size_t
> Dimension>
> int DoIt(int argc, char *argv[]){
>
> typedef InputPixelType OutputPixelType;
>
> typedef itk::Image<InputPixelType, Dimension> InputImageType;
> typedef itk::Image<OutputPixelType, Dimension> OutputImageType;
>
>
> typedef itk::ImageFileReader<InputImageType> ReaderType;
> typename ReaderType::Pointer reader = ReaderType::New();
>
> reader->SetFileName(argv[1]);
> FilterWatcher watcherI(reader, "reading");
> watcherI.QuietOn();
> watcherI.ReportTimeOn();
> try
> {
> reader->Update();
> }
> catch (itk::ExceptionObject &ex)
> {
> if (!strcmp(ex.GetDescription(), "Filter does not have progress.")){
> std::cerr << ex.GetDescription() << std::endl;
> return EXIT_FAILURE;
> }
> }
>
>
> typedef double CoordType;
>
> //typedef itk::Mesh<CoordType> MeshType;
> typedef itk::QuadEdgeMesh<CoordType, Dimension> InputMeshType;
>
> const InputPixelType isoValue= static_cast<InputPixelType>(atof(argv[2]));
>
> typedef itk::BinaryMask3DMeshSource<InputImageType, InputMeshType>
> MeshSourceType;
> typename MeshSourceType::Pointer meshSource = MeshSourceType::New();
>
> meshSource->SetInput(reader->GetOutput());
> meshSource->SetObjectValue(isoValue);
>
> FilterWatcher watcher(meshSource, "BinaryMask3DMeshSource");
> watcher.QuietOn();
> watcher.ReportTimeOn();
> try
> {
> meshSource->Update();
> }
> catch (itk::ExceptionObject &ex)
> {
> if (!strcmp(ex.GetDescription(), "Filter does not have progress.")){
> std::cerr << ex.GetDescription() << std::endl;
> return EXIT_FAILURE;
> }
> }
>
>
> std::cout << "Nodes = " << meshSource->GetNumberOfNodes() << std::endl;
> std::cout << "Cells = " << meshSource->GetNumberOfCells() << std::endl;
>
> typename InputMeshType::Pointer mesh= meshSource->GetOutput();
>
> //// check if meshSource can be represented as an itkQuadEdgeMesh
> //// if BinaryMask3DMeshSource caps boundary surfaces this check is not
> necessary!
>
> typedef itk::QuadEdgeMeshTopologyChecker<InputMeshType> TQEchecker;
> typename TQEchecker::Pointer QEchecker = TQEchecker::New();
> QEchecker->SetMesh(mesh);
> std::cerr << "Executing
> QuadEdgeMeshTopologyChecker->ValidateEulerCharacteristic()";
> if(!QEchecker->ValidateEulerCharacteristic()){
> std::cerr << "Not a valid QuadEdgeMesh! Aborting." << std::endl;
> return EXIT_FAILURE;
> }
> std::cerr << " done." << std::endl;
>
> //// calc Euler characteristics \xi
> //// For closed smooth manifolds, the Euler characteristic coincides with
> the Euler number:
> http://en.wikipedia.org/wiki/Euler_characteristic#Relations_to_other_invariants
> //// general formula see:
> http://en.wikipedia.org/wiki/Genus_%28mathematics%29#Orientable_surface
> ////
> ////
> ////
> //// good QuadEdge explanation:
> http://www.cs.cmu.edu/afs/andrew/scs/cs/15-463/2001/pub/src/a2/quadedge.html
> /// from itkQuadEdgeMeshTopologyChecker.hxx
>
> typedef itk::QuadEdgeMeshBoundaryEdgesMeshFunction<InputMeshType>
> BoundaryEdges;
>
> typename BoundaryEdges::Pointer boundaryEdges = BoundaryEdges::New();
>
> // Number of USED points
> typedef typename InputMeshType::PointIdentifier
> PointIdentifier;
> std::cerr << "ComputeNumberOfPoints";
> PointIdentifier numPoints = mesh->ComputeNumberOfPoints();
> std::cerr << " done." << std::endl;
> // Number of USED edges
> typedef typename InputMeshType::CellIdentifier
> CellIdentifier;
> std::cerr << "ComputeNumberOfEdges";
> CellIdentifier numEdges = mesh->ComputeNumberOfEdges();
> std::cerr << " done." << std::endl;
> // Number of USED faces
> std::cerr << "ComputeNumberOfFaces";
> CellIdentifier numFaces = mesh->ComputeNumberOfFaces();
> std::cerr << " done." << std::endl;
> // Number of Boundaries
> std::cerr << "ComputeNumberOfBoundaryEdges";
> typename BoundaryEdges::OutputType
> listOfBoundaries = boundaryEdges->Evaluate((*mesh));
> std::cerr << " done." << std::endl;
> CellIdentifier numBounds = listOfBoundaries->size();
> delete listOfBoundaries;
>
> /**
> * Number of points
> *
> * There are two methods to get the number of points.
> * 1. itk::QuadEdgeMesh::ComputeNumberOfPoints()
> * 2. itk::Mesh::GetNumberOfPoints()
> *
> * As an itk::QuadEdgeMesh is an itk::Mesh by inheritance, the user
> * can use both. 1. will returned the number of points actually
> * used by at least one edge, while 2. will give you the number
> * of points in the container. Number of unused points can be found
> * by making the difference between the two values.
> */
> if(mesh->GetNumberOfPoints() != numPoints){
> std::cerr << "There are isolated vertices! Aborting." << std::endl;
> return EXIT_FAILURE;
> }
>
> // The euler formula states:
> // numFaces - numEdges + numPoints == 2 - 2 * genus - numBounds == \xi
> // hence 2 * genus= 2 - numBounds - numFaces + numEdges - numPoints must
> be an oddeven number.
>
> typedef ::itk::OffsetValueType OffsetValueType;
>
> OffsetValueType xi=
> + OffsetValueType(numFaces)
> - OffsetValueType(numEdges)
> + OffsetValueType(numPoints);
>
> std::cout << "The Euler characteristic Xi is: " << xi << std::endl;
>
> OffsetValueType twiceGenus=
> + OffsetValueType(2) - OffsetValueType(numBounds)
> - xi;
>
> std::cout << "The genus is: " << twiceGenus/OffsetValueType(2) <<
> std::endl;
>
> return EXIT_SUCCESS;
>
> }
>
>
> void dispatch_cT(itk::ImageIOBase::IOComponentType componentType,
> itk::ImageIOBase::IOPixelType pixelType, size_t dimensionType, int argc,
> char *argv[]){
>
>
> //http://www.itk.org/Doxygen45/html/classitk_1_1ImageIOBase.html#a8dc783055a0af6f0a5a26cb080feb178
> //http://www.itk.org/Doxygen45/html/itkImageIOBase_8h_source.html#l00107
> //IOComponentType: UNKNOWNCOMPONENTTYPE, UCHAR, CHAR, USHORT, SHORT, UINT,
> INT, ULONG, LONG, FLOAT, DOUBLE
>
> switch (componentType){
> case itk::ImageIOBase::UCHAR:{
> typedef unsigned char InputComponentType;
> dispatch_pT<InputComponentType>(pixelType, dimensionType, argc, argv);
> }
> break;
> case itk::ImageIOBase::CHAR:{
> typedef char InputComponentType;
> dispatch_pT<InputComponentType>(pixelType, dimensionType, argc, argv);
> }
> break;
> case itk::ImageIOBase::USHORT:{
> typedef unsigned short InputComponentType;
> dispatch_pT<InputComponentType>(pixelType, dimensionType, argc, argv);
> }
> break;
> case itk::ImageIOBase::SHORT:{
> typedef short InputComponentType;
> dispatch_pT<InputComponentType>(pixelType, dimensionType, argc, argv);
> }
> break;
> case itk::ImageIOBase::UINT:{
> typedef unsigned int InputComponentType;
> dispatch_pT<InputComponentType>(pixelType, dimensionType, argc, argv);
> }
> break;
> case itk::ImageIOBase::INT:{
> typedef int InputComponentType;
> dispatch_pT<InputComponentType>(pixelType, dimensionType, argc, argv);
> }
> break;
> // case itk::ImageIOBase::ULONG:{
> // typedef unsigned long InputComponentType;
> // dispatch_pT<InputComponentType>(pixelType, dimensionType, argc,
> argv);
> // }
> // break;
> // case itk::ImageIOBase::LONG:{
> // typedef long InputComponentType;
> // dispatch_pT<InputComponentType>(pixelType, dimensionType, argc,
> argv);
> // }
> // break;
> // case itk::ImageIOBase::FLOAT:{
> // typedef float InputComponentType;
> // dispatch_pT<InputComponentType>(pixelType, dimensionType, argc,
> argv);
> // }
> // break;
> // case itk::ImageIOBase::DOUBLE:{
> // typedef double InputComponentType;
> // dispatch_pT<InputComponentType>(pixelType, dimensionType, argc,
> argv);
> // }
> // break;
> case itk::ImageIOBase::UNKNOWNCOMPONENTTYPE:
> default:
> std::cout << "unknown component type" << std::endl;
> break;
> }//switch
> }
>
> template<typename InputComponentType>
> void dispatch_pT(itk::ImageIOBase::IOPixelType pixelType, size_t
> dimensionType, int argc, char *argv[]){
>
>
> //http://www.itk.org/Doxygen45/html/classitk_1_1ImageIOBase.html#abd189f096c2a1b3ea559bc3e4849f658
> //http://www.itk.org/Doxygen45/html/itkImageIOBase_8h_source.html#l00099
> //IOPixelType:: UNKNOWNPIXELTYPE, SCALAR, RGB, RGBA, OFFSET, VECTOR,
> POINT, COVARIANTVECTOR, SYMMETRICSECONDRANKTENSOR, DIFFUSIONTENSOR3D,
> COMPLEX, FIXEDARRAY, MATRIX
>
> switch (pixelType){
> case itk::ImageIOBase::SCALAR:{
> typedef InputComponentType InputPixelType;
> //typedef itk::SCALARPixel<InputComponentType> InputPixelType;
> dispatch_D<InputComponentType, InputPixelType>(dimensionType, argc,
> argv);
> }
> break;
> // case itk::ImageIOBase::RGB:{
> // typedef itk::RGBPixel<InputComponentType> InputPixelType;
> // dispatch_D<InputComponentType, InputPixelType>(dimensionType, argc,
> argv);
> // }
> // break;
> // case itk::ImageIOBase::RGBA:{
> // typedef itk::RGBAPixel<InputComponentType> InputPixelType;
> // dispatch_D<InputComponentType, InputPixelType>(dimensionType, argc,
> argv);
> // }
> // break;
> // case itk::ImageIOBase::VECTOR:{
> // //typedef itk::VECTORPixel<InputComponentType> InputPixelType; //does
> not work!
> // typedef itk::VariableLengthVector<InputComponentType> InputPixelType;
> // dispatch_D<InputComponentType, InputPixelType>(dimensionType, argc,
> argv);
> // }
> // break;
> case itk::ImageIOBase::UNKNOWNPIXELTYPE:
> default:
> std::cout << "unknown pixel type" << std::endl;
> break;
> }//switch
> }
>
>
> template<typename InputComponentType, typename InputPixelType>
> void dispatch_D(size_t dimensionType, int argc, char *argv[]){
> switch (dimensionType){
> // case 1:
> // DoIt<InputComponentType, InputPixelType, 1>(argc, argv);
> // break;
> // case 2:
> // DoIt<InputComponentType, InputPixelType, 2>(argc, argv);
> // break;
> case 3:
> DoIt<InputComponentType, InputPixelType, 3>(argc, argv);
> break;
> case 4:
> DoIt<InputComponentType, InputPixelType, 3>(argc, argv);
> break;
> case 5:
> DoIt<InputComponentType, InputPixelType, 3>(argc, argv);
> break;
> default:
> std::cout << "Images of dimension " << dimensionType << " are not
> supported!" << std::endl;
> break;
> }
> }
>
>
>
> ////from
> http://itk-users.7.n7.nabble.com/Pad-image-with-0-but-keep-its-type-what-ever-it-is-td27442.html
> //namespace itk{
> // Description:
> // Get the PixelType and ComponentType from fileName
>
> void GetImageType (std::string fileName,
> itk::ImageIOBase::IOPixelType &pixelType,
> itk::ImageIOBase::IOComponentType &componentType,
> size_t &dimensionType
> //ImageIOBase::IODimensionType &dimensionType
> ){
> typedef itk::Image<unsigned char, 3> ImageType;
> itk::ImageFileReader<ImageType>::Pointer imageReader=
> itk::ImageFileReader<ImageType>::New();
> imageReader->SetFileName(fileName.c_str());
> imageReader->UpdateOutputInformation();
>
> pixelType = imageReader->GetImageIO()->GetPixelType();
> componentType = imageReader->GetImageIO()->GetComponentType();
> dimensionType= imageReader->GetImageIO()->GetNumberOfDimensions();
>
> //std::cerr << "Pixel Type is " <<
> imageReader->GetImageIO()->GetComponentTypeAsString(pixelType) << std::endl;
> std::cerr << "component type is: " <<
> imageReader->GetImageIO()->GetComponentTypeAsString(componentType) <<
> std::endl;
> std::cerr << "numDimensions: " << dimensionType << std::endl;
> std::cerr << "component size: " <<
> imageReader->GetImageIO()->GetComponentSize() << std::endl;
> std::cerr << "pixel type (string): " <<
> imageReader->GetImageIO()->GetPixelTypeAsString(imageReader->GetImageIO()->GetPixelType())
> << std::endl;
> std::cerr << "pixel type: " << pixelType << std::endl <<
> "_________________________" << std::endl << std::endl;
>
>
> }
>
>
>
> int main(int argc, char *argv[]){
> if ( argc != 3 )
> {
> std::cerr << "Missing Parameters: "
> << argv[0]
> << " Input_Image"
> << " isoValue"
> << std::endl;
>
> return EXIT_FAILURE;
> }
>
> std::string ifn = argv[1];
> // std::string ofn = argv[2];
> // int compress= atoi(argv[3]);
>
> itk::ImageIOBase::IOPixelType pixelType;
> typename itk::ImageIOBase::IOComponentType componentType;
> //itk::ImageIOBase::IOComponentType componentType1;
> //itk::ImageIOBase::IODimensionType dimensionType1;
> size_t dimensionType;
>
>
> try
> {
> GetImageType(argv[1], pixelType, componentType, dimensionType);
>
> dispatch_cT(componentType, pixelType, dimensionType, argc, argv);
> }//try
> catch( itk::ExceptionObject &excep)
> {
> std::cerr << argv[0] << ": exception caught !" << std::endl;
> std::cerr << excep << std::endl;
> return EXIT_FAILURE;
> }
>
> return EXIT_SUCCESS;
>
> }
>
>
>
>
>
>
>
>
>
>
> On 08/06/14 15:51, Luis Ibanez wrote:
>>
>>
>> Roman,
>>
>>
>> Could you please share with the list an example
>> of a program where you are observing this behavior ?
>>
>> In particular, we would like to see if this is occurring with
>> Image filters, or rather with other types of ITK objects.
>>
>>
>>
>> The itk::Image does listen to the (potential) memory allocation messages
>> (that is, the bad_alloc exception).
>>
>>
>>
>> The itk::Image class uses a helper class to store the pixel data.
>>
>>
>>
>> This helper is the:
>>
>> itk::ImportImageContainer
>>
>>
>> https://github.com/InsightSoftwareConsortium/ITK/blob/master/Modules/Core/Common/include/itkImportImageContainer.h
>>
>>
>>
>>
>> and its memory allocation is typically done in the AllocateElements()
>> function
>> in line:
>>
>>
>> https://github.com/InsightSoftwareConsortium/ITK/blob/master/Modules/Core/Common/include/itkImportImageContainer.hxx#L173
>>
>> the code looks like
>>
>> TElement *data;
>>
>> try
>> {
>> if ( UseDefaultConstructor )
>> {
>> data = new TElement[size](); //POD types initialized to 0, others
>> use
>> default constructor.
>> }
>> else
>> {
>> data = new TElement[size]; //Faster but uninitialized
>> }
>> }
>> catch ( ... )
>> {
>> data = ITK_NULLPTR;
>> }
>> if ( !data )
>> {
>> // We cannot construct an error string here because we may be out
>> // of memory. Do not use the exception macro.
>> throw MemoryAllocationError(__FILE__, __LINE__,
>> "Failed to allocate memory for image.",
>> ITK_LOCATION);
>> }
>> return data;
>>
>>
>>
>> As you can see, the case of exceptions is managed in that function.
>>
>>
>> Here is an example of a test program demonstrating that ITK image
>> will throw exceptions when allocating memory beyond limits.
>>
>>
>> #include "itkImage.h"
>> #include <iostream>
>>
>> int main(int argc, char * argv[] )
>> {
>> if( argc < 2 )
>> {
>> std::cerr << "Missing arguments" << std::endl;
>> return 1;
>> }
>>
>> typedef itk::Image< unsigned short, 3 > ImageType;
>>
>> ImageType::Pointer image = ImageType::New();
>>
>> ImageType::RegionType region;
>> ImageType::SizeType size;
>>
>> size_t side = atoi( argv[1] );
>>
>> size[0] = side;
>> size[1] = side;
>> size[2] = side;
>>
>> region.SetSize( size );
>>
>> image->SetRegions( region );
>>
>> try
>> {
>> image->Allocate();
>> }
>> catch( std::exception & excp )
>> {
>> std::cerr << excp.what() << std::endl;
>> return 1;
>> }
>>
>> std::cout << "ITK Hello World !" << std::endl;
>>
>> return 0;
>> }
>>
>>
>> and what happens when executing in an Ubuntu Linux machine with 32GB or
>> RAM
>>
>>
>> $ ./HelloWorld 1000
>> ITK Hello World !
>>
>>
>> $ ./HelloWorld 1000000
>>
>> /home/ibanez/src/ITK/Modules/Core/Common/include/itkImportImageContainer.hxx:199:
>> Failed to allocate memory for image.
>>
>>
>>
>>
>> Please share with the list a minimal example of the problem you are
>> observing,
>> and in this way we will be able to help track the issue.
>>
>>
>> Thanks
>>
>>
>> Luis
>>
>>
>
> --
> Dr. Roman Grothausmann
>
> Tomographie und Digitale Bildverarbeitung
> Tomography and Digital Image Analysis
>
> Institut für Funktionelle und Angewandte Anatomie, OE 4120
> Medizinische Hochschule Hannover
> Carl-Neuberg-Str. 1
> D-30625 Hannover
>
> Tel. +49 511 532-9574
>
> _____________________________________
> Powered by www.kitware.com
>
> Visit other Kitware open-source projects at
> http://www.kitware.com/opensource/opensource.html
>
> Kitware offers ITK Training Courses, for more information visit:
> http://www.kitware.com/products/protraining.php
>
> Please keep messages on-topic and check the ITK FAQ at:
> http://www.itk.org/Wiki/ITK_FAQ
>
> Follow this link to subscribe/unsubscribe:
> http://public.kitware.com/mailman/listinfo/insight-users
>
_____________________________________
Powered by www.kitware.com
Visit other Kitware open-source projects at
http://www.kitware.com/opensource/opensource.html
Kitware offers ITK Training Courses, for more information visit:
http://www.kitware.com/products/protraining.php
Please keep messages on-topic and check the ITK FAQ at:
http://www.itk.org/Wiki/ITK_FAQ
Follow this link to subscribe/unsubscribe:
http://public.kitware.com/mailman/listinfo/insight-users
More information about the Community
mailing list