[5d12a0]: / src / labelOverlapMeasures.cxx

Download this file

111 lines (89 with data), 3.9 kB

  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
#include <nanobind/nanobind.h>
#include <nanobind/stl/vector.h>
#include <nanobind/stl/string.h>
#include <nanobind/stl/tuple.h>
#include <nanobind/stl/list.h>
#include <nanobind/ndarray.h>
#include <nanobind/stl/shared_ptr.h>
#include <exception>
#include <vector>
#include <string>
#include "itkImage.h"
#include "itkLabelOverlapMeasuresImageFilter.h"
#include "antsImage.h"
namespace nb = nanobind;
using namespace nb::literals;
template<class PrecisionType, unsigned int ImageDimension>
nb::dict labelOverlapMeasures( AntsImage<itk::Image<PrecisionType, ImageDimension>> & antsSourceImage,
AntsImage<itk::Image<PrecisionType, ImageDimension>> & antsTargetImage )
{
using ImageType = itk::Image<PrecisionType, ImageDimension>;
using ImagePointerType = typename ImageType::Pointer;
typename ImageType::Pointer itkSourceImage = antsSourceImage.ptr;
typename ImageType::Pointer itkTargetImage = antsTargetImage.ptr;
using FilterType = itk::LabelOverlapMeasuresImageFilter<ImageType>;
typename FilterType::Pointer filter = FilterType::New();
filter->SetSourceImage( itkSourceImage );
filter->SetTargetImage( itkTargetImage );
filter->Update();
typename FilterType::MapType labelMap = filter->GetLabelSetMeasures();
// Sort the labels
std::vector<PrecisionType> allLabels;
allLabels.clear();
for( typename FilterType::MapType::const_iterator it = labelMap.begin();
it != labelMap.end(); ++it )
{
if( (*it).first == 0 )
{
continue;
}
const int label = (*it).first;
allLabels.push_back( label );
}
std::sort( allLabels.begin(), allLabels.end() );
// Now put the results in an Rcpp data frame
unsigned int vectorLength = 1 + allLabels.size();
std::vector<PrecisionType> labels( vectorLength );
std::vector<double> totalOrTargetOverlap( vectorLength );
std::vector<double> unionOverlap( vectorLength );
std::vector<double> meanOverlap( vectorLength );
std::vector<double> volumeSimilarity( vectorLength );
std::vector<double> falseNegativeError( vectorLength );
std::vector<double> falsePositiveError( vectorLength );
// We'll replace label '0' with "All" in the R wrapper.
labels[0] = itk::NumericTraits<PrecisionType>::Zero;
totalOrTargetOverlap[0] = filter->GetTotalOverlap();
unionOverlap[0] = filter->GetUnionOverlap();
meanOverlap[0] = filter->GetMeanOverlap();
volumeSimilarity[0] = filter->GetVolumeSimilarity();
falseNegativeError[0] = filter->GetFalseNegativeError();
falsePositiveError[0] = filter->GetFalsePositiveError();
unsigned int i = 1;
typename std::vector<PrecisionType>::const_iterator itL = allLabels.begin();
for( itL = allLabels.begin(); itL != allLabels.end(); ++itL )
{
labels[i] = *itL;
totalOrTargetOverlap[i] = filter->GetTargetOverlap( *itL );
unionOverlap[i] = filter->GetUnionOverlap( *itL );
meanOverlap[i] = filter->GetMeanOverlap( *itL );
volumeSimilarity[i] = filter->GetVolumeSimilarity( *itL );
falseNegativeError[i] = filter->GetFalseNegativeError( *itL );
falsePositiveError[i] = filter->GetFalsePositiveError( *itL );
i++;
}
nb::dict labelOverlapMeasures;
labelOverlapMeasures["Label"] = labels;
labelOverlapMeasures["TotalOrTargetOverlap"] = totalOrTargetOverlap;
labelOverlapMeasures["UnionOverlap"] = unionOverlap;
labelOverlapMeasures["MeanOverlap"] = meanOverlap;
labelOverlapMeasures["VolumeSimilarity"] = volumeSimilarity;
labelOverlapMeasures["FalseNegativeError"] = falseNegativeError;
labelOverlapMeasures["FalsePositiveError"] = falsePositiveError;
return labelOverlapMeasures;
}
void local_labelOverlapMeasures(nb::module_ &m)
{
m.def("labelOverlapMeasures2D", &labelOverlapMeasures<unsigned int, 2>);
m.def("labelOverlapMeasures3D", &labelOverlapMeasures<unsigned int, 3>);
m.def("labelOverlapMeasures4D", &labelOverlapMeasures<unsigned int, 4>);
}