Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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_code
8/// \macro_output
9///
10/// \date July 2008
11/// \author Wouter Verkerke
12
13#include "RooRealVar.h"
14#include "RooDataSet.h"
15#include "RooGaussian.h"
16#include "RooPolynomial.h"
17#include "RooAddPdf.h"
18#include "RooProdPdf.h"
19#include "TCanvas.h"
20#include "TAxis.h"
21#include "RooPlot.h"
22using namespace RooFit;
23
24void rf603_multicpu()
25{
26
27 // C r e a t e 3 D p d f a n d d a t a
28 // -------------------------------------------
29
30 // Create observables
31 RooRealVar x("x", "x", -5, 5);
32 RooRealVar y("y", "y", -5, 5);
33 RooRealVar z("z", "z", -5, 5);
34
35 // Create signal pdf gauss(x)*gauss(y)*gauss(z)
36 RooGaussian gx("gx", "gx", x, 0.0, 1.0);
37 RooGaussian gy("gy", "gy", y, 0.0, 1.0);
38 RooGaussian gz("gz", "gz", z, 0.0, 1.0);
39 RooProdPdf sig("sig", "sig", RooArgSet(gx, gy, gz));
40
41 // Create background pdf poly(x)*poly(y)*poly(z)
42 RooPolynomial px("px", "px", x, RooArgSet(-0.1, 0.004));
43 RooPolynomial py("py", "py", y, RooArgSet(0.1, -0.004));
44 RooPolynomial pz("pz", "pz", z);
45 RooProdPdf bkg("bkg", "bkg", RooArgSet(px, py, pz));
46
47 // Create composite pdf sig+bkg
48 RooRealVar fsig("fsig", "signal fraction", 0.1, 0., 1.);
49 RooAddPdf model("model", "model", RooArgList(sig, bkg), fsig);
50
51 // Generate large dataset
52 std::unique_ptr<RooDataSet> data{model.generate({x, y, z}, 200000)};
53
54 // P a r a l l e l f i t t i n g
55 // -------------------------------
56
57 // In parallel mode the likelihood calculation is split in N pieces,
58 // that are calculated in parallel and added a posteriori before passing
59 // it back to MINUIT.
60
61 // Use four processes and time results both in wall time and CPU time
62 model.fitTo(*data, NumCPU(4), Timer(true), PrintLevel(-1));
63
64 // P a r a l l e l M C p r o j e c t i o n s
65 // ----------------------------------------------
66
67 // Construct signal, total likelihood projection on (y,z) observables and likelihood ratio
68 RooAbsPdf *sigyz = sig.createProjection(x);
69 RooAbsPdf *totyz = model.createProjection(x);
70 RooFormulaVar llratio_func("llratio", "log10(@0)-log10(@1)", RooArgList(*sigyz, *totyz));
71
72 // Calculate likelihood ratio for each event, define subset of events with high signal likelihood
73 data->addColumn(llratio_func);
74 RooDataSet *dataSel = (RooDataSet *)data->reduce(Cut("llratio>0.7"));
75
76 // Make plot frame and plot data
77 RooPlot *frame = x.frame(Title("Projection on X with LLratio(y,z)>0.7"), Bins(40));
78 dataSel->plotOn(frame);
79
80 // Perform parallel projection using MC integration of pdf using given input dataSet.
81 // In this mode the data-weighted average of the pdf is calculated by splitting the
82 // input dataset in N equal pieces and calculating in parallel the weighted average
83 // one each subset. The N results of those calculations are then weighted into the
84 // final result
85
86 // Use four processes
87 model.plotOn(frame, ProjWData(*dataSel), NumCPU(4));
88
89 new TCanvas("rf603_multicpu", "rf603_multicpu", 600, 600);
90 gPad->SetLeftMargin(0.15);
91 frame->GetYaxis()->SetTitleOffset(1.6);
92 frame->Draw();
93}
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
#define gPad
virtual RooPlot * plotOn(RooPlot *frame, const RooCmdArg &arg1={}, const RooCmdArg &arg2={}, const RooCmdArg &arg3={}, const RooCmdArg &arg4={}, const RooCmdArg &arg5={}, const RooCmdArg &arg6={}, const RooCmdArg &arg7={}, const RooCmdArg &arg8={}) const
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
virtual RooAbsPdf * createProjection(const RooArgSet &iset)
Return a p.d.f that represent a projection of this p.d.f integrated over given observables.
Efficient implementation of a sum of PDFs of the form.
Definition RooAddPdf.h:33
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:55
RooDataSet is a container class to hold unbinned data.
Definition RooDataSet.h:57
A RooFormulaVar is a generic implementation of a real-valued object, which takes a RooArgList of serv...
Plain Gaussian p.d.f.
Definition RooGaussian.h:24
A RooPlot is a plot frame and a container for graphics objects within that frame.
Definition RooPlot.h:43
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:239
TAxis * GetYaxis() const
Definition RooPlot.cxx:1279
void Draw(Option_t *options=nullptr) override
Draw this plot and all of the elements it contains.
Definition RooPlot.cxx:652
RooPolynomial implements a polynomial p.d.f of the form.
Efficient implementation of a product of PDFs of the form.
Definition RooProdPdf.h:33
RooRealVar represents a 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
The Canvas class.
Definition TCanvas.h:23
RooCmdArg Bins(Int_t nbin)
RooCmdArg PrintLevel(Int_t code)
RooCmdArg NumCPU(Int_t nCPU, Int_t interleave=0)
RooCmdArg ProjWData(const RooAbsData &projData, bool binData=false)
RooCmdArg Cut(const char *cutSpec)
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
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