'MULTIDIMENSIONAL MODELS' RooFit tutorial macro #309
Projecting p.d.f and data slices in discrete observables
Processing /mnt/build/workspace/root-makedoc-v614/rootspi/rdoc/src/v6-14-00-patches/tutorials/roofit/rf310_sliceplot.C...
[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:Plotting -- RooAbsReal::plotOn(bmix) plot on dt integrates over variables (tagFlav,mixState)
[#1] INFO:Plotting -- RooTreeData::plotOn: plotting 3787 events out of 20000 total events
[#1] INFO:Plotting -- RooAbsReal::plotOn(bmix) plot on dt represents a slice in (mixState)
[#1] INFO:Plotting -- RooAbsReal::plotOn(bmix) plot on dt integrates over variables (tagFlav)
[#1] INFO:Plotting -- RooTreeData::plotOn: plotting 16213 events out of 20000 total events
[#1] INFO:Plotting -- RooAbsReal::plotOn(bmix) plot on dt represents a slice in (mixState)
[#1] INFO:Plotting -- RooAbsReal::plotOn(bmix) plot on dt integrates over variables (tagFlav)
void rf310_sliceplot()
{
RooCategory mixState(
"mixState",
"B0/B0bar mixing state") ;
RooCategory tagFlav(
"tagFlav",
"Flavour of the tagged B0") ;
mixState.defineType("mixed",-1) ;
mixState.defineType("unmixed",1) ;
tagFlav.defineType("B0",1) ;
tagFlav.defineType("B0bar",-1) ;
RooRealVar w(
"w",
"Flavor Mistag rate",0.03,0.0,1.0) ;
RooRealVar dw(
"dw",
"Flavor Mistag rate difference between B0 and B0bar",0.01) ;
RooBMixDecay bmix_gm1(
"bmix",
"decay",dt,mixState,tagFlav,tau,dm,w,dw,gm1,
RooBMixDecay::DoubleSided) ;
RooPlot* frame = dt.frame(
Title(
"Inclusive decay distribution")) ;
bmix_gm1.plotOn(frame) ;
RooPlot* frame2 = dt.frame(
Title(
"Decay distribution of mixed events")) ;
data->
plotOn(frame2,
Cut(
"mixState==mixState::mixed")) ;
bmix_gm1.plotOn(frame2,
Slice(mixState,
"mixed")) ;
RooPlot* frame3 = dt.frame(
Title(
"Decay distribution of unmixed events")) ;
data->
plotOn(frame3,
Cut(
"mixState==mixState::unmixed")) ;
bmix_gm1.plotOn(frame3,
Slice(mixState,
"unmixed")) ;
TCanvas* c =
new TCanvas(
"rf310_sliceplot",
"rf310_sliceplot",1200,400) ;
}
- Author
- 07/2008 - Wouter Verkerke
Definition in file rf310_sliceplot.C.