Skip to content

Commit f1473cd

Browse files
qogksdnrguitargeek
authored andcommitted
[RF] Fixing issues in RooSimultaneous with ConditionalObservables
`RooSimultaneous` is a PDF class that performs the simultaneous fitting of multiple pairs of a PDF and a dataset. `RooFit::ConditionalObservables(RooAbsArgs...)` is a command argument to the `fitTo()` method of the `RooAbsPdf` classes, which specifies variables that should be excluded from the normalization set of a PDF when fitting. If I use the command argument for fitting using the class, the class does not deal with the conditional variable well, and the conditional variables become floating parameters. This should be fixed. Closes #19166
1 parent 99f8c33 commit f1473cd

File tree

3 files changed

+66
-1
lines changed

3 files changed

+66
-1
lines changed

roofit/roofitcore/src/FitHelpers.cxx

Lines changed: 11 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -745,7 +745,17 @@ std::unique_ptr<RooAbsReal> createNLL(RooAbsPdf &pdf, RooAbsData &data, const Ro
745745

746746
RooArgSet normSet;
747747
pdf.getObservables(data.get(), normSet);
748-
normSet.remove(projDeps, true, true);
748+
749+
if (dynamic_cast<RooSimultaneous const*>(&pdf)) {
750+
for (auto i : projDeps) {
751+
auto res = normSet.find(i->GetName());
752+
if (res != nullptr) {
753+
res->setAttribute("__conditional__");
754+
}
755+
}
756+
} else {
757+
normSet.remove(projDeps);
758+
}
749759

750760
pdf.setAttribute("SplitRange", splitRange);
751761
pdf.setStringAttribute("RangeName", rangeName);

roofit/roofitcore/src/RooSimultaneous.cxx

Lines changed: 4 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -1217,6 +1217,10 @@ RooSimultaneous::compileForNormSet(RooArgSet const &normSet, RooFit::Detail::Com
12171217

12181218
std::unique_ptr<RooArgSet> pdfNormSet{
12191219
std::unique_ptr<RooArgSet>(pdfClone->getVariables())->selectByAttrib("__obs__", true)};
1220+
std::unique_ptr<RooArgSet> condVarSet{
1221+
std::unique_ptr<RooArgSet>(pdfClone->getVariables())->selectByAttrib("__conditional__", true)};
1222+
1223+
pdfNormSet->remove(*condVarSet, true, true);
12201224

12211225
if (rangeName) {
12221226
pdfClone->setNormRange(RooHelpers::getRangeNameForSimComponent(rangeName, splitRange, catName).c_str());

roofit/roofitcore/test/testRooSimultaneous.cxx

Lines changed: 51 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -456,3 +456,54 @@ TEST(RooSimultaneous, DuplicateExtendedPdfs)
456456
EXPECT_FLOAT_EQ(simPdfVal, 0.05);
457457
EXPECT_DOUBLE_EQ(simPdfVal, simPdfRef.getVal(normSet));
458458
}
459+
460+
/// GitHub issue #19166.
461+
/// RooSimultaneous::fitTo() should work well with ConditionalOvservables.
462+
/// We test this by performing the similar test to a test case in
463+
/// RooSimultaneousSingleChannelCrossCheck.
464+
TEST_P(TestStatisticTest, RooSimultaneousSingleChannelCrossCheckWithCondVar)
465+
{
466+
using namespace RooFit;
467+
468+
// silence log output
469+
RooHelpers::LocalChangeMsgLevel changeMsgLvl(RooFit::WARNING);
470+
471+
RooWorkspace ws;
472+
ws.factory("Uniform::uniform(width[1, 3])");
473+
ws.factory("expr::sigma(\"@0*@1\", width, width_scale[0.8, 0.5, 1.5])");
474+
ws.factory("Gaussian::model(x[-20, 20], mean[3, -20., 20.], sigma)");
475+
476+
RooRealVar &width = *ws.var("width");
477+
RooAbsPdf &widthModel = *ws.pdf("uniform");
478+
std::unique_ptr<RooDataSet> protoData{widthModel.generate(width, 100)};
479+
480+
RooRealVar &x = *ws.var("x");
481+
RooAbsPdf &model = *ws.pdf("model");
482+
483+
std::unique_ptr<RooDataSet> data{model.generate(x, *protoData)};
484+
485+
RooCategory cat("cat", "cat");
486+
cat.defineType("physics");
487+
488+
RooArgSet params;
489+
RooArgSet initialParams;
490+
model.getParameters(data->get(), params);
491+
params.snapshot(initialParams);
492+
493+
RooSimultaneous modelSim("modelSim", "modelSim", {{"physics", &model}}, cat);
494+
495+
RooDataSet combData("combData", "combData", {x, width}, Index(cat), Import({{"physics", data.get()}}));
496+
497+
using namespace RooFit;
498+
499+
params.assign(initialParams);
500+
std::unique_ptr<RooFitResult> resDirect{
501+
model.fitTo(*data, ConditionalObservables(width), Save(), PrintLevel(-1), _evalBackend)};
502+
503+
params.assign(initialParams);
504+
std::unique_ptr<RooFitResult> resSimWrapped{
505+
modelSim.fitTo(combData, ConditionalObservables(width), Save(), PrintLevel(-1), _evalBackend)};
506+
507+
EXPECT_TRUE(resSimWrapped->isIdentical(*resDirect))
508+
<< "Inconsistency in RooSimultaneous wrapping with ConditionalObservables";
509+
}

0 commit comments

Comments
 (0)