Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf619_discrete_profiling.C
Go to the documentation of this file.
1///\file
2///\ingroup tutorial_roofit_main
3/// \notebook -js
4/// Basic functionality: demonstrate fitting multiple models using RooMultiPdf and selecting the best one via Discrete
5/// Profiling method.
6///
7/// \macro_image
8/// \macro_code
9/// \macro_output
10///
11/// \date July 2025
12/// \author Galin Bistrev
13
14using namespace RooFit;
15
16void rf619_discrete_profiling()
17{
18 RooRealVar x("x", "Observable", 0, 50);
19
20 // Category 0 Pdf-s: Exponential + Chebyshev.
21 RooRealVar lambda1("lambda1", "slope1", -0.025, -0.1, -0.02);
22 RooExponential expo1("expo1", "Exponential 1", x, lambda1);
23
24 RooRealVar c0("c0", "Cheby coeff 0", -1.0, -1.0, 1.0);
25 RooRealVar c1("c1", "Cheby coeff 1", 0.4, 0.05, 0.5);
27 RooChebychev cheb("cheb", "Chebyshev PDF", x, chebCoeffs);
28
29 RooCategory pdfIndex0("pdfIndex0", "pdf index 0");
30 RooMultiPdf multiPdf0("multiPdf0", "multiPdf0", pdfIndex0, RooArgList(expo1, cheb));
31
32 // Adding complexity to the model via introdcing exponential pdf.
33 RooRealVar lambdaExtra("lambdaExtra", "extra slope", -0.05, -1.0, -0.01);
34 RooExponential expoExtra("expoExtra", "extra exponential", x, lambdaExtra);
35
36 // Category 1 Pdf-s: Gaussian + Landau.
37 RooRealVar mean("mean", "shared mean", 25, 0, 50);
38 RooRealVar sigmaG("sigmaG", "Gaussian width", 2.0, 0.0, 5.0);
39 RooRealVar sigmaL("sigmaL", "Landau width", 3.0, 1.0, 8.0);
40
41 RooGaussian gauss1("gauss1", "Gaussian", x, mean, sigmaG);
42 RooLandau landau1("landau1", "Landau", x, mean, sigmaL);
43
44 RooCategory pdfIndex1("pdfIndex1", "pdf index 1");
45 RooMultiPdf multiPdf1("multiPdf1", "multiPdf1", pdfIndex1, RooArgList(gauss1, landau1));
46
47 // Adding complexity to the model via introdcing gaussExtra pdf. Profile scan in this model
48 // will also be done over a shared parameter from both models (mean).
49 RooRealVar sigmaExtra("sigmaExtra", "extra Gaussian width", 3.0, 1.0, 6.0);
50 RooGaussian gaussExtra("gaussExtra", "extra Gaussian", x, mean, sigmaExtra);
51
52 // Creation of AddPdf objects.Add multiple PDFs into a single model for each category.
53 RooRealVar frac0("frac0", "fraction for cat0", 0.7, 0.0, 1.0);
54 RooAddPdf addPdf0("addPdf0", "multiPdf0 + extra expo", RooArgList(multiPdf0, gaussExtra), frac0);
55
56 RooRealVar frac1("frac1", "fraction for cat1", 0.5, 0.0, 1.0);
57 RooAddPdf addPdf1("addPdf1", "multiPdf1 + extra gauss", RooArgList(multiPdf1, expoExtra), frac1);
58
59 // Simultaneous Pdf across categories.
60 RooCategory catIndex("catIndex", "Category");
61 catIndex.defineType("cat0", 0);
62 catIndex.defineType("cat1", 1);
63
64 RooSimultaneous simPdf("simPdf", "simultaneous model", catIndex);
65 simPdf.addPdf(addPdf0, "cat0");
66 simPdf.addPdf(addPdf1, "cat1");
67
68 // Generate toy data for each AddPdf.
69 RooDataSet *data0 = addPdf0.generate(x, 800);
70 RooDataSet *data1 = addPdf1.generate(x, 1000);
71
72 // Ploting individual Pdf-s
73
74 RooPlot *frame0 = x.frame();
75 data0->plotOn(frame0);
76 addPdf0.plotOn(frame0);
77 pdfIndex0.setIndex(1);
78 addPdf0.plotOn(frame0, LineColor(kRed));
79 frame0->SetTitle("");
80 frame0->GetXaxis()->SetTitle("Observable");
81 frame0->GetYaxis()->SetTitle("Events");
82
83 TLegend *leg0 = new TLegend(0.6, 0.7, 0.9, 0.9);
84 leg0->AddEntry(frame0->getObject(0), "Data", "lep"); // Data points
85 leg0->AddEntry(frame0->getObject(1), "Expo ", "l");
86 leg0->AddEntry(frame0->getObject(2), "Poly", "l");
87
88 // Ploting individual Pdf-s - Category 1
89
90 RooPlot *frame1 = x.frame();
91 data1->plotOn(frame1);
92 addPdf1.plotOn(frame1);
93 pdfIndex1.setIndex(1);
94 addPdf1.plotOn(frame1, LineColor(kRed));
95 frame1->SetTitle("");
96 frame1->GetXaxis()->SetTitle("Observable");
97 frame1->GetYaxis()->SetTitle("Events");
98
99 // Create legend for Category 1
100 TLegend *leg1 = new TLegend(0.6, 0.7, 0.9, 0.9);
101 leg1->AddEntry(frame1->getObject(0), "Data", "lep"); // Data points
102 leg1->AddEntry(frame1->getObject(1), "Gauss", "l");
103 leg1->AddEntry(frame1->getObject(2), "Landau", "l");
104
105 // Combine datasets for simultaneous fit.
106 RooDataSet *data = new RooDataSet("data", "combined", RooArgSet(x, catIndex));
107
108 RooArgSet vars(x, catIndex);
109
110 for (int i = 0; i < data0->numEntries(); ++i) {
111 x.setVal(data0->get(i)->getRealValue("x"));
112 catIndex.setLabel("cat0");
113 data->add(vars);
114 }
115 for (int i = 0; i < data1->numEntries(); ++i) {
116 x.setVal(data1->get(i)->getRealValue("x"));
117 catIndex.setLabel("cat1");
118 data->add(vars);
119 }
120
121 // Create NLL with codegen and minimize it via the discrete profiling method.
122 std::unique_ptr<RooAbsReal> nll1(simPdf.createNLL(*data, EvalBackend("codegen")));
124
125 minim.setStrategy(1);
126 minim.setEps(1e-7);
127 minim.setPrintLevel(-1);
128
129 // Setup profiling of 'mean' over discrete combinations.
130 const int nMeanPoints = 40;
131 const double meanMin = 17;
132 const double meanMax = 33;
133
134 // Generate all discrete combinations of PDF indices for profiling.
135 std::vector<std::vector<int>> combosToPlot;
136 for (int i = 0; i < pdfIndex0.size(); ++i) {
137 for (int j = 0; j < pdfIndex1.size(); ++j) {
138 combosToPlot.push_back({i, j});
139 }
140 }
141
142 // Creates a canvas where all NLL vs mean will be drawn.
143 TCanvas *c = new TCanvas("c_rf619", "NLL vs Mean for Different Discrete Combinations", 1200, 400);
144 c->Divide(3, 1);
145
146 // Define arrays of ROOT colors and marker styles.
147 int colors[] = {kRed, kBlue, kGreen + 2, kMagenta, kOrange + 7};
148 int markers[] = {20, 21, 22, 23, 33};
149
150 // Container to store all TGraph objects:
151 std::vector<TGraph *> graphs;
152
153 // Create one TGraph for each discrete combination of indices.
154 // Each graph will hold the NLL vs mean curve for that combination.
155 // Colors and marker styles are cycled through the predefined arrays.
156 for (size_t idx = 0; idx < combosToPlot.size(); ++idx) {
157 TGraph *g = new TGraph(nMeanPoints);
158 g->SetLineColor(colors[idx % 5]);
159 g->SetMarkerColor(colors[idx % 5]);
160 g->SetMarkerStyle(markers[idx % 5]);
161 g->SetTitle(Form("Combo [%d,%d]", combosToPlot[idx][0], combosToPlot[idx][1]));
162 graphs.push_back(g);
163 }
164
165 // Create a special TGraph for the profiled NLL curve across all combinations.
166 // Drawn in bold black to stand out.
168 profileGraph->SetLineColor(kBlack);
169 profileGraph->SetLineWidth(4);
170 profileGraph->SetMarkerColor(kBlack);
171 profileGraph->SetMarkerStyle(22);
172 profileGraph->SetTitle("Profile");
173 graphs.push_back(profileGraph);
174
175 // Loop over mean values, compute NLL and profile likelihood.
176 for (int i = 0; i < nMeanPoints; ++i) {
177
178 // Fix the "mean" parameter at this scan point
179 double meanVal = meanMin + i * (meanMax - meanMin) / (nMeanPoints - 1);
180 mean.setVal(meanVal);
181
182 // Loop over all discrete PDF index combinations
183 for (size_t comboIdx = 0; comboIdx < combosToPlot.size(); ++comboIdx) {
184 const auto &combo = combosToPlot[comboIdx];
185
186 // Fix both category indices to this combination
187 pdfIndex0.setIndex(combo[0]);
188 pdfIndex1.setIndex(combo[1]);
189
190 // Freeze discrete indices and the mean parameter
191 // so the minimizer only optimizes continuous parameters
192 pdfIndex0.setConstant(true);
193 pdfIndex1.setConstant(true);
194 mean.setConstant(true);
195
196 // Minimize over the continuous parameters
197 // Evaluate NLL at this configuration
198 minim.minimize("Minuit2", "Migrad");
199 double nllVal = nll1->getVal();
200
201 // Store NLL vs. mean value in the corresponding graph
202 graphs[comboIdx]->SetPoint(i, meanVal, nllVal);
203
204 // Unfreeze categories and mean so the next loop iteration can change thems
205 pdfIndex0.setConstant(false);
206 pdfIndex1.setConstant(false);
207 mean.setConstant(false);
208 }
209
210 // Freeze mean (profiling at this value)
211 mean.setConstant(true);
212
213 // Scan all discrete combinations internally
214 minim.minimize("Minuit2", "Migrad");
215
216 // Store profiled NLL and plot it
217 double profNLL = nll1->getVal();
218 graphs.back()->SetPoint(i, meanVal, profNLL);
219 }
220 // Panel 1: Category 0
221 c->cd(1);
222 gPad->SetLeftMargin(0.15);
223 frame0->GetYaxis()->SetTitleOffset(1.4);
224 frame0->Draw();
225 if (leg0)
226 leg0->Draw();
227
228 // Panel 2: Category 1
229 c->cd(2);
230 gPad->SetLeftMargin(0.15);
231 frame1->GetYaxis()->SetTitleOffset(1.4);
232 frame1->Draw();
233 if (leg1)
234 leg1->Draw();
235
236 c->cd(3);
237 gPad->SetLeftMargin(0.15);
238
239 TMultiGraph *mg = new TMultiGraph();
240 for (auto &g : graphs) {
241 mg->Add(g, "PL");
242 mg->GetYaxis()->SetTitleOffset(1.8);
243 }
244
245 mg->Draw("APL");
246 mg->GetXaxis()->SetTitle("Mean");
247 mg->GetYaxis()->SetTitle("NLL");
248
249 gPad->BuildLegend();
250}
#define c(i)
Definition RSha256.hxx:101
#define g(i)
Definition RSha256.hxx:105
#define e(i)
Definition RSha256.hxx:103
@ kRed
Definition Rtypes.h:67
@ kOrange
Definition Rtypes.h:68
@ kBlack
Definition Rtypes.h:66
@ kGreen
Definition Rtypes.h:67
@ kMagenta
Definition Rtypes.h:67
@ kBlue
Definition Rtypes.h:67
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
char * Form(const char *fmt,...)
Formats a string in a circular formatting buffer.
Definition TString.cxx:2495
#define gPad
Color * colors
Definition X3DBuffer.c:21
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:24
Object to represent discrete states.
Definition RooCategory.h:28
Chebychev polynomial p.d.f.
Container class to hold unbinned data.
Definition RooDataSet.h:32
Exponential PDF.
Plain Gaussian p.d.f.
Definition RooGaussian.h:24
Landau distribution p.d.f.
Definition RooLandau.h:24
Wrapper class around ROOT::Math::Minimizer that provides a seamless interface between the minimizer f...
The class RooMultiPdf allows for the creation of a RooMultiPdf object, which can switch between previ...
Definition RooMultiPdf.h:9
Plot frame and a container for graphics objects within that frame.
Definition RooPlot.h:43
Variable that can be changed from the outside.
Definition RooRealVar.h:37
Facilitates simultaneous fitting of multiple PDFs to subsets of a given dataset.
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title.
Definition TAttAxis.cxx:279
The Canvas class.
Definition TCanvas.h:23
A TGraph is an object made of two arrays X and Y with npoints each.
Definition TGraph.h:41
This class displays a legend box (TPaveText) containing several legend entries.
Definition TLegend.h:23
A TMultiGraph is a collection of TGraph (or derived) objects.
Definition TMultiGraph.h:34
virtual void Add(TGraph *graph, Option_t *chopt="")
Add a new graph to the list of graphs.
void Draw(Option_t *chopt="") override
Draw this multigraph with its current attributes.
TAxis * GetYaxis()
Get y axis of the graph.
TAxis * GetXaxis()
Get x axis of the graph.
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition TNamed.cxx:173
RooCmdArg LineColor(TColorNumber color)
return c1
Definition legend1.C:41
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 CodegenImpl.h:67