Rf 2 0 4B_Extended Likelihood_Ranged Fit

This macro demonstrates how to set up a fit in two ranges such that it does not only fit the shapes in each region, but also takes into account the relative normalization of the two.

1. Shape fits (plain likelihood)

If you perform a fit in two ranges in RooFit, e.g. pdf->fitTo(data,Range("Range1,Range2")) it will construct a simple simultaneous fit of the two regions.

In case the pdf is not extended, i.e., a shape fit, it will only fit the shapes in the two selected ranges, and not take into account the relative predicted yields in those ranges.

In certain models (like exponential decays) and configurations (e.g. narrow ranges that are far apart), the relative normalization of the ranges may carry much more information about the function parameters than the shape of the distribution inside those ranges. Therefore, it is important to take that into account.

This is particularly important for cases where the 2-range fit is meant to be representative of a full-range fit, but with a blinded signal region inside it.

2. Shape+rate fits (extended likelihood)

Also if your pdf is already 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 interpretatin 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.

Author: Stephan Hageboeck, Wouter Verkerke
This notebook tutorial was automatically generated with ROOTBOOK-izer from the macro found in the ROOT repository on Saturday, November 28, 2020 at 10:40 AM.

In [1]:
%%cpp -d
#include "RooRealVar.h"
#include "RooExponential.h"
#include "RooGaussian.h"
#include "RooAddPdf.h"
#include "RooDataSet.h"
#include "RooPlot.h"
#include "RooExtendPdf.h"
#include "TCanvas.h"
In [2]:
%%cpp -d
// This is a workaround to make sure the namespace is used inside functions
using namespace RooFit;
In [3]:
// 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);

 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, o the 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]
 // that will introduce additional sensitivity of the likelihood to the slope parameter alpha of the exponential model through the `frac_LEFT` and `frac_RIGHT` integrals.
 //
 // In the extreme case of an exponential function and a fit in narrow LEFT and RIGHT ranges, this sensitivity may actually be larger
 // than from the shapes.
 //
 // This is also nicely demonstrated in the example below where the uncertainty on alpha is almost 5x smaller if the extended term is included.


 TCanvas* c = new TCanvas("c", "c", 2100, 700);
 c->Divide(3);
 c->cd(1);

 RooFitResult* r = model.fitTo(*data, Range("LEFT,RIGHT"), Save());

 RooPlot* frame = x.frame();
 data->plotOn(frame);
 model.plotOn(frame, VisualizeError(*r));
 model.plotOn(frame);
 model.paramOn(frame, Label("Bkg fit. Large errors since\nnormalisation ignored"));
 frame->Draw();

 c->cd(2);

 RooFitResult* r2 = extmodel.fitTo(*data, Range("LEFT,RIGHT"), Save());
 RooPlot* frame2 = x.frame();
 data->plotOn(frame2);
 extmodel.plotOn(frame2);
 extmodel.plotOn(frame2, VisualizeError(*r2));
 extmodel.paramOn(frame2, Label("Bkg fit. Normalisation\nincluded"), 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", RooArgSet(sig, model), RooArgSet(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
 // ```

 RooFitResult* r3 = modelsum.fitTo(*data, Range("LEFT,RIGHT"), Save());

 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();
RooFit v3.60 -- Developed by Wouter Verkerke and David Kirkby 
                Copyright (C) 2000-2013 NIKHEF, University of California & Stanford University
                All rights reserved, please read http://roofit.sourceforge.net/license.txt

[#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:Minization -- p.d.f. provides expected number of events, including extended term in likelihood.
[#1] INFO:Fitting -- RooAbsOptTestStatistic::ctor(nll_extmodel_modelData) constructing test statistic for sub-range named LEFT
[#1] INFO:Eval -- RooRealVar::setRange(x) new range named 'NormalizationRangeForLEFT' created with bounds [10,100]
[#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) fixing interpretation of coefficients of any RooAddPdf to full domain of observables 
[#1] INFO:Minization -- RooMinimizer::optimizeConst: activating const optimization
[#1] INFO:Minization -- RooMinimizer::optimizeConst: deactivating const optimization
RooRealVar::N = 10232.9 +/- 1125.64  L(0 - 20000) 
[#1] INFO:Fitting -- RooAbsOptTestStatistic::ctor(nll_model_modelData_LEFT) constructing test statistic for sub-range named LEFT
[#1] INFO:Eval -- RooRealVar::setRange(x) new range named 'fit_nll_model_modelData_LEFT' created with bounds [10,20]
[#1] INFO:Fitting -- RooAbsOptTestStatistic::ctor(nll_model_modelData_LEFT) fixing interpretation of coefficients of any RooAddPdf to full domain of observables 
[#1] INFO:Fitting -- RooAbsOptTestStatistic::ctor(nll_model_modelData_RIGHT) constructing test statistic for sub-range named RIGHT
[#1] INFO:Eval -- RooRealVar::setRange(x) new range named 'NormalizationRangeForRIGHT' created with bounds [10,100]
[#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_RIGHT) fixing interpretation of coefficients of any RooAddPdf to full domain of observables 
[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_model_modelData) Summation contains a RooNLLVar, using its error level
[#1] INFO:Minization -- RooMinimizer::optimizeConst: activating const optimization
 **********
 **    1 **SET PRINT           1
 **********
 **********
 **    2 **SET NOGRAD
 **********
 PARAMETER DEFINITIONS:
    NO.   NAME         VALUE      STEP SIZE      LIMITS
     1 alpha       -3.88140e-02  5.94722e-03   -1.00000e-01 -0.00000e+00
 **********
 **    3 **SET ERR         0.5
 **********
 **********
 **    4 **SET PRINT           1
 **********
 **********
 **    5 **SET STR           1
 **********
 NOW USING STRATEGY  1: TRY TO BALANCE SPEED AGAINST RELIABILITY
 **********
 **    6 **MIGRAD         500           1
 **********
 FIRST CALL TO USER FUNCTION AT NEW START POINT, WITH IFLAG=4.
 START MIGRAD MINIMIZATION.  STRATEGY  1.  CONVERGENCE WHEN EDM .LT. 1.00e-03
 FCN=11765.3 FROM MIGRAD    STATUS=INITIATE        4 CALLS           5 TOTAL
                     EDM= unknown      STRATEGY= 1      NO ERROR MATRIX       
  EXT PARAMETER               CURRENT GUESS       STEP         FIRST   
  NO.   NAME      VALUE            ERROR          SIZE      DERIVATIVE 
   1  alpha       -3.88140e-02   5.94722e-03   1.22392e-01  -2.41993e+01
                               ERR DEF= 0.5
 MIGRAD MINIMIZATION HAS CONVERGED.
 MIGRAD WILL VERIFY CONVERGENCE AND ERROR MATRIX.
 COVARIANCE MATRIX CALCULATED SUCCESSFULLY
 FCN=11764.5 FROM MIGRAD    STATUS=CONVERGED      17 CALLS          18 TOTAL
                     EDM=1.00065e-07    STRATEGY= 1      ERROR MATRIX ACCURATE 
  EXT PARAMETER                                   STEP         FIRST   
  NO.   NAME      VALUE            ERROR          SIZE      DERIVATIVE 
   1  alpha       -3.57051e-02   2.49242e-03   3.89675e-03  -6.07822e-03
                               ERR DEF= 0.5
 EXTERNAL ERROR MATRIX.    NDIM=  25    NPAR=  1    ERR DEF=0.5
  6.218e-06 
 **********
 **    7 **SET ERR         0.5
 **********
 **********
 **    8 **SET PRINT           1
 **********
 **********
 **    9 **HESSE         500
 **********
 COVARIANCE MATRIX CALCULATED SUCCESSFULLY
 FCN=11764.5 FROM HESSE     STATUS=OK              5 CALLS          23 TOTAL
                     EDM=1.00032e-07    STRATEGY= 1      ERROR MATRIX ACCURATE 
  EXT PARAMETER                                INTERNAL      INTERNAL  
  NO.   NAME      VALUE            ERROR       STEP SIZE       VALUE   
   1  alpha       -3.57051e-02   2.49241e-03   1.55870e-04   2.89944e-01
                               ERR DEF= 0.5
 EXTERNAL ERROR MATRIX.    NDIM=  25    NPAR=  1    ERR DEF=0.5
  6.218e-06 
[#1] INFO:Minization -- RooMinimizer::optimizeConst: deactivating const optimization
[#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>.setStringAttribute("fitrange", nullptr);
	- 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>.setStringAttribute("fitrange", nullptr);
	- 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:Minization -- p.d.f. provides expected number of events, including extended term in likelihood.
[#1] INFO:Fitting -- RooAbsOptTestStatistic::ctor(nll_extmodel_modelData_LEFT) constructing test statistic for sub-range named LEFT
[#1] INFO:Eval -- RooRealVar::setRange(x) new range named 'fit_nll_extmodel_modelData_LEFT' created with bounds [10,20]
[#1] INFO:Fitting -- RooAbsOptTestStatistic::ctor(nll_extmodel_modelData_LEFT) fixing interpretation of coefficients of any RooAddPdf to full domain of observables 
[#1] INFO:Fitting -- RooAbsOptTestStatistic::ctor(nll_extmodel_modelData_RIGHT) constructing test statistic for sub-range named RIGHT
[#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_RIGHT) fixing interpretation of coefficients of any RooAddPdf to full domain of observables 
[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_extmodel_modelData) Summation contains a RooNLLVar, using its error level
[#1] INFO:Minization -- RooMinimizer::optimizeConst: activating const optimization
 **********
 **    1 **SET PRINT           1
 **********
 **********
 **    2 **SET NOGRAD
 **********
 PARAMETER DEFINITIONS:
    NO.   NAME         VALUE      STEP SIZE      LIMITS
     1 N            1.02329e+04  1.12564e+03    0.00000e+00  2.00000e+04
     2 alpha       -3.57051e-02  2.49241e-03   -1.00000e-01 -0.00000e+00
 **********
 **    3 **SET ERR         0.5
 **********
 **********
 **    4 **SET PRINT           1
 **********
 **********
 **    5 **SET STR           1
 **********
 NOW USING STRATEGY  1: TRY TO BALANCE SPEED AGAINST RELIABILITY
 **********
 **    6 **MIGRAD        1000           1
 **********
 FIRST CALL TO USER FUNCTION AT NEW START POINT, WITH IFLAG=4.
 START MIGRAD MINIMIZATION.  STRATEGY  1.  CONVERGENCE WHEN EDM .LT. 1.00e-03
 FCN=-19005.9 FROM MIGRAD    STATUS=INITIATE        8 CALLS           9 TOTAL
                     EDM= unknown      STRATEGY= 1      NO ERROR MATRIX       
  EXT PARAMETER               CURRENT GUESS       STEP         FIRST   
  NO.   NAME      VALUE            ERROR          SIZE      DERIVATIVE 
   1  N            1.02329e+04   1.12564e+03   1.12835e-01   6.55035e+01
   2  alpha       -3.57051e-02   2.49241e-03   5.20493e-02   6.95485e+02
                               ERR DEF= 0.5
 MIGRAD MINIMIZATION HAS CONVERGED.
 MIGRAD WILL VERIFY CONVERGENCE AND ERROR MATRIX.
 COVARIANCE MATRIX CALCULATED SUCCESSFULLY
 FCN=-19036.9 FROM MIGRAD    STATUS=CONVERGED      37 CALLS          38 TOTAL
                     EDM=8.39226e-08    STRATEGY= 1      ERROR MATRIX ACCURATE 
  EXT PARAMETER                                   STEP         FIRST   
  NO.   NAME      VALUE            ERROR          SIZE      DERIVATIVE 
   1  N            9.98775e+03   1.49879e+02   1.41908e-03   1.65276e-02
   2  alpha       -3.99744e-02   5.59846e-04   1.08236e-03   1.09845e-02
                               ERR DEF= 0.5
 EXTERNAL ERROR MATRIX.    NDIM=  25    NPAR=  2    ERR DEF=0.5
  2.247e+04  9.166e-03 
  9.166e-03  3.134e-07 
 PARAMETER  CORRELATION COEFFICIENTS  
       NO.  GLOBAL      1      2
        1  0.10923   1.000  0.109
        2  0.10923   0.109  1.000
 **********
 **    7 **SET ERR         0.5
 **********
 **********
 **    8 **SET PRINT           1
 **********
 **********
 **    9 **HESSE        1000
 **********
 COVARIANCE MATRIX CALCULATED SUCCESSFULLY
 FCN=-19036.9 FROM HESSE     STATUS=OK             10 CALLS          48 TOTAL
                     EDM=8.39789e-08    STRATEGY= 1      ERROR MATRIX ACCURATE 
  EXT PARAMETER                                INTERNAL      INTERNAL  
  NO.   NAME      VALUE            ERROR       STEP SIZE       VALUE   
   1  N            9.98775e+03   1.49884e+02   5.67631e-05  -1.22528e-03
   2  alpha       -3.99744e-02   5.59867e-04   4.32944e-05   2.01880e-01
                               ERR DEF= 0.5
 EXTERNAL ERROR MATRIX.    NDIM=  25    NPAR=  2    ERR DEF=0.5
  2.247e+04  9.194e-03 
  9.194e-03  3.135e-07 
 PARAMETER  CORRELATION COEFFICIENTS  
       NO.  GLOBAL      1      2
        1  0.10955   1.000  0.110
        2  0.10955   0.110  1.000
[#1] INFO:Minization -- RooMinimizer::optimizeConst: deactivating const optimization
[#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>.setStringAttribute("fitrange", nullptr);
	- 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>.setStringAttribute("fitrange", nullptr);
	- 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 [-1e+30, 1e+30] of the RooGaussian 'sig' exceeds the safe range of (0, inf). Advise to limit its range.
[#1] INFO:Minization -- p.d.f. provides expected number of events, including extended term in likelihood.
[#1] INFO:Fitting -- RooAbsOptTestStatistic::ctor(nll_modelsum_modelData_LEFT) constructing test statistic for sub-range named LEFT
[#1] INFO:Eval -- RooRealVar::setRange(x) new range named 'fit_nll_modelsum_modelData_LEFT' created with bounds [10,20]
[#1] INFO:Fitting -- RooAbsOptTestStatistic::ctor(nll_modelsum_modelData_LEFT) fixing interpretation of coefficients of any RooAddPdf to full domain of observables 
[#1] INFO:Fitting -- RooAbsOptTestStatistic::ctor(nll_modelsum_modelData_RIGHT) constructing test statistic for sub-range named RIGHT
[#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_RIGHT) fixing interpretation of coefficients of any RooAddPdf to full domain of observables 
[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_modelsum_modelData) Summation contains a RooNLLVar, using its error level
[#1] INFO:Minization -- RooMinimizer::optimizeConst: activating const optimization
[#1] INFO:Minization --  The following expressions have been identified as constant and will be precalculated and cached: (sig)
[#1] INFO:Minization --  The following expressions will be evaluated in cache-and-track mode: (model)
[#1] INFO:Minization --  The following expressions have been identified as constant and will be precalculated and cached: (sig)
[#1] INFO:Minization --  The following expressions will be evaluated in cache-and-track mode: (model)
 **********
 **    1 **SET PRINT           1
 **********
 **********
 **    2 **SET NOGRAD
 **********
 PARAMETER DEFINITIONS:
    NO.   NAME         VALUE      STEP SIZE      LIMITS
     1 Nbkg         1.00000e+04  2.00000e+03    0.00000e+00  2.00000e+04
     2 Nsig         1.00000e+03  2.00000e+02    0.00000e+00  2.00000e+03
     3 alpha       -3.99744e-02  5.59867e-04   -1.00000e-01 -0.00000e+00
 **********
 **    3 **SET ERR         0.5
 **********
 **********
 **    4 **SET PRINT           1
 **********
 **********
 **    5 **SET STR           1
 **********
 NOW USING STRATEGY  1: TRY TO BALANCE SPEED AGAINST RELIABILITY
 **********
 **    6 **MIGRAD        1500           1
 **********
 FIRST CALL TO USER FUNCTION AT NEW START POINT, WITH IFLAG=4.
 START MIGRAD MINIMIZATION.  STRATEGY  1.  CONVERGENCE WHEN EDM .LT. 1.00e-03
 FCN=-19036.9 FROM MIGRAD    STATUS=INITIATE       10 CALLS          11 TOTAL
                     EDM= unknown      STRATEGY= 1      NO ERROR MATRIX       
  EXT PARAMETER               CURRENT GUESS       STEP         FIRST   
  NO.   NAME      VALUE            ERROR          SIZE      DERIVATIVE 
   1  Nbkg         1.00000e+04   2.00000e+03   2.01358e-01   5.58605e+00
   2  Nsig         1.00000e+03   2.00000e+02   2.01358e-01   3.95359e-03
   3  alpha       -3.99744e-02   5.59867e-04   1.14297e-02  -7.64048e-01
                               ERR DEF= 0.5
 MIGRAD MINIMIZATION HAS CONVERGED.
 MIGRAD WILL VERIFY CONVERGENCE AND ERROR MATRIX.
 COVARIANCE MATRIX CALCULATED SUCCESSFULLY
 FCN=-19036.9 FROM MIGRAD    STATUS=CONVERGED      87 CALLS          88 TOTAL
                     EDM=1.83501e-09    STRATEGY= 1      ERROR MATRIX ACCURATE 
  EXT PARAMETER                                   STEP         FIRST   
  NO.   NAME      VALUE            ERROR          SIZE      DERIVATIVE 
   1  Nbkg         9.98770e+03   1.49879e+02   1.41948e-03  -1.86415e-03
   2  Nsig         8.62282e-08   1.93761e+03   5.00000e-01** at limit **
   3  alpha       -3.99745e-02   5.59848e-04   1.08244e-03  -2.58583e-03
                               ERR DEF= 0.5
 EXTERNAL ERROR MATRIX.    NDIM=  25    NPAR=  3    ERR DEF=0.5
  2.247e+04 -1.156e-01  9.166e-03 
 -1.156e-01  4.533e-02 -1.401e-07 
  9.166e-03 -1.401e-07  3.134e-07 
 PARAMETER  CORRELATION COEFFICIENTS  
       NO.  GLOBAL      1      2      3
        1  0.10929   1.000 -0.004  0.109
        2  0.00371  -0.004  1.000 -0.001
        3  0.10924   0.109 -0.001  1.000
 **********
 **    7 **SET ERR         0.5
 **********
 **********
 **    8 **SET PRINT           1
 **********
 **********
 **    9 **HESSE        1500
 **********
 COVARIANCE MATRIX CALCULATED SUCCESSFULLY
 FCN=-19036.9 FROM HESSE     STATUS=OK             16 CALLS         104 TOTAL
                     EDM=1.83305e-09    STRATEGY= 1      ERROR MATRIX ACCURATE 
  EXT PARAMETER                                INTERNAL      INTERNAL  
  NO.   NAME      VALUE            ERROR       STEP SIZE       VALUE   
   1  Nbkg         9.98770e+03   1.49884e+02   5.67792e-05  -1.22967e-03
   2  Nsig         8.62282e-08   1.97068e+03   1.00000e-01  -2.47715e+03
                                 WARNING -   - ABOVE PARAMETER IS AT LIMIT.
   3  alpha       -3.99745e-02   5.59868e-04   4.32977e-05   2.01878e-01
                               ERR DEF= 0.5
 EXTERNAL ERROR MATRIX.    NDIM=  25    NPAR=  3    ERR DEF=0.5
  2.247e+04 -2.316e-02  9.194e-03 
 -2.316e-02  4.444e-02 -2.815e-08 
  9.194e-03 -2.815e-08  3.135e-07 
 PARAMETER  CORRELATION COEFFICIENTS  
       NO.  GLOBAL      1      2      3
        1  0.10955   1.000 -0.001  0.110
        2  0.00075  -0.001  1.000 -0.000
        3  0.10955   0.110 -0.000  1.000
[#1] INFO:Minization -- RooMinimizer::optimizeConst: deactivating const optimization
[#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>.setStringAttribute("fitrange", nullptr);
	- 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>.setStringAttribute("fitrange", nullptr);
	- 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'