Skip to content

[RF] Plotting shifted RooHistPdfs doesn't go well #13030

Open
@guitargeek

Description

@guitargeek

If you plot a shifted RooHistPdf, it doesn't look good. Probably because it is sampled at the wrong points, the numeric bin integrator doesn't work in this case, or both.

This is a problem that comes up relatively often in the ROOT forum, the last time here:
https://root-forum.cern.ch/t/roofit-in-root-6-28-04/55267

By looking for the keyword "shift" in the RooFit forum, one can see that this problem also came up in the past already:
https://root-forum.cern.ch/search?expanded=true&q=shift%20%23roofit-and-roostats
For example here:
https://root-forum.cern.ch/t/how-to-implement-a-horizontal-shift-for-roohistpdf/16787/6

I think I can approximately reproduce the workflow of the user that posted on the forum with this code.

void repro()
{
   using namespace RooFit;

   RooRealVar x("x", "", 1000, 1500);
   x.setBins(50);

   RooRealVar shift("shift", "", 10.0, -100, 100);

   RooFormulaVar xShifted("x_shifted", "x - shift", {x, shift});
   // Doesn't work better either...
   // RooLinearVar xShifted("x_shifted", "", x, RooConst(1.0), shift);

   std::vector<std::unique_ptr<RooDataHist>> templateHists;
   RooArgSet pdfs;
   RooArgSet yields;

   // Fill the templates
   for (std::size_t i = 0; i < 2; ++i) {
      auto suffix = std::to_string(i);

      TF1 pdf("pdf", "gaus", x.getMin(), x.getMax());
      pdf.SetParameters(1, 1200 + 100 * i, 50);
      pdf.Print();
      TH1D h("h", "", x.numBins(), x.getMin(), x.getMax());
      h.FillRandom("pdf", 10000);

      templateHists.emplace_back(std::make_unique<RooDataHist>(("template_data_hist_" + suffix).c_str(), "", x, &h));
      pdfs.addOwned(
         *new RooHistPdf(("template_hist_pdf_" + suffix).c_str(), "", xShifted, x, *templateHists.back(), 0));
      yields.addOwned(*new RooRealVar(("n_" + suffix).c_str(), "", 1000, 100, 100000));
   }

   // Construct final model
   RooAddPdf model("model", "", pdfs, yields);

   std::unique_ptr<RooDataHist> data{model.generateBinned(x)};

   // If the range is larger than the variable range, the problem is even amplified
   model.fitTo(*data, Range(0.0, 5000));

   TCanvas c1;

   auto frame = x.frame();
   data->plotOn(frame);
   model.plotOn(frame, Range(0.0, 2000.0));
   model.plotOn(frame, Components(*pdfs[0]), LineStyle(kDashed), LineColor(kRed));
   model.plotOn(frame, Components(*pdfs[1]), LineStyle(kDashed), LineColor(kBlue));
   frame->Draw();

   c1.SaveAs("plot.png");
}

The plot looks like this:
plot

It would be nice if this plotting would work better, sampling the RooHistPdf at the correct points.

Metadata

Metadata

Assignees

Type

No type

Projects

No projects

Milestone

No milestone

Relationships

None yet

Development

No branches or pull requests

Issue actions