Rf 7 0 1_Efficiencyfit

Special pdf's: unbinned maximum likelihood fit of an efficiency eff(x) function

to a dataset D(x,cut), where cut is a category encoding a selection, of which the efficiency as function of x should be described by eff(x)

Author: 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 11:06 AM.

In [1]:
%%cpp -d
#include "RooRealVar.h"
#include "RooDataSet.h"
#include "RooGaussian.h"
#include "RooConstVar.h"
#include "RooFormulaVar.h"
#include "RooProdPdf.h"
#include "RooEfficiency.h"
#include "RooPolynomial.h"
#include "RooCategory.h"
#include "TCanvas.h"
#include "TAxis.h"
#include "RooPlot.h"
In [2]:
%%cpp -d
// This is a workaround to make sure the namespace is used inside functions
using namespace RooFit;

Construct efficiency function e(x)

Declare variables x,mean,sigma with associated name, title, initial value and allowed range

In [3]:
RooRealVar x("x", "x", -10, 10);
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

Efficiency function eff(x;a,b)

In [4]:
RooRealVar a("a", "a", 0.4, 0, 1);
RooRealVar b("b", "b", 5);
RooRealVar c("c", "c", -1, -10, 10);
RooFormulaVar effFunc("effFunc", "(1-a)+a*cos((x-c)/b)", RooArgList(a, b, c, x));

Construct conditional efficiency pdf e(cut|x)

Acceptance state cut (1 or 0)

In [5]:
RooCategory cut("cut", "cutr", { {"accept", 1}, {"reject", 0} });

Construct efficiency pdf eff(cut|x)

In [6]:
RooEfficiency effPdf("effPdf", "effPdf", effFunc, cut, "accept");

Generate data (x, cut) from a toy model

Construct global shape pdf shape(x) and product model(x,cut) = eff(cut|x)*shape(x) (These are only needed to generate some toy MC here to be used later)

In [7]:
RooPolynomial shapePdf("shapePdf", "shapePdf", x, RooConst(-0.095));
RooProdPdf model("model", "model", shapePdf, Conditional(effPdf, cut));

Generate some toy data from model

In [8]:
RooDataSet *data = model.generate(RooArgSet(x, cut), 10000);

Fit conditional efficiency pdf to data

Fit conditional efficiency pdf to data

In [9]:
effPdf.fitTo(*data, ConditionalObservables(x));
[#1] INFO:Minization -- RooMinimizer::optimizeConst: activating const optimization
 **********
 **    1 **SET PRINT           1
 **********
 **********
 **    2 **SET NOGRAD
 **********
 PARAMETER DEFINITIONS:
    NO.   NAME         VALUE      STEP SIZE      LIMITS
     1 a            4.00000e-01  1.00000e-01    0.00000e+00  1.00000e+00
     2 c           -1.00000e+00  2.00000e+00   -1.00000e+01  1.00000e+01
 **********
 **    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=3887.11 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  a            4.00000e-01   1.00000e-01   2.05758e-01   8.51804e+01
   2  c           -1.00000e+00   2.00000e+00   2.02430e-01   7.00873e+01
                               ERR DEF= 0.5
 MIGRAD MINIMIZATION HAS CONVERGED.
 MIGRAD WILL VERIFY CONVERGENCE AND ERROR MATRIX.
 COVARIANCE MATRIX CALCULATED SUCCESSFULLY
 FCN=3886.26 FROM MIGRAD    STATUS=CONVERGED      30 CALLS          31 TOTAL
                     EDM=2.277e-05    STRATEGY= 1      ERROR MATRIX ACCURATE 
  EXT PARAMETER                                   STEP         FIRST   
  NO.   NAME      VALUE            ERROR          SIZE      DERIVATIVE 
   1  a            3.89817e-01   8.14984e-03   6.58343e-04  -2.65715e-01
   2  c           -9.91750e-01   6.56208e-02   2.58326e-04  -6.57933e-01
                               ERR DEF= 0.5
 EXTERNAL ERROR MATRIX.    NDIM=  25    NPAR=  2    ERR DEF=0.5
  6.643e-05 -2.189e-04 
 -2.189e-04  4.306e-03 
 PARAMETER  CORRELATION COEFFICIENTS  
       NO.  GLOBAL      1      2
        1  0.40934   1.000 -0.409
        2  0.40934  -0.409  1.000
 **********
 **    7 **SET ERR         0.5
 **********
 **********
 **    8 **SET PRINT           1
 **********
 **********
 **    9 **HESSE        1000
 **********
 COVARIANCE MATRIX CALCULATED SUCCESSFULLY
 FCN=3886.26 FROM HESSE     STATUS=OK             10 CALLS          41 TOTAL
                     EDM=2.27764e-05    STRATEGY= 1      ERROR MATRIX ACCURATE 
  EXT PARAMETER                                INTERNAL      INTERNAL  
  NO.   NAME      VALUE            ERROR       STEP SIZE       VALUE   
   1  a            3.89817e-01   8.14814e-03   1.31669e-04  -2.22191e-01
   2  c           -9.91750e-01   6.56071e-02   5.16653e-05  -9.93383e-02
                               ERR DEF= 0.5
 EXTERNAL ERROR MATRIX.    NDIM=  25    NPAR=  2    ERR DEF=0.5
  6.640e-05 -2.186e-04 
 -2.186e-04  4.304e-03 
 PARAMETER  CORRELATION COEFFICIENTS  
       NO.  GLOBAL      1      2
        1  0.40891   1.000 -0.409
        2  0.40891  -0.409  1.000
[#1] INFO:Minization -- RooMinimizer::optimizeConst: deactivating const optimization

Plot fitted, data efficiency

Plot distribution of all events and accepted fraction of events on frame

In [10]:
RooPlot *frame1 = x.frame(Bins(20), Title("Data (all, accepted)"));
data->plotOn(frame1);
data->plotOn(frame1, Cut("cut==cut::accept"), MarkerColor(kRed), LineColor(kRed));
[#1] INFO:Plotting -- RooTreeData::plotOn: plotting 8176 events out of 10000 total events

Plot accept/reject efficiency on data overlay fitted efficiency curve

In [11]:
RooPlot *frame2 = x.frame(Bins(20), Title("Fitted efficiency"));
data->plotOn(frame2, Efficiency(cut)); // needs ROOT version >= 5.21
effFunc.plotOn(frame2, LineColor(kRed));

Draw all frames on a canvas

In [12]:
TCanvas *ca = new TCanvas("rf701_efficiency", "rf701_efficiency", 800, 400);
ca->Divide(2);
ca->cd(1);
gPad->SetLeftMargin(0.15);
frame1->GetYaxis()->SetTitleOffset(1.6);
frame1->Draw();
ca->cd(2);
gPad->SetLeftMargin(0.15);
frame2->GetYaxis()->SetTitleOffset(1.4);
frame2->Draw();

Draw all canvases

In [13]:
%jsroot on
gROOT->GetListOfCanvases()->Draw()