Logo ROOT   6.12/07
Reference Guide
rf603_multicpu.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// \notebook -js
4 /// 'LIKELIHOOD AND MINIMIZATION' RooFit tutorial macro #603
5 ///
6 /// Setting up a multi-core parallelized unbinned maximum likelihood fit
7 ///
8 /// \macro_image
9 /// \macro_output
10 /// \macro_code
11 /// \author 07/2008 - Wouter Verkerke
12 
13 
14 #include "RooRealVar.h"
15 #include "RooDataSet.h"
16 #include "RooGaussian.h"
17 #include "RooConstVar.h"
18 #include "RooPolynomial.h"
19 #include "RooAddPdf.h"
20 #include "RooProdPdf.h"
21 #include "TCanvas.h"
22 #include "TAxis.h"
23 #include "RooPlot.h"
24 using namespace RooFit ;
25 
26 
27 void rf603_multicpu()
28 {
29 
30  // C r e a t e 3 D p d f a n d d a t a
31  // -------------------------------------------
32 
33  // Create observables
34  RooRealVar x("x","x",-5,5) ;
35  RooRealVar y("y","y",-5,5) ;
36  RooRealVar z("z","z",-5,5) ;
37 
38  // Create signal pdf gauss(x)*gauss(y)*gauss(z)
39  RooGaussian gx("gx","gx",x,RooConst(0),RooConst(1)) ;
40  RooGaussian gy("gy","gy",y,RooConst(0),RooConst(1)) ;
41  RooGaussian gz("gz","gz",z,RooConst(0),RooConst(1)) ;
42  RooProdPdf sig("sig","sig",RooArgSet(gx,gy,gz)) ;
43 
44  // Create background pdf poly(x)*poly(y)*poly(z)
45  RooPolynomial px("px","px",x,RooArgSet(RooConst(-0.1),RooConst(0.004))) ;
46  RooPolynomial py("py","py",y,RooArgSet(RooConst(0.1),RooConst(-0.004))) ;
47  RooPolynomial pz("pz","pz",z) ;
48  RooProdPdf bkg("bkg","bkg",RooArgSet(px,py,pz)) ;
49 
50  // Create composite pdf sig+bkg
51  RooRealVar fsig("fsig","signal fraction",0.1,0.,1.) ;
52  RooAddPdf model("model","model",RooArgList(sig,bkg),fsig) ;
53 
54  // Generate large dataset
55  RooDataSet* data = model.generate(RooArgSet(x,y,z),200000) ;
56 
57 
58 
59  // P a r a l l e l f i t t i n g
60  // -------------------------------
61 
62  // In parallel mode the likelihood calculation is split in N pieces,
63  // that are calculated in parallel and added a posteriori before passing
64  // it back to MINUIT.
65 
66  // Use four processes and time results both in wall time and CPU time
67  model.fitTo(*data,NumCPU(4),Timer(kTRUE)) ;
68 
69 
70 
71  // P a r a l l e l M C p r o j e c t i o n s
72  // ----------------------------------------------
73 
74  // Construct signal, total likelihood projection on (y,z) observables and likelihood ratio
75  RooAbsPdf* sigyz = sig.createProjection(x) ;
76  RooAbsPdf* totyz = model.createProjection(x) ;
77  RooFormulaVar llratio_func("llratio","log10(@0)-log10(@1)",RooArgList(*sigyz,*totyz)) ;
78 
79  // Calculate likelihood ratio for each event, define subset of events with high signal likelihood
80  data->addColumn(llratio_func) ;
81  RooDataSet* dataSel = (RooDataSet*) data->reduce(Cut("llratio>0.7")) ;
82 
83  // Make plot frame and plot data
84  RooPlot* frame = x.frame(Title("Projection on X with LLratio(y,z)>0.7"),Bins(40)) ;
85  dataSel->plotOn(frame) ;
86 
87  // Perform parallel projection using MC integration of pdf using given input dataSet.
88  // In this mode the data-weighted average of the pdf is calculated by splitting the
89  // input dataset in N equal pieces and calculating in parallel the weighted average
90  // one each subset. The N results of those calculations are then weighted into the
91  // final result
92 
93  // Use four processes
94  model.plotOn(frame,ProjWData(*dataSel),NumCPU(4)) ;
95 
96 
97  new TCanvas("rf603_multicpu","rf603_multicpu",600,600) ;
98  gPad->SetLeftMargin(0.15) ; frame->GetYaxis()->SetTitleOffset(1.6) ; frame->Draw() ;
99 
100 }
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
RooAddPdf is an efficient implementation of a sum of PDFs of the form.
Definition: RooAddPdf.h:29
TAxis * GetYaxis() const
Definition: RooPlot.cxx:1118
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
Plot dataset on specified frame.
Definition: RooAbsData.cxx:568
RooCmdArg Cut(const char *cutSpec)
RooProdPdf is an efficient implementation of a product of PDFs of the form.
Definition: RooProdPdf.h:31
RooAbsData * reduce(const RooCmdArg &arg1, const RooCmdArg &arg2=RooCmdArg(), const RooCmdArg &arg3=RooCmdArg(), const RooCmdArg &arg4=RooCmdArg(), const RooCmdArg &arg5=RooCmdArg(), const RooCmdArg &arg6=RooCmdArg(), const RooCmdArg &arg7=RooCmdArg(), const RooCmdArg &arg8=RooCmdArg())
Create a reduced copy of this dataset.
Definition: RooAbsData.cxx:359
RooCmdArg ProjWData(const RooAbsData &projData, Bool_t binData=kFALSE)
RooCmdArg Title(const char *name)
Double_t x[n]
Definition: legend1.C:17
virtual RooAbsArg * addColumn(RooAbsArg &var, Bool_t adjustRange=kTRUE)
Add a column with the values of the given (function) argument to this dataset.
Plain Gaussian p.d.f.
Definition: RooGaussian.h:25
RooCmdArg NumCPU(Int_t nCPU, Int_t interleave=0)
virtual RooAbsPdf * createProjection(const RooArgSet &iset)
Return a p.d.f that represent a projection of this p.d.f integrated over given observables.
Definition: RooAbsPdf.cxx:3012
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:36
RooCmdArg Timer(Bool_t flag=kTRUE)
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:29
A RooPlot is a plot frame and a container for graphics objects within that frame. ...
Definition: RooPlot.h:41
The Canvas class.
Definition: TCanvas.h:31
Double_t y[n]
Definition: legend1.C:17
you should not use this method at all Int_t Int_t z
Definition: TRolke.cxx:630
RooConstVar & RooConst(Double_t val)
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
RooCmdArg Bins(Int_t nbin)
#define gPad
Definition: TVirtualPad.h:285
RooPolynomial implements a polynomial p.d.f of the form By default coefficient a_0 is chosen to be 1...
Definition: RooPolynomial.h:28
const Bool_t kTRUE
Definition: RtypesCore.h:87
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition: RooPlot.cxx:559