Skip to content

Commit 25f8bbd

Browse files
committed
[RF] Don't normalize over non-dependents when compiling RooFormulaVar
Don't normalize pdfs over non-dependents in `RooFormulaVar::compileForNormSet()`. This fixes use cases where pdfs that don't depend on any observables are used for their unnormalized `RooAbsPdf::evaluate()` shape as functions. Even thought technically not allowed because evaluating pdfs without a normalization set is undefined, this pattern is used a lot in the wild, so we need to support it also in the new vectorizing evaluation backend. A unit test that covers the original user-reported problem is also added.
1 parent 202843a commit 25f8bbd

File tree

4 files changed

+64
-3
lines changed

4 files changed

+64
-3
lines changed

roofit/roofitcore/inc/RooFormulaVar.h

Lines changed: 3 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -77,6 +77,9 @@ class RooFormulaVar : public RooAbsReal {
7777

7878
std::string getUniqueFuncName() const;
7979

80+
std::unique_ptr<RooAbsArg>
81+
compileForNormSet(RooArgSet const &normSet, RooFit::Detail::CompileContext &ctx) const override;
82+
8083
protected:
8184
// Post-processing of server redirection
8285
bool redirectServersHook(const RooAbsCollection& newServerList, bool mustReplaceAll, bool nameChange, bool isRecursive) override ;

roofit/roofitcore/src/RooFormulaVar.cxx

Lines changed: 19 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -316,3 +316,22 @@ std::string RooFormulaVar::getUniqueFuncName() const
316316
{
317317
return getFormula().getTFormula()->GetUniqueFuncName().Data();
318318
}
319+
320+
std::unique_ptr<RooAbsArg>
321+
RooFormulaVar::compileForNormSet(RooArgSet const &normSet, RooFit::Detail::CompileContext &ctx) const
322+
{
323+
// Some users exploit unnormalized RooAbsPdfs as inputs for RooFormulaVars,
324+
// relying on what the pdf returns from RooAbsPdf::evaluate(). This is in
325+
// principle not allowed because every pdf needs to be evaluated with a
326+
// normalization set, but it's so common in user code that we need to
327+
// support it. To make this work, we need to make sure that the no
328+
// normalization over non-dependents is happening at this point, reducing
329+
// the normalization set to the subset of actual dependents.
330+
// See also the "PdfAsFunctionInFormulaVar" test in testRooAbsPdf.
331+
RooArgSet depList;
332+
getObservables(&normSet, depList);
333+
auto newArg = std::unique_ptr<RooAbsArg>{static_cast<RooAbsArg *>(Clone())};
334+
ctx.markAsCompiled(*newArg);
335+
ctx.compileServers(*newArg, depList);
336+
return newArg;
337+
}

roofit/roofitcore/test/CMakeLists.txt

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -18,7 +18,7 @@ ROOT_ADD_GTEST(testRooDataHist testRooDataHist.cxx LIBRARIES RooFitCore
1818
ROOT_ADD_GTEST(testRooBinSamplingPdf testRooBinSamplingPdf.cxx LIBRARIES RooFitCore)
1919
ROOT_ADD_GTEST(testRooWrapperPdf testRooWrapperPdf.cxx LIBRARIES Gpad RooFitCore)
2020
ROOT_ADD_GTEST(testGenericPdf testGenericPdf.cxx LIBRARIES RooFitCore)
21-
ROOT_ADD_GTEST(testRooAbsPdf testRooAbsPdf.cxx LIBRARIES RooFitCore)
21+
ROOT_ADD_GTEST(testRooAbsPdf testRooAbsPdf.cxx LIBRARIES RooFitCore RooFit)
2222
ROOT_ADD_GTEST(testRooAbsCollection testRooAbsCollection.cxx LIBRARIES RooFitCore)
2323

2424
ROOT_ADD_GTEST(testRooDataSet testRooDataSet.cxx LIBRARIES Tree RooFitCore

roofit/roofitcore/test/testRooAbsPdf.cxx

Lines changed: 41 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -2,23 +2,26 @@
22
// Authors: Stephan Hageboeck, CERN 04/2020
33
// Jonas Rembser, CERN 04/2021
44

5-
#include <RooAddition.h>
65
#include <RooAddPdf.h>
6+
#include <RooAddition.h>
77
#include <RooCategory.h>
88
#include <RooConstVar.h>
99
#include <RooDataHist.h>
1010
#include <RooDataSet.h>
1111
#include <RooFitResult.h>
1212
#include <RooFormulaVar.h>
13+
#include <RooGaussian.h>
1314
#include <RooGenericPdf.h>
1415
#include <RooHelpers.h>
16+
#include <RooParametricStepFunction.h>
1517
#include <RooProdPdf.h>
1618
#include <RooProduct.h>
19+
#include <RooRandom.h>
1720
#include <RooRealVar.h>
1821
#include <RooSimultaneous.h>
1922
#include <RooWorkspace.h>
20-
#include <RooRandom.h>
2123

24+
#include <TArrayD.h>
2225
#include <TClass.h>
2326
#include <TRandom.h>
2427

@@ -375,6 +378,42 @@ TEST_P(FitTest, ProblemsWith2DSimultaneousFit)
375378
simPdf.fitTo(*data, PrintLevel(-1), _evalBackend);
376379
}
377380

381+
// This test covers a usecase by an ATLAS collaborator. The unnormalized shape
382+
// of a RooFit pdf is used as a function inside a RooFormulaVar, which is used
383+
// for the bins of a RooParametricStep function. This case is potentially
384+
// fragile, because it requires that the top-level normalization set is ignored
385+
// for the inner pdfs that don't depend on the observable, as normalizing over
386+
// a non-dependent is a corner case where the variable should be dropped.
387+
TEST_P(FitTest, PdfAsFunctionInFormulaVar)
388+
{
389+
using namespace RooFit;
390+
391+
RooRealVar x{"x", "x", 0., 1.};
392+
393+
const double arg = std::sqrt(-8 * std::log(0.5));
394+
RooGaussian gauss1{"gauss1", "", 10 - arg, 10, 2};
395+
396+
RooFormulaVar func1{"func1", "x[0]", {gauss1}};
397+
398+
// The parametric step function doesn't make any attempt at self
399+
// normalization, so the pdf value in the first bin is just the value of the
400+
// unnormalized Gaussian. And the first bin is the almost the whole domain.
401+
TArrayD limits(3);
402+
limits[0] = 0.;
403+
limits[1] = 1 - 1e-9;
404+
limits[2] = 1.;
405+
RooParametricStepFunction pdf("pdf", "pdf", x, {func1}, limits, limits.size() - 1);
406+
407+
int nEvents = 10000;
408+
std::unique_ptr<RooAbsData> data{pdf.generateBinned(x, nEvents)};
409+
410+
std::unique_ptr<RooAbsReal> nll{pdf.createNLL(*data, _evalBackend)};
411+
412+
// The test is designed to have an analytical reference value
413+
double ref = nEvents * -std::log(0.5);
414+
EXPECT_FLOAT_EQ(nll->getVal(), ref);
415+
}
416+
378417
// Verifies that a server pdf gets correctly reevaluated when the normalization
379418
// set is changed.
380419
TEST(RooAbsPdf, NormSetChange)

0 commit comments

Comments
 (0)