Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf603_multicpu.C File Reference

Detailed Description

View in nbviewer Open in SWAN
Likelihood and minimization: setting up a multi-core parallelized unbinned maximum likelihood fit

#include "RooRealVar.h"
#include "RooDataSet.h"
#include "RooGaussian.h"
#include "RooConstVar.h"
#include "RooPolynomial.h"
#include "RooAddPdf.h"
#include "RooProdPdf.h"
#include "TCanvas.h"
#include "TAxis.h"
#include "RooPlot.h"
using namespace RooFit;
{
// C r e a t e 3 D p d f a n d d a t a
// -------------------------------------------
// Create observables
RooRealVar x("x", "x", -5, 5);
RooRealVar y("y", "y", -5, 5);
RooRealVar z("z", "z", -5, 5);
// Create signal pdf gauss(x)*gauss(y)*gauss(z)
RooGaussian gx("gx", "gx", x, RooConst(0), RooConst(1));
RooGaussian gy("gy", "gy", y, RooConst(0), RooConst(1));
RooGaussian gz("gz", "gz", z, RooConst(0), RooConst(1));
RooProdPdf sig("sig", "sig", RooArgSet(gx, gy, gz));
// Create background pdf poly(x)*poly(y)*poly(z)
RooPolynomial px("px", "px", x, RooArgSet(-0.1, 0.004));
RooPolynomial py("py", "py", y, RooArgSet(0.1, -0.004));
RooPolynomial pz("pz", "pz", z);
RooProdPdf bkg("bkg", "bkg", RooArgSet(px, py, pz));
// Create composite pdf sig+bkg
RooRealVar fsig("fsig", "signal fraction", 0.1, 0., 1.);
RooAddPdf model("model", "model", RooArgList(sig, bkg), fsig);
// Generate large dataset
RooDataSet *data = model.generate(RooArgSet(x, y, z), 200000);
// P a r a l l e l f i t t i n g
// -------------------------------
// In parallel mode the likelihood calculation is split in N pieces,
// that are calculated in parallel and added a posteriori before passing
// it back to MINUIT.
// Use four processes and time results both in wall time and CPU time
model.fitTo(*data, NumCPU(4), Timer(true), PrintLevel(-1));
// P a r a l l e l M C p r o j e c t i o n s
// ----------------------------------------------
// Construct signal, total likelihood projection on (y,z) observables and likelihood ratio
RooAbsPdf *sigyz = sig.createProjection(x);
RooAbsPdf *totyz = model.createProjection(x);
RooFormulaVar llratio_func("llratio", "log10(@0)-log10(@1)", RooArgList(*sigyz, *totyz));
// Calculate likelihood ratio for each event, define subset of events with high signal likelihood
data->addColumn(llratio_func);
RooDataSet *dataSel = (RooDataSet *)data->reduce(Cut("llratio>0.7"));
// Make plot frame and plot data
RooPlot *frame = x.frame(Title("Projection on X with LLratio(y,z)>0.7"), Bins(40));
dataSel->plotOn(frame);
// Perform parallel projection using MC integration of pdf using given input dataSet.
// In this mode the data-weighted average of the pdf is calculated by splitting the
// input dataset in N equal pieces and calculating in parallel the weighted average
// one each subset. The N results of those calculations are then weighted into the
// final result
// Use four processes
model.plotOn(frame, ProjWData(*dataSel), NumCPU(4));
new TCanvas("rf603_multicpu", "rf603_multicpu", 600, 600);
gPad->SetLeftMargin(0.15);
frame->GetYaxis()->SetTitleOffset(1.6);
frame->Draw();
}
PrintLevel
Definition RooMinuit.h:6
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=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
virtual RooAbsPdf * createProjection(const RooArgSet &iset)
Return a p.d.f that represent a projection of this p.d.f integrated over given observables.
RooAddPdf is an efficient implementation of a sum of PDFs of the form.
Definition RooAddPdf.h:34
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:1276
void Draw(Option_t *options=nullptr) override
Draw this plot and all of the elements it contains.
Definition RooPlot.cxx:649
RooPolynomial implements a polynomial p.d.f of the form.
RooProdPdf is an 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:40
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
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 Common.h:18
[#1] INFO:Eval -- RooAbsTestStatistic::initMPMode: started 4 remote server process.
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:Minimization -- The following expressions have been identified as constant and will be precalculated and cached: (sig,bkg)
[#1] INFO:Minimization -- The following expressions have been identified as constant and will be precalculated and cached: (sig,bkg)
[#1] INFO:Minimization -- The following expressions have been identified as constant and will be precalculated and cached: (sig,bkg)
[#1] INFO:Minimization -- The following expressions have been identified as constant and will be precalculated and cached: (sig,bkg)
[#1] INFO:Minimization -- Command timer: Real time 0:00:00, CP time 0.080
[#1] INFO:Minimization -- Session timer: Real time 0:00:00, CP time 0.080
[#1] INFO:Minimization -- Command timer: Real time 0:00:00, CP time 0.010
[#1] INFO:Minimization -- Session timer: Real time 0:00:00, CP time 0.090, 2 slices
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
[#1] INFO:InputArguments -- The formula llratio>0.7 claims to use the variables (x,y,z,llratio) but only (llratio) seem to be in use.
inputs: llratio>0.7
[#1] INFO:Plotting -- RooAbsReal::plotOn(model) plot on x averages using data variables (y,z)
[#1] INFO:Plotting -- RooAbsReal::plotOn(model) only the following components of the projection data will be used: (y,z)
[#1] INFO:Plotting -- RooDataWeightedAverage::ctor(modelDataWgtAvg) constructing data weighted average of function model_Norm[x] over 14497 data points of (y,z) with a total weight of 14497
[#1] INFO:Eval -- RooAbsTestStatistic::initMPMode: started 4 remote server process.
[#1] INFO:Minimization -- The following expressions have been identified as constant and will be precalculated and cached: (gy,gz,py,pz)
[#1] INFO:Minimization -- The following expressions have been identified as constant and will be precalculated and cached: (gy,gz,py,pz)
[#1] INFO:Minimization -- The following expressions have been identified as constant and will be precalculated and cached: (gy,gz,py,pz)
[#1] INFO:Minimization -- The following expressions have been identified as constant and will be precalculated and cached: (gy,gz,py,pz)
Date
July 2008
Author
Wouter Verkerke

Definition in file rf603_multicpu.C.