View Issue Details [ Jump to Notes ] | [ Print ] | ||||||||
ID | Project | Category | View Status | Date Submitted | Last Update | ||||
0009139 | ITK | public | 2009-06-11 12:08 | 2009-12-16 08:25 | |||||
Reporter | Hans Johnson | ||||||||
Assigned To | Hans Johnson | ||||||||
Priority | high | Severity | crash | Reproducibility | always | ||||
Status | closed | Resolution | fixed | ||||||
Platform | OS | OS Version | |||||||
Product Version | ITK-3-14 | ||||||||
Target Version | Fixed in Version | ||||||||
Summary | 0009139: itkOptImageToImageMetric segmentation faults | ||||||||
Description | Much of the code in ::MultiThreadingInitialize has a dependency on the variable m_NumberOfFixedImageSamples for determining the length of the weights arrays. This is OK if that variable were only set from the class API before computing the registration. Unfortunately, after initialization it is possible that the FixedImageRegion, or the number of voxels available in the FixedImage mask can change this value after initialization to a smaller number. If that is the case, then the size of the weights arrays also need to be re-adjusted in size. Ultimately the segmentation fault arises in a for loop over several arrays that should all be the same length (i.e. The new smaller value of m_NumberOfFixedImageSamples), but the loop index is specified as one of the arrays that has not been adjusted to the smaller size. | ||||||||
Additional Information | This problem has been difficult to reproduce in a small artificial test. In the real world application it occurs when using FixedImage Mask to restrict the area to be smaller than the number of requested samples (i.e. request 1 million sampling points from a FixedImageMask with only 900K voxels). This would induce a "race" condition where cached bspline weights were initialized to be of length 1M, and then the arrays of samples were reduced to 900K without re-initializing the weights cache. The solution was to make sure that re-initialization is done EVERY TIME the number of samples is changed. | ||||||||
Tags | No tags attached. | ||||||||
Resolution Date | |||||||||
Sprint | |||||||||
Sprint Status | |||||||||
Attached Files | itkOptImageToImageMetric.h [^] (26,501 bytes) 2009-06-11 12:14 itkOptImageToImageMetric.txx [^] (45,204 bytes) 2009-06-11 12:14 | ||||||||
Relationships | |
Relationships |
Notes | |
(0016708) Hans Johnson (developer) 2009-06-11 12:11 |
BRAINSFit run that shows was a first clue to what the problem was. [hjohnson@pandora ~]$ valgrind /IPLlinux/ipldev/scratch/hjohnson/src/BRAINS-COMPILE/Linux/DEBUG_64-build/BRAINS3-build/src/bin/BRAINSFit --fixedVolume /nca/structural/MR6/9503953/35249099/10_AUTO.v020Global/35249099_class_csf_clipped.nii.gz --movingVolume /ipldev/scratch/hjohnson/src/BRAINS3/src/iplProg/itkEMS/BRAINS3Atlas/template_class_csf_clipped.nii.gz --initializeTransformMode CenterOfHead --maskProcessingMode ROIAUTO --outputVolume /hjohnson/HDNI/PREDICT_TRAINING/regina_ann/TrainingModels/2009_Apr_Training/Deformations/Final_Refined_Registration/AtlasToSub/35249099.mat_BSpline_output.nii.gz --outputVolumePixelType short --outputTransform /hjohnson/HDNI/PREDICT_TRAINING/regina_ann/TrainingModels/2009_Apr_Training/Deformations/Final_Refined_Registration/AtlasToSub/35249099.mat --transformType Rigid,ScaleVersor3D,BSpline --numberOfIterations 1500 --spatialScale 1000.0 --reproportionScale 1.0 --minimumStepSize 0.005 --numberOfSamples 1000000 --skewScale 1.0 --splineGridSize 14,10,12 --numberOfHistogramBins 50 --numberOfMatchPoints 10 --histogramMatch |tee logger 2>&1 ==27809== Memcheck, a memory error detector. ==27809== Copyright (C) 2002-2005, and GNU GPL'd, by Julian Seward et al. ==27809== Using LibVEX rev 1575, a library for dynamic binary translation. ==27809== Copyright (C) 2004-2005, and GNU GPL'd, by OpenWorks LLP. ==27809== Using valgrind-3.1.1, a dynamic binary instrumentation framework. ==27809== Copyright (C) 2000-2005, and GNU GPL'd, by Julian Seward et al. ==27809== For more details, rerun with: -v ==27809== Original Fixed image origin[0, 255, 0, 0] Original moving volume Direction matrix: 1 0 0 0 0 1 0 0 0 0 1 0 0 0 0 1 extracted moving volume Direction matrix: 1 0 0 0 1 0 0 0 1 Original fixed volume Direction matrix: 1 0 0 0 0 -0 -1 0 0 1 0 0 0 0 0 1 extracted fixed volume Direction matrix: 1 0 0 0 -0 -1 0 1 0 Original moving volume Origin: [-127.5, -127.5, -127.5, 0] extracted moving volume Origin: [-127.5, -127.5, -127.5] Original fixed volume Origin: [0, 255, 0, 0] extracted fixed volume Origin: [0, 255, 0] Computing otsu thresholds with percentile: 0.01. LowHigh Thresholds: [0.0280731,89.25,237.74] -------- Start ConnectedComponentImageFilter "" Progress Quiet Filter took 61.7658 seconds. -------- End ConnectedComponentImageFilter "" Removed 282 smaller objects. ^[[C Computing otsu thresholds with percentile: 0.01. LowHigh Thresholds: [0.0283293,89.25,240.209] -------- Start ConnectedComponentImageFilter "" Progress Quiet Filter took 61.0111 seconds. -------- End ConnectedComponentImageFilter "" Removed 36 smaller objects. =============================== Starting Transform Estimations for Rigid=============================== Computing otsu thresholds with percentile: 0.01. LowHigh Thresholds: [0,89,237] -------- Start ConnectedComponentImageFilter "" Progress Quiet Filter took 61.2844 seconds. -------- End ConnectedComponentImageFilter "" Removed 34 smaller objects. extremum = 127.5 voxelSize = 0.001 cubic CC MIN MAX 37 185 Histo values 37 ... 185 -> 148 headSizeLimit = 1000 CCs zero bin count to be skipped = 72 ForegroundLevel = 185 topOfHeadDistFromExtremeSI 45 Y_Location_from_Top_of_Head: = 8.05938 Computing otsu thresholds with percentile: 0.01. LowHigh Thresholds: [0,89,236] -------- Start ConnectedComponentImageFilter "" Progress Quiet Filter took 61.2616 seconds. -------- End ConnectedComponentImageFilter "" Removed 282 smaller objects. extremum = 255 voxelSize = 0.001 cubic CC MIN MAX 49 186 Histo values 49 ... 186 -> 137 headSizeLimit = 1000 CCs zero bin count to be skipped = 320 ForegroundLevel = 186 topOfHeadDistFromExtremeSI 53 Y_Location_from_Top_of_Head: = 129.9 Initializing transform with CenterOfHead to VersorRigid3DTransform (0x8107fd8) RTTI typeinfo: itk::VersorRigid3DTransform<double> Reference Count: 2 Modified Time: 2580 Debug: Off Observers: none Matrix: 1 0 0 0 1 0 0 0 1 Offset: [-128.845, -125.117, -121.841] Center: [128.005, 148.272, 129.9] Translation: [-128.845, -125.117, -121.841] Inverse: 1 0 0 0 1 0 0 0 1 Singular: 0 Versor: [ 0, 0, 0, 1 ] =============================================== Initializer, optimizerScales: [1, 1, 1, 0.001, 0.001, 0.001]. 0 -0.128821 [0.00315007, -0.000612427, -0.00804151, -128.759, -125.213, -121.993] 1 -0.137211 [0.00684012, -0.00233695, -0.0132159, -128.738, -125.232, -122.19] 2 -0.144138 [0.00854006, -0.00337992, -0.0150535, -128.698, -125.229, -122.386] 3 -0.147919 [0.010394, -0.00390493, -0.0164256, -128.655, -125.233, -122.581] 4 -0.151599 [0.0123935, -0.00454459, -0.0173688, -128.614, -125.241, -122.777] 5 -0.155332 [0.0142511, -0.00445197, -0.0179263, -128.563, -125.26, -122.969] 6 -0.159096 [0.0162956, -0.00438678, -0.0185906, -128.512, -125.282, -123.161] 7 -0.162657 [0.0186019, -0.0047017, -0.0190266, -128.469, -125.316, -123.354] 8 -0.166057 [0.0213642, -0.00512257, -0.0191276, -128.429, -125.352, -123.546] 9 -0.169642 [0.0237647, -0.00569085, -0.0194515, -128.396, -125.392, -123.739] 10 -0.172895 [0.0263812, -0.00625109, -0.0194389, -128.366, -125.453, -123.927] 11 -0.175823 [0.0286768, -0.00639747, -0.0193174, -128.328, -125.525, -124.11] 12 -0.178194 [0.0296593, -0.00707959, -0.0189574, -128.292, -125.615, -124.285] 13 -0.180085 [0.0303526, -0.00753382, -0.0195449, -128.254, -125.713, -124.455] 14 -0.181912 [0.0303147, -0.00719864, -0.0190568, -128.199, -125.823, -124.613] 15 -0.183418 [0.0305764, -0.00763744, -0.0200438, -128.168, -125.938, -124.774] 16 -0.184737 [0.0298366, -0.007493, -0.0191762, -128.117, -126.056, -124.927] 17 -0.186005 [0.0288667, -0.00921316, -0.0208045, -128.101, -126.194, -125.07] 18 -0.186791 [0.0282575, -0.00758654, -0.0191426, -128.006, -126.315, -125.199] 19 -0.187676 [0.0283775, -0.00929761, -0.0219562, -128.043, -126.428, -125.359] 20 -0.188122 [0.0276027, -0.0079668, -0.0180641, -127.953, -126.586, -125.443] 21 -0.188671 [0.028641, -0.0127787, -0.0240802, -128.065, -126.711, -125.55] 22 -0.187407 [0.0271454, -0.00838592, -0.0174496, -127.972, -126.74, -125.567] 23 -0.188977 [0.0272845, -0.00997659, -0.0196024, -128.007, -126.765, -125.591] 24 -0.189344 [0.02723, -0.00953713, -0.0195371, -127.995, -126.801, -125.623] 25 -0.189366 [0.0272346, -0.00990227, -0.0199078, -127.993, -126.832, -125.662] 26 -0.189335 [0.0271179, -0.00916726, -0.0199218, -127.984, -126.861, -125.702] 27 -0.189402 [0.0267066, -0.0103017, -0.0202926, -128.002, -126.888, -125.741] 28 -0.189398 [0.0277576, -0.00833556, -0.0191297, -127.97, -126.922, -125.758] 29 -0.189286 [0.0274334, -0.00978334, -0.0201393, -127.991, -126.933, -125.765] 30 -0.18945 [0.0277115, -0.00919098, -0.0196752, -127.986, -126.956, -125.771] 31 -0.189496 [0.0276979, -0.00966595, -0.0200138, -127.997, -126.977, -125.78] 32 -0.189476 [0.0277195, -0.00913262, -0.0197375, -127.988, -127, -125.784] 33 -0.189502 [0.0270427, -0.0098896, -0.0204035, -128.005, -127.016, -125.792] 34 -0.189453 [0.0276844, -0.00886729, -0.0196072, -127.987, -127.033, -125.799] 35 -0.189522 [0.027356, -0.00958007, -0.0202772, -128, -127.033, -125.8] 36 -0.189498 [0.0274112, -0.00932893, -0.0201334, -127.995, -127.037, -125.801] 37 -0.189546 [0.0273703, -0.00913528, -0.0200941, -127.992, -127.042, -125.801] 38 -0.18956 [0.0272891, -0.0094375, -0.0204258, -127.995, -127.047, -125.803] 39 -0.189524 [0.0273275, -0.00922639, -0.0201913, -127.989, -127.05, -125.805] 40 -0.189594 [0.0272435, -0.00919906, -0.0204024, -127.99, -127.056, -125.804] 41 -0.189531 [0.0273188, -0.00924714, -0.0201718, -127.986, -127.058, -125.808] Matrix = 0.999015 0.0398133 -0.019585 -0.0408238 0.997694 -0.0542307 0.0173807 0.0549769 0.998336 Offset = [-131.219, -114.446, -135.968] =============================== Starting Transform Estimations for ScaleVersor3D=============================== Initializer, optimizerScales: [1, 1, 1, 0.001, 0.001, 0.001, 1, 1, 1]. 0 -0.188366 [0.0199496, -0.013814, -0.0183863, -127.939, -127.124, -125.642, 1.05784, 0.979216, 1.04581] 1 -0.201766 [0.0216539, -0.0114085, -0.0203277, -127.98, -127.127, -125.733, 1.0573, 0.978173, 1.04672] 2 -0.203417 [0.0222002, -0.0101834, -0.0205865, -128.031, -127.136, -125.818, 1.0568, 0.977333, 1.04771] 3 -0.204092 [0.0226688, -0.00855159, -0.0207398, -128.073, -127.144, -125.909, 1.05654, 0.976574, 1.0491] 4 -0.204609 [0.0223696, -0.00716034, -0.0207012, -128.111, -127.166, -125.998, 1.05664, 0.975971, 1.05071] 5 -0.205031 [0.0208482, -0.00776376, -0.0200765, -128.159, -127.185, -126.084, 1.05662, 0.975246, 1.05214] 6 -0.205327 [0.019965, -0.00710902, -0.0209706, -128.187, -127.213, -126.176, 1.0574, 0.974252, 1.05369] 7 -0.205633 [0.0206648, -0.00798805, -0.0209407, -128.228, -127.264, -126.252, 1.05894, 0.972706, 1.05511] 8 -0.205691 [0.0184655, -0.00440914, -0.0191889, -128.177, -127.347, -126.271, 1.05964, 0.972092, 1.05768] 9 -0.204998 [0.0189885, -0.00691607, -0.0228477, -128.214, -127.339, -126.302, 1.0595, 0.971684, 1.05761] 10 -0.205561 [0.0204702, -0.00743821, -0.0185095, -128.218, -127.362, -126.345, 1.05954, 0.970108, 1.05791] 11 -0.205528 [0.0194043, -0.00765185, -0.0282718, -128.209, -127.404, -126.33, 1.05875, 0.971582, 1.05823] 12 -0.203735 [0.0192003, -0.00491258, -0.0197742, -128.165, -127.419, -126.33, 1.0587, 0.971089, 1.05918] 13 -0.205231 [0.0190462, -0.00608451, -0.0210873, -128.186, -127.409, -126.337, 1.05866, 0.9708, 1.05926] 14 -0.205674 [0.0189665, -0.00712029, -0.0212367, -128.204, -127.4, -126.352, 1.05875, 0.970468, 1.05946] 15 -0.205685 [0.0192533, -0.00840991, -0.0209266, -128.224, -127.4, -126.367, 1.05934, 0.969964, 1.05997] 16 -0.205664 [0.0195052, -0.00763662, -0.0209263, -128.213, -127.401, -126.363, 1.05921, 0.970032, 1.0602] 17 -0.205742 [0.0197027, -0.00779215, -0.02067, -128.206, -127.399, -126.373, 1.0599, 0.969881, 1.06028] 18 -0.205735 [0.0193924, -0.0076042, -0.0209455, -128.201, -127.399, -126.369, 1.05979, 0.970011, 1.06024] Matrix = 1.0588 0.0415777 -0.0160141 -0.0421676 0.968382 -0.0384494 0.0143894 0.0390865 1.05938 Offset = [-139.811, -112.319, -141.719] =============================== Starting Transform Estimations for BSpline=============================== Starting Registration ==27809== Warning: set address range perms: large range 479927808, a 0, v 1 ==27809== Warning: set address range perms: large range 479927808, a 0, v 1 ==27809== Thread 8: ==27809== Conditional jump or move depends on uninitialised value(s) ==27809== at 0x4D2B96A: itk::ImageToImageMetric<itk::Image<float, 3>, itk::Image<float, 3> >::TransformPointWithDerivatives(unsigned, itk::Point<double, 3>&, bool&, double&, itk::CovariantVector<double, 3>&, unsigned) const (itkOptImageToImageMetric.txx:890) ==27809== by 0x4D8F830: itk::ImageToImageMetric<itk::Image<float, 3>, itk::Image<float, 3> >::GetValueAndDerivativeThread(unsigned) const (itkOptImageToImageMetric.txx:1316) ==27809== by 0x4D71115: itk::ImageToImageMetric<itk::Image<float, 3>, itk::Image<float, 3> >::GetValueAndDerivativeMultiThreaded(void*) (itkOptImageToImageMetric.txx:1260) ==27809== by 0x6FD3B13: itk::MultiThreader::SingleMethodProxy(void*) (itkMultiThreader.cxx:853) ==27809== by 0x30AA506136: start_thread (in /lib64/tls/libpthread-2.3.4.so) ==27809== by 0x30A9EC9F02: clone (in /lib64/tls/libc-2.3.4.so) ==27809== ==27809== Conditional jump or move depends on uninitialised value(s) ==27809== at 0x4D2BBB6: itk::ImageToImageMetric<itk::Image<float, 3>, itk::Image<float, 3> >::TransformPointWithDerivatives(unsigned, itk::Point<double, 3>&, bool&, double&, itk::CovariantVector<double, 3>&, unsigned) const (itkOptImageToImageMetric.txx:938) ==27809== by 0x4D8F830: itk::ImageToImageMetric<itk::Image<float, 3>, itk::Image<float, 3> >::GetValueAndDerivativeThread(unsigned) const (itkOptImageToImageMetric.txx:1316) ==27809== by 0x4D71115: itk::ImageToImageMetric<itk::Image<float, 3>, itk::Image<float, 3> >::GetValueAndDerivativeMultiThreaded(void*) (itkOptImageToImageMetric.txx:1260) ==27809== by 0x6FD3B13: itk::MultiThreader::SingleMethodProxy(void*) (itkMultiThreader.cxx:853) ==27809== by 0x30AA506136: start_thread (in /lib64/tls/libpthread-2.3.4.so) ==27809== by 0x30A9EC9F02: clone (in /lib64/tls/libc-2.3.4.so) ==27809== ==27809== Conditional jump or move depends on uninitialised value(s) ==27809== at 0x4D8F835: itk::ImageToImageMetric<itk::Image<float, 3>, itk::Image<float, 3> >::GetValueAndDerivativeThread(unsigned) const (itkOptImageToImageMetric.txx:1320) ==27809== by 0x4D71115: itk::ImageToImageMetric<itk::Image<float, 3>, itk::Image<float, 3> >::GetValueAndDerivativeMultiThreaded(void*) (itkOptImageToImageMetric.txx:1260) ==27809== by 0x6FD3B13: itk::MultiThreader::SingleMethodProxy(void*) (itkMultiThreader.cxx:853) ==27809== by 0x30AA506136: start_thread (in /lib64/tls/libpthread-2.3.4.so) ==27809== by 0x30A9EC9F02: clone (in /lib64/tls/libc-2.3.4.so) ==27809== ==27809== Invalid read of size 8 ==27809== at 0x4D5A2EF: std::_Bit_reference::operator bool() const (stl_bvector.h:79) ==27809== by 0x4D5A2D2: std::_Bit_const_iterator::operator*() const (stl_bvector.h:279) ==27809== by 0x4D5A296: std::vector<bool, std::allocator<bool> >::operator[](unsigned long) const (stl_bvector.h:576) ==27809== by 0x4D2B95A: itk::ImageToImageMetric<itk::Image<float, 3>, itk::Image<float, 3> >::TransformPointWithDerivatives(unsigned, itk::Point<double, 3>&, bool&, double&, itk::CovariantVector<double, 3>&, unsigned) const (itkOptImageToImageMetric.txx:888) ==27809== by 0x4D8F830: itk::ImageToImageMetric<itk::Image<float, 3>, itk::Image<float, 3> >::GetValueAndDerivativeThread(unsigned) const (itkOptImageToImageMetric.txx:1316) ==27809== by 0x4D71115: itk::ImageToImageMetric<itk::Image<float, 3>, itk::Image<float, 3> >::GetValueAndDerivativeMultiThreaded(void*) (itkOptImageToImageMetric.txx:1260) ==27809== by 0x6FD3B13: itk::MultiThreader::SingleMethodProxy(void*) (itkMultiThreader.cxx:853) ==27809== by 0x30AA506136: start_thread (in /lib64/tls/libpthread-2.3.4.so) ==27809== by 0x30A9EC9F02: clone (in /lib64/tls/libc-2.3.4.so) ==27809== Address 0x82210B0 is 0 bytes after a block of size 117,176 alloc'd ==27809== at 0x4904DB5: operator new(unsigned long) (vg_replace_malloc.c:168) ==27809== by 0x4D91920: __gnu_cxx::new_allocator<unsigned long>::allocate(unsigned long, void const*) (new_allocator.h:81) ==27809== by 0x4D95EB0: std::_Bvector_base<std::allocator<bool> >::_M_allocate(unsigned long) (stl_bvector.h:379) ==27809== by 0x4D905E6: std::vector<bool, std::allocator<bool> >::_M_fill_insert(std::_Bit_iterator, unsigned long, bool) (stl_bvector.h:824) ==27809== by 0x4D71EB9: std::vector<bool, std::allocator<bool> >::insert(std::_Bit_iterator, unsigned long, bool) (stl_bvector.h:837) ==27809== by 0x4D59B6F: std::vector<bool, std::allocator<bool> >::resize(unsigned long, bool) (stl_bvector.h:861) ==27809== by 0x4D29405: itk::ImageToImageMetric<itk::Image<float, 3>, itk::Image<float, 3> >::MultiThreadingInitialize() (itkOptImageToImageMetric.txx:383) ==27809== by 0x4D265F1: itk::MattesMutualInformationImageToImageMetric<itk::Image<float, 3>, itk::Image<float, 3> >::Initialize() (itkOptMattesMutualInformationImageToImageMetric.txx:180) ==27809== by 0x4C9ECC3: itk::ImageRegistrationMethod<itk::Image<float, 3>, itk::Image<float, 3> >::Initialize() (itkImageRegistrationMethod.txx:218) ==27809== by 0x4C3F26B: itk::ImageRegistrationMethod<itk::Image<float, 3>, itk::Image<float, 3> >::StartRegistration() (itkImageRegistrationMethod.txx:268) ==27809== by 0x4D2192A: itk::ImageRegistrationMethod<itk::Image<float, 3>, itk::Image<float, 3> >::GenerateData() (itkImageRegistrationMethod.txx:342) ==27809== by 0x6FF1AE5: itk::ProcessObject::UpdateOutputData(itk::DataObject*) (itkProcessObject.cxx:987) ==27809== ==27809== Invalid read of size 8 ==27809== at 0x4D5A388: vnl_matrix<double>::operator[](unsigned) const (vnl_matrix.h:171) ==27809== by 0x4D2B982: itk::ImageToImageMetric<itk::Image<float, 3>, itk::Image<float, 3> >::TransformPointWithDerivatives(unsigned, itk::Point<double, 3>&, bool&, double&, itk::CovariantVector<double, 3>&, unsigned) const (itkOptImageToImageMetric.txx:895) ==27809== by 0x4D8F830: itk::ImageToImageMetric<itk::Image<float, 3>, itk::Image<float, 3> >::GetValueAndDerivativeThread(unsigned) const (itkOptImageToImageMetric.txx:1316) ==27809== by 0x4D71115: itk::ImageToImageMetric<itk::Image<float, 3>, itk::Image<float, 3> >::GetValueAndDerivativeMultiThreaded(void*) (itkOptImageToImageMetric.txx:1260) ==27809== by 0x6FD3B13: itk::MultiThreader::SingleMethodProxy(void*) (itkMultiThreader.cxx:853) ==27809== by 0x30AA506136: start_thread (in /lib64/tls/libpthread-2.3.4.so) ==27809== by 0x30A9EC9F02: clone (in /lib64/tls/libc-2.3.4.so) ==27809== Address 0x10B073B8 is 352 bytes inside a block of size 79,560 free'd ==27809== at 0x4905B14: operator delete[](void*) (vg_replace_malloc.c:256) ==27809== by 0x4C4DD49: vnl_sse_dealloc(void*, unsigned, unsigned) (vnl_sse.h:101) ==27809== by 0x4C4DD24: vnl_c_vector_dealloc(void*, int, int) (vnl_c_vector.txx:349) ==27809== by 0x4C4DD01: vnl_c_vector<double>::deallocate(double*, int) (vnl_c_vector.txx:418) ==27809== by 0x4C410BC: vnl_vector<double>::destroy() (vnl_vector.txx:277) ==27809== by 0x4C35B59: vnl_vector<double>::~vnl_vector() (vnl_vector.txx:268) ==27809== by 0x4C35B36: itk::Array<double>::~Array() (itkArray.txx:71) ==27809== by 0x4D2B2DD: itk::ImageToImageMetric<itk::Image<float, 3>, itk::Image<float, 3> >::PreComputeTransformValues() (itkOptImageToImageMetric.txx:714) ==27809== by 0x4D2941B: itk::ImageToImageMetric<itk::Image<float, 3>, itk::Image<float, 3> >::MultiThreadingInitialize() (itkOptImageToImageMetric.txx:385) ==27809== by 0x4D265F1: itk::MattesMutualInformationImageToImageMetric<itk::Image<float, 3>, itk::Image<float, 3> >::Initialize() (itkOptMattesMutualInformationImageToImageMetric.txx:180) ==27809== by 0x4C9ECC3: itk::ImageRegistrationMethod<itk::Image<float, 3>, itk::Image<float, 3> >::Initialize() (itkImageRegistrationMethod.txx:218) ==27809== by 0x4C3F26B: itk::ImageRegistrationMethod<itk::Image<float, 3>, itk::Image<float, 3> >::StartRegistration() (itkImageRegistrationMethod.txx:268) ==27809== ==27809== Invalid read of size 8 ==27809== at 0x4D2B9EA: itk::ImageToImageMetric<itk::Image<float, 3>, itk::Image<float, 3> >::TransformPointWithDerivatives(unsigned, itk::Point<double, 3>&, bool&, double&, itk::CovariantVector<double, 3>&, unsigned) const (itkOptImageToImageMetric.txx:901) ==27809== by 0x4D8F830: itk::ImageToImageMetric<itk::Image<float, 3>, itk::Image<float, 3> >::GetValueAndDerivativeThread(unsigned) const (itkOptImageToImageMetric.txx:1316) ==27809== by 0x4D71115: itk::ImageToImageMetric<itk::Image<float, 3>, itk::Image<float, 3> >::GetValueAndDerivativeMultiThreaded(void*) (itkOptImageToImageMetric.txx:1260) ==27809== by 0x6FD3B13: itk::MultiThreader::SingleMethodProxy(void*) (itkMultiThreader.cxx:853) ==27809== by 0x30AA506136: start_thread (in /lib64/tls/libpthread-2.3.4.so) ==27809== by 0x30A9EC9F02: clone (in /lib64/tls/libc-2.3.4.so) ==27809== Address 0x3D46CA48 is 1,152 bytes inside a block of size 79,560 free'd ==27809== at 0x4905B14: operator delete[](void*) (vg_replace_malloc.c:256) ==27809== by 0x4C4DD49: vnl_sse_dealloc(void*, unsigned, unsigned) (vnl_sse.h:101) ==27809== by 0x4C4DD24: vnl_c_vector_dealloc(void*, int, int) (vnl_c_vector.txx:349) ==27809== by 0x4C4DD01: vnl_c_vector<double>::deallocate(double*, int) (vnl_c_vector.txx:418) ==27809== by 0x4C410BC: vnl_vector<double>::destroy() (vnl_vector.txx:277) ==27809== by 0x4C35B59: vnl_vector<double>::~vnl_vector() (vnl_vector.txx:268) ==27809== by 0x4C35B36: itk::Array<double>::~Array() (itkArray.txx:71) ==27809== by 0x4D28434: itk::MattesMutualInformationImageToImageMetric<itk::Image<float, 3>, itk::Image<float, 3> >::Initialize() (itkOptMattesMutualInformationImageToImageMetric.txx:444) ==27809== by 0x4C9ECC3: itk::ImageRegistrationMethod<itk::Image<float, 3>, itk::Image<float, 3> >::Initialize() (itkImageRegistrationMethod.txx:218) ==27809== by 0x4C3F26B: itk::ImageRegistrationMethod<itk::Image<float, 3>, itk::Image<float, 3> >::StartRegistration() (itkImageRegistrationMethod.txx:268) ==27809== by 0x4D2192A: itk::ImageRegistrationMethod<itk::Image<float, 3>, itk::Image<float, 3> >::GenerateData() (itkImageRegistrationMethod.txx:342) ==27809== by 0x6FF1AE5: itk::ProcessObject::UpdateOutputData(itk::DataObject*) (itkProcessObject.cxx:987) ==27809== ==27809== Invalid read of size 8 ==27809== at 0x4D2BA6F: itk::ImageToImageMetric<itk::Image<float, 3>, itk::Image<float, 3> >::TransformPointWithDerivatives(unsigned, itk::Point<double, 3>&, bool&, double&, itk::CovariantVector<double, 3>&, unsigned) const (itkOptImageToImageMetric.txx:908) ==27809== by 0x4D8F830: itk::ImageToImageMetric<itk::Image<float, 3>, itk::Image<float, 3> >::GetValueAndDerivativeThread(unsigned) const (itkOptImageToImageMetric.txx:1316) ==27809== by 0x4D71115: itk::ImageToImageMetric<itk::Image<float, 3>, itk::Image<float, 3> >::GetValueAndDerivativeMultiThreaded(void*) (itkOptImageToImageMetric.txx:1260) ==27809== by 0x6FD3B13: itk::MultiThreader::SingleMethodProxy(void*) (itkMultiThreader.cxx:853) ==27809== by 0x30AA506136: start_thread (in /lib64/tls/libpthread-2.3.4.so) ==27809== by 0x30A9EC9F02: clone (in /lib64/tls/libc-2.3.4.so) ==27809== Address 0x0 is not stack'd, malloc'd or (recently) free'd ==27809== ==27809== Process terminating with default action of signal 11 (SIGSEGV) ==27809== Access not within mapped region at address 0x0 ==27809== at 0x4D2BA6F: itk::ImageToImageMetric<itk::Image<float, 3>, itk::Image<float, 3> >::TransformPointWithDerivatives(unsigned, itk::Point<double, 3>&, bool&, double&, itk::CovariantVector<double, 3>&, unsigned) const (itkOptImageToImageMetric.txx:908) ==27809== by 0x4D8F830: itk::ImageToImageMetric<itk::Image<float, 3>, itk::Image<float, 3> >::GetValueAndDerivativeThread(unsigned) const (itkOptImageToImageMetric.txx:1316) ==27809== by 0x4D71115: itk::ImageToImageMetric<itk::Image<float, 3>, itk::Image<float, 3> >::GetValueAndDerivativeMultiThreaded(void*) (itkOptImageToImageMetric.txx:1260) ==27809== by 0x6FD3B13: itk::MultiThreader::SingleMethodProxy(void*) (itkMultiThreader.cxx:853) ==27809== by 0x30AA506136: start_thread (in /lib64/tls/libpthread-2.3.4.so) ==27809== by 0x30A9EC9F02: clone (in /lib64/tls/libc-2.3.4.so) ==27809== ==27809== ERROR SUMMARY: 154 errors from 7 contexts (suppressed: 5 from 2) ==27809== malloc/free: in use at exit: 1,413,843,410 bytes in 14,677 blocks. ==27809== malloc/free: 3,677,748 allocs, 3,663,071 frees, 4,781,492,805 bytes allocated. ==27809== For counts of detected errors, rerun with: -v ==27809== searching for pointers to 14,677 not-freed blocks. ==27809== checked 1,477,967,448 bytes. ==27809== ==27809== LEAK SUMMARY: ==27809== definitely lost: 0 bytes in 0 blocks. ==27809== possibly lost: 392,990 bytes in 11,096 blocks. ==27809== still reachable: 1,413,450,420 bytes in 3,581 blocks. ==27809== suppressed: 0 bytes in 0 blocks. ==27809== Reachable blocks (those to which a pointer was found) are not shown. ==27809== To see them, rerun with: --show-reachable=yes |
(0016711) Hans Johnson (developer) 2009-06-11 16:38 |
Committed fix to previously reported bug. /cvsroot/Insight/Insight/Code/Review/itkOptImageToImageMetric.h,v <-- itkOptImageToImageMetric.h new revision: 1.23; previous revision: 1.22 /cvsroot/Insight/Insight/Code/Review/itkOptImageToImageMetric.txx,v <-- itkOptImageToImageMetric.txx new revision: 1.30; previous revision: 1.29 -- Hans J. Johnson, Ph.D. Associate Faculty hans-johnson@uiowa.edu 319 353 8587 ? From: Hans Johnson <hans-johnson@uiowa.edu> Date: Thu, 11 Jun 2009 14:15:14 -0500 To: Hans Johnson <hans-johnson@uiowa.edu>, Luis Ibanez <luis.ibanez@kitware.com>, Kent Williams <norman-k-williams@uiowa.edu> Cc: ITK <insight-developers@itk.org> Subject: Re: [Insight-developers] Bug in Code/Review/itkOptImageToImageMetric.txx? Luis, I have committed a test that now exercises the situation described below. It should start failing on the dashboard for those builds that use ITK_USE_OPTIMIZED_REGISTRATION_METHODS:BOOL=ON. cvs commit -m"BUG:0009139: Added a test that fails because the number of samples is internally reduced when the fixed image mask is smaller than the number of samples requested. This leads to an inconsistent initialization of several variables. The solution listed in bug 0009139 ensures that all internal arrays are re-initialized to the correct length if this internal reduction in samples occurs." Testing/Code/Algorithms/itkMattesMutualInformationImageToImageMetricTest.cxx Committer: Hans Johnson <hjohnson@mail.psychiatry.uiowa.edu> /cvsroot/Insight/Insight/Testing/Code/Algorithms/itkMattesMutualInformationImageToImageMetricTest.cxx,v <-- Testing/Code/Algorithms/itkMattesMutualInformationImageToImageMetricTest.cxx new revision: 1.21; previous revision: 1.20 ==================== I’ve not yet committed the fix to the code that is attached to the bug report. I’ll do that in about 1 hour. Hans -- Hans J. Johnson, Ph.D. Associate Faculty hans-johnson@uiowa.edu 319 353 8587 ? From: Hans Johnson <hans-johnson@uiowa.edu> Date: Thu, 11 Jun 2009 11:15:26 -0500 To: Luis Ibanez <luis.ibanez@kitware.com>, Kent Williams <norman-k-williams@uiowa.edu> Cc: ITK <insight-developers@itk.org> Subject: Re: [Insight-developers] Bug in Code/Review/itkOptImageToImageMetric.txx? Luis, I think I’ve finally tracked down the exact problem. Much of the code in ::MultiThreadingInitialize has a dependency on the variable m_NumberOfFixedImageSamples for determining the length of the weights arrays. This is OK if that variable were only set from the class API before computing the registration. Unfortunately, after initialization it is possible that the FixedImageRegion, or the number of voxels available in the FixedImage mask can change this value after initialization to a smaller number. If that is the case, then the size of the weights arrays also need to be re-adjusted in size. Ultimately the segmentation fault arises in a for loop over several arrays that should all be the same length (i.e. The new smaller value of m_NumberOfFixedImageSamples), but the loop index is specified as one of the arrays that has not been adjusted to the smaller size. Currently the code changes I have pass all the old regression test, and produce a good result in the “real world” code I have (where it was previously segmentation faulting). I do not yet have a regression test that fails with the old code and succeeds with the new code. I’ve spent more time trying to make a small test case than it took to actually fix the bug that was discovered by running a real test case. Currently when attempting to make a test case that only causes this problem, other error checking is stepping in front of this and not even getting to cause this problem. I’ve added bug number “0009139” with a proposed patch that still needs testing. Hans -- Hans J. Johnson, Ph.D. Associate Faculty hans-johnson@uiowa.edu 319 353 8587 ? From: Luis Ibanez <luis.ibanez@kitware.com> Date: Wed, 10 Jun 2009 12:57:56 -0400 To: Kent Williams <norman-k-williams@uiowa.edu> Cc: ITK <insight-developers@itk.org> Subject: Re: [Insight-developers] Bug in Code/Review/itkOptImageToImageMetric.txx? Hi Kent, The number of fixed samples is updated according the number of pixels inside the (potential) fixed image mask in lines: 454-513, more specifically in 467-472 if ( count > maxcount || randIter.IsAtEnd() ) { m_NumberOfFixedImageSamples = samples_found; samples.resize(samples_found); break; } but... of course... it is always possible that there is a bug and we missed some combination of circumstances. Any chance that we could reduce this to a minimal case ? E.g. I would think in taking one of the examples in Insight/Examples/Registration and adding a Mask to it, so that we can reproduce this in a small test case. I'll be happy to help setting this up. BTW: In your case are us using AllPixels or are you letting the Metric subsample the image ? Depending on that, the code will take two different paths. Thanks Luis -------------------------------------------------------------------------------------------------- On Wed, Jun 10, 2009 at 12:27 PM, kent williams <norman-k-williams@uiowa.edu> wrote: We ran into a problem that Hans could probably explain better than I: We are getting core dumps in itkOptImageToImageMetric.txx. The problem is that there is a member array m_BSplineTransformWeightsArray, that gets sized based on m_NumberOfFixedImageSamples and m_NumBSplineWeights. It's getting allocated (in our case to 10,000 samples), but when it gets filled in (in PreComputeTransformValues()), every element is not set -- some are left with a value of zero or some random value. I think the problem is that the iteration over weights doesn't take into account the reduced actual sample count due to using masks. if(sampleOk) { // If the transform is BSplineDeformable, we can use the precomputed // weights and indices to obtained the mapped position const WeightsValueType * weights = if(sampleOk) { // If the transform is BSplineDeformable, we can use the precomputed // weights and indices to obtained the mapped position const WeightsValueType * weights = m_BSplineTransformWeightsArray[sampleNumber]; const IndexValueType * indices = m_BSplineTransformIndicesArray[sampleNumber]; for( unsigned int j = 0; j < FixedImageDimension; j++ ) { mappedPoint[j] = m_BSplinePreTransformPointsArray[sampleNumber][j]; } for ( unsigned int k = 0; k < m_NumBSplineWeights; k++ ) { for ( unsigned int j = 0; j < FixedImageDimension; j++ ) { //CRASHES HERE line 908 mappedPoint[j] += weights[k] * m_Parameters[ indices[k] + m_BSplineParametersOffset[j] ]; } } } } m_BSplineTransformWeightsArray[sampleNumber]; const IndexValueType * indices = m_BSplineTransformIndicesArray[sampleNumber]; for( unsigned int j = 0; j < FixedImageDimension; j++ ) { mappedPoint[j] = m_BSplinePreTransformPointsArray[sampleNumber][j]; } for ( unsigned int k = 0; k < m_NumBSplineWeights; k++ ) { for ( unsigned int j = 0; j < FixedImageDimension; j++ ) { mappedPoint[j] += weights[k] * m_Parameters[ indices[k] + m_BSplineParametersOffset[j] ]; } } } } |
(0018913) Hans Johnson (developer) 2009-12-16 08:25 |
I believe that this was addressed and fixed in the summer of 2009. |
Notes |
Issue History | |||
Date Modified | Username | Field | Change |
2009-06-11 12:08 | Hans Johnson | New Issue | |
2009-06-11 12:11 | Hans Johnson | Note Added: 0016708 | |
2009-06-11 12:14 | Hans Johnson | Assigned To | => Hans Johnson |
2009-06-11 12:14 | Hans Johnson | Status | new => assigned |
2009-06-11 12:14 | Hans Johnson | Description Updated | |
2009-06-11 12:14 | Hans Johnson | File Added: itkOptImageToImageMetric.h | |
2009-06-11 12:14 | Hans Johnson | File Added: itkOptImageToImageMetric.txx | |
2009-06-11 16:38 | Hans Johnson | Note Added: 0016711 | |
2009-12-16 08:25 | Hans Johnson | Note Added: 0018913 | |
2009-12-16 08:25 | Hans Johnson | Status | assigned => closed |
2009-12-16 08:25 | Hans Johnson | Resolution | open => fixed |
Issue History |
Copyright © 2000 - 2018 MantisBT Team |