File: itkDanielssonDistanceMapImageFilterTest.cxx

package info (click to toggle)
insighttoolkit5 5.4.3-5
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid, trixie
  • size: 704,384 kB
  • sloc: cpp: 783,592; ansic: 628,724; xml: 44,704; fortran: 34,250; python: 22,874; sh: 4,078; pascal: 2,636; lisp: 2,158; makefile: 464; yacc: 328; asm: 205; perl: 203; lex: 146; tcl: 132; javascript: 98; csh: 81
file content (226 lines) | stat: -rw-r--r-- 7,448 bytes parent folder | download
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
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
/*=========================================================================
 *
 *  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 "itkDanielssonDistanceMapImageFilter.h"
#include "itkStdStreamStateSave.h"
#include "itkTestingMacros.h"

int
itkDanielssonDistanceMapImageFilterTest(int, char *[])
{
  // Save the format stream variables for std::cout
  // They will be restored when coutState goes out of scope
  itk::StdStreamStateSave coutState(std::cout);

  std::cout << "Test ITK Danielsson Distance Map" << std::endl << std::endl;
  std::cout << "Compute the distance map of a 9x9 image" << std::endl;
  std::cout << "with a point at (4,4) (value=1)" << std::endl << std::endl;
  std::cout << "with a point at (1,6) (value=2)" << std::endl << std::endl;


  using myImageType2D1 = itk::Image<unsigned char, 2>;
  using myImageType2D2 = itk::Image<float, 2>;

  /* Allocate the 2D image */
  myImageType2D1::SizeType   size2D = { { 9, 9 } };
  myImageType2D1::IndexType  index2D = { { 0, 0 } };
  myImageType2D1::RegionType region2D{ index2D, size2D };

  auto inputImage2D = myImageType2D1::New();
  inputImage2D->SetRegions(region2D);
  inputImage2D->AllocateInitialized();

  // 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.

  index2D[0] = 4;
  index2D[1] = 4;
  inputImage2D->SetPixel(index2D, 1);
  index2D[0] = 1;
  index2D[1] = 6;
  inputImage2D->SetPixel(index2D, 2);

  // Create Danielsson Distance Map filter
  using myFilterType2D = itk::DanielssonDistanceMapImageFilter<myImageType2D1, myImageType2D2>;

  auto filter2D = myFilterType2D::New();

  ITK_EXERCISE_BASIC_OBJECT_METHODS(filter2D, DanielssonDistanceMapImageFilter, ImageToImageFilter);


  filter2D->SetInput(inputImage2D);

  myImageType2D2::Pointer outputDistance2D = filter2D->GetOutput();

  using VoronoiImageType = myFilterType2D::VoronoiImageType;

  VoronoiImageType::Pointer outputVoronoi2D = filter2D->GetVoronoiMap();

  myFilterType2D::VectorImagePointer outputComponents = filter2D->GetVectorDistanceMap();

  ITK_TRY_EXPECT_NO_EXCEPTION(filter2D->Update());


  ShowDistanceMap(outputDistance2D);
  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 Squared Distance functionality
  myImageType2D2::IndexType index;
  index[0] = 0;
  index[1] = 0;
  const double distance1 = outputDistance2D->GetPixel(index);

  bool squaredDistance = true;
  ITK_TEST_SET_GET_BOOLEAN(filter2D, SquaredDistance, squaredDistance);

  ITK_TRY_EXPECT_NO_EXCEPTION(filter2D->Update());


  const double                    distance2 = outputDistance2D->GetPixel(index);
  const myImageType2D2::PixelType epsilon = 1e-5;
  if (itk::Math::abs(distance2 - distance1 * distance1) > epsilon)
  {
    std::cerr << "Error in use of the SetSquaredDistance() method" << std::endl;
    return EXIT_FAILURE;
  }

  std::cout << "Squared Distance Map " << std::endl;
  ShowDistanceMap(outputDistance2D);


  // Test for images with anisotropic spacing
  myImageType2D1::SpacingType anisotropicSpacing;

  anisotropicSpacing[0] = 1.0;
  anisotropicSpacing[1] = 5.0;

  inputImage2D->SetSpacing(anisotropicSpacing);

  inputImage2D->FillBuffer(0);

  index2D[0] = 4;
  index2D[1] = 4;
  inputImage2D->SetPixel(index2D, 1);

  filter2D->SetInput(inputImage2D);

  bool inputIsBinary = true;
  ITK_TEST_SET_GET_BOOLEAN(filter2D, InputIsBinary, inputIsBinary);

  bool useImageSpacing = true;
  ITK_TEST_SET_GET_BOOLEAN(filter2D, UseImageSpacing, useImageSpacing);

  ITK_TRY_EXPECT_NO_EXCEPTION(filter2D->Update());

  index2D[1] = 5;
  auto expectedValue = static_cast<myImageType2D2::PixelType>(anisotropicSpacing[1]);
  expectedValue *= expectedValue;
  myImageType2D2::PixelType pixelValue = filter2D->GetOutput()->GetPixel(index2D);
  if (itk::Math::abs(expectedValue - pixelValue) > epsilon)
  {
    std::cerr << "Error when image spacing is anisotropic." << std::endl;
    std::cerr << "Pixel value was " << pixelValue << ", expected " << expectedValue << std::endl;
    return EXIT_FAILURE;
  }

  ShowDistanceMap(outputDistance2D);

  // Create a large 3D image with a small foreground object.  The foreground is denoted by a pixel value of 0,
  // and the background by a non-zero pixel value.  This will test speedups to the code that ignore all background
  // pixels in the computation since those pixels do not influence the result.

  // Allocate the 3D image
  using ImageType3D = itk::Image<float, 3>;
  ImageType3D::SizeType   size3D = { { 200, 200, 200 } };
  ImageType3D::IndexType  index3D = { { 0, 0 } };
  ImageType3D::RegionType region3D{ index3D, size3D };
  auto                    inputImage3D = ImageType3D::New();
  inputImage3D->SetRegions(region3D);
  inputImage3D->Allocate();
  inputImage3D->FillBuffer(1);

  // Set a few pixels in the middle of the image to 0.  These are the foreground pixels for which the distance will be
  // solved.
  ImageType3D::IndexType foregroundIndex;
  ImageType3D::SizeType  foregroundSize;
  for (unsigned int i = 0; i < 3; ++i)
  {
    foregroundSize[i] = 5;
    foregroundIndex[i] = (size3D[i] / 2) - foregroundSize[i] / 2;
  }
  ImageType3D::RegionType foregroundRegion{ foregroundIndex, foregroundSize };

  itk::ImageRegionIteratorWithIndex<ImageType3D> it3D(inputImage3D, foregroundRegion);
  for (it3D.GoToBegin(); !it3D.IsAtEnd(); ++it3D)
  {
    it3D.Set(0);
  }

  // Create Danielsson Distance Map filter
  using myFilterType3D = itk::DanielssonDistanceMapImageFilter<ImageType3D, ImageType3D>;

  auto filter3D = myFilterType3D::New();

  filter3D->SetInput(inputImage3D);

  ITK_TRY_EXPECT_NO_EXCEPTION(filter3D->Update());


  std::cout << "Test finished." << std::endl;
  return EXIT_SUCCESS;
}