Skip to content

Commit 10801b6

Browse files
N-Dekkerhjmjohnson
authored andcommitted
BUG: Add ImageAdaptor::ComputeOffset, to fix ImageRegionIterator support
Fixed issue #5870, "ImageAdaptor should have its own ComputeOffset, ImageRegionIterator support broken" `ImageRegionConstIterator<AdapterType>::Increment()` depends on `AdapterType::ComputeOffset`.
1 parent 22f6d27 commit 10801b6

File tree

2 files changed

+77
-2
lines changed

2 files changed

+77
-2
lines changed

Modules/Core/ImageAdaptors/include/itkImageAdaptor.h

Lines changed: 7 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -234,6 +234,13 @@ class ITK_TEMPLATE_EXPORT ImageAdaptor : public ImageBase<TImage::ImageDimension
234234
IndexType
235235
ComputeIndex(OffsetValueType offset) const;
236236

237+
/** Computes the offset from the beginning of the buffer. */
238+
OffsetValueType
239+
ComputeOffset(const IndexType & index) const
240+
{
241+
return m_Image->ComputeOffset(index);
242+
}
243+
237244
/** PixelContainer type alias support Used to construct a container for
238245
* the pixel data. */
239246
using PixelContainer = typename TImage::PixelContainer;

Modules/Core/ImageAdaptors/test/itkImageAdaptorGTest.cxx

Lines changed: 70 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -17,11 +17,12 @@
1717
*=========================================================================*/
1818

1919
// First include the header file to be tested:
20-
#include "itkImageBase.h"
20+
#include "itkImageAdaptor.h"
2121

2222
#include "itkCovariantVector.h"
2323
#include "itkImage.h"
24-
#include "itkImageAdaptor.h"
24+
#include "itkImageRegionConstIteratorWithIndex.h"
25+
#include "itkImageRegionIterator.h"
2526
#include "itkVector.h"
2627

2728
#include <gtest/gtest.h>
@@ -185,3 +186,70 @@ TEST(ImageAdapter, TransformMemberFunctionsReturnSameAsForImage)
185186
const auto image3D = itk::Image<unsigned char, 3>::New();
186187
Expect_Transform_member_functions_return_the_same_for_an_ImageAdapter_as_for_its_image(*image3D);
187188
}
189+
190+
191+
// Tests that ComputeOffset(index) returns the same value for an adaptor as for its internal image.
192+
TEST(ImageAdaptor, ComputeOffset)
193+
{
194+
constexpr unsigned int dimension{ 2 };
195+
using PixelType = int;
196+
using ImageType = itk::Image<PixelType, dimension>;
197+
198+
const auto image = ImageType::New();
199+
const auto adaptor = itk::ImageAdaptor<ImageType, DummyPixelAccessor<PixelType>>::New();
200+
201+
// Use a very small image region, to speed up the unit test.
202+
const itk::ImageRegion imageRegion{ itk::Index<dimension>::Filled(1), itk::Size<dimension>::Filled(2) };
203+
204+
adaptor->SetImage(image);
205+
adaptor->SetRegions(imageRegion);
206+
adaptor->AllocateInitialized();
207+
208+
itk::ImageRegionConstIteratorWithIndex<ImageType> it(image, imageRegion);
209+
for (it.GoToBegin(); !it.IsAtEnd(); ++it)
210+
{
211+
EXPECT_EQ(adaptor->ComputeOffset(it.GetIndex()), image->ComputeOffset(it.GetIndex()));
212+
}
213+
}
214+
215+
216+
// Tests that ImageAdaptor supports iteration by means of ImageRegionIterator.
217+
TEST(ImageAdaptor, SupportsRegionIterator)
218+
{
219+
constexpr unsigned int dimension{ 2 };
220+
using PixelType = int;
221+
using ImageType = itk::Image<PixelType, dimension>;
222+
223+
const auto image = ImageType::New();
224+
const auto adaptor = itk::ImageAdaptor<ImageType, DummyPixelAccessor<PixelType>>::New();
225+
226+
adaptor->SetImage(image);
227+
adaptor->SetRegions(itk::Size<dimension>::Filled(4));
228+
adaptor->AllocateInitialized();
229+
230+
constexpr auto regionIndex = itk::Index<dimension>::Filled(1);
231+
constexpr auto regionSize = itk::Size<dimension>::Filled(2);
232+
const itk::ImageRegion region{ regionIndex, regionSize };
233+
234+
// Set the pixel values by means of the adaptor, for the specified region.
235+
PixelType pixelValue{};
236+
constexpr PixelType maxPixelValue{ regionSize.CalculateProductOfElements() };
237+
238+
using AdaptorType = itk::ImageAdaptor<ImageType, DummyPixelAccessor<PixelType>>;
239+
for (itk::ImageRegionIterator<AdaptorType> iterator(adaptor, region); !iterator.IsAtEnd(); ++iterator)
240+
{
241+
++pixelValue;
242+
ASSERT_LE(pixelValue, maxPixelValue);
243+
iterator.Set(pixelValue);
244+
EXPECT_EQ(iterator.Get(), pixelValue);
245+
}
246+
247+
// Now check if the intenal image has got the expected pixel values, in the specified region.
248+
PixelType expectedPixelValue{};
249+
250+
for (itk::ImageRegionConstIterator<ImageType> iterator(image, region); !iterator.IsAtEnd(); ++iterator)
251+
{
252+
++expectedPixelValue;
253+
EXPECT_EQ(iterator.Get(), expectedPixelValue);
254+
}
255+
}

0 commit comments

Comments
 (0)