Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf204b_extendedLikelihood_rangedFit.C File Reference

Detailed Description

View in nbviewer Open in SWAN
This macro demonstrates how to set up a fit in two ranges for plain likelihoods and extended likelihoods.

1. Shape fits (plain likelihood)

If you fit a non-extended pdf in two ranges, e.g. pdf->fitTo(data,Range("Range1,Range2")), it will fit the shapes in the two selected ranges and also take into account the relative predicted yields in those ranges.

This is useful for example to represent a full-range fit, but with a blinded signal region inside it.

2. Shape+rate fits (extended likelihood)

If your pdf is extended, i.e. measuring both the distribution in the observable as well as the event count in the fitted region, some intervention is needed to make fits in ranges work in a way that corresponds to intuition.

If an extended fit is performed in a sub-range, the observed yield is only that of the subrange, hence the expected event count will converge to a number that is smaller than what's visible in a plot. In such cases, it is often preferred to interpret the extended term with respect to the full range that's plotted, i.e., apply a correction to the extended likelihood term in such a way that the interpretation of the expected event count remains that of the full range. This can be done by applying a correcion factor (equal to the fraction of the pdf that is contained in the fitted range) in the Poisson term that represents the extended likelihood term.

If an extended likelihood fit is performed over two sub-ranges, this correction is even more important: without it, each component likelihood would have a different interpretation of the expected event count (each corresponding to the count in its own region), and a joint fit of these regions with different interpretations of the same model parameter results in a number that is not easily interpreted.

If both regions correct their interpretations such that N_expected refers to the full range, it is interpreted easily, and consistent in both regions.

This requires that the likelihood model is extended using RooAddPdf in the form SumPdf = Nsig * sigPdf + Nbkg * bkgPdf.

#include "RooRealVar.h"
#include "RooExponential.h"
#include "RooGaussian.h"
#include "RooAddPdf.h"
#include "RooDataSet.h"
#include "RooPlot.h"
#include "RooExtendPdf.h"
#include "RooFitResult.h"
#include "TCanvas.h"
{
using namespace RooFit;
// PART 1: Background-only fits
// ----------------------------
// Build plain exponential model
RooRealVar x("x", "x", 10, 100);
RooRealVar alpha("alpha", "alpha", -0.04, -0.1, -0.0);
RooExponential model("model", "Exponential model", x, alpha);
// Define side band regions and full range
x.setRange("LEFT",10,20);
x.setRange("RIGHT",60,100);
x.setRange("FULL",10,100);
std::unique_ptr<RooDataSet> data{model.generate(x, 10000)};
// Construct an extended pdf, which measures the event count N **on the full range**.
// If the actual domain of x that is fitted is identical to FULL, this has no affect.
//
// If the fitted domain is a subset of `FULL`, though, the expected event count is divided by
// \f[
// \mathrm{frac} = \frac{
// \int_{\mathrm{Fit range}} \mathrm{model}(x) \; \mathrm{d}x }{
// \int_{\mathrm{Full range}} \mathrm{model}(x) \; \mathrm{d}x }.
// \f]
// `N` will therefore return the count extrapolated to the full range instead of the fit range.
//
// **Note**: When using a RooAddPdf for extending the likelihood, the same effect can be achieved with
// [RooAddPdf::fixCoefRange()](https://root.cern.ch/doc/master/classRooAddPdf.html#ab631caf4b59e4c4221f8967aecbf2a65),
RooRealVar N("N", "Extended term", 0, 20000);
RooExtendPdf extmodel("extmodel", "Extended model", model, N, "FULL");
// It can be instructive to fit the above model to either the LEFT or RIGHT range. `N` should approximately converge to the expected number
// of events in the full range. One may try to leave out `"FULL"` in the constructor, or the interpretation of `N` changes.
extmodel.fitTo(*data, Range("LEFT"), PrintLevel(-1));
N.Print();
// If we now do a simultaneous fit to the extended model, instead of the original model, the LEFT and RIGHT range will each correct their local
// `N` such that it refers to the `FULL` range.
//
// This joint fit of the extmodel will include (w.r.t. the plain model fit) a product of extended terms
// \f[
// L_\mathrm{ext} = L
// \cdot \mathrm{Poisson} \left( N_\mathrm{obs}^\mathrm{LEFT} | N_\mathrm{exp} / \mathrm{frac LEFT} \right)
// \cdot \mathrm{Poisson} \left( N_\mathrm{obs}^\mathrm{RIGHT} | N_\mathrm{exp} / \mathrm{frac RIGHT} \right)
// \f]
TCanvas* c = new TCanvas("c", "c", 2100, 700);
c->Divide(3);
c->cd(1);
std::unique_ptr<RooFitResult> r{model.fitTo(*data, Range("LEFT,RIGHT"), PrintLevel(-1), Save())};
r->Print();
RooPlot* frame = x.frame();
data->plotOn(frame);
model.plotOn(frame, VisualizeError(*r));
model.plotOn(frame);
model.paramOn(frame, Label("Non-extended fit"));
frame->Draw();
c->cd(2);
std::unique_ptr<RooFitResult> r2{extmodel.fitTo(*data, Range("LEFT,RIGHT"), PrintLevel(-1), Save())};
r2->Print();
RooPlot* frame2 = x.frame();
data->plotOn(frame2);
extmodel.plotOn(frame2);
extmodel.plotOn(frame2, VisualizeError(*r2));
extmodel.paramOn(frame2, Label("Extended fit"), Layout(0.4,0.95));
frame2->Draw();
// PART 2: Extending with RooAddPdf
// --------------------------------
//
// Now we repeat the above exercise, but instead of explicitly adding an extended term to a single shape pdf (RooExponential),
// we assume that we already have an extended likelihood model in the form of a RooAddPdf constructed in the form `Nsig * sigPdf + Nbkg * bkgPdf`.
//
// We add a Gaussian to the previously defined exponential background.
// The signal shape parameters are chosen constant, since the signal region is entirely in the blinded region, i.e., the fit has no sensitivity.
c->cd(3);
RooRealVar Nsig("Nsig", "Number of signal events", 1000, 0, 2000);
RooRealVar Nbkg("Nbkg", "Number of background events", 10000, 0, 20000);
RooRealVar mean("mean", "Mean of signal model", 40.);
RooRealVar width("width", "Width of signal model", 5.);
RooGaussian sig("sig", "Signal model", x, mean, width);
RooAddPdf modelsum("modelsum", "NSig*signal + NBkg*background", {sig, model}, {Nsig, Nbkg});
// This model will automatically insert the correction factor for the reinterpretation of Nsig and Nnbkg in the full ranges.
//
// When this happens, it reports this with lines like the following:
// ```
// [#1] INFO:Fitting -- RooAbsOptTestStatistic::ctor(nll_modelsum_modelsumData_LEFT) fixing interpretation of coefficients of any RooAddPdf to full domain of observables
// [#1] INFO:Fitting -- RooAbsOptTestStatistic::ctor(nll_modelsum_modelsumData_RIGHT) fixing interpretation of coefficients of any RooAddPdf to full domain of observables
// ```
std::unique_ptr<RooFitResult> r3{modelsum.fitTo(*data, Range("LEFT,RIGHT"), PrintLevel(-1), Save())};
r3->Print();
RooPlot* frame3 = x.frame();
data->plotOn(frame3);
modelsum.plotOn(frame3);
modelsum.plotOn(frame3, VisualizeError(*r3));
modelsum.paramOn(frame3, Label("S+B fit with RooAddPdf"), Layout(0.3,0.95));
frame3->Draw();
c->Draw();
}
#define c(i)
Definition RSha256.hxx:101
#define N
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t r
Option_t Option_t width
RooFit::OwningPtr< RooFitResult > fitTo(RooAbsData &data, CmdArgs_t const &... cmdArgs)
Fit PDF to given dataset.
Definition RooAbsPdf.h:156
Efficient implementation of a sum of PDFs of the form.
Definition RooAddPdf.h:33
Exponential PDF.
RooExtendPdf is a wrapper around an existing PDF that adds a parameteric extended likelihood term to ...
Plain Gaussian p.d.f.
Definition RooGaussian.h:24
A RooPlot is a plot frame and a container for graphics objects within that frame.
Definition RooPlot.h:43
static RooPlot * frame(const RooAbsRealLValue &var, double xmin, double xmax, Int_t nBins)
Create a new frame for a given variable in x.
Definition RooPlot.cxx:239
void Draw(Option_t *options=nullptr) override
Draw this plot and all of the elements it contains.
Definition RooPlot.cxx:652
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:37
The Canvas class.
Definition TCanvas.h:23
RooCmdArg Layout(double xmin, double xmax=0.99, double ymin=0.95)
RooCmdArg Save(bool flag=true)
RooCmdArg PrintLevel(Int_t code)
RooCmdArg VisualizeError(const RooDataSet &paramData, double Z=1)
Double_t x[n]
Definition legend1.C:17
const Rcpp::internal::NamedPlaceHolder & Label
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition JSONIO.h:26
Ta Range(0, 0, 1, 1)
[#1] INFO:Eval -- RooRealVar::setRange(x) new range named 'LEFT' created with bounds [10,20]
[#1] INFO:Eval -- RooRealVar::setRange(x) new range named 'RIGHT' created with bounds [60,100]
[#1] INFO:Eval -- RooRealVar::setRange(x) new range named 'FULL' created with bounds [10,100]
[#1] INFO:Minimization -- p.d.f. provides expected number of events, including extended term in likelihood.
[#1] INFO:Eval -- RooRealVar::setRange(x) new range named 'fit_nll_extmodel_modelData' created with bounds [10,20]
[#1] INFO:Fitting -- RooAbsOptTestStatistic::ctor(nll_extmodel_modelData) constructing test statistic for sub-range named LEFT
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
RooRealVar::N = 3395.66 +/- 58.2751 L(0 - 20000)
[#1] INFO:Eval -- RooRealVar::setRange(x) new range named 'fit_nll_model_modelData_LEFT' created with bounds [10,20]
[#1] INFO:Eval -- RooRealVar::setRange(x) new range named 'fit_nll_model_modelData_RIGHT' created with bounds [60,100]
[#1] INFO:Fitting -- RooAbsOptTestStatistic::ctor(nll_model_modelData) constructing test statistic for sub-range named LEFT,RIGHT
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
RooFitResult: minimized FCN value: 14265.9, estimated distance to minimum: 1.48012e-07
covariance matrix quality: Full, accurate covariance matrix
Status : MINIMIZE=0 HESSE=0
Floating Parameter FinalValue +/- Error
-------------------- --------------------------
alpha -3.9974e-02 +/- 5.60e-04
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f was fitted in a subrange and no explicit Range() and NormRange() was specified. Plotting / normalising in fit range. To override, do one of the following
- Clear the automatic fit range attribute: <pdf>.removeStringAttribute("fitrange");
- Explicitly specify the plotting range: Range("<rangeName>").
- Explicitly specify where to compute the normalisation: NormRange("<rangeName>").
The default (full) range can be denoted with Range("") / NormRange("").
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) only plotting range 'fit_nll_model_modelData_LEFT,fit_nll_model_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_model_modelData_LEFT,fit_nll_model_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) only plotting range 'fit_nll_model_modelData_LEFT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_model_modelData_LEFT,fit_nll_model_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) only plotting range 'fit_nll_model_modelData_LEFT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_model_modelData_LEFT,fit_nll_model_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) only plotting range 'fit_nll_model_modelData_LEFT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_model_modelData_LEFT,fit_nll_model_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) only plotting range 'fit_nll_model_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_model_modelData_LEFT,fit_nll_model_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) only plotting range 'fit_nll_model_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_model_modelData_LEFT,fit_nll_model_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) only plotting range 'fit_nll_model_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_model_modelData_LEFT,fit_nll_model_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f was fitted in a subrange and no explicit Range() and NormRange() was specified. Plotting / normalising in fit range. To override, do one of the following
- Clear the automatic fit range attribute: <pdf>.removeStringAttribute("fitrange");
- Explicitly specify the plotting range: Range("<rangeName>").
- Explicitly specify where to compute the normalisation: NormRange("<rangeName>").
The default (full) range can be denoted with Range("") / NormRange("").
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) only plotting range 'fit_nll_model_modelData_LEFT,fit_nll_model_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_model_modelData_LEFT,fit_nll_model_modelData_RIGHT'
[#1] INFO:Minimization -- p.d.f. provides expected number of events, including extended term in likelihood.
[#1] INFO:Eval -- RooRealVar::setRange(x) new range named 'fit_nll_extmodel_modelData_LEFT' created with bounds [10,20]
[#1] INFO:Eval -- RooRealVar::setRange(x) new range named 'fit_nll_extmodel_modelData_RIGHT' created with bounds [60,100]
[#1] INFO:Fitting -- RooAbsOptTestStatistic::ctor(nll_extmodel_modelData) constructing test statistic for sub-range named LEFT,RIGHT
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
RooFitResult: minimized FCN value: -15500.6, estimated distance to minimum: 0.000427967
covariance matrix quality: Full, accurate covariance matrix
Status : MINIMIZE=0 HESSE=0
Floating Parameter FinalValue +/- Error
-------------------- --------------------------
N 4.4939e+03 +/- 6.70e+01
alpha -4.8258e-02 +/- 8.32e-04
[#1] INFO:Plotting -- RooAbsPdf::plotOn(extmodel) p.d.f was fitted in a subrange and no explicit Range() and NormRange() was specified. Plotting / normalising in fit range. To override, do one of the following
- Clear the automatic fit range attribute: <pdf>.removeStringAttribute("fitrange");
- Explicitly specify the plotting range: Range("<rangeName>").
- Explicitly specify where to compute the normalisation: NormRange("<rangeName>").
The default (full) range can be denoted with Range("") / NormRange("").
[#1] INFO:Plotting -- RooAbsPdf::plotOn(extmodel) only plotting range 'fit_nll_extmodel_modelData_LEFT,fit_nll_extmodel_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(extmodel) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_extmodel_modelData_LEFT,fit_nll_extmodel_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(extmodel) p.d.f was fitted in a subrange and no explicit Range() and NormRange() was specified. Plotting / normalising in fit range. To override, do one of the following
- Clear the automatic fit range attribute: <pdf>.removeStringAttribute("fitrange");
- Explicitly specify the plotting range: Range("<rangeName>").
- Explicitly specify where to compute the normalisation: NormRange("<rangeName>").
The default (full) range can be denoted with Range("") / NormRange("").
[#1] INFO:Plotting -- RooAbsPdf::plotOn(extmodel) only plotting range 'fit_nll_extmodel_modelData_LEFT,fit_nll_extmodel_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(extmodel) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_extmodel_modelData_LEFT,fit_nll_extmodel_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(extmodel) only plotting range 'fit_nll_extmodel_modelData_LEFT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(extmodel) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_extmodel_modelData_LEFT,fit_nll_extmodel_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(extmodel) only plotting range 'fit_nll_extmodel_modelData_LEFT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(extmodel) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_extmodel_modelData_LEFT,fit_nll_extmodel_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(extmodel) only plotting range 'fit_nll_extmodel_modelData_LEFT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(extmodel) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_extmodel_modelData_LEFT,fit_nll_extmodel_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(extmodel) only plotting range 'fit_nll_extmodel_modelData_LEFT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(extmodel) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_extmodel_modelData_LEFT,fit_nll_extmodel_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(extmodel) only plotting range 'fit_nll_extmodel_modelData_LEFT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(extmodel) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_extmodel_modelData_LEFT,fit_nll_extmodel_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(extmodel) only plotting range 'fit_nll_extmodel_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(extmodel) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_extmodel_modelData_LEFT,fit_nll_extmodel_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(extmodel) only plotting range 'fit_nll_extmodel_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(extmodel) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_extmodel_modelData_LEFT,fit_nll_extmodel_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(extmodel) only plotting range 'fit_nll_extmodel_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(extmodel) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_extmodel_modelData_LEFT,fit_nll_extmodel_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(extmodel) only plotting range 'fit_nll_extmodel_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(extmodel) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_extmodel_modelData_LEFT,fit_nll_extmodel_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(extmodel) only plotting range 'fit_nll_extmodel_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(extmodel) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_extmodel_modelData_LEFT,fit_nll_extmodel_modelData_RIGHT'
[#0] WARNING:InputArguments -- The parameter 'width' with range [-inf, inf] of the RooGaussian 'sig' exceeds the safe range of (0, inf). Advise to limit its range.
[#1] INFO:Minimization -- p.d.f. provides expected number of events, including extended term in likelihood.
[#1] INFO:Eval -- RooRealVar::setRange(x) new range named 'fit_nll_modelsum_modelData_LEFT' created with bounds [10,20]
[#1] INFO:Eval -- RooRealVar::setRange(x) new range named 'fit_nll_modelsum_modelData_RIGHT' created with bounds [60,100]
[#1] INFO:Fitting -- RooAbsOptTestStatistic::ctor(nll_modelsum_modelData) constructing test statistic for sub-range named LEFT,RIGHT
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:Minimization -- The following expressions have been identified as constant and will be precalculated and cached: (sig)
[#1] INFO:Minimization -- The following expressions will be evaluated in cache-and-track mode: (model)
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
RooFitResult: minimized FCN value: -19036.9, estimated distance to minimum: 1.041e-09
covariance matrix quality: Full, accurate covariance matrix
Status : MINIMIZE=0 HESSE=0
Floating Parameter FinalValue +/- Error
-------------------- --------------------------
Nbkg 9.9877e+03 +/- 1.50e+02
Nsig 4.4502e-06 +/- 1.97e+03
alpha -3.9975e-02 +/- 5.60e-04
[#1] INFO:Plotting -- RooAbsPdf::plotOn(modelsum) p.d.f was fitted in a subrange and no explicit Range() and NormRange() was specified. Plotting / normalising in fit range. To override, do one of the following
- Clear the automatic fit range attribute: <pdf>.removeStringAttribute("fitrange");
- Explicitly specify the plotting range: Range("<rangeName>").
- Explicitly specify where to compute the normalisation: NormRange("<rangeName>").
The default (full) range can be denoted with Range("") / NormRange("").
[#1] INFO:Plotting -- RooAbsPdf::plotOn(modelsum) only plotting range 'fit_nll_modelsum_modelData_LEFT,fit_nll_modelsum_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(modelsum) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_modelsum_modelData_LEFT,fit_nll_modelsum_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(modelsum) p.d.f was fitted in a subrange and no explicit Range() and NormRange() was specified. Plotting / normalising in fit range. To override, do one of the following
- Clear the automatic fit range attribute: <pdf>.removeStringAttribute("fitrange");
- Explicitly specify the plotting range: Range("<rangeName>").
- Explicitly specify where to compute the normalisation: NormRange("<rangeName>").
The default (full) range can be denoted with Range("") / NormRange("").
[#1] INFO:Plotting -- RooAbsPdf::plotOn(modelsum) only plotting range 'fit_nll_modelsum_modelData_LEFT,fit_nll_modelsum_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(modelsum) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_modelsum_modelData_LEFT,fit_nll_modelsum_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(modelsum) only plotting range 'fit_nll_modelsum_modelData_LEFT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(modelsum) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_modelsum_modelData_LEFT,fit_nll_modelsum_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(modelsum) only plotting range 'fit_nll_modelsum_modelData_LEFT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(modelsum) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_modelsum_modelData_LEFT,fit_nll_modelsum_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(modelsum) only plotting range 'fit_nll_modelsum_modelData_LEFT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(modelsum) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_modelsum_modelData_LEFT,fit_nll_modelsum_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(modelsum) only plotting range 'fit_nll_modelsum_modelData_LEFT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(modelsum) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_modelsum_modelData_LEFT,fit_nll_modelsum_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(modelsum) only plotting range 'fit_nll_modelsum_modelData_LEFT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(modelsum) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_modelsum_modelData_LEFT,fit_nll_modelsum_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(modelsum) only plotting range 'fit_nll_modelsum_modelData_LEFT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(modelsum) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_modelsum_modelData_LEFT,fit_nll_modelsum_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(modelsum) only plotting range 'fit_nll_modelsum_modelData_LEFT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(modelsum) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_modelsum_modelData_LEFT,fit_nll_modelsum_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(modelsum) only plotting range 'fit_nll_modelsum_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(modelsum) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_modelsum_modelData_LEFT,fit_nll_modelsum_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(modelsum) only plotting range 'fit_nll_modelsum_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(modelsum) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_modelsum_modelData_LEFT,fit_nll_modelsum_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(modelsum) only plotting range 'fit_nll_modelsum_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(modelsum) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_modelsum_modelData_LEFT,fit_nll_modelsum_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(modelsum) only plotting range 'fit_nll_modelsum_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(modelsum) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_modelsum_modelData_LEFT,fit_nll_modelsum_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(modelsum) only plotting range 'fit_nll_modelsum_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(modelsum) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_modelsum_modelData_LEFT,fit_nll_modelsum_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(modelsum) only plotting range 'fit_nll_modelsum_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(modelsum) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_modelsum_modelData_LEFT,fit_nll_modelsum_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(modelsum) only plotting range 'fit_nll_modelsum_modelData_RIGHT'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(modelsum) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_modelsum_modelData_LEFT,fit_nll_modelsum_modelData_RIGHT'
Authors
Stephan Hageboeck, Wouter Verkerke

Definition in file rf204b_extendedLikelihood_rangedFit.C.