Skip to content

Commit ff1ffe0

Browse files
author
Simon Rit
committed
BUG: Remove const_cast from VariableLengthVector constructor
The `Get` function of `VectorImage` iterators calls the `VariableLengthVector` constructor parameterized by a const pointer and a vector length (see `DefaultVectorPixelAccessor::Get`). The `m_Data` pointer of the constructed VariableLengthVector which was returned was then pointing to the returned pixel data. This led to an unexpected behavior in `CastImageFilter::DynamicThreadedGenerateDataDispatched` where modifications to the local variable ``` OutputPixelType value{ outputIt.Get() }; ``` of `CastImageFilter::DynamicThreadedGenerateDataDispatched` were actually modifying the pixel of the image at the iterator position. The problem might occur elsewhere so it is safer to copy the data when a VariableLengthVector is constructed from a const pointer despite the potential loss in performances. A `CastImageFilter` test has been added to illustrate the issue: when a `VectorImage` was casted to an `Image< Vector <`, the first pixel of the region was taking the value of the last pixel of the region. The test illustrated this by creating an image filled with 0 except the first pixel. The casted image was completely full of 0 before bug correction.
1 parent a8dd05b commit ff1ffe0

File tree

2 files changed

+58
-2
lines changed

2 files changed

+58
-2
lines changed

Modules/Core/Common/include/itkVariableLengthVector.hxx

Lines changed: 13 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -45,9 +45,20 @@ VariableLengthVector<TValue>::VariableLengthVector(ValueType * datain, unsigned
4545
template <typename TValue>
4646
VariableLengthVector<TValue>::VariableLengthVector(const ValueType * datain, unsigned int sz, bool LetArrayManageMemory)
4747
: m_LetArrayManageMemory(LetArrayManageMemory)
48-
, m_Data(const_cast<ValueType *>(datain))
4948
, m_NumElements(sz)
50-
{}
49+
{
50+
if (m_NumElements != 0)
51+
{
52+
m_Data = this->AllocateElements(m_NumElements);
53+
itkAssertInDebugAndIgnoreInReleaseMacro(m_Data != nullptr);
54+
itkAssertInDebugAndIgnoreInReleaseMacro(datain != nullptr);
55+
std::copy_n(datain, m_NumElements, &this->m_Data[0]);
56+
}
57+
else
58+
{
59+
m_Data = nullptr;
60+
}
61+
}
5162

5263
template <typename TValue>
5364
VariableLengthVector<TValue>::VariableLengthVector(const VariableLengthVector<TValue> & v)

Modules/Filtering/ImageFilterBase/test/itkCastImageFilterTest.cxx

Lines changed: 45 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -400,6 +400,50 @@ TestVectorImageCast2()
400400
}
401401

402402

403+
bool
404+
TestVectorImageCast3()
405+
{
406+
// This function casts an Image<Vector<float, 3>, 3>
407+
// to a VectorImage<float, 3>
408+
std::cout << "Casting from an Image<Vector<float,3>, 3> to VectorImage<float, 3> ... ";
409+
410+
using VectorFloatImageType = itk::Image<itk::Vector<float, 3>, 3>;
411+
using FloatVectorImageType = itk::VectorImage<float, 3>;
412+
413+
// Create a 1x3 image of 2D vectors
414+
auto image = VectorFloatImageType::New();
415+
416+
constexpr itk::Size<3> size{ { 1, 1, 4 } };
417+
const itk::ImageRegion<3> region{ size };
418+
image->SetRegions(region);
419+
image->AllocateInitialized();
420+
VectorFloatImageType::PixelType vec{ { 1.3, 5.3, 6.4 } };
421+
// Only the first pixel will be the vector (1.3, 5.3, 6.4)
422+
image->SetPixel(itk::Index<3>(), vec);
423+
424+
using CastImageFilterType = itk::CastImageFilter<VectorFloatImageType, FloatVectorImageType>;
425+
auto castImageFilter = CastImageFilterType::New();
426+
castImageFilter->SetInput(image);
427+
castImageFilter->SetNumberOfWorkUnits(1);
428+
castImageFilter->Update();
429+
430+
bool success = true;
431+
if (std::memcmp(castImageFilter->GetOutput()->GetBufferPointer(),
432+
image->GetBufferPointer(),
433+
sizeof(float) * size.CalculateProductOfElements()) == 0)
434+
{
435+
std::cout << "[PASSED]" << std::endl;
436+
}
437+
else
438+
{
439+
std::cout << "[FAILED]" << std::endl;
440+
success = false;
441+
}
442+
443+
return success;
444+
}
445+
446+
403447
int
404448
itkCastImageFilterTest(int, char *[])
405449
{
@@ -431,6 +475,7 @@ itkCastImageFilterTest(int, char *[])
431475
success &= TestCastFrom<double>();
432476
success &= TestVectorImageCast1();
433477
success &= TestVectorImageCast2();
478+
success &= TestVectorImageCast3();
434479

435480
std::cout << std::endl;
436481
if (!success)

0 commit comments

Comments
 (0)