Logo ROOT   6.18/05
Reference Guide
rf306_condpereventerrors.C File Reference

Detailed Description

View in nbviewer Open in SWAN Multidimensional models: conditional p.d.f.

with per-event errors

␛[1mRooFit v3.60 -- Developed by Wouter Verkerke and David Kirkby␛[0m
Copyright (C) 2000-2013 NIKHEF, University of California & Stanford University
All rights reserved, please read http://roofit.sourceforge.net/license.txt
[#1] INFO:Minization -- RooMinimizer::optimizeConst: activating const optimization
**********
** 1 **SET PRINT 1
**********
**********
** 2 **SET NOGRAD
**********
PARAMETER DEFINITIONS:
NO. NAME VALUE STEP SIZE LIMITS
1 bias 0.00000e+00 2.00000e+00 -1.00000e+01 1.00000e+01
2 sigma 1.00000e+00 4.50000e-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 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=23876.4 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 bias 0.00000e+00 2.00000e+00 2.01358e-01 -1.70342e+02
2 sigma 1.00000e+00 4.50000e-01 1.63378e-01 8.62474e+01
ERR DEF= 0.5
MIGRAD MINIMIZATION HAS CONVERGED.
MIGRAD WILL VERIFY CONVERGENCE AND ERROR MATRIX.
COVARIANCE MATRIX CALCULATED SUCCESSFULLY
FCN=23876.2 FROM MIGRAD STATUS=CONVERGED 30 CALLS 31 TOTAL
EDM=3.03467e-08 STRATEGY= 1 ERROR MATRIX ACCURATE
EXT PARAMETER STEP FIRST
NO. NAME VALUE ERROR SIZE DERIVATIVE
1 bias 5.26463e-03 1.72103e-02 1.83632e-04 1.00797e-01
2 sigma 9.87130e-01 2.04183e-02 7.70435e-04 2.52620e-03
ERR DEF= 0.5
EXTERNAL ERROR MATRIX. NDIM= 25 NPAR= 2 ERR DEF=0.5
2.962e-04 -4.402e-06
-4.402e-06 4.169e-04
PARAMETER CORRELATION COEFFICIENTS
NO. GLOBAL 1 2
1 0.01253 1.000 -0.013
2 0.01253 -0.013 1.000
**********
** 7 **SET ERR 0.5
**********
**********
** 8 **SET PRINT 1
**********
**********
** 9 **HESSE 1000
**********
COVARIANCE MATRIX CALCULATED SUCCESSFULLY
FCN=23876.2 FROM HESSE STATUS=OK 12 CALLS 43 TOTAL
EDM=3.12688e-08 STRATEGY= 1 ERROR MATRIX ACCURATE
EXT PARAMETER INTERNAL INTERNAL
NO. NAME VALUE ERROR STEP SIZE VALUE
1 bias 5.26463e-03 1.72105e-02 3.67263e-05 5.26463e-04
2 sigma 9.87130e-01 2.04260e-02 3.08174e-05 -9.62778e-01
ERR DEF= 0.5
EXTERNAL ERROR MATRIX. NDIM= 25 NPAR= 2 ERR DEF=0.5
2.962e-04 -4.786e-06
-4.786e-06 4.172e-04
PARAMETER CORRELATION COEFFICIENTS
NO. GLOBAL 1 2
1 0.01361 1.000 -0.014
2 0.01361 -0.014 1.000
[#1] INFO:Minization -- RooMinimizer::optimizeConst: deactivating const optimization
[#1] INFO:NumericIntegration -- RooRealIntegral::init(gm1_conv_exp(-abs(@0)/@1)_dt_tau_[decay_gm]_Int[dt,dterr]) using numeric integrator RooIntegrator1D to calculate Int(dterr)
[#1] INFO:Plotting -- RooAbsReal::plotOn(decay_gm) plot on dt averages using data variables (dterr)
[#1] INFO:Plotting -- RooDataWeightedAverage::ctor(decay_gmDataWgtAvg) constructing data weighted average of function decay_gm_Norm[dt] over 100 data points of (dterr) with a total weight of 10000
.........................................................................................................................................................................................................................
#include "RooRealVar.h"
#include "RooDataSet.h"
#include "RooGaussian.h"
#include "RooConstVar.h"
#include "RooGaussModel.h"
#include "RooDecay.h"
#include "RooLandau.h"
#include "RooPlot.h"
#include "TCanvas.h"
#include "TAxis.h"
#include "TH2D.h"
using namespace RooFit;
{
// B - p h y s i c s p d f w i t h p e r - e v e n t G a u s s i a n r e s o l u t i o n
// ----------------------------------------------------------------------------------------------
// Observables
RooRealVar dt("dt", "dt", -10, 10);
RooRealVar dterr("dterr", "per-event error on dt", 0.01, 10);
// Build a gaussian resolution model scaled by the per-event error = gauss(dt,bias,sigma*dterr)
RooRealVar bias("bias", "bias", 0, -10, 10);
RooRealVar sigma("sigma", "per-event error scale factor", 1, 0.1, 10);
RooGaussModel gm("gm1", "gauss model scaled bt per-event error", dt, bias, sigma, dterr);
// Construct decay(dt) (x) gauss1(dt|dterr)
RooRealVar tau("tau", "tau", 1.548);
RooDecay decay_gm("decay_gm", "decay", dt, tau, gm, RooDecay::DoubleSided);
// C o n s t r u c t f a k e ' e x t e r n a l ' d a t a w i t h p e r - e v e n t e r r o r
// ------------------------------------------------------------------------------------------------------
// Use landau p.d.f to get somewhat realistic distribution with long tail
RooLandau pdfDtErr("pdfDtErr", "pdfDtErr", dterr, RooConst(1), RooConst(0.25));
RooDataSet *expDataDterr = pdfDtErr.generate(dterr, 10000);
// S a m p l e d a t a f r o m c o n d i t i o n a l d e c a y _ g m ( d t | d t e r r )
// ---------------------------------------------------------------------------------------------
// Specify external dataset with dterr values to use decay_dm as conditional p.d.f.
RooDataSet *data = decay_gm.generate(dt, ProtoData(*expDataDterr));
// F i t c o n d i t i o n a l d e c a y _ d m ( d t | d t e r r )
// ---------------------------------------------------------------------
// Specify dterr as conditional observable
decay_gm.fitTo(*data, ConditionalObservables(dterr));
// P l o t c o n d i t i o n a l d e c a y _ d m ( d t | d t e r r )
// ---------------------------------------------------------------------
// Make two-dimensional plot of conditional p.d.f in (dt,dterr)
TH1 *hh_decay = decay_gm.createHistogram("hh_decay", dt, Binning(50), YVar(dterr, Binning(50)));
hh_decay->SetLineColor(kBlue);
// Plot decay_gm(dt|dterr) at various values of dterr
RooPlot *frame = dt.frame(Title("Slices of decay(dt|dterr) at various dterr"));
for (Int_t ibin = 0; ibin < 100; ibin += 20) {
dterr.setBin(ibin);
decay_gm.plotOn(frame, Normalization(5.));
}
// Make projection of data an dt
RooPlot *frame2 = dt.frame(Title("Projection of decay(dt|dterr) on dt"));
data->plotOn(frame2);
// Make projection of decay(dt|dterr) on dt.
//
// Instead of integrating out dterr, make a weighted average of curves
// at values dterr_i as given in the external dataset.
// (The kTRUE argument bins the data before projection to speed up the process)
decay_gm.plotOn(frame2, ProjWData(*expDataDterr, kTRUE));
// Draw all frames on canvas
TCanvas *c = new TCanvas("rf306_condpereventerrors", "rf306_condperventerrors", 1200, 400);
c->Divide(3);
c->cd(1);
gPad->SetLeftMargin(0.20);
hh_decay->GetZaxis()->SetTitleOffset(2.5);
hh_decay->Draw("surf");
c->cd(2);
gPad->SetLeftMargin(0.15);
frame->GetYaxis()->SetTitleOffset(1.6);
frame->Draw();
c->cd(3);
gPad->SetLeftMargin(0.15);
frame2->GetYaxis()->SetTitleOffset(1.6);
frame2->Draw();
}
#define c(i)
Definition: RSha256.hxx:101
int Int_t
Definition: RtypesCore.h:41
const Bool_t kTRUE
Definition: RtypesCore.h:87
@ kBlue
Definition: Rtypes.h:64
#define gPad
Definition: TVirtualPad.h:286
virtual RooPlot * plotOn(RooPlot *frame, const RooCmdArg &arg1=RooCmdArg::none(), const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none(), const RooCmdArg &arg6=RooCmdArg::none(), const RooCmdArg &arg7=RooCmdArg::none(), const RooCmdArg &arg8=RooCmdArg::none()) const
Calls RooPlot* plotOn(RooPlot* frame, const RooLinkedList& cmdList) const ;.
Definition: RooAbsData.cxx:552
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:31
Single or double sided decay function that can be analytically convolved with any RooResolutionModel ...
Definition: RooDecay.h:22
@ DoubleSided
Definition: RooDecay.h:25
Class RooGaussModel implements a RooResolutionModel that models a Gaussian distribution.
Definition: RooGaussModel.h:26
Landau distribution p.d.f.
Definition: RooLandau.h:24
A RooPlot is a plot frame and a container for graphics objects within that frame.
Definition: RooPlot.h:41
TAxis * GetYaxis() const
Definition: RooPlot.cxx:1123
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition: RooPlot.cxx:558
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:36
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title Offset is a correction factor with respect to the "s...
Definition: TAttAxis.cxx:294
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition: TAttLine.h:40
The Canvas class.
Definition: TCanvas.h:31
The TH1 histogram class.
Definition: TH1.h:56
TAxis * GetZaxis()
Definition: TH1.h:318
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2981
const Double_t sigma
Template specialisation used in RooAbsArg:
RooCmdArg ProjWData(const RooAbsData &projData, Bool_t binData=kFALSE)
RooCmdArg ProtoData(const RooDataSet &protoData, Bool_t randomizeOrder=kFALSE, Bool_t resample=kFALSE)
RooCmdArg Binning(const RooAbsBinning &binning)
RooCmdArg YVar(const RooAbsRealLValue &var, const RooCmdArg &arg=RooCmdArg::none())
RooConstVar & RooConst(Double_t val)
RooCmdArg ConditionalObservables(const RooArgSet &set)
RooCmdArg Normalization(Double_t scaleFactor)
const char * Title
Definition: TXMLSetup.cxx:67
Author
07/2008 - Wouter Verkerke

Definition in file rf306_condpereventerrors.C.