In [1]:
from ROOTaaS.iPyROOT import ROOT
ROOT.toCpp()


Welcome to ROOTaas Beta
Notebook is in Cpp mode

In [2]:
using namespace RooFit ;

  RooRealVar dt("dt","dt",-20,20) ;

  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 dm("dm","delta m(B)",0.472,0.,1.0) ;
  RooRealVar tau("tau","B0 decay time",1.547,1.0,2.0) ;
  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) ;

  RooRealVar bias1("bias1","bias1",0) ;
  RooRealVar sigma1("sigma1","sigma1",0.01) ;  
  RooGaussModel gm1("gm1","gauss model 1",dt,bias1,sigma1) ;

  RooBMixDecay bmix_gm1("bmix","decay",dt,mixState,tagFlav,tau,dm,w,dw,gm1,RooBMixDecay::DoubleSided) ;
  
  RooDataSet *data = bmix_gm1.generate(RooArgSet(dt,tagFlav,mixState),20000) ;

  RooPlot* frame = dt.frame(Title("Inclusive decay distribution")) ;
  data->plotOn(frame) ;
  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) ;
  c->Divide(3) ;
  c->cd(1) ; gPad->SetLeftMargin(0.15) ; frame->GetYaxis()->SetTitleOffset(1.4) ; gPad->SetLogy() ; frame->Draw() ;
  c->cd(2) ; gPad->SetLeftMargin(0.15) ; frame2->GetYaxis()->SetTitleOffset(1.4) ; gPad->SetLogy() ; frame2->Draw() ;
  c->cd(3) ; gPad->SetLeftMargin(0.15) ; frame3->GetYaxis()->SetTitleOffset(1.4) ; gPad->SetLogy() ; frame3->Draw() ;


Out[2]:
0
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

[#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)

In [ ]: