diff --git a/Modules/ImageStatistics/Testing/CMakeLists.txt b/Modules/ImageStatistics/Testing/CMakeLists.txt index ec552a7be1..45387426df 100644 --- a/Modules/ImageStatistics/Testing/CMakeLists.txt +++ b/Modules/ImageStatistics/Testing/CMakeLists.txt @@ -1,5 +1,34 @@ MITK_CREATE_MODULE_TESTS(EXTRA_DEPENDS DICOMReader) -mitkAddCustomModuleTest(mitkImageStatisticsHotspotTest_Case1 mitkImageStatisticsHotspotTest ${CMAKE_CURRENT_SOURCE_DIR}/Data/Hotspot_Case1.xml) - mitkAddCustomModuleTest(mitkRoiMeasurementsTests mitkRoiMeasurementsTest ${MITK_DATA_DIR}/ImageStatisticsTestData/) + +file(GLOB allHotSpotTests RELATIVE "${CMAKE_CURRENT_SOURCE_DIR}/Data" "${CMAKE_CURRENT_SOURCE_DIR}/Data/Hotspot_Case*.xml") +foreach(testcase ${allHotSpotTests}) + string(REGEX REPLACE "[^a-zA-Z0-9_]" "_" testcaseName ${testcase}) + mitkAddCustomModuleTest(mitkImageStatisticsHotspotTest_${testcaseName} mitkImageStatisticsHotspotTest ${CMAKE_CURRENT_SOURCE_DIR}/Data/${testcase}) +endforeach() + +# +# The following lines may be activated to generate new test cases for mitkImageStatisticsHotspotTest. +# Test cases are generated by mitkMultiGaussianTest. All .xml files in Data/TestGeneration/Input will +# be processed and transformed into new .xml files containing statistics in Data/TestGeneration/Output. +# +if (false) + set(testInputDir ${CMAKE_CURRENT_SOURCE_DIR}/Data/TestGeneration/Input) + set(testOutputDir ${CMAKE_CURRENT_SOURCE_DIR}/Data/TestGeneration/Output) + file(GLOB testcasesToGenerate RELATIVE "${testInputDir}" "${testInputDir}/*.xml") + + if (NOT EXISTS ${testOutputDir}) + file(MAKE_DIRECTORY ${testOutputDir}) + endif() + + foreach(testinput ${testcasesToGenerate}) + string(REGEX REPLACE "[^a-zA-Z0-9_]\\+" "_" testcaseName ${testinput}) + string(REGEX REPLACE "\\.xml" "" testoutput ${testinput}) + message("Generate hotspot test case '${testinput}'. Output in '${testoutput}.xml' and '${testoutput}.nrrd'") + mitkAddCustomModuleTest(mitkMultiGaussianTest_${testcaseName} + mitkMultiGaussianTest + ${testOutputDir}/${testoutput} + ${testInputDir}/${testinput}) + endforeach() +endif() diff --git a/Modules/ImageStatistics/Testing/files.cmake b/Modules/ImageStatistics/Testing/files.cmake index e583626672..def2ec655c 100644 --- a/Modules/ImageStatistics/Testing/files.cmake +++ b/Modules/ImageStatistics/Testing/files.cmake @@ -1,11 +1,11 @@ set(MODULE_TESTS mitkImageStatisticsCalculatorTest.cpp mitkPointSetStatisticsCalculatorTest.cpp mitkPointSetDifferenceStatisticsCalculatorTest.cpp - #mitkMultiGaussianTest.cpp # TODO: activate test if bug 15821 has been fixed - ) +) set(MODULE_CUSTOM_TESTS mitkRoiMeasurementsTest.cpp mitkImageStatisticsHotspotTest.cpp +# mitkMultiGaussianTest.cpp # TODO: activate test to generate new test cases for mitkImageStatisticsHotspotTest ) diff --git a/Modules/ImageStatistics/Testing/mitkImageStatisticsHotspotTest.cpp b/Modules/ImageStatistics/Testing/mitkImageStatisticsHotspotTest.cpp index 799a7926c1..9c54b2bdcf 100644 --- a/Modules/ImageStatistics/Testing/mitkImageStatisticsHotspotTest.cpp +++ b/Modules/ImageStatistics/Testing/mitkImageStatisticsHotspotTest.cpp @@ -1,592 +1,618 @@ /*=================================================================== The Medical Imaging Interaction Toolkit (MITK) Copyright (c) German Cancer Research Center, Division of Medical and Biological Informatics. All rights reserved. This software is distributed WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See LICENSE.txt or http://www.mitk.org for details. ===================================================================*/ #include "mitkImageStatisticsCalculator.h" #include "itkMultiGaussianImageSource.h" #include "mitkTestingMacros.h" #include #include #include #include /** \section hotspotCalculationTestCases Testcases To see the different Hotspot-Testcases have a look at the \ref hotspottestdoc. */ struct mitkImageStatisticsHotspotTestClass { /** \brief Test parameters for one test case. Describes all aspects of a single test case: - parameters to generate a test image - parameters of a ROI that describes where to calculate statistics - expected statistics results */ struct Parameters { public: // XML-Tag /** \brief XML-Tag "image-rows": size of x-dimension */ int m_ImageRows; /** \brief XML-Tag "image-columns": size of y-dimension */ int m_ImageColumns; /** \brief XML-Tag "image-slices": size of z-dimension */ int m_ImageSlices; /** \brief XML-Tag "numberOfGaussians": number of used gauss-functions */ int m_NumberOfGaussian; /** \brief XML-Tags "spacingX", "spacingY", "spacingZ": spacing of image in every direction */ float m_Spacing[3]; /** \brief XML-Tag "entireHotSpotInImage" */ unsigned int m_EntireHotspotInImage; // XML-Tag /** \brief XML-Tag "centerIndexX: gaussian parameter*/ std::vector m_CenterX; /** \brief XML-Tag "centerIndexY: gaussian parameter */ std::vector m_CenterY; /** \brief XML-Tag "centerIndexZ: gaussian parameter */ std::vector m_CenterZ; /** \brief XML-Tag "deviationX: gaussian parameter */ std::vector m_SigmaX; /** \brief XML-Tag "deviationY: gaussian parameter */ std::vector m_SigmaY; /** \brief XML-Tag "deviationZ: gaussian parameter */ std::vector m_SigmaZ; /** \brief XML-Tag "altitude: gaussian parameter */ std::vector m_Altitude; // XML-Tag /** \brief XML-Tag "numberOfLabels": number of different labels which appear in the mask */ unsigned int m_NumberOfLabels; /** \brief XML-Tag "hotspotRadiusInMM": radius of hotspot */ float m_HotspotRadiusInMM; // XML-Tag /** \brief XML-Tag "maximumSizeX": maximum position of ROI in x-dimension */ vnl_vector m_MaxSizeX; /** \brief XML-Tag "minimumSizeX": minimum position of ROI in x-dimension */ vnl_vector m_MinSizeX; /** \brief XML-Tag "maximumSizeX": maximum position of ROI in y-dimension */ vnl_vector m_MaxSizeY; /** \brief XML-Tag "minimumSizeX": minimum position of ROI in y-dimension */ vnl_vector m_MinSizeY; /** \brief XML-Tag "maximumSizeX": maximum position of ROI in z-dimension */ vnl_vector m_MaxSizeZ; /** \brief XML-Tag "minimumSizeX": minimum position of ROI in z-dimension */ vnl_vector m_MinSizeZ; /** \brief XML-Tag "label": value of label */ vnl_vector m_Label; //XML-Tag /** \brief XML-Tag "minimum": minimum inside hotspot */ vnl_vector m_HotspotMin; /** \brief XML-Tag "maximum": maximum inside hotspot */ vnl_vector m_HotspotMax; /** \brief XML-Tag "mean": mean value of hotspot */ vnl_vector m_HotspotMean; /** \brief XML-Tag "maximumIndexX": x-coordinate of maximum-location inside hotspot */ vnl_vector m_HotspotMaxIndexX; /** \brief XML-Tag "maximumIndexX": y-coordinate of maximum-location inside hotspot */ vnl_vector m_HotspotMaxIndexY; /** \brief XML-Tag "maximumIndexX": z-coordinate of maximum-location inside hotspot */ vnl_vector m_HotspotMaxIndexZ; /** \brief XML-Tag "maximumIndexX": x-coordinate of maximum-location inside hotspot */ vnl_vector m_HotspotMinIndexX; /** \brief XML-Tag "maximumIndexX": y-coordinate of maximum-location inside hotspot */ vnl_vector m_HotspotMinIndexY; /** \brief XML-Tag "maximumIndexX": z-coordinate of maximum-location inside hotspot */ vnl_vector m_HotspotMinIndexZ; /** \brief XML-Tag "maximumIndexX": x-coordinate of hotspot-location */ vnl_vector m_HotspotIndexX; /** \brief XML-Tag "maximumIndexX": y-coordinate of hotspot-location */ vnl_vector m_HotspotIndexY; /** \brief XML-Tag "maximumIndexX": z-coordinate of hotspot-location */ vnl_vector m_HotspotIndexZ; }; /** \brief Find/Convert integer attribute in itk::DOMNode. */ static int GetIntegerAttribute(itk::DOMNode* domNode, const std::string& tag) { assert(domNode); MITK_TEST_CONDITION_REQUIRED( domNode->HasAttribute(tag), "Tag '" << tag << "' is defined in test parameters" ); std::string attributeValue = domNode->GetAttribute(tag); int resultValue; try { //MITK_TEST_OUTPUT( << "Converting tag value '" << attributeValue << "' for tag '" << tag << "' to integer"); std::stringstream(attributeValue) >> resultValue; return resultValue; } catch(std::exception& e) { MITK_TEST_CONDITION_REQUIRED(false, "Convert tag value '" << attributeValue << "' for tag '" << tag << "' to integer"); return 0; // just to satisfy compiler } } /** \brief Find/Convert double attribute in itk::DOMNode. */ static double GetDoubleAttribute(itk::DOMNode* domNode, const std::string& tag) { assert(domNode); MITK_TEST_CONDITION_REQUIRED( domNode->HasAttribute(tag), "Tag '" << tag << "' is defined in test parameters" ); std::string attributeValue = domNode->GetAttribute(tag); double resultValue; try { //MITK_TEST_OUTPUT( << "Converting tag value '" << attributeValue << "' for tag '" << tag << "' to double"); std::stringstream(attributeValue) >> resultValue; return resultValue; } catch(std::exception& e) { MITK_TEST_CONDITION_REQUIRED(false, "Convert tag value '" << attributeValue << "' for tag '" << tag << "' to double"); return 0.0; // just to satisfy compiler } } /** \brief Read XML file describing the test parameters. Reads XML file given in first commandline parameter in order to construct a Parameters structure. The XML file should be structurs as the following example, i.e. we describe the three test aspects of Parameters in four different tags, with all the details described as tag attributes. */ /** \verbatim \endverbatim */ static Parameters ParseParameters(int argc, char* argv[]) { MITK_TEST_CONDITION_REQUIRED(argc == 2, "Test is invoked with exactly 1 parameter (XML parameters file)"); MITK_INFO << "Reading parameters from file '" << argv[1] << "'"; std::string filename = argv[1]; Parameters result; itk::DOMNodeXMLReader::Pointer xmlReader = itk::DOMNodeXMLReader::New(); xmlReader->SetFileName( filename ); try { xmlReader->Update(); itk::DOMNode::Pointer domRoot = xmlReader->GetOutput(); typedef std::vector NodeList; NodeList testimages; domRoot->GetChildren("testimage", testimages); MITK_TEST_CONDITION_REQUIRED( testimages.size() == 1, "One test image defined" ) itk::DOMNode* testimage = testimages[0]; result.m_ImageRows = GetIntegerAttribute( testimage, "image-rows" ); result.m_ImageColumns = GetIntegerAttribute( testimage, "image-columns" ); result.m_ImageSlices = GetIntegerAttribute( testimage, "image-slices" ); result.m_NumberOfGaussian = GetIntegerAttribute( testimage, "numberOfGaussians" ); result.m_Spacing[0] = GetDoubleAttribute(testimage, "spacingX"); result.m_Spacing[1] = GetDoubleAttribute(testimage, "spacingY"); result.m_Spacing[2] = GetDoubleAttribute(testimage, "spacingZ"); result.m_EntireHotspotInImage = GetIntegerAttribute( testimage, "entireHotSpotInImage" ); MITK_TEST_OUTPUT( << "Read size parameters (x,y,z): " << result.m_ImageRows << "," << result.m_ImageColumns << "," << result.m_ImageSlices); MITK_TEST_OUTPUT( << "Read spacing parameters (x,y,z): " << result.m_Spacing[0] << "," << result.m_Spacing[1] << "," << result.m_Spacing[2]); NodeList gaussians; testimage->GetChildren("gaussian", gaussians); MITK_TEST_CONDITION_REQUIRED( gaussians.size() >= 1, "At least one gaussian is defined" ) result.m_CenterX.resize(result.m_NumberOfGaussian); result.m_CenterY.resize(result.m_NumberOfGaussian); result.m_CenterZ.resize(result.m_NumberOfGaussian); result.m_SigmaX.resize(result.m_NumberOfGaussian); result.m_SigmaY.resize(result.m_NumberOfGaussian); result.m_SigmaZ.resize(result.m_NumberOfGaussian); result.m_Altitude.resize(result.m_NumberOfGaussian); for(int i = 0; i < result.m_NumberOfGaussian ; ++i) { itk::DOMNode* gaussian = gaussians[i]; result.m_CenterX[i] = GetIntegerAttribute(gaussian, "centerIndexX"); result.m_CenterY[i] = GetIntegerAttribute(gaussian, "centerIndexY"); result.m_CenterZ[i] = GetIntegerAttribute(gaussian, "centerIndexZ"); result.m_SigmaX[i] = GetIntegerAttribute(gaussian, "deviationX"); result.m_SigmaY[i] = GetIntegerAttribute(gaussian, "deviationY"); result.m_SigmaZ[i] = GetIntegerAttribute(gaussian, "deviationZ"); result.m_Altitude[i] = GetIntegerAttribute(gaussian, "altitude"); result.m_CenterX[i] *= result.m_Spacing[0]; result.m_CenterY[i] *= result.m_Spacing[1]; result.m_CenterZ[i] *= result.m_Spacing[2]; result.m_SigmaX[i] *= result.m_Spacing[0]; result.m_SigmaY[i] *= result.m_Spacing[1]; result.m_SigmaZ[i] *= result.m_Spacing[2]; } NodeList segmentations; domRoot->GetChildren("segmentation", segmentations); MITK_TEST_CONDITION_REQUIRED( segmentations.size() == 1, "One segmentation defined"); itk::DOMNode* segmentation = segmentations[0]; result.m_NumberOfLabels = GetIntegerAttribute(segmentation, "numberOfLabels"); result.m_HotspotRadiusInMM = GetDoubleAttribute(segmentation, "hotspotRadiusInMM"); // read ROI parameters, fill result structure NodeList rois; segmentation->GetChildren("roi", rois); MITK_TEST_CONDITION_REQUIRED( rois.size() >= 1, "At least one ROI defined" ) result.m_MaxSizeX.set_size(result.m_NumberOfLabels); result.m_MinSizeX.set_size(result.m_NumberOfLabels); result.m_MaxSizeY.set_size(result.m_NumberOfLabels); result.m_MinSizeY.set_size(result.m_NumberOfLabels); result.m_MaxSizeZ.set_size(result.m_NumberOfLabels); result.m_MinSizeZ.set_size(result.m_NumberOfLabels); result.m_Label.set_size(result.m_NumberOfLabels); for(int i = 0; i < rois.size(); ++i) { result.m_MaxSizeX[i] = GetIntegerAttribute(rois[i], "maximumSizeX"); result.m_MinSizeX[i] = GetIntegerAttribute(rois[i], "minimumSizeX"); result.m_MaxSizeY[i] = GetIntegerAttribute(rois[i], "maximumSizeY"); result.m_MinSizeY[i] = GetIntegerAttribute(rois[i], "minimumSizeY"); result.m_MaxSizeZ[i] = GetIntegerAttribute(rois[i], "maximumSizeZ"); result.m_MinSizeZ[i] = GetIntegerAttribute(rois[i], "minimumSizeZ"); result.m_Label[i] = GetIntegerAttribute(rois[i], "label"); } // read statistic parameters, fill result structure NodeList statistics; domRoot->GetChildren("statistic", statistics); MITK_TEST_CONDITION_REQUIRED( statistics.size() >= 1 , "At least one statistic defined" ) MITK_TEST_CONDITION_REQUIRED( statistics.size() == rois.size(), "Same number of rois and corresponding statistics defined"); result.m_HotspotMin.set_size(statistics.size()); result.m_HotspotMax.set_size(statistics.size()); result.m_HotspotMean.set_size(statistics.size()); result.m_HotspotMinIndexX.set_size(statistics.size()); result.m_HotspotMinIndexY.set_size(statistics.size()); result.m_HotspotMinIndexZ.set_size(statistics.size()); result.m_HotspotMaxIndexX.set_size(statistics.size()); result.m_HotspotMaxIndexY.set_size(statistics.size()); result.m_HotspotMaxIndexZ.set_size(statistics.size()); result.m_HotspotIndexX.set_size(statistics.size()); result.m_HotspotIndexY.set_size(statistics.size()); result.m_HotspotIndexZ.set_size(statistics.size()); for(int i = 0; i < statistics.size(); ++i) { result.m_HotspotMin[i] = GetDoubleAttribute(statistics[i], "minimum"); result.m_HotspotMax[i] = GetDoubleAttribute(statistics[i], "maximum"); result.m_HotspotMean[i] = GetDoubleAttribute(statistics[i], "mean"); result.m_HotspotMinIndexX[i] = GetIntegerAttribute(statistics[i], "minimumIndexX"); result.m_HotspotMinIndexY[i] = GetIntegerAttribute(statistics[i], "minimumIndexY"); result.m_HotspotMinIndexZ[i] = GetIntegerAttribute(statistics[i], "minimumIndexZ"); result.m_HotspotMaxIndexX[i] = GetIntegerAttribute(statistics[i], "maximumIndexX"); result.m_HotspotMaxIndexY[i] = GetIntegerAttribute(statistics[i], "maximumIndexY"); result.m_HotspotMaxIndexZ[i] = GetIntegerAttribute(statistics[i], "maximumIndexZ"); result.m_HotspotIndexX[i] = GetIntegerAttribute(statistics[i], "hotspotIndexX"); result.m_HotspotIndexY[i] = GetIntegerAttribute(statistics[i], "hotspotIndexY"); result.m_HotspotIndexZ[i] = GetIntegerAttribute(statistics[i], "hotspotIndexZ"); } return result; } catch (std::exception& e) { MITK_TEST_CONDITION_REQUIRED(false, "Reading test parameters from XML file. Error message: " << e.what()); } } /** \brief Generate an image that contains a couple of 3D gaussian distributions. Uses the given parameters to produce a test image using class MultiGaussianImageSource. */ static mitk::Image::Pointer BuildTestImage(const Parameters& testParameters) { mitk::Image::Pointer result; typedef double PixelType; const unsigned int Dimension = 3; typedef itk::Image ImageType; ImageType::Pointer image = ImageType::New(); typedef itk::MultiGaussianImageSource< ImageType > MultiGaussianImageSource; MultiGaussianImageSource::Pointer gaussianGenerator = MultiGaussianImageSource::New(); ImageType::SizeValueType size[3]; size[0] = testParameters.m_ImageColumns; size[1] = testParameters.m_ImageRows; size[2] = testParameters.m_ImageSlices; itk::MultiGaussianImageSource::VectorType centerXVec, centerYVec, centerZVec, sigmaXVec, sigmaYVec, sigmaZVec, altitudeVec; for(int i = 0; i < testParameters.m_NumberOfGaussian; ++i) { centerXVec.push_back(testParameters.m_CenterX[i]); centerYVec.push_back(testParameters.m_CenterY[i]); centerZVec.push_back(testParameters.m_CenterZ[i]); sigmaXVec.push_back(testParameters.m_SigmaX[i]); sigmaYVec.push_back(testParameters.m_SigmaY[i]); sigmaZVec.push_back(testParameters.m_SigmaZ[i]); altitudeVec.push_back(testParameters.m_Altitude[i]); } ImageType::SpacingType spacing; for(int i = 0; i < Dimension; ++i) spacing[i] = testParameters.m_Spacing[i]; gaussianGenerator->SetSize( size ); gaussianGenerator->SetSpacing( spacing ); gaussianGenerator->SetRadius(testParameters.m_HotspotRadiusInMM); gaussianGenerator->SetNumberOfGausssians(testParameters.m_NumberOfGaussian); gaussianGenerator->AddGaussian(centerXVec, centerYVec, centerZVec, sigmaXVec, sigmaYVec, sigmaZVec, altitudeVec); gaussianGenerator->Update(); image = gaussianGenerator->GetOutput(); mitk::CastToMitkImage(image, result); return result; } /** \brief Calculates hotspot statistics for given test image and ROI parameters. Uses ImageStatisticsCalculator to find a hotspot in a defined ROI within the given image. */ static mitk::ImageStatisticsCalculator::Statistics CalculateStatistics(mitk::Image* image, const Parameters& testParameters, unsigned int label) { mitk::ImageStatisticsCalculator::Statistics result; const unsigned int Dimension = 3; typedef itk::Image MaskImageType; MaskImageType::Pointer mask = MaskImageType::New(); MaskImageType::SizeType size; MaskImageType::SpacingType spacing; MaskImageType::IndexType start; mitk::ImageStatisticsCalculator::Pointer statisticsCalculator = mitk::ImageStatisticsCalculator::New(); statisticsCalculator->SetImage(image); mitk::Image::Pointer mitkMaskImage; if((testParameters.m_MaxSizeX[label] > testParameters.m_MinSizeX[label] && testParameters.m_MinSizeX[label] >= 0) && (testParameters.m_MaxSizeY[label] > testParameters.m_MinSizeY[label] && testParameters.m_MinSizeY[label] >= 0) && (testParameters.m_MaxSizeZ[label] > testParameters.m_MinSizeZ[label] && testParameters.m_MinSizeZ[label] >= 0)) { for(int i = 0; i < Dimension; ++i) { start[i] = 0; spacing[i] = testParameters.m_Spacing[i]; } size[0] = testParameters.m_ImageColumns; size[1] = testParameters.m_ImageRows; size[2] = testParameters.m_ImageSlices; MaskImageType::RegionType region; region.SetIndex(start); region.SetSize(size); mask->SetSpacing(spacing); mask->SetRegions(region); mask->Allocate(); typedef itk::ImageRegionIteratorWithIndex MaskImageIteratorType; MaskImageIteratorType maskIt(mask, region); for(maskIt.GoToBegin(); !maskIt.IsAtEnd(); ++maskIt) { maskIt.Set(0); } for(int i = 0; i < testParameters.m_NumberOfLabels; ++i) { for(maskIt.GoToBegin(); !maskIt.IsAtEnd(); ++maskIt) { MaskImageType::IndexType index = maskIt.GetIndex(); if((index[0] >= testParameters.m_MinSizeX[i] && index[0] < testParameters.m_MaxSizeX[i] ) && (index[1] >= testParameters.m_MinSizeY[i] && index[1] < testParameters.m_MaxSizeY[i] ) && (index[2] >= testParameters.m_MinSizeZ[i] && index[2] < testParameters.m_MaxSizeZ[i] )) { maskIt.Set(testParameters.m_Label[i]); } } } MITK_DEBUG << "Masking mode has set to image"; mitk::CastToMitkImage(mask, mitkMaskImage); statisticsCalculator->SetImageMask(mitkMaskImage); statisticsCalculator->SetMaskingModeToImage(); } else { MITK_DEBUG << "Masking mode has set to none"; statisticsCalculator->SetMaskingModeToNone(); } statisticsCalculator->SetHotspotRadiusInMM(testParameters.m_HotspotRadiusInMM); statisticsCalculator->SetCalculateHotspot(true); if(testParameters.m_EntireHotspotInImage == 1) { MITK_INFO << "Hotspot must be completly inside image"; statisticsCalculator->SetHotspotMustBeCompletlyInsideImage(true); } else { MITK_INFO << "Hotspot must not be completly inside image"; statisticsCalculator->SetHotspotMustBeCompletlyInsideImage(false); } statisticsCalculator->ComputeStatistics(); result = statisticsCalculator->GetStatistics(0, label); return result; } + static void ValidateStatisticsItem(const std::string& label, double testvalue, double reference, double tolerance) + { + double diff = ::fabs(reference - testvalue); + MITK_TEST_CONDITION( diff < tolerance, "'" << label << "' value close enough to reference value " + "(value=" << testvalue << + ", reference=" << reference << + ", diff=" << diff << ")" ); + } + + static void ValidateStatisticsItem(const std::string& label, const vnl_vector& testvalue, const vnl_vector& reference) + { + double diffX = testvalue[0] - reference[0]; + double diffY = testvalue[1] - reference[1]; + double diffZ = testvalue[2] - reference[2]; + + std::stringstream testPosition; + testPosition << testvalue[0] << "," << testvalue[1] << "," << testvalue[2]; + std::stringstream referencePosition; + referencePosition << reference[0] << "," << reference[1] << "," << reference[2]; + MITK_TEST_CONDITION( diffX < mitk::eps && diffY < mitk::eps && diffZ < mitk::eps, + "'" << label << "' close enough to reference value " << + "(value=[" << testPosition.str() << "]," << + " reference=[" << referencePosition.str() << "]"); + } + + /** \brief Compares calculated against actual statistics values. Checks validness of all statistics aspects. Lets test fail if any aspect is not sufficiently equal. */ static void ValidateStatistics(const mitk::ImageStatisticsCalculator::Statistics& statistics, const Parameters& testParameters, unsigned int label) { // check all expected test result against actual results double eps = 1.6; - MITK_TEST_CONDITION( ::fabs(testParameters.m_HotspotMean[label] - statistics.GetHotspotStatistics().GetMean() ) < eps, "Mean value of hotspot in XML-File: " << testParameters.m_HotspotMean[label] << " (Mean value of hotspot calculated in mitkImageStatisticsCalculator: " << statistics.GetHotspotStatistics().GetMean() << ")" ); - MITK_TEST_CONDITION( ::fabs(testParameters.m_HotspotMax[label]- statistics.GetHotspotStatistics().GetMax() ) < eps, "Maximum of hotspot in XML-File: " << testParameters.m_HotspotMax[label] << " (Maximum of hotspot calculated in mitkImageStatisticsCalculator: " << statistics.GetHotspotStatistics().GetMax() << ")" ); - MITK_TEST_CONDITION( ::fabs(testParameters.m_HotspotMin[label] - statistics.GetHotspotStatistics().GetMin() ) < eps, "Minimum of hotspot in XML-File: " << testParameters.m_HotspotMin[label] << " (Minimum of hotspot calculated in mitkImageStatisticsCalculator: " << statistics.GetHotspotStatistics().GetMin() << ")" ); - - MITK_TEST_CONDITION( statistics.GetHotspotStatistics().GetHotspotIndex()[0] == testParameters.m_HotspotIndexX[label] && - statistics.GetHotspotStatistics().GetHotspotIndex()[1] == testParameters.m_HotspotIndexY[label] && - statistics.GetHotspotStatistics().GetHotspotIndex()[2] == testParameters.m_HotspotIndexZ[label] , - "Index of hotspot in XML-File: " << testParameters.m_HotspotIndexX[label] << " " << testParameters.m_HotspotIndexY[label] << " " << testParameters.m_HotspotIndexZ[label] - << " (Index of hotspot calculated in mitkImageStatisticsCalculator: " << statistics.GetHotspotStatistics().GetHotspotIndex() << ")" ); - - MITK_TEST_CONDITION( statistics.GetHotspotStatistics().GetMaxIndex()[0] == testParameters.m_HotspotMaxIndexX[label] && - statistics.GetHotspotStatistics().GetMaxIndex()[1] == testParameters.m_HotspotMaxIndexY[label] && - statistics.GetHotspotStatistics().GetMaxIndex()[2] == testParameters.m_HotspotMaxIndexZ[label] , - "Index of hotspot-maximum in XML-File: " << testParameters.m_HotspotMaxIndexX[label] << " " << testParameters.m_HotspotMaxIndexY[label] << " " << testParameters.m_HotspotMaxIndexZ[label] - << " (Index of hotspot-maximum calculated in mitkImageStatisticsCalculator: " << statistics.GetHotspotStatistics().GetMaxIndex() << ")" ); - - MITK_TEST_CONDITION( statistics.GetHotspotStatistics().GetMinIndex()[0] == testParameters.m_HotspotMinIndexX[label] && - statistics.GetHotspotStatistics().GetMinIndex()[1] == testParameters.m_HotspotMinIndexY[label] && - statistics.GetHotspotStatistics().GetMinIndex()[2] == testParameters.m_HotspotMinIndexZ[label] , - "Index of hotspot-minimum in XML-File: " << testParameters.m_HotspotMinIndexX[label] << " " << testParameters.m_HotspotMinIndexY[label] << " " << testParameters.m_HotspotMinIndexZ[label] - << " (Index of hotspot-minimum calculated in mitkImageStatisticsCalculator: " << statistics.GetHotspotStatistics().GetMinIndex() << ")" ); + ValidateStatisticsItem("Hotspot mean", statistics.GetHotspotStatistics().GetMean(), testParameters.m_HotspotMean[label], eps); + ValidateStatisticsItem("Hotspot maximum", statistics.GetHotspotStatistics().GetMax(), testParameters.m_HotspotMax[label], eps); + ValidateStatisticsItem("Hotspot minimum", statistics.GetHotspotStatistics().GetMin(), testParameters.m_HotspotMin[label], eps); + + vnl_vector referenceHotspotCenterIndex; referenceHotspotCenterIndex.set_size(3); + referenceHotspotCenterIndex[0] = testParameters.m_HotspotIndexX[label]; + referenceHotspotCenterIndex[1] = testParameters.m_HotspotIndexY[label]; + referenceHotspotCenterIndex[2] = testParameters.m_HotspotIndexZ[label]; + ValidateStatisticsItem("Hotspot center position", statistics.GetHotspotStatistics().GetHotspotIndex(), referenceHotspotCenterIndex); + + vnl_vector referenceHotspotMaxIndex; referenceHotspotMaxIndex.set_size(3); + referenceHotspotMaxIndex[0] = testParameters.m_HotspotMaxIndexX[label]; + referenceHotspotMaxIndex[1] = testParameters.m_HotspotMaxIndexY[label]; + referenceHotspotMaxIndex[2] = testParameters.m_HotspotMaxIndexZ[label]; + ValidateStatisticsItem("Hotspot maximum position", statistics.GetHotspotStatistics().GetHotspotIndex(), referenceHotspotMaxIndex); + + vnl_vector referenceHotspotMinIndex; referenceHotspotMinIndex.set_size(3); + referenceHotspotMinIndex[0] = testParameters.m_HotspotMinIndexX[label]; + referenceHotspotMinIndex[1] = testParameters.m_HotspotMinIndexY[label]; + referenceHotspotMinIndex[2] = testParameters.m_HotspotMinIndexZ[label]; + ValidateStatisticsItem("Hotspot minimum position", statistics.GetHotspotStatistics().GetHotspotIndex(), referenceHotspotMinIndex); } }; /** \brief Verifies that hotspot statistics part of ImageStatisticsCalculator. The test reads parameters from an XML-file to generate a test-image, calculates the hotspot statistics of the image and checks if the calculated statistics are the same as the specified values of the XML-file. */ int mitkImageStatisticsHotspotTest(int argc, char* argv[]) { MITK_TEST_BEGIN("mitkImageStatisticsHotspotTest") try { mitkImageStatisticsHotspotTestClass::Parameters parameters = mitkImageStatisticsHotspotTestClass::ParseParameters(argc,argv); mitk::Image::Pointer image = mitkImageStatisticsHotspotTestClass::BuildTestImage(parameters); MITK_TEST_CONDITION_REQUIRED( image.IsNotNull(), "Generate test image" ); for(int label = 0; label < parameters.m_NumberOfLabels; ++label) { mitk::ImageStatisticsCalculator::Statistics statistics = mitkImageStatisticsHotspotTestClass::CalculateStatistics(image, parameters, label); mitkImageStatisticsHotspotTestClass::ValidateStatistics(statistics, parameters, label); std::cout << std::endl; } } catch (std::exception& e) { std::cout << "Error: " << e.what() << std::endl; } MITK_TEST_END() } diff --git a/Modules/ImageStatistics/itkMultiGaussianImageSource.h b/Modules/ImageStatistics/itkMultiGaussianImageSource.h index ba4d250dd2..7569c3c0f2 100644 --- a/Modules/ImageStatistics/itkMultiGaussianImageSource.h +++ b/Modules/ImageStatistics/itkMultiGaussianImageSource.h @@ -1,375 +1,374 @@ /*========================================================================= * * Copyright Insight Software Consortium * * 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 * * http://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. * *=========================================================================*/ /*========================================================================= * * Portions of this file are subject to the VTK Toolkit Version 3 copyright. * * Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen * * For complete copyright, license and disclaimer of warranty information * please refer to the NOTICE file at the top of the ITK source tree. * *=========================================================================*/ #ifndef __itkMultiGaussianImageSource_h #define __itkMultiGaussianImageSource_h #include "itkImageSource.h" #include "itkNumericTraits.h" #include "itkImageRegionIteratorWithIndex.h" #include "itkImageFileWriter.h" #include namespace itk { /** \class MultiGaussianImageSource \brief Generate an 3-dimensional multigaussian image. -This class defines an 3-dimensional Image, in which the value at one voxel equals the value of a multigaussian function evaluated at the voxel's coordinates. The multigaussian function is build as a sum of N gaussian function. This is defined by the following parameters: +This class defines an 3-dimensional Image, in which the value at one voxel equals the value of a multigaussian function evaluated at the voxel's coordinates. The multigaussian function is built as a sum of N gaussian function and is defined by the following parameters (\ref Generation-of-a-multigauss-image): 1. CenterX, CenterY, CenterZ - vectors of the size of N determining the expectancy value at the x-, y- and the z-axis. That means: The i-th gaussian bell curve takes its maximal value at the voxel with index [CenterX(i); CenterY(i); Centerz(i)]. 2. SigmaX, SigmaY, SigmaZ - vectors of the size of N determining the deviation at the x-, y- and the z-axis. That means: The width of the i-th gaussian bell curve is determined by the deviation in the x-axis, which is SigmaX(i), in the y-axis is SigmaY(i) and in the z-axis is SigmaZ(i). 3. Altitude - vector of the size of N determining the altitude: the i-th gaussian bell curve has a height of Altitude(i). This class allows by the method CalculateMidpointAndMeanValue() to find a sphere with a specified radius that has a maximal mean value over all sphere with that radius with midpoint inside or on the boundary of the image. Furthermore it can calculate the maximal und minimal pixel intensities and whose indices in the founded sphere. +To serve as a test tool for ImageStatisticsCalculator, esp. the "hotspot search" feature of this class, MultiGaussianImageSource is also able to calculate the position of a sphere that maximizes the mean value of the voxels within the sphere (\ref Algorithm-for-calculating-statistic-in-a-sphere). -\section Multigaussian-image-source-generator Multigaussian image source generator +\section Generation-of-a-multigauss-image Generation of a multigauss image -\subsection Abstract Abstract - -The class name is MultiGaussianImageSource and can be found in the directory MINT\Mint_Bin\CMakeExternals\Source\MITK\Modules\ImageStatistics. This class generates an image, which pixel intensities are the values of a multigauss function. The class allows calculating a spherical region of defined size which has a maximal mean value and returns the position and the mean value. Furthermore it can calculate the maximal und minimal pixel intensities and whose indices in the founded sphere. - -\subsection Generation-of-a-multigauss-image Generation of a multigauss image - -A multigauss function consists of the sum of \f$ N \f$ gauss function. The \f$ i \f$-th \linebreak (\f$0 \leq i \leq N \f$) gaussian is described with the following seven parameters: +A multigauss function consists of the sum of \f$ N \f$ gauss function. The \f$ i \f$-th \linebreak (\f$0 \leq i \leq N \f$) gaussian is described with the following seven parameters (see above): - \f$ x_0^{(i)} \f$ is the expectancy value in the \f$ x \f$-Axis - \f$ y_0^{(i)} \f$ is the expectancy value in the \f$ y \f$-Axis - \f$ z_0^{(i)} \f$ is the expectancy value in the \f$ z \f$-Axis - \f$ \sigma_x^{(i)} \f$ is the deviation in the \f$ x \f$-Axis - \f$ \sigma_y^{(i)} \f$ is the deviation in the \f$ y \f$-Axis - \f$ \sigma_z^{(i)} \f$ is the deviation in the \f$ z \f$-Axis - \f$ a^{(i)} \f$ is the altitude of the gaussian. A gauss function has the following form: \f{eqnarray}{ \nonumber f^{(i)}(x,y,z) = a^{(i)} exp \left[ - \left( \frac{(x - x_0^{(i)})^2}{2 (\sigma_x^{(i)})^2} + \frac{(y - y_0^{(i)})^2}{2 (\sigma_y^{(i)})^2} + \frac{(z - z_0^{(i)})^2}{2 (\sigma_z^{(i)})^2} \right) \right]. \f} A multigauss function has then the form: \f{align*}{ f_N(x,y,z) =& \sum_{i=0}^{N}f^{(i)}(x,y,z)\\ =&\sum_{0}^{N} a^{(i)} exp \left[ - \left( \frac{(x - x_0^{(i)})^2}{2 (\sigma_x^{(i)})^2} + \frac{(y - y_0^{(i)})^2}{2 (\sigma_y^{(i)})^2} + \frac{(z - z_0^{(i)})^2}{2 (\sigma_z^{(i)})^2} \right) \right]. \f} The multigauss function \f$f_N\f$ will be evaluated at each voxel coordinate to become the voxel intensity. -\subsection Algorithm-for-calculating-statistic-in-a-sphere Algorithm for calculating statistic in a sphere +\section Algorithm-for-calculating-statistic-in-a-sphere Algorithm for calculating statistic in a sphere This section explains how we can find a sphere region which has a maximal mean value over all sphere regions with a fixed radius. Furthermore we want to calculate the maximal and minimal value in the wanted sphere. To calculate the mean value in a sphere we integrate the gaussians over the whole sphere. The antiderivative is unknown as an explicit function, but we have a value table for the distribution function of the normal distribution \f$ \Phi(x) \f$ for \f$ x \f$ between \f$ -3.99 \f$ and \f$ 3.99 \f$ with step size \f$ 0.01 \f$. The only problem is that we cannot integrate over a spherical region, because we have an 3-dim integral and therefore are the integral limits dependent from each other and we cannot evaluate \f$ \Phi \f$. So we approximate the sphere with cuboids inside the sphere and prisms on the boundary of the sphere. We calculate these cuboids with the octree recursive method: We start by subdividing the wrapping box of the sphere in eight cuboids. Further we subdivide each cuboid in eight cuboids and check for each of them, whether it is inside or outside the sphere or whether it intersects the sphere surface. We save those of them, which are inside the sphere and continue to subdivide the cuboids that intersect the sphere until the recursion breaks. In the last step we take the half of the cuboids on the boundary and this are the prisms. Now we can calculate and sum the integrals over the cuboids and divide through the volume of the body to obtain the mean value. For each cuboid \f$ Q = [a_1, b_1]\times[a_2, b_2]\times[a_3, b_3] \f$ we apply Fubini's theorem for integrable functions and become for the integral the following: \f{align*}{ m_k =& \sum_{i=0}^{N} \int_{Q} f^{(i)}(x,y,z)dx\\ =&\sum_{i=0}^{N} a^{(i)} \int_{Q} exp \left[ - \left( \frac{(x - x_0^{(i)})^2}{2 (\sigma_x^{(i)})^2} + \frac{(y - y_0^{(i)})^2}{2 (\sigma_y^{(i)})^2} + \frac{(z - z_0^{(i)})^2}{2 (\sigma_z^{(i)})^2} \right) \right] dx \\ =& \sum_{i=0}^{N} a^{(i)} \int_{a_1}^{b_1} exp \left[ - \frac{(x - x_0^{(i)})^2}{2 (\sigma_x^{(i)})^2} \right] dx \int_{a_2}^{b_2}exp \left[ -\frac{(y - y_0^{(i)})^2}{2 (\sigma_y^{(i)})^2} \right]dx \int_{a_3}^{b_3}exp \left[ -\frac{(z - z_0^{(i)})^2}{2 (\sigma_z^{(i)})^2} \right]dx. \f} So we calculate three one dimensional integrals: \f{align*}{ \int_{a}^{b} & exp \left[ - \frac{(x - x_0^{(i)})^2}{2 (\sigma_x^{(i)})^2} \right] dx \\ =&\int_{-\infty}^{b} exp \left[ - \frac{(x - x_0^{(i)})^2}{2 (\sigma_x^{(i)})^2} \right] dx - \int_{-\infty}^{a} exp \left[ - \frac{(x - x_0^{(i)})^2}{2 (\sigma_x^{(i)})^2} \right] dx \\ =& \sigma_x^{(i)} \left[\int_{-\infty}^{(a - x_0^{(i)})/ \sigma_x^{(i)}} e^{-\frac{t^2}{2}} dt - \int_{-\infty}^{(b - x_0^{(i)})/ \sigma_x^{(i)}} e^{-\frac{t^2}{2}}dt \right] \\ =&\sigma_x^{(i)} \sqrt{(\pi)} \left[ \Phi \left( \frac{(a - x_0^{(i)})^2}{\sigma_x^{(i)}} \right) - \Phi \left ( \frac{(b - x_0^{(i)})^2}{\sigma_x^{(i)}} \right) \right]. \f} and become for the integral over \f$ Q \f$: \f{align*}{ m_k =& \sum_{i=0}^{N} \sigma_x^{(i)} \sigma_y^{(i)} \sigma_z^{(i)} \pi^{1.5} \left[ \Phi \left( \frac{(a_1 - x_0^{(i)})^2}{\sigma_x^{(i)}} \right) - \Phi \left ( \frac{(b_1 - x_0^{(i)})^2}{\sigma_x^{(i)}} \right) \right]\times \\ &\left[ \Phi \left( \frac{(a_2 - y_0^{(i)})^2}{\sigma_y^{(i)}} \right) - \Phi \left ( \frac{(b_2 - y_0^{(i)})^2}{\sigma_y^{(i)}} \right) \right]\times \left[ \Phi \left( \frac{(a_3 - z_0^{(i)})^2}{\sigma_z^{(i)}} \right) - \Phi \left ( \frac{(b_3 - z_0^{(i)})^2}{\sigma_z^{(i)}} \right) \right]. \f} For the integral over the prism we take the half of the integral over the corresponding cuboid. -Altogether we become for the mean value in the sphere: +Altogether we find the mean value in the sphere as: \f{align*}{ \left( \sum_{Q_k \text{ Cuboid}} m_k + \sum_{P_l \text{ Prism}} 0.5 m_l \right )/Volume(B), \f} where Volume(B) is the volume of the body that approximate the sphere. Now we know how to calculate the mean value in a sphere for given midpoint and radius. So we assume each voxel in the given image to be the sphere midpoint and we calculate the mean value as described above. If the new mean value is greater than the "maximum-until-now", we take the new value to be the "maximum-until-now". Then we go to the next voxel and make the same calculation and so on. At the same time we save the coordinates of the midpoint voxel. After we found the midpoint and the maximal mean value, we can calculate the maximum and the minimum in the sphere: we just traverse all the voxels in the region and take the maximum and minimum value and the respective coordinates. -\subsection Input-and-output Input and output +\section Input-and-output Input and output -An example for an input in the command-line is: mitkMultiGaussianTest C:/temp/outputFile C:/temp/inputFile.xml -Here is outputFile the name of the gaussian image with extension .nrrd and at the same time the name of the output file with extension .xml, which is the same as the inputFile, only added the calculated mean value, max and min and the corresponding indexes in the statistic tag. Here we see an example for the input and output .xml file: +An example for an input in the command-line is: +\verbatim + mitkMultiGaussianTest C:/temp/outputFile C:/temp/inputFile.xml +\endverbatim ---------------------- INPUT --------------------------------------------------- +Here is outputFile the name of the gaussian image with extension .nrrd and at the same time the name of the output file with extension .xml, which is the same as the inputFile, only added the calculated mean value, max and min and the corresponding indexes in the statistic tag. Here we see an example for the input and output .xml file: \verbatim +INPUT: + \endverbatim ---------------------- OUTPUT -------------------------------------------------- - \verbatim +OUTPUT: + \endverbatim \subsection Parameter-for-the-input Parameter for the input In the tag \a testimage we describe the image that we generate. Image rows/columns/slices gives the number of rows/columns/slices of the image; \a numberOfGaussians is the number of gauss functions (\f$ N \f$); spacing defines the extend of one voxel for each direction. The parameter \a entireHotSpotInImage determines whether the whole sphere is in the image included (\f$ = 1 \f$) or only the midpoint of the sphere is inside the image. NOTE: When the \a entireHotSpotInImage \f$ = 0 \f$ it is possible that we find the midpoint of the sphere on the border of the image. In this case we cut the approximation of the sphere, so that we become a body, which is completely inside the image, but not a "sphere" anymore. To that effect is the volume of the body decreased and that could lead to unexpected results. In the subtag \a gaussian we describe each gauss function as mentioned in the second section. In the tag \a segmentation we define the radius of the wanted sphere in mm (\a hotspotRadiusInMM ). We can also set the number of labels (\a numberOfLabels ) to be an positive number and this determines the number of regions of interest (ROI). In each ROI we find the sphere with the wanted properties and midpoint inside the ROI, but not necessarily the whole sphere. In the subtag \a roi we set label number and the index coordinates for the borders of the roi. */ template< typename TOutputImage > class ITK_EXPORT MultiGaussianImageSource:public ImageSource< TOutputImage > { public: /** Standard class typedefs. */ typedef MultiGaussianImageSource Self; typedef ImageSource< TOutputImage > Superclass; typedef SmartPointer< Self > Pointer; typedef SmartPointer< const Self > ConstPointer; /** Typedef for the output image PixelType. */ typedef typename TOutputImage::PixelType OutputImagePixelType; /** Typedef to describe the output image region type. */ typedef typename TOutputImage::RegionType OutputImageRegionType; /** Method for creation through the object factory. */ itkNewMacro(Self); /** Basic types from the OutputImageType */ typedef typename TOutputImage::SizeType SizeType; typedef typename TOutputImage::IndexType IndexType; typedef typename TOutputImage::SpacingType SpacingType; typedef typename TOutputImage::PointType PointType; typedef typename SizeType::SizeValueType SizeValueType; typedef SizeValueType SizeValueArrayType[TOutputImage::ImageDimension]; typedef typename TOutputImage::SpacingValueType SpacingValueType; typedef SpacingValueType SpacingValueArrayType[TOutputImage::ImageDimension]; typedef typename TOutputImage::PointValueType PointValueType; typedef PointValueType PointValueArrayType[TOutputImage::ImageDimension]; typedef typename itk::ImageRegion<3> ::SizeValueType SizeRegionType; /** Typedef to describe the sphere radius type. */ typedef double RadiusType; /** Typedef to describe the standard vector type. */ typedef std::vector VectorType; /** Typedef to describe the itk vector type. */ typedef Vector ItkVectorType; /** Typedef to describe the ImageRegionIteratorWithIndex type. */ typedef ImageRegionIteratorWithIndex IteratorType; /** Typedef to describe the Poiner type at the output image. */ typedef typename TOutputImage::Pointer ImageType; typedef MapContainer MapContainerPoints; typedef MapContainer MapContainerRadius; /** Set/Get size of the output image. */ itkSetMacro(Size, SizeType); virtual void SetSize(SizeValueArrayType sizeArray); virtual const SizeValueType * GetSize() const; /** Set/Get spacing of the output image. */ itkSetMacro(Spacing, SpacingType); virtual void SetSpacing(SpacingValueArrayType spacingArray); virtual const SpacingValueType * GetSpacing() const; /** Set/Get origin of the output image. This programm works proper only with origin [0.0, 0.0, 0.0] */ itkSetMacro(Origin, PointType); virtual void SetOrigin(PointValueArrayType originArray); virtual const PointValueType * GetOrigin() const; /** Get the number of gaussian functions in the output image. */ virtual unsigned int GetNumberOfGaussians() const; /** Set the number of gaussian function. */ virtual void SetNumberOfGausssians( unsigned int ); /** Set/Get the radius of the sphere. */ virtual const RadiusType GetRadius() const; virtual void SetRadius( RadiusType radius ); /** Get the maximal mean value in a sphere over all possible spheres with midpoint in the image. */ virtual const OutputImagePixelType GetMaxMeanValue() const; /** Get the index of the midpoint of a sphere with the maximal mean value.*/ virtual const IndexType GetSphereMidpoint() const; /** Calculates the value of the multigaussian function at a Point given by its coordinates [x, y, z]. */ virtual const double MultiGaussianFunctionValueAtPoint(double , double, double); /** Adds a multigaussian defined by the parameter: CenterX, CenterY, CenterZ, SigmaX, SigmaY, SigmaZ, Altitude. All parameters should have the same size, which determinates the number of the gaussian added. */ - virtual void AddGaussian( VectorType, VectorType, VectorType, VectorType, VectorType, VectorType, VectorType); + virtual void AddGaussian( VectorType centerX, VectorType centerY, VectorType centerZ, VectorType sigmaX, VectorType sigmaY, VectorType sigmaZ, VectorType altitude); /** Calculates and set the index of the midpoint of the sphere with the maximal mean value as well as the mean value. */ virtual void CalculateTheMidpointAndTheMeanValueWithOctree(); /** Calculates and set the index an the value of maximulm and minimum in the wanted sphere. */ virtual void CalculateMaxAndMinInSphere(); /** Get the index in the sphere with maximal value. */ virtual const IndexType GetMaxValueIndexInSphere() const; /** Get the maximal value in the sphere. */ virtual const OutputImagePixelType GetMaxValueInSphere() const; /** Get the index in the sphere with minimal value. */ virtual const IndexType GetMinValueIndexInSphere() const; /** Get the minimal value in the sphere. */ virtual const OutputImagePixelType GetMinValueInSphere() const; /** Set the region of interest. */ virtual void SetRegionOfInterest(ItkVectorType, ItkVectorType); /** Write a .mps file to visualise the point in the sphere. */ virtual void WriteXMLToTestTheCuboidInsideTheSphere(); /**This recursive method realise the octree method. It subdivide a cuboid in eight cuboids, when this cuboid crosses the boundary of sphere. If the cuboid is inside the sphere, it calculates the integral. */ virtual void CalculateEdgesInSphere( PointType globalCoordinateMidpointCuboid, PointType globalCoordinateMidpointSphere, double cuboidRadius, int level); /**Calculate and return value of the integral of the gaussian in a cuboid region with the dimension 3: in the x-axis between xMin and xMax and in the y-axis between yMin and yMax and in the z-axis also between zMin and zMax. */ virtual double MultiGaussianFunctionValueAtCuboid(double xMin, double xMax, double yMin, double yMax, double zMin, double zMax); /** Inseret the midpoints of cuboid in a vector m_Midpoints, so that we can visualise it. */ virtual void InsertPoints( PointType globalCoordinateMidpointCuboid, double cuboidRadius); /** Start the octree recursion in eigth directions for the sphere with midpoint globalCoordinateMidpointSphere. */ virtual void GenerateCuboidSegmentationInSphere( PointType globalCoordinateMidpointSphere ); /** Get the the values of the cumulative distribution function of the normal distribution. */ virtual double FunctionPhi(double value); /** Check if a cuboid with midpoint globalCoordinateMidpointCuboid and side length sideLength intersect the sphere with midpoint globalCoordinateMidpointSphere boundary. */ virtual unsigned int IntersectTheSphere( PointType globalCoordinateMidpointCuboid, PointType globalCoordinateMidpointSphere, double sideLength); /** Set the tabel values of the distribution function of the normal distribution. */ void SetNormalDistributionValues(); /** Set the minimum possible pixel value. By default, it is * NumericTraits::min(). */ itkSetClampMacro( Min, OutputImagePixelType, NumericTraits< OutputImagePixelType >::NonpositiveMin(), NumericTraits< OutputImagePixelType >::max() ); /** Check if a index is inside the image*/ bool IsInImage(IndexType index); /** Get the minimum possible pixel value. */ itkGetConstMacro(Min, OutputImagePixelType); /** Set the maximum possible pixel value. By default, it is * NumericTraits::max(). */ itkSetClampMacro( Max, OutputImagePixelType, NumericTraits< OutputImagePixelType >::NonpositiveMin(), NumericTraits< OutputImagePixelType >::max() ); /** Get the maximum possible pixel value. */ itkGetConstMacro(Max, OutputImagePixelType); protected: MultiGaussianImageSource(); ~MultiGaussianImageSource(); void PrintSelf(std::ostream & os, Indent indent) const; virtual void GenerateData(); virtual void GenerateOutputInformation(); private: MultiGaussianImageSource(const MultiGaussianImageSource &); //purposely not implemented void operator=(const MultiGaussianImageSource &); //purposely not implemented SizeType m_Size; //size of the output image SpacingType m_Spacing; //spacing PointType m_Origin; //origin OutputImagePixelType m_MaxValueInSphere; //maximal value in the wanted sphere IndexType m_MaxValueIndexInSphere; //index of the maximal value in the wanted sphere OutputImagePixelType m_MinValueInSphere; //minimal value in the wanted sphere IndexType m_MinValueIndexInSphere; //index of the minimal value in the wanted sphere unsigned int m_NumberOfGaussians; //number of Gaussians RadiusType m_Radius; //radius of the sphere unsigned int m_RadiusStepNumber; //number of steps to traverse the sphere radius OutputImagePixelType m_MeanValue; //mean value in the wanted sphere OutputImagePixelType m_ValueAtMidpoint; //value at the midpoint of the wanted sphere IndexType m_SphereMidpoint; //midpoint of the wanted sphere VectorType m_SigmaX; //deviation in the x-axis VectorType m_SigmaY; //deviation in the y-axis VectorType m_SigmaZ; //deviation in the z-axis VectorType m_CenterX; //x-coordinate of the mean value of Gaussians VectorType m_CenterY; //y-coordinate of the mean value of Gaussians VectorType m_CenterZ; //z-coordinate of the mean value of Gaussians VectorType m_Altitude; //amplitude ItkVectorType m_RegionOfInterestMax; //maximal values for the coordinates in the region of interest ItkVectorType m_RegionOfInterestMin; //minimal values for the coordinates in the region of interest typename TOutputImage::PixelType m_Min; //minimum possible value typename TOutputImage::PixelType m_Max; //maximum possible value PointType m_GlobalCoordinate; //physical coordiante of the sphere midpoint bool m_WriteMPS; //1 = write a MPS File to visualise the cuboid midpoints of one approximation of the sphere MapContainerPoints m_Midpoints; //the midpoints of the cuboids MapContainerRadius m_RadiusCuboid; //the radius ( = 0.5 * side length) of the cuboids (in the same order as the midpoints in m_Midpoints) double m_Volume; //the volume of the body, that approximize the sphere double m_NormalDistValues [410];//normal distribution values double m_meanValueTemp; //= m_Volume * meanValue in each sphere // The following variables are deprecated, and provided here just for // backward compatibility. It use is discouraged. mutable PointValueArrayType m_OriginArray; mutable SpacingValueArrayType m_SpacingArray; }; } // end namespace itk #ifndef ITK_MANUAL_INSTANTIATION #include "itkMultiGaussianImageSource.hxx" #endif #endif