Skip to content

Commit 66e6d8b

Browse files
N-Dekkerdzenanz
authored andcommitted
ENH: Add GoogleTest unit test ImageRegionIterator.IsConstructedAtBegin
Checks that an iterator, constructed by `IteratorType(image, region)` is at the begin. Checks 22 different iterator types from Core/Common.
1 parent ffe60a3 commit 66e6d8b

File tree

2 files changed

+137
-0
lines changed

2 files changed

+137
-0
lines changed

Modules/Core/Common/test/CMakeLists.txt

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -1735,6 +1735,7 @@ set(ITKCommonGTests
17351735
itkImageIORegionGTest.cxx
17361736
itkImageRandomConstIteratorWithIndexGTest.cxx
17371737
itkImageRegionGTest.cxx
1738+
itkImageRegionIteratorGTest.cxx
17381739
itkIndexGTest.cxx
17391740
itkIndexRangeGTest.cxx
17401741
itkLightObjectGTest.cxx
Lines changed: 136 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,136 @@
1+
/*=========================================================================
2+
*
3+
* Copyright NumFOCUS
4+
*
5+
* Licensed under the Apache License, Version 2.0 (the "License");
6+
* you may not use this file except in compliance with the License.
7+
* You may obtain a copy of the License at
8+
*
9+
* https://www.apache.org/licenses/LICENSE-2.0.txt
10+
*
11+
* Unless required by applicable law or agreed to in writing, software
12+
* distributed under the License is distributed on an "AS IS" BASIS,
13+
* WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14+
* See the License for the specific language governing permissions and
15+
* limitations under the License.
16+
*
17+
*=========================================================================*/
18+
19+
// First include the header files to be tested:
20+
#include "itkImageConstIterator.h"
21+
#include "itkImageConstIteratorWithIndex.h"
22+
#include "itkImageConstIteratorWithOnlyIndex.h"
23+
#include "itkImageIterator.h"
24+
#include "itkImageIteratorWithIndex.h"
25+
#include "itkImageLinearConstIteratorWithIndex.h"
26+
#include "itkImageLinearIteratorWithIndex.h"
27+
#include "itkImageRegionConstIterator.h"
28+
#include "itkImageRegionConstIteratorWithIndex.h"
29+
#include "itkImageRegionConstIteratorWithOnlyIndex.h"
30+
#include "itkImageRegionExclusionConstIteratorWithIndex.h"
31+
#include "itkImageRegionExclusionIteratorWithIndex.h"
32+
#include "itkImageRegionIterator.h"
33+
#include "itkImageRegionIteratorWithIndex.h"
34+
#include "itkImageRegionReverseConstIterator.h"
35+
#include "itkImageRegionReverseIterator.h"
36+
#include "itkImageReverseConstIterator.h"
37+
#include "itkImageReverseIterator.h"
38+
#include "itkImageScanlineConstIterator.h"
39+
#include "itkImageScanlineIterator.h"
40+
#include "itkImageSliceConstIteratorWithIndex.h"
41+
#include "itkImageSliceIteratorWithIndex.h"
42+
43+
#include "itkImage.h"
44+
#include <gtest/gtest.h>
45+
#include <type_traits> // For false_type and true_type.
46+
47+
namespace
48+
{
49+
// Has_IsAtBegin tells whether the specified iterator type has an `IsAtBegin()` member function, returning a `bool`. It
50+
// is inspired by "How to detect whether there is a specific member variable in class?", answered by Andy Prowl, Jan 25,
51+
// 2013 at https://stackoverflow.com/a/14523787 (the `has_id` example).
52+
template <typename TIterator, typename = bool>
53+
struct Has_IsAtBegin : std::false_type
54+
{};
55+
56+
template <typename TIterator>
57+
struct Has_IsAtBegin<TIterator, decltype(std::declval<TIterator>().IsAtBegin())> : std::true_type
58+
{};
59+
60+
61+
template <typename TIterator>
62+
void
63+
CheckConstructedAtBegin()
64+
{
65+
using ImageType = typename TIterator::ImageType;
66+
using IndexType = typename TIterator::IndexType;
67+
using SizeType = typename TIterator::SizeType;
68+
using RegionType = typename TIterator::RegionType;
69+
70+
const auto image = ImageType::New();
71+
72+
// Use a small image size, so that the unit test won't take a long time.
73+
static constexpr itk::SizeValueType imageSizeValue{ 4 };
74+
75+
image->SetRegions(SizeType::Filled(imageSizeValue));
76+
image->Allocate();
77+
78+
// Check various regions, specified by the following `indexValue` and `sizeValue` combinations:
79+
for (const itk::IndexValueType indexValue : { 0, 1 })
80+
{
81+
for (const auto sizeValue : { itk::SizeValueType{ 1 }, imageSizeValue - 1 })
82+
{
83+
const RegionType imageRegion(IndexType::Filled(indexValue), SizeType::Filled(sizeValue));
84+
85+
const TIterator iterator(image, imageRegion);
86+
TIterator iteratorThatGoesToBegin = iterator;
87+
iteratorThatGoesToBegin.GoToBegin();
88+
EXPECT_EQ(iterator, iteratorThatGoesToBegin);
89+
90+
if constexpr (Has_IsAtBegin<TIterator>::value)
91+
{
92+
// Extra check, using IsAtBegin(), if the iterator has that member function. (Some iterators
93+
// do not have an IsAtBegin() member function, for example, ImageRegionConstIteratorWithIndex.)
94+
EXPECT_TRUE(TIterator(image, imageRegion).IsAtBegin());
95+
}
96+
}
97+
}
98+
}
99+
100+
101+
template <template <typename> typename... TIteratorTemplate>
102+
void
103+
CheckIteratorsConstructedAtBegin()
104+
{
105+
(CheckConstructedAtBegin<TIteratorTemplate<itk::Image<int>>>(), ...);
106+
(CheckConstructedAtBegin<TIteratorTemplate<itk::Image<double, 3>>>(), ...);
107+
}
108+
} // namespace
109+
110+
111+
// Checks that an iterator that is just constructed by `IteratorType(image, region)` is at the begin.
112+
TEST(ImageRegionIterator, IsConstructedAtBegin)
113+
{
114+
CheckIteratorsConstructedAtBegin<itk::ImageConstIterator,
115+
itk::ImageConstIteratorWithIndex,
116+
itk::ImageConstIteratorWithOnlyIndex,
117+
itk::ImageIterator,
118+
itk::ImageIteratorWithIndex,
119+
itk::ImageLinearConstIteratorWithIndex,
120+
itk::ImageLinearIteratorWithIndex,
121+
itk::ImageRegionConstIterator,
122+
itk::ImageRegionConstIteratorWithIndex,
123+
itk::ImageRegionConstIteratorWithOnlyIndex,
124+
itk::ImageRegionExclusionConstIteratorWithIndex,
125+
itk::ImageRegionExclusionIteratorWithIndex,
126+
itk::ImageRegionIterator,
127+
itk::ImageRegionIteratorWithIndex,
128+
itk::ImageRegionReverseConstIterator,
129+
itk::ImageRegionReverseIterator,
130+
itk::ImageReverseConstIterator,
131+
itk::ImageReverseIterator,
132+
itk::ImageScanlineConstIterator,
133+
itk::ImageScanlineIterator,
134+
itk::ImageSliceConstIteratorWithIndex,
135+
itk::ImageSliceIteratorWithIndex>();
136+
}

0 commit comments

Comments
 (0)