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

Detailed Description

View in nbviewer Open in SWAN
Multidimensional models: full pdf with per-event errors

#include "RooRealVar.h"
#include "RooDataSet.h"
#include "RooGaussian.h"
#include "RooGaussModel.h"
#include "RooDecay.h"
#include "RooLandau.h"
#include "RooProdPdf.h"
#include "RooHistPdf.h"
#include "RooPlot.h"
#include "TCanvas.h"
#include "TAxis.h"
#include "TH1.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 e m p i r i c a l p d f f o r p e r - e v e n t e r r o r
// -----------------------------------------------------------------
// Use landau pdf to get empirical distribution with long tail
RooLandau pdfDtErr("pdfDtErr", "pdfDtErr", dterr, 1.0, 0.25);
std::unique_ptr<RooDataSet> expDataDterr{pdfDtErr.generate(dterr, 10000)};
// Construct a histogram pdf to describe the shape of the dtErr distribution
std::unique_ptr<RooDataHist> expHistDterr{expDataDterr->binnedClone()};
RooHistPdf pdfErr("pdfErr", "pdfErr", dterr, *expHistDterr);
// C o n s t r u c t c o n d i t i o n a l p r o d u c t d e c a y _ d m ( d t | d t e r r ) * p d f ( d t e r
// r )
// ----------------------------------------------------------------------------------------------------------------------
// Construct production of conditional decay_dm(dt|dterr) with empirical pdfErr(dterr)
RooProdPdf model("model", "model", pdfErr, Conditional(decay_gm, dt));
// (Alternatively you could also use the landau shape pdfDtErr)
// RooProdPdf model("model","model",pdfDtErr,Conditional(decay_gm,dt)) ;
// S a m p l e, f i t a n d p l o t p r o d u c t m o d e l
// ------------------------------------------------------------------
// Specify external dataset with dterr values to use model_dm as conditional pdf
std::unique_ptr<RooDataSet> data{model.generate({dt, dterr}, 10000)};
// 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
model.fitTo(*data, PrintLevel(-1));
// 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 pdf in (dt,dterr)
TH1 *hh_model = model.createHistogram("hh_model", dt, Binning(50), YVar(dterr, Binning(50)));
hh_model->SetLineColor(kBlue);
// Make projection of data an dt
RooPlot *frame = dt.frame(Title("Projection of model(dt|dterr) on dt"));
data->plotOn(frame);
model.plotOn(frame);
// Draw all frames on canvas
TCanvas *c = new TCanvas("rf307_fullpereventerrors", "rf307_fullperventerrors", 800, 400);
c->Divide(2);
c->cd(1);
gPad->SetLeftMargin(0.20);
hh_model->GetZaxis()->SetTitleOffset(2.5);
hh_model->Draw("surf");
c->cd(2);
gPad->SetLeftMargin(0.15);
frame->GetYaxis()->SetTitleOffset(1.6);
frame->Draw();
}
#define c(i)
Definition RSha256.hxx:101
@ kBlue
Definition Rtypes.h:66
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
#define gPad
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.
A propability density function sampled from a multidimensional histogram.
Definition RooHistPdf.h:30
Landau distribution p.d.f.
Definition RooLandau.h:24
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:185
TAxis * GetYaxis() const
Definition RooPlot.cxx:1224
void Draw(Option_t *options=nullptr) override
Draw this plot and all of the elements it contains.
Definition RooPlot.cxx:597
Efficient implementation of a product of PDFs of the form.
Definition RooProdPdf.h:39
Variable that can be changed from the outside.
Definition RooRealVar.h:37
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title.
Definition TAttAxis.cxx:298
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition TAttLine.h:40
The Canvas class.
Definition TCanvas.h:23
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:59
TAxis * GetZaxis()
Definition TH1.h:327
void Draw(Option_t *option="") override
Draw this histogram with options.
Definition TH1.cxx:3068
RooCmdArg YVar(const RooAbsRealLValue &var, const RooCmdArg &arg={})
RooCmdArg Conditional(const RooArgSet &pdfSet, const RooArgSet &depSet, bool depsAreCond=false)
RooCmdArg PrintLevel(Int_t code)
RooCmdArg Binning(const RooAbsBinning &binning)
const Double_t sigma
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition CodegenImpl.h:64
const char * Title
Definition TXMLSetup.cxx:68
Date
July 2008
Author
Wouter Verkerke

Definition in file rf307_fullpereventerrors.C.