|
| 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 | +// The header file to be tested: |
| 20 | +#include "itkComposeImageFilter.h" |
| 21 | +#include "itkVectorIndexSelectionCastImageFilter.h" |
| 22 | +#include "itkMinimumMaximumImageCalculator.h" |
| 23 | +#include "itkExtractImageFilter.h" |
| 24 | +#include "itkPasteImageFilter.h" |
| 25 | + |
| 26 | + |
| 27 | +#include "itkImage.h" |
| 28 | + |
| 29 | +// Google Test header file: |
| 30 | +#include <gtest/gtest.h> |
| 31 | + |
| 32 | +// Standard C++ header files: |
| 33 | +#include <limits> |
| 34 | +#include <random> |
| 35 | + |
| 36 | + |
| 37 | +TEST(ComposeImageFilter, simpleitk_2220) |
| 38 | +{ |
| 39 | + // Test case for SimpleITK issue #2220 |
| 40 | + // Reported data corruption with 4G compose on Windows system |
| 41 | + |
| 42 | + |
| 43 | + const unsigned int size = 400; |
| 44 | + const unsigned int nchannels = 100; |
| 45 | + using ImageType = itk::Image<unsigned char, 3>; |
| 46 | + using VectorImageType = itk::VectorImage<unsigned char, 3>; |
| 47 | + std::vector<ImageType::Pointer> images; |
| 48 | + |
| 49 | + for (unsigned int i = 0; i < nchannels; ++i) |
| 50 | + { |
| 51 | + auto image = ImageType::New(); |
| 52 | + ImageType::RegionType region; |
| 53 | + ImageType::SizeType imageSize = { { size, size, size } }; |
| 54 | + region.SetSize(imageSize); |
| 55 | + image->SetRegions(region); |
| 56 | + image->Allocate(); |
| 57 | + image->FillBuffer(0); |
| 58 | + |
| 59 | + using PasteFilterType = itk::PasteImageFilter<ImageType>; |
| 60 | + PasteFilterType::Pointer paster = PasteFilterType::New(); |
| 61 | + ImageType::IndexType destinationIndex = { { 0, 0, static_cast<int>(i) } }; |
| 62 | + ImageType::SizeType sourceSize = { { size, size, 1 } }; |
| 63 | + paster->SetConstant(i%250); |
| 64 | + paster->SetSourceRegion(sourceSize); |
| 65 | + paster->SetDestinationImage(image); |
| 66 | + paster->SetDestinationIndex(destinationIndex); |
| 67 | + paster->Update(); |
| 68 | + image = paster->GetOutput(); |
| 69 | + images.push_back(image); |
| 70 | + } |
| 71 | + |
| 72 | + using ComposeFilterType = itk::ComposeImageFilter<ImageType>; |
| 73 | + ComposeFilterType::Pointer composeFilter = ComposeFilterType::New(); |
| 74 | + for (unsigned int i = 0; i < nchannels; ++i) |
| 75 | + { |
| 76 | + composeFilter->SetInput(i, images[i]); |
| 77 | + } |
| 78 | + |
| 79 | + composeFilter->Update(); |
| 80 | + VectorImageType::Pointer img = composeFilter->GetOutput(); |
| 81 | + std::cout << "Compose filter executed." << std::endl; |
| 82 | + |
| 83 | + for (unsigned int i = 0; i < nchannels; ++i) |
| 84 | + { |
| 85 | + itk::IndexValueType z = i; |
| 86 | + std::set<unsigned int> values; |
| 87 | + |
| 88 | + for (itk::IndexValueType y = 0; y < size; ++y) |
| 89 | + { |
| 90 | + for (itk::IndexValueType x = 0; x < size; ++x) |
| 91 | + { |
| 92 | + itk::Index<3> idx = { { x, y, z } }; |
| 93 | + |
| 94 | + auto pixel = img->GetPixel(idx); |
| 95 | + if (pixel[i] != i%250) |
| 96 | + { |
| 97 | + values.insert(pixel[i]); |
| 98 | + } |
| 99 | + } |
| 100 | + } |
| 101 | + |
| 102 | + EXPECT_TRUE(values.empty()) << "channel " << i << " has incorrect values;" << std::endl; |
| 103 | + |
| 104 | + if (!values.empty()) |
| 105 | + { |
| 106 | + std::cout << "values: "; |
| 107 | + for (auto v : values) |
| 108 | + { |
| 109 | + std::cout << v << " "; |
| 110 | + } |
| 111 | + std::cout << std::endl; |
| 112 | + } |
| 113 | + |
| 114 | + } |
| 115 | +} |
0 commit comments