-
Notifications
You must be signed in to change notification settings - Fork 1.2k
New issue
Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.
By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.
Already on GitHub? Sign in to your account
[RF] Make plotting of histograms work for all orders of magnitude #13316
Conversation
847e319
to
87e3f43
Compare
87e3f43
to
47be0ec
Compare
Starting build on |
Build failed on ROOT-ubuntu2004/python3. Failing tests: |
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
LGTM!
It makes sense to use a relative epsilon value instead of an absolute one
The RooFit plotting did not work well for variables of small orders of magnitude, because there was a hardcoded check with an absolute value to see if two points are considered the same: `x2 - x1 < 1e-20`. If this condition is true, the point at x2 was not plotted. This meant that some points that were actually different were skipped. The fix is to change this condition to a relative one, considering xmin and xmax. The code to generate sampling hints for histograms was also moved inside RooCurve, such that it's easier to syncronize the relevant magic constants. Here is a demo of a plot that didn't work before but now works: ``` void repro() { int nBins = 4; double scale = 1e-31; double xmin = 0. * scale; double xmax = 1. * scale; double binWidth = (xmax - xmin) / nBins; // Fill the bin boundaries std::vector<double> binBoundaries(nBins + 1); for (int i = 0; i <= nBins; ++i) { binBoundaries[i] = i * binWidth; } // The RooParametricStepFunction needs a TArrayD TArrayD step_edges{int(binBoundaries.size()), binBoundaries.data()}; RooRealVar Gamma("Gamma", "Gamma", xmin, xmax); RooArgList step_values; step_values.addOwned(std::make_unique<RooRealVar>("coef1", "coef1", 0.0 / scale, 0.0, 1.0 / scale)); step_values.addOwned(std::make_unique<RooRealVar>("coef2", "coef2", 0.1 / scale, 0.0, 1.0 / scale)); step_values.addOwned(std::make_unique<RooRealVar>("coef3", "coef3", 0.1 / scale, 0.0, 1.0 / scale)); RooParametricStepFunction prior("prior", "Prior for decay rate", Gamma, step_values, step_edges, nBins); RooPlot *xframe = Gamma.frame(xmin, xmax); prior.plotOn(xframe); auto c1 = new TCanvas(); xframe->Draw(); c1->SaveAs("plot.png"); } ``` This bugfix was inspired by the following forum post: https://root-forum.cern.ch/t/rooparametricstepfunction/55820
47be0ec
to
976a55d
Compare
Starting build on |
Build failed on mac11/noimt. Failing tests: |
Build failed on ROOT-debian10-i386/soversion. Errors:
|
The RooFit plotting did not work well for variables of small orders of
magnitude, because there was a hardcoded check with an absolute value to
see if two points are considered the same:
x2 - x1 < 1e-20
. If thiscondition is true, the point at x2 was not plotted. This meant that some
points that were actually different were skipped.
The fix is to change this condition to a relative one, considering xmin
and xmax. The code to generate sampling hints for histograms was also
moved inside RooCurve, such that it's easier to syncronize the relevant
magic constants.
Here is a demo of a plot that didn't work before but now works:
This bugfix was inspired by the following forum post:
https://root-forum.cern.ch/t/rooparametricstepfunction/55820