Skip to content

Commit 63be8d3

Browse files
committed
[RF] Correctly normalize RooSimultaneous when projecting over data
Since a8ef8b0, any RooAbsPdf that makes a prediction for the expected number of events will use this prediction for the normalization when plotting. This made is more meaningful to compare models with data because the model normalization was not automatically scaled to match the data. However, it also broke plots with projections of RooSimultaneous pdfs over data, because the temporary RooAddPdf that is created in this case did not use coefficients that were scaled to the actual dataset, but to unity. Correctly normalizing to the sum of entries in the dataset fixes this problem. A unit test is also implemented. Closes #20383. (cherry picked from commit 61735d8)
1 parent 7e82797 commit 63be8d3

File tree

3 files changed

+37
-3
lines changed

3 files changed

+37
-3
lines changed

roofit/roofitcore/src/RooSimultaneous.cxx

Lines changed: 2 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -810,8 +810,8 @@ RooPlot* RooSimultaneous::plotOn(RooPlot *frame, RooLinkedList& cmdList) const
810810
}
811811
if (skip) continue ;
812812

813-
// Instantiate a RRV holding this pdfs weight fraction
814-
wgtCompList.addOwned(std::make_unique<RooRealVar>(proxy->name(),"coef",wTable->getFrac(proxy->name())));
813+
// Instantiate a RRV holding this pdfs weight
814+
wgtCompList.addOwned(std::make_unique<RooRealVar>(proxy->name(),"coef",wTable->get(proxy->name())));
815815
sumWeight += wTable->getFrac(proxy->name()) ;
816816

817817
// Add the PDF to list list

roofit/roofitcore/test/CMakeLists.txt

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -65,7 +65,7 @@ ROOT_ADD_GTEST(testRooPolyFunc testRooPolyFunc.cxx LIBRARIES Gpad RooFitCore)
6565
ROOT_ADD_GTEST(testRooRealL TestStatistics/testRooRealL.cxx LIBRARIES RooFitCore)
6666
ROOT_ADD_GTEST(testRooRombergIntegrator testRooRombergIntegrator.cxx LIBRARIES MathCore RooFitCore)
6767
ROOT_ADD_GTEST(testRooSTLRefCountList testRooSTLRefCountList.cxx LIBRARIES RooFitCore)
68-
ROOT_ADD_GTEST(testRooSimultaneous testRooSimultaneous.cxx LIBRARIES RooFitCore)
68+
ROOT_ADD_GTEST(testRooSimultaneous testRooSimultaneous.cxx LIBRARIES RooFitCore RooFit)
6969
ROOT_ADD_GTEST(testRooTruthModel testRooTruthModel.cxx LIBRARIES RooFitCore RooFit)
7070
ROOT_ADD_GTEST(testSumW2Error testSumW2Error.cxx LIBRARIES Gpad RooFitCore)
7171
ROOT_ADD_GTEST(testTestStatistics testTestStatistics.cxx LIBRARIES RooFitCore)

roofit/roofitcore/test/testRooSimultaneous.cxx

Lines changed: 34 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -9,11 +9,13 @@
99
#include <RooGenericPdf.h>
1010
#include <RooHelpers.h>
1111
#include <RooMinimizer.h>
12+
#include <RooPlot.h>
1213
#include <RooProdPdf.h>
1314
#include <RooRandom.h>
1415
#include <RooRealVar.h>
1516
#include <RooSimultaneous.h>
1617
#include <RooThresholdCategory.h>
18+
#include <RooUniform.h>
1719
#include <RooWorkspace.h>
1820

1921
#include "gtest_wrapper.h"
@@ -507,3 +509,35 @@ TEST_P(TestStatisticTest, RooSimultaneousSingleChannelCrossCheckWithCondVar)
507509
EXPECT_TRUE(resSimWrapped->isIdentical(*resDirect))
508510
<< "Inconsistency in RooSimultaneous wrapping with ConditionalObservables";
509511
}
512+
513+
/// GitHub issue #20383.
514+
/// Check that the the simultaneous pdf is normalized correctly when plotting
515+
/// with a projection dataset.
516+
TEST(RooSimultaneous, PlotProjWData)
517+
{
518+
RooRealVar x("x", "x", -8, 8);
519+
x.setBins(1);
520+
521+
RooUniform model{"model", "", x};
522+
RooUniform model_ctl{"model_ctl", "", x};
523+
524+
RooCategory sample("sample", "sample", {{"physics", 0}, {"control", 1}});
525+
526+
RooArgSet vars{x, sample};
527+
RooDataHist combData{"combData", "", vars};
528+
sample.setLabel("physics");
529+
combData.add(vars, 1000);
530+
sample.setLabel("control");
531+
combData.add(vars, 2000);
532+
533+
RooSimultaneous simPdf("simPdf", "simultaneous pdf", {{"physics", &model}, {"control", &model_ctl}}, sample);
534+
535+
RooPlot *frame = x.frame();
536+
combData.plotOn(frame);
537+
simPdf.plotOn(frame, RooFit::ProjWData(sample, combData));
538+
539+
// The pdf should be normalized to match the data. In this test, we plot a
540+
// single bin and the model is uniform, to the curve should be equal to the
541+
// sum of data entries in the center.
542+
EXPECT_DOUBLE_EQ(frame->getCurve()->interpolate(0.), combData.sumEntries());
543+
}

0 commit comments

Comments
 (0)