Rf 3 1 4_Paramfitrange

Multidimensional models: working with parametrized ranges in a fit. This an example of a fit with an acceptance that changes per-event

pdf = exp(-t/tau) with t[tmin,5]

where t and tmin are both observables in the dataset

Author: Wouter Verkerke
This notebook tutorial was automatically generated with ROOTBOOK-izer from the macro found in the ROOT repository on Monday, January 17, 2022 at 10:00 AM.

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

Define observables and decay pdf

Declare observables

In [3]:
RooRealVar t("t", "t", 0, 5);
RooRealVar tmin("tmin", "tmin", 0, 0, 5);
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

Make parametrized range in t : [tmin,5]

In [4]:
t.setRange(tmin, RooConst(t.getMax()));

Make pdf

In [5]:
RooRealVar tau("tau", "tau", -1.54, -10, -0.1);
RooExponential model("model", "model", t, tau);

Create input data

Generate complete dataset without acceptance cuts (for reference)

In [6]:
RooDataSet *dall = model.generate(t, 10000);

Generate a (fake) prototype dataset for acceptance limit values

In [7]:
RooDataSet *tmp = RooGaussian("gmin", "gmin", tmin, RooConst(0), RooConst(0.5)).generate(tmin, 5000);

Generate dataset with t values that observe (t>tmin)

In [8]:
RooDataSet *dacc = model.generate(t, ProtoData(*tmp));

Fit pdf to data in acceptance region

In [9]:
RooFitResult *r = model.fitTo(*dacc, Save());
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
 **********
 **    1 **SET PRINT           1
 **********
 **********
 **    2 **SET NOGRAD
 **********
 PARAMETER DEFINITIONS:
    NO.   NAME         VALUE      STEP SIZE      LIMITS
     1 tau         -1.54000e+00  7.20000e-01   -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         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=2824.02 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  tau         -1.54000e+00   7.20000e-01   2.12947e-01  -4.56069e+01
                               ERR DEF= 0.5
 MIGRAD MINIMIZATION HAS CONVERGED.
 MIGRAD WILL VERIFY CONVERGENCE AND ERROR MATRIX.
 COVARIANCE MATRIX CALCULATED SUCCESSFULLY
 FCN=2823.97 FROM MIGRAD    STATUS=CONVERGED      12 CALLS          13 TOTAL
                     EDM=6.75021e-08    STRATEGY= 1      ERROR MATRIX ACCURATE 
  EXT PARAMETER                                   STEP         FIRST   
  NO.   NAME      VALUE            ERROR          SIZE      DERIVATIVE 
   1  tau         -1.53353e+00   2.21980e-02   2.34054e-04   4.07752e-02
                               ERR DEF= 0.5
 EXTERNAL ERROR MATRIX.    NDIM=  25    NPAR=  1    ERR DEF=0.5
  4.928e-04 
 **********
 **    7 **SET ERR         0.5
 **********
 **********
 **    8 **SET PRINT           1
 **********
 **********
 **    9 **HESSE         500
 **********
 COVARIANCE MATRIX CALCULATED SUCCESSFULLY
 FCN=2823.97 FROM HESSE     STATUS=OK              5 CALLS          18 TOTAL
                     EDM=6.74739e-08    STRATEGY= 1      ERROR MATRIX ACCURATE 
  EXT PARAMETER                                INTERNAL      INTERNAL  
  NO.   NAME      VALUE            ERROR       STEP SIZE       VALUE   
   1  tau         -1.53353e+00   2.21980e-02   4.68108e-05   7.90063e-01
                               ERR DEF= 0.5
 EXTERNAL ERROR MATRIX.    NDIM=  25    NPAR=  1    ERR DEF=0.5
  4.928e-04 
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization

Plot fitted pdf on full and accepted data

Make plot frame, add datasets and overlay model

In [10]:
RooPlot *frame = t.frame(Title("Fit to data with per-event acceptance"));
dall->plotOn(frame, MarkerColor(kRed), LineColor(kRed));
model.plotOn(frame);
dacc->plotOn(frame);
[#0] WARNING:Plotting -- Cannot apply a bin width correction and use Poisson errors. Not correcting for bin width.
[#0] WARNING:Plotting -- Cannot apply a bin width correction and use Poisson errors. Not correcting for bin width.
[#1] INFO:Plotting -- RooPlot::updateFitRangeNorm: New event count of 5000 will supercede previous event count of 10000 for normalization of PDF projections

Print fit results to demonstrate absence of bias

In [11]:
r->Print("v");

new TCanvas("rf314_paramranges", "rf314_paramranges", 600, 600);
gPad->SetLeftMargin(0.15);
frame->GetYaxis()->SetTitleOffset(1.6);
frame->Draw();

return;
  RooFitResult: minimized FCN value: 2823.97, estimated distance to minimum: 6.74739e-08
                covariance matrix quality: Full, accurate covariance matrix
                Status : MINIMIZE=0 HESSE=0 

    Floating Parameter  InitialValue    FinalValue +/-  Error     GblCorr.
  --------------------  ------------  --------------------------  --------
                   tau   -1.5400e+00   -1.5335e+00 +/-  2.22e-02  <none>

Draw all canvases

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