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

Detailed Description

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

#include "RooRealVar.h"
#include "RooDataSet.h"
#include "RooGaussian.h"
#include "RooGaussModel.h"
#include "RooDecay.h"
#include "RooLandau.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 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 pdf to get somewhat realistic distribution with long tail
RooLandau pdfDtErr("pdfDtErr", "pdfDtErr", dterr, 1.0, 0.25);
std::unique_ptr<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 pdf
std::unique_ptr<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), 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_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 true argument bins the data before projection to speed up the process)
decay_gm.plotOn(frame2, ProjWData(*expDataDterr, true));
// 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:45
@ 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.
Landau distribution p.d.f.
Definition RooLandau.h:24
Plot frame and a container for graphics objects within that frame.
Definition RooPlot.h:45
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:225
TAxis * GetYaxis() const
Definition RooPlot.cxx:1264
void Draw(Option_t *options=nullptr) override
Draw this plot and all of the elements it contains.
Definition RooPlot.cxx:637
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:326
void Draw(Option_t *option="") override
Draw this histogram with options.
Definition TH1.cxx:3066
RooCmdArg YVar(const RooAbsRealLValue &var, const RooCmdArg &arg={})
RooCmdArg PrintLevel(Int_t code)
RooCmdArg ConditionalObservables(Args_t &&... argsOrArgSet)
Create a RooCmdArg to declare conditional observables.
RooCmdArg Binning(const RooAbsBinning &binning)
RooCmdArg ProjWData(const RooAbsData &projData, bool binData=false)
RooCmdArg Normalization(double scaleFactor)
const Double_t sigma
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition JSONIO.h:26
const char * Title
Definition TXMLSetup.cxx:68
[#1] INFO:Fitting -- RooAbsPdf::fitTo(gm1_conv_exp(-abs(@0)/@1)_dt_tau_[decay_gm]_over_gm1_conv_exp(-abs(@0)/@1)_dt_tau_[decay_gm]_Int[dt]) fixing normalization set for coefficient determination to observables in data
[#1] INFO:Fitting -- using CPU computation library compiled with -mavx2
[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_gm1_conv_exp(-abs(@0)/@1)_dt_tau_[decay_gm]_over_gm1_conv_exp(-abs(@0)/@1)_dt_tau_[decay_gm]_Int[dt]_decay_gmData) Summation contains a RooNLLVar, using its error level
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: 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)
Date
July 2008
Author
Wouter Verkerke

Definition in file rf306_condpereventerrors.C.