Skip to content

[RF] Buggy range overlap check in createNLL when SplitRange option is used #11396

@AlkaidCheng

Description

@AlkaidCheng
  • Checked for duplicates

Describe the bug

Since ROOT 6.26, there is a new range overlap check when calling the createNLL from a RooAbsPdf instance with a Range argument that contains multiple range names. However, it ignores the fact that when the SplitRange option is used, the name of the range to be checked should be appended by the appropriate category label from each category of the simultaneous pdf. This causes an exception to be falsely raised even if the ranges do not overlap. This is because all the named ranges will return the full observable range since it fetches the range from a range name that does not exist.

To Reproduce

In PyROOT, the issue can be reproduced using

import ROOT

ws_cat1 = ROOT.RooWorkspace("ws_cat1")
ws_cat1.factory("Gaussian::pdf_cat1(x_cat1[0,10],mu_cat1[4,0,10],sigma_cat1[1.0,0.1,10.0])")
pdf_cat1 = ws_cat1.pdf("pdf_cat1")
x_cat1 = ws_cat1.var("x_cat1")
x_cat1.setRange("SideBandLo_cat1", 0, 2)
x_cat1.setRange("SideBandHi_cat1", 6, 10)
ds_cat1 = pdf_cat1.generate(ROOT.RooArgSet(x_cat1), 11000)

ws_cat2 = ROOT.RooWorkspace("ws_cat2")
ws_cat2.factory("Gaussian::pdf_cat2(x_cat2[0,10],mu_cat2[6,0,10],sigma_cat2[1.0,0.1,10.0])")
pdf_cat2 = ws_cat2.pdf("pdf_cat2")
x_cat2 = ws_cat2.var("x_cat2")
x_cat2.setRange("SideBandLo_cat2", 0, 4)
x_cat2.setRange("SideBandHi_cat2", 8, 10)
ds_cat2 = pdf_cat2.generate(ROOT.RooArgSet(x_cat2), 11000)

index_cat = ROOT.RooCategory("cat", "cat")
index_cat.defineType("cat1")
index_cat.defineType("cat2")

sim_pdf = ROOT.RooSimultaneous("sim_pdf", "", index_cat)
sim_pdf.addPdf(pdf_cat1, "cat1")
sim_pdf.addPdf(pdf_cat2, "cat2")

ROOT.gInterpreter.GenerateDictionary("std::map<std::string, RooDataSet*>", "map;string;RooDataSet.h")
ROOT.gInterpreter.GenerateDictionary("std::pair<std::string, RooDataSet*>", "map;string;RooDataSet.h")
dsmap = ROOT.std.map('string, RooDataSet*')()
dsmap.keepalive = list()
dsmap.keepalive.append(ds_cat1)
dsmap.keepalive.append(ds_cat2)
dsmap.insert(dsmap.begin(), ROOT.std.pair("const string, RooDataSet*")("cat1", ds_cat1))
dsmap.insert(dsmap.begin(), ROOT.std.pair("const string, RooDataSet*")("cat2", ds_cat2))
combData = ROOT.RooDataSet("combData", "", ROOT.RooArgSet(x_cat1, x_cat2),
                           ROOT.RooFit.Index(index_cat),
                           ROOT.RooFit.Import(dsmap))

nll = sim_pdf.createNLL(combData, ROOT.RooFit.Range("SideBandLo,SideBandHi"), ROOT.RooFit.SplitRange())

The last step raises the error "runtime_error: Error in RooAbsPdf::createNLL! The ranges SideBandLo,SideBandHi are overlapping!" when using ROOT 6.26+.

Metadata

Metadata

Assignees

Type

No type

Projects

No projects

Milestone

Relationships

None yet

Development

No branches or pull requests

Issue actions