Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
rf603_multicpu.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// \notebook -js
4 /// Likelihood and minimization: setting up a multi-core parallelized unbinned maximum likelihood fit
5 ///
6 /// \macro_image
7 /// \macro_output
8 /// \macro_code
9 /// \author 07/2008 - Wouter Verkerke
10 
11 #include "RooRealVar.h"
12 #include "RooDataSet.h"
13 #include "RooGaussian.h"
14 #include "RooConstVar.h"
15 #include "RooPolynomial.h"
16 #include "RooAddPdf.h"
17 #include "RooProdPdf.h"
18 #include "TCanvas.h"
19 #include "TAxis.h"
20 #include "RooPlot.h"
21 using namespace RooFit;
22 
23 void rf603_multicpu()
24 {
25 
26  // C r e a t e 3 D p d f a n d d a t a
27  // -------------------------------------------
28 
29  // Create observables
30  RooRealVar x("x", "x", -5, 5);
31  RooRealVar y("y", "y", -5, 5);
32  RooRealVar z("z", "z", -5, 5);
33 
34  // Create signal pdf gauss(x)*gauss(y)*gauss(z)
35  RooGaussian gx("gx", "gx", x, RooConst(0), RooConst(1));
36  RooGaussian gy("gy", "gy", y, RooConst(0), RooConst(1));
37  RooGaussian gz("gz", "gz", z, RooConst(0), RooConst(1));
38  RooProdPdf sig("sig", "sig", RooArgSet(gx, gy, gz));
39 
40  // Create background pdf poly(x)*poly(y)*poly(z)
41  RooPolynomial px("px", "px", x, RooArgSet(RooConst(-0.1), RooConst(0.004)));
42  RooPolynomial py("py", "py", y, RooArgSet(RooConst(0.1), RooConst(-0.004)));
43  RooPolynomial pz("pz", "pz", z);
44  RooProdPdf bkg("bkg", "bkg", RooArgSet(px, py, pz));
45 
46  // Create composite pdf sig+bkg
47  RooRealVar fsig("fsig", "signal fraction", 0.1, 0., 1.);
48  RooAddPdf model("model", "model", RooArgList(sig, bkg), fsig);
49 
50  // Generate large dataset
51  RooDataSet *data = model.generate(RooArgSet(x, y, z), 200000);
52 
53  // P a r a l l e l f i t t i n g
54  // -------------------------------
55 
56  // In parallel mode the likelihood calculation is split in N pieces,
57  // that are calculated in parallel and added a posteriori before passing
58  // it back to MINUIT.
59 
60  // Use four processes and time results both in wall time and CPU time
61  model.fitTo(*data, NumCPU(4), Timer(kTRUE));
62 
63  // P a r a l l e l M C p r o j e c t i o n s
64  // ----------------------------------------------
65 
66  // Construct signal, total likelihood projection on (y,z) observables and likelihood ratio
67  RooAbsPdf *sigyz = sig.createProjection(x);
68  RooAbsPdf *totyz = model.createProjection(x);
69  RooFormulaVar llratio_func("llratio", "log10(@0)-log10(@1)", RooArgList(*sigyz, *totyz));
70 
71  // Calculate likelihood ratio for each event, define subset of events with high signal likelihood
72  data->addColumn(llratio_func);
73  RooDataSet *dataSel = (RooDataSet *)data->reduce(Cut("llratio>0.7"));
74 
75  // Make plot frame and plot data
76  RooPlot *frame = x.frame(Title("Projection on X with LLratio(y,z)>0.7"), Bins(40));
77  dataSel->plotOn(frame);
78 
79  // Perform parallel projection using MC integration of pdf using given input dataSet.
80  // In this mode the data-weighted average of the pdf is calculated by splitting the
81  // input dataset in N equal pieces and calculating in parallel the weighted average
82  // one each subset. The N results of those calculations are then weighted into the
83  // final result
84 
85  // Use four processes
86  model.plotOn(frame, ProjWData(*dataSel), NumCPU(4));
87 
88  new TCanvas("rf603_multicpu", "rf603_multicpu", 600, 600);
89  gPad->SetLeftMargin(0.15);
90  frame->GetYaxis()->SetTitleOffset(1.6);
91  frame->Draw();
92 }