1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166
|
/*=========================================================================
*
* Copyright NumFOCUS
*
* Licensed under the Apache License, Version 2.0 (the "License");
* you may not use this file except in compliance with the License.
* You may obtain a copy of the License at
*
* https://www.apache.org/licenses/LICENSE-2.0.txt
*
* Unless required by applicable law or agreed to in writing, software
* distributed under the License is distributed on an "AS IS" BASIS,
* WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
* See the License for the specific language governing permissions and
* limitations under the License.
*
*=========================================================================*/
#include "itkShowDistanceMap.h"
#include "itkSignedDanielssonDistanceMapImageFilter.h"
#include "itkStdStreamStateSave.h"
void
test(int);
int
itkSignedDanielssonDistanceMapImageFilterTest(int, char *[])
{
test(1); // Test with a 9x9 square, with a 5x5 subsquare in the middle ON
test(0); // Test with 2 points.. same test and results as
// DanielssonDistanceMap
return EXIT_SUCCESS;
}
void
test(int testIdx)
{
// Save the format stream variables for std::cout
// They will be restored when coutState goes out of scope
// scope.
itk::StdStreamStateSave coutState(std::cout);
constexpr unsigned int Dimension = 2;
using PixelType = float;
using myImageType2D1 = itk::Image<unsigned char, Dimension>;
using myImageType2D2 = itk::Image<PixelType, Dimension>;
/* TEST 1: For a point image, SignedDaniessonDistanceMapImageFilter should
* give the same output as DaniessonDistanceMapImageFilter */
/* Allocate the 2D image */
myImageType2D1::SizeType size2D;
size2D.Fill(9);
myImageType2D1::IndexType index2D;
index2D.Fill(0);
myImageType2D1::RegionType region2D{ index2D, size2D };
auto inputImage2D = myImageType2D1::New();
inputImage2D->SetRegions(region2D);
inputImage2D->AllocateInitialized();
if (!testIdx)
{
std::cout << "Compute with a 9x9 image, pixels (4,4) and (1,6) set to ON."
<< " This subtest is the same as the DanielssonDistanceMapTest "
<< "and should yield the same results." << std::endl
<< std::endl;
/* Set pixel (4,4) with the value 1
* and pixel (1,6) with the value 2
* The Danielsson Distance is performed for each pixel with a value > 0
* The ClosestPoints computation is based on the value of the pixel.
* Test 1 contains two isolated points and should return the same result
* as if the filter SignedDanielssonDistanceMapImageFilter was changed
* to DanielssonDistanceMapImageFilter
*/
index2D[0] = 4;
index2D[1] = 4;
inputImage2D->SetPixel(index2D, 1);
index2D[0] = 1;
index2D[1] = 6;
inputImage2D->SetPixel(index2D, 2);
}
else
{
std::cout << "Compute with a 9x9 image, with a 5x5 square at the center set to ON." << std::endl << std::endl;
// Test the signed Danielsson Output for the a 5x5 square in a 9x9 image
int i, j;
for (i = 2; i <= 6; ++i)
{
for (j = 2; j <= 6; ++j)
{
index2D[0] = i;
index2D[1] = j;
inputImage2D->SetPixel(index2D, 1);
}
}
}
/* Create Danielsson Distance Map filter */
using myFilterType2D = itk::SignedDanielssonDistanceMapImageFilter<myImageType2D1, myImageType2D2>;
auto filter2D = myFilterType2D::New();
filter2D->SetInput(inputImage2D);
myImageType2D2::Pointer outputDistance2D = filter2D->GetOutput();
using VoronoiImageType = myFilterType2D::VoronoiImageType;
VoronoiImageType::Pointer outputVoronoi2D = filter2D->GetVoronoiMap();
myFilterType2D::VectorImageType::Pointer outputComponents = filter2D->GetVectorDistanceMap();
filter2D->Update();
ShowDistanceMap(outputDistance2D);
std::cout << std::endl << std::endl;
std::cout << "Voronoi Map Image 2D" << std::endl << std::endl;
ShowDistanceMap(outputVoronoi2D);
/* Show VectorsComponents Points map */
std::cout << std::endl << std::endl;
std::cout << "Components Map Image 2D" << std::endl << std::endl;
itk::ImageSliceConstIteratorWithIndex<myFilterType2D::VectorImageType> it2D4(outputComponents,
outputComponents->GetRequestedRegion());
it2D4.SetFirstDirection(0);
it2D4.SetSecondDirection(1);
while (!it2D4.IsAtEnd())
{
while (!it2D4.IsAtEndOfSlice())
{
while (!it2D4.IsAtEndOfLine())
{
std::cout << '[';
for (unsigned int i = 0; i < 2; ++i)
{
std::cout << it2D4.Get()[i];
if (i == 0)
{
std::cout << ',';
}
}
std::cout << ']';
std::cout << '\t';
++it2D4;
}
std::cout << std::endl;
it2D4.NextLine();
}
it2D4.NextSlice();
}
/* Test Distance functionality */
filter2D->Update();
std::cout << "Distance Map " << std::endl;
ShowDistanceMap(outputDistance2D);
}
|