Skip to content
Snippets Groups Projects
Commit 7d39c076 authored by Jonas Rembser's avatar Jonas Rembser
Browse files

[RF] Enable analytic integration of RooHistPdfs with RooLinearVars

The RooHistPdf and RooHistFunc should be able to do analytic integration
if the input is a linear transformation of a variable using
RooLinearVar.

This makes the fits faster where one wants to shift a template on the
x-axis, which is for example talked about in this forum post:
https://root-forum.cern.ch/t/roofit-pdf-normalization-integration/53905
parent fce73f05
No related branches found
No related tags found
No related merge requests found
......@@ -74,6 +74,8 @@ public:
Int_t getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName=nullptr) const override ;
double analyticalIntegral(Int_t code, const char* rangeName=nullptr) const override ;
bool forceAnalyticalInt(const RooAbsArg& dep) const override;
/// Set use of special boundary conditions for c.d.f.s
void setCdfBoundaries(bool flag) {
_cdfBoundaries = flag ;
......
......@@ -65,6 +65,8 @@ public:
Int_t getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName=nullptr) const override ;
double analyticalIntegral(Int_t code, const char* rangeName=nullptr) const override ;
bool forceAnalyticalInt(const RooAbsArg& dep) const override;
void setCdfBoundaries(bool flag) {
// Set use of special boundary conditions for c.d.f.s
_cdfBoundaries = flag ;
......@@ -119,17 +121,19 @@ private:
friend class RooHistFunc;
static bool forceAnalyticalInt(RooArgSet const& pdfObsList, RooAbsArg const& dep);
static Int_t getAnalyticalIntegral(RooArgSet& allVars,
RooArgSet& analVars,
const char* rangeName,
RooArgSet const& histObsList,
RooSetProxy const& pdfObsList,
RooArgSet const& pdfObsList,
Int_t intOrder) ;
static double analyticalIntegral(Int_t code,
const char* rangeName,
RooArgSet const& histObsList,
RooSetProxy const& pdfObsList,
RooArgSet const& pdfObsList,
RooDataHist& dataHist,
bool histFuncMode) ;
......
......@@ -314,6 +314,11 @@ double RooHistFunc::analyticalIntegral(Int_t code, const char* rangeName) const
return RooHistPdf::analyticalIntegral(code, rangeName, _histObsList, _depList, *_dataHist, true);
}
bool RooHistFunc::forceAnalyticalInt(const RooAbsArg& dep) const
{
return RooHistPdf::forceAnalyticalInt(_depList, dep);
}
////////////////////////////////////////////////////////////////////////////////
/// Return sampling hint for making curves of (projections) of this function
......
......@@ -261,6 +261,7 @@ double RooHistPdf::totVolume() const
}
namespace {
bool fullRange(const RooAbsArg& x, const RooAbsArg& y ,const char* range)
{
const RooAbsRealLValue *_x = dynamic_cast<const RooAbsRealLValue*>(&x);
......@@ -276,15 +277,34 @@ bool fullRange(const RooAbsArg& x, const RooAbsArg& y ,const char* range)
}
return (_x->getMin(range) == _y->getMin() && _x->getMax(range) == _y->getMax());
}
bool okayForAnalytical(RooAbsArg const& obs, RooArgSet const& allVars)
{
auto lobs = dynamic_cast<RooAbsRealLValue const*>(&obs);
if(lobs == nullptr) return false;
bool isOkayForAnalyticalInt = false;
for(RooAbsArg *var : allVars) {
if(obs.dependsOn(*var)) {
if(!lobs->isJacobianOK(*var)) return false;
isOkayForAnalyticalInt = true;
}
}
return isOkayForAnalyticalInt;
}
} // namespace
Int_t RooHistPdf::getAnalyticalIntegral(RooArgSet& allVars,
RooArgSet& analVars,
const char* rangeName,
RooArgSet const& histObsList,
RooSetProxy const& pdfObsList,
Int_t intOrder) {
RooArgSet const& pdfObsList,
Int_t intOrder)
{
// First make list of pdf observables to histogram observables
// and select only those for which the integral is over the full range
......@@ -294,7 +314,7 @@ Int_t RooHistPdf::getAnalyticalIntegral(RooArgSet& allVars,
const auto pa = pdfObsList[n];
const auto ha = histObsList[n];
if (allVars.find(*pa)) {
if (okayForAnalytical(*pa, allVars)) {
code |= 2 << n;
analVars.add(*pa);
if (fullRange(*pa, *ha, rangeName)) {
......@@ -323,7 +343,7 @@ Int_t RooHistPdf::getAnalyticalIntegral(RooArgSet& allVars,
double RooHistPdf::analyticalIntegral(Int_t code,
const char* rangeName,
RooArgSet const& histObsList,
RooSetProxy const& pdfObsList,
RooArgSet const& pdfObsList,
RooDataHist& dataHist,
bool histFuncMode) {
// Simplest scenario, full-range integration over all dependents
......@@ -383,6 +403,30 @@ double RooHistPdf::analyticalIntegral(Int_t code, const char* rangeName) const
}
bool RooHistPdf::forceAnalyticalInt(RooArgSet const& pdfObsList, const RooAbsArg& dep)
{
bool isOkayForAnalyticalInt = false;
for (RooAbsArg * obs : pdfObsList) {
if(obs->dependsOn(dep)) {
// If the observable doesn't depend linearly on the integration
// variable we will not do analytical integration.
auto lvalue = dynamic_cast<RooAbsRealLValue const*>(obs);
if(!(lvalue && lvalue->isJacobianOK(dep))) return false;
isOkayForAnalyticalInt = true;
}
}
return isOkayForAnalyticalInt;
}
bool RooHistPdf::forceAnalyticalInt(const RooAbsArg& dep) const
{
return forceAnalyticalInt(_pdfObsList, dep);
}
////////////////////////////////////////////////////////////////////////////////
/// Return sampling hint for making curves of (projections) of this function
/// as the recursive division strategy of RooCurve cannot deal efficiently
......@@ -611,4 +655,3 @@ void RooHistPdf::Streamer(TBuffer &R__b)
R__b.WriteClassBuffer(RooHistPdf::Class(),this);
}
}
......@@ -61,6 +61,7 @@ ROOT_ADD_GTEST(testGlobalObservables testGlobalObservables.cxx LIBRARIES RooFitC
ROOT_ADD_GTEST(testRooPolyFunc testRooPolyFunc.cxx LIBRARIES Gpad RooFitCore)
ROOT_ADD_GTEST(testSumW2Error testSumW2Error.cxx LIBRARIES Gpad RooFitCore)
ROOT_ADD_GTEST(testRooHist testRooHist.cxx LIBRARIES RooFitCore)
ROOT_ADD_GTEST(testRooHistPdf testRooHistPdf.cxx LIBRARIES RooFitCore)
if (roofit_multiprocess)
ROOT_ADD_GTEST(testTestStatisticsPlot TestStatistics/testPlot.cxx LIBRARIES RooFitMultiProcess RooFitCore RooFit
COPY_TO_BUILDDIR ${CMAKE_CURRENT_SOURCE_DIR}/TestStatistics/TestStatistics_ref.root)
......
// Tests for the RooHistPdf
// Authors: Jonas Rembser, CERN 03/2023
#include <RooDataHist.h>
#include <RooHistPdf.h>
#include <RooLinearVar.h>
#include <RooRealIntegral.h>
#include <RooRealVar.h>
#include <gtest/gtest.h>
// Verify that RooFit correctly uses analytic integration when having a
// RooLinearVar as the observable of a RooHistPdf.
TEST(RooHistPdf, AnalyticIntWithRooLinearVar)
{
RooRealVar x{"x", "x", 0, -10, 10};
x.setBins(10);
RooDataHist dataHist("dataHist", "dataHist", x);
for(int i = 0; i < x.numBins(); ++i) {
dataHist.set(i, 10.0, 0.0);
}
RooRealVar shift{"shift", "shift", 2.0, -10, 10};
RooRealVar slope{"slope", "slope", 1.0};
RooLinearVar xShifted{"x_shifted", "x_shifted", x, slope, shift};
RooHistPdf pdf{"pdf", "pdf", xShifted, x, dataHist};
RooRealIntegral integ{"integ", "integ", pdf, x};
EXPECT_DOUBLE_EQ(integ.getVal(), 90.);
EXPECT_EQ(integ.anaIntVars().size(), 1);
}
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment