Skip to content

Commit

Permalink
[RF] Forbid calling getVal() with r-value references to norm set
Browse files Browse the repository at this point in the history
Calling RooAbsReal::getVal() with r-value references to the
normalization set (i.e. with RooArgSets that are created in place) is a
bad idea, because it breaks RooFits caching logic and potentially
introduces significant overhead.

There was even a dedicated Pythonization to deal with the related memory
problems.

Since to avoid the problems and the need for a workaround in Python,
this commit simply suggests to forbid calling getVal() with r-value
references at runtime by printing a clear exception.
  • Loading branch information
guitargeek committed Nov 27, 2024
1 parent 1adb93b commit 337b56a
Show file tree
Hide file tree
Showing 7 changed files with 52 additions and 27 deletions.
Original file line number Diff line number Diff line change
Expand Up @@ -116,22 +116,6 @@ def chi2FitTo(self, *args, **kwargs):
args, kwargs = _kwargs_to_roocmdargs(*args, **kwargs)
return self._chi2FitTo(*args, **kwargs)

def getVal(self, normalizationSet=None):
# We do the conversion to RooArgSet now, such that we can keep alive
# the normalization set by setting it as an attribute of this
# RooAbsReal.
if isinstance(normalizationSet, (set, list, tuple)):
import ROOT

normalizationSet = ROOT.RooArgSet(normalizationSet)
# With the pythonizations, we have the opportunity to use the Python
# reference counting to make sure the last normalization set doesn't
# get deleted under our feet (RooFit tries to use it by pointer when
# you call getVal() without any normalization set the next time).
if normalizationSet:
self._getVal_normSet = normalizationSet
return self._getVal(normalizationSet) if normalizationSet else self._getVal()

def setEvalErrorLoggingMode(m):

import ROOT
Expand Down
2 changes: 2 additions & 0 deletions roofit/roofitcore/inc/RooAbsReal.h
Original file line number Diff line number Diff line change
Expand Up @@ -131,6 +131,8 @@ class RooAbsReal : public RooAbsArg {
return _fast ? _value : getValV(normalisationSet.empty() ? nullptr : &normalisationSet) ;
}

double getVal(RooArgSet &&) const;

virtual double getValV(const RooArgSet* normalisationSet = nullptr) const ;

double getPropagatedError(const RooFitResult &fr, const RooArgSet &nset = {}) const;
Expand Down
26 changes: 26 additions & 0 deletions roofit/roofitcore/src/RooAbsReal.cxx
Original file line number Diff line number Diff line change
Expand Up @@ -4475,3 +4475,29 @@ void RooAbsReal::enableOffsetting(bool flag)


RooAbsReal::Ref::Ref(double val) : _ref{RooFit::RooConst(val)} {}

////////////////////////////////////////////////////////////////////////////////
/// Calling RooAbsReal::getVal() with an r-value reference is a common
/// performance trap, because this usually happens when implicitly constructing
/// the RooArgSet to be used as the parameter (for example, in calls like
/// `pdf.getVal(x)`).
///
/// Creating the RooArgSet can cause quite some overhead, especially when the
/// evaluated object is just a simple variable. Even worse, many RooFit objects
/// internally cache information using the uniqueId() of the normalization set
/// as the key. So by constructing normalization sets in place, RooFits caching
/// logic is broken.
///
/// To avoid these kind of problems, getVal() will just throw an error when
/// it's called with an r-value reference. This also catches the cases where
/// one uses it in Python, implicitly creating the normalization set from a
/// Python list or set.
double RooAbsReal::getVal(RooArgSet &&) const
{
std::stringstream errMsg;
errMsg << "calling RooAbsReal::getVal() with r-value references to the normalization set is not allowed, because "
"it breaks RooFits caching logic and potentially introduces significant overhead. Please explicitly "
"create the RooArgSet outside the call to getVal().";
coutF(Eval) << errMsg.str() << std::endl;
throw std::runtime_error(errMsg.str());
}
4 changes: 3 additions & 1 deletion roofit/roofitcore/test/testRooAbsPdf.cxx
Original file line number Diff line number Diff line change
Expand Up @@ -389,8 +389,10 @@ TEST(RooAbsPdf, NormSetChange)

RooAddition add("add", "add", {gauss});

RooArgSet normSet{x};

double v1 = add.getVal();
double v2 = add.getVal(x);
double v2 = add.getVal(normSet);

// The change of normalization set should trigger a recomputation of the
// value, so val2 should be different from val1. }
Expand Down
3 changes: 2 additions & 1 deletion roofit/roofitcore/test/testRooDataHist.cxx
Original file line number Diff line number Diff line change
Expand Up @@ -740,9 +740,10 @@ TEST_P(WeightsTest, VectorizedWeights)
}

std::vector<double> weightsGetVal(nVals);
RooArgSet normSet{x};
for (std::size_t i = 0; i < nVals; ++i) {
x.setVal(xVals[i]);
weightsGetVal[i] = absReal->getVal(x);
weightsGetVal[i] = absReal->getVal(normSet);
}
x.setVal(0.0);

Expand Down
25 changes: 17 additions & 8 deletions roofit/roofitcore/test/testRooProdPdf.cxx
Original file line number Diff line number Diff line change
Expand Up @@ -99,14 +99,23 @@ TEST(RooProdPdf, TestGetPartIntList)
// Product of all the pdfs.
auto &prod = static_cast<RooProdPdf &>(*ws.factory("PROD::prod(pdf1, pdf2)"));

EXPECT_DOUBLE_EQ(prod.getVal({}), 1.0);
EXPECT_DOUBLE_EQ(prod.getVal({x}), 1. / a);
EXPECT_DOUBLE_EQ(prod.getVal({y}), 1. / b);
EXPECT_DOUBLE_EQ(prod.getVal({z}), 1. / c);
EXPECT_DOUBLE_EQ(prod.getVal({x, y}), 1. / a / b);
EXPECT_DOUBLE_EQ(prod.getVal({x, z}), 1. / a / c);
EXPECT_DOUBLE_EQ(prod.getVal({y, z}), 1. / b / c);
EXPECT_DOUBLE_EQ(prod.getVal({x, y, z}), 1. / a / b / c);
RooArgSet normSetNada{};
RooArgSet normSetX{x};
RooArgSet normSetY{y};
RooArgSet normSetZ{z};
RooArgSet normSetXY{x, y};
RooArgSet normSetXZ{x, z};
RooArgSet normSetYZ{y, z};
RooArgSet normSetXYZ{x, y, z};

EXPECT_DOUBLE_EQ(prod.getVal(normSetNada), 1.0);
EXPECT_DOUBLE_EQ(prod.getVal(normSetX), 1. / a);
EXPECT_DOUBLE_EQ(prod.getVal(normSetY), 1. / b);
EXPECT_DOUBLE_EQ(prod.getVal(normSetZ), 1. / c);
EXPECT_DOUBLE_EQ(prod.getVal(normSetXY), 1. / a / b);
EXPECT_DOUBLE_EQ(prod.getVal(normSetXZ), 1. / a / c);
EXPECT_DOUBLE_EQ(prod.getVal(normSetYZ), 1. / b / c);
EXPECT_DOUBLE_EQ(prod.getVal(normSetXYZ), 1. / a / b / c);
}

TEST(RooProdPdf, TestDepsAreCond)
Expand Down
3 changes: 2 additions & 1 deletion roofit/roofitcore/test/testRooWrapperPdf.cxx
Original file line number Diff line number Diff line change
Expand Up @@ -29,7 +29,8 @@ TEST(RooWrapperPdf, Basics)

RooWrapperPdf polPdf("polPdf", "polynomial PDF", pol);

EXPECT_GT(pol.getVal(x)*1.05, polPdf.getVal(x)) << "Wrapper pdf normalises.";
RooArgSet normSet{x};
EXPECT_GT(pol.getVal(normSet)*1.05, polPdf.getVal(normSet)) << "Wrapper pdf normalises.";

RooArgSet intSet(x);
RooArgSet numSet;
Expand Down

0 comments on commit 337b56a

Please sign in to comment.