Skip to content

[hist] Using TKDE::Fill works with empty tkde #14740

New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Open
wants to merge 6 commits into
base: master
Choose a base branch
from
Open
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
13 changes: 12 additions & 1 deletion hist/hist/src/TKDE.cxx
Original file line number Diff line number Diff line change
Expand Up @@ -518,7 +518,15 @@ void TKDE::InitFromNewData() {
if (fUseMirroring) {
SetMirroredEvents();
}
else {
// to set fBinCount and fSumOfCounts
SetBinCountData();
}
// in case of I/O reset the kernel
if (!fKernelFunction) SetKernelFunction(nullptr);

SetKernel();

}

///////////////////////////////////////////////////////////////////////////////////////////////////////
Expand Down Expand Up @@ -876,6 +884,9 @@ void TKDE::SetBinCountData() {
}
fBinCount.clear();
}
if (fSumOfCounts == 0) {
Warning("SetBinCountData()", "Empty sum of counts, might lead to nan/inf errors when normalizing.");
}
}

////////////////////////////////////////////////////////////////////////////////
Expand Down Expand Up @@ -1015,7 +1026,7 @@ Double_t TKDE::TKernel::operator()(Double_t x) const {
Bool_t useCount = (fKDE->fBinCount.size() == n);
// also in case of unbinned unweighted data fSumOfCounts is sum of events in range
// events outside range should be used to normalize the TKDE ??
Double_t nSum = fKDE->fSumOfCounts; //(useBins) ? fKDE->fSumOfCounts : fKDE->fNEvents;
Double_t nSum = fKDE->fSumOfCounts; //(useCount) ? fKDE->fSumOfCounts : fKDE->fNEvents;
//if (!useCount) nSum = fKDE->fNEvents;
// in case of non-adaptive fWeights is a vector of size 1
Bool_t hasAdaptiveWeights = (fWeights.size() == n);
Expand Down
10 changes: 9 additions & 1 deletion hist/hist/test/test_tkde.cxx
Original file line number Diff line number Diff line change
Expand Up @@ -496,4 +496,12 @@ TEST(TKDE, tkde_io_mirror)
for (size_t i = 0; i < t.xtest.size(); ++i) {
EXPECT_NEAR(t.values1[i], t.values2[i], delta);
}
}
}

TEST(TKDE, tkde_fill_empty) // ROOT #7808
{
TRandom3 r;
auto kde = new TKDE(0, nullptr, 0, 5, "KernelType:Gaussian;Iteration:Adaptive;Mirror:noMirror;Binning:RelaxedBinning", 1);
for (unsigned int i = 0; i < 100; i++) { kde->Fill(r.Gaus(2,1)); }
EXPECT_NEAR(kde->GetValue(2), 0.487581, 1e-1); // Worked on ROOT <= v6.16, failed from 6.18 until this test with result "nan"
}
Loading