CMS 3D CMS Logo

FitWithRooFit.cc
Go to the documentation of this file.
2 
4 
5 // Import TH1 histogram into a RooDataHist
6 std::unique_ptr<RooDataHist> FitWithRooFit::importTH1(TH1* histo, double xMin, double xMax) {
7  if ((xMin == xMax) && xMin == 0) {
8  xMin = histo->GetXaxis()->GetXmin();
9  xMax = histo->GetXaxis()->GetXmax();
10  }
11  // Declare observable x
12  RooRealVar x("InvMass", "di-muon mass M(#mu^{+}#mu^{-}) [GeV]", xMin, xMax);
13  // Create a binned dataset that imports contents of TH1 and associates its contents to observable 'x'
14  return std::make_unique<RooDataHist>("dh", "dh", x, RooFit::Import(*histo));
15 }
16 
17 // Plot and fit a RooDataHist fitting signal and background
19  TH1* histo, const TString signalType, const TString backgroundType, double xMin, double xMax, bool sumW2Error) {
21 
22  std::unique_ptr<RooDataHist> dh = importTH1(histo, xMin, xMax);
23  RooRealVar x(*static_cast<RooRealVar*>(dh->get()->find("x")));
24 
25  // Make plot of binned dataset showing Poisson error bars (RooFit default)
26  RooPlot* frame = x.frame(RooFit::Title("di-muon mass M(#mu^{+}#mu^{-}) [GeV]"));
27  frame->SetName(TString(histo->GetName()) + "_frame");
28  dh->plotOn(frame);
29 
30  // Build the composite model
31  std::unique_ptr<RooAbsPdf> model = buildModel(&x, signalType, backgroundType);
32 
33  std::unique_ptr<RooAbsReal> chi2{model->createChi2(*dh, RooFit::DataError(RooAbsData::SumW2))};
34 
35  // Fit the composite model
36  // -----------------------
37  // Fit with likelihood
38  if (!useChi2_) {
39  model->fitTo(*dh, RooFit::SumW2Error(sumW2Error));
40  }
41  // Fit with chi^2
42  else {
43  edm::LogWarning("FitWithRooFit") << "FITTING WITH CHI^2";
44  RooMinimizer m(*chi2);
45  m.migrad();
46  m.hesse();
47  // RooFitResult* r_chi2_wgt = m.save();
48  }
49  model->plotOn(frame, RooFit::LineColor(kRed));
50  model->plotOn(frame, RooFit::Components(backgroundType), RooFit::LineStyle(kDashed));
51  model->paramOn(frame,
52  RooFit::Label("fit result"),
53  RooFit::Layout(0.65, 0.90, 0.90),
54  RooFit::Format("NEU", RooFit::AutoPrecision(2)));
55 
56  // TODO: fix next lines to get the prob(chi2) (ndof should be dynamically set according to the choosen pdf)
57  // double chi2 = xframe.chiSquare("model","data",ndof);
58  // double ndoff = xframeGetNbinsX();
59  // double chi2prob = TMath::Prob(chi2,ndoff);
60 
61  // P l o t a n d f i t a R o o D a t a H i s t w i t h i n t e r n a l e r r o r s
62  // ---------------------------------------------------------------------------------------------
63 
64  // If histogram has custom error (i.e. its contents is does not originate from a Poisson process
65  // but e.g. is a sum of weighted events) you can data with symmetric 'sum-of-weights' error instead
66  // (same error bars as shown by ROOT)
67  RooPlot* frame2 = x.frame(RooFit::Title("di-muon mass M(#mu^{+}#mu^{-}) [GeV]"));
68  dh->plotOn(frame2, RooFit::DataError(RooAbsData::SumW2));
69  model->plotOn(frame2, RooFit::LineColor(kRed));
70  model->plotOn(frame2, RooFit::Components(backgroundType), RooFit::LineStyle(kDashed));
71  model->paramOn(frame2,
72  RooFit::Label("fit result"),
73  RooFit::Layout(0.65, 0.90, 0.90),
74  RooFit::Format("NEU", RooFit::AutoPrecision(2)));
75 
76  // Please note that error bars shown (Poisson or SumW2) are for visualization only, the are NOT used
77  // in a maximum likelihood fit
78  //
79  // A (binned) ML fit will ALWAYS assume the Poisson error interpretation of data (the mathematical definition
80  // of likelihood does not take any external definition of errors). Data with non-unit weights can only be correctly
81  // fitted with a chi^2 fit (see rf602_chi2fit.C)
82 
83  // Draw all frames on a canvas
84  // if( canvas == 0 ) {
85  // canvas = new TCanvas("rf102_dataimport","rf102_dataimport",800,800);
86  // canvas->cd(1);
87  // }
88  // canvas->Divide(2,1);
89  // canvas->cd(1) ; frame->Draw();
90  // canvas->cd(2) ; frame2->Draw();
91 
92  frame2->Draw();
93 }
94 
95 // Initialization methods for all the parameters
96 void FitWithRooFit::initMean(double value, double min, double max, const TString& name, const TString& title) {
97  mean_ = std::make_unique<RooRealVar>(name, title, value, min, max);
99 }
100 
101 void FitWithRooFit::initMean2(double value, double min, double max, const TString& name, const TString& title) {
102  mean2_ = std::make_unique<RooRealVar>(name, title, value, min, max);
104 }
105 
106 void FitWithRooFit::initMean3(double value, double min, double max, const TString& name, const TString& title) {
107  mean3_ = std::make_unique<RooRealVar>(name, title, value, min, max);
109 }
110 
111 void FitWithRooFit::initSigma(double value, double min, double max, const TString& name, const TString& title) {
112  sigma_ = std::make_unique<RooRealVar>(name, title, value, min, max);
114 }
115 
116 void FitWithRooFit::initSigma2(double value, double min, double max, const TString& name, const TString& title) {
117  sigma2_ = std::make_unique<RooRealVar>(name, title, value, min, max);
119 }
120 
121 void FitWithRooFit::initSigma3(double value, double min, double max, const TString& name, const TString& title) {
122  sigma3_ = std::make_unique<RooRealVar>(name, title, value, min, max);
124 }
125 
126 void FitWithRooFit::initGamma(double value, double min, double max, const TString& name, const TString& title) {
127  gamma_ = std::make_unique<RooRealVar>(name, title, value, min, max);
129 }
130 
131 void FitWithRooFit::initGaussFrac(double value, double min, double max, const TString& name, const TString& title) {
132  gaussFrac_ = std::make_unique<RooRealVar>(name, title, value, min, max);
134 }
135 
136 void FitWithRooFit::initGaussFrac2(double value, double min, double max, const TString& name, const TString& title) {
137  gaussFrac2_ = std::make_unique<RooRealVar>(name, title, value, min, max);
139 }
140 
141 void FitWithRooFit::initExpCoeffA0(double value, double min, double max, const TString& name, const TString& title) {
142  expCoeffa0_ = std::make_unique<RooRealVar>(name, title, value, min, max);
144 }
145 
146 void FitWithRooFit::initExpCoeffA1(double value, double min, double max, const TString& name, const TString& title) {
147  expCoeffa1_ = std::make_unique<RooRealVar>(name, title, value, min, max);
149 }
150 
151 void FitWithRooFit::initExpCoeffA2(double value, double min, double max, const TString& name, const TString& title) {
152  expCoeffa2_ = std::make_unique<RooRealVar>(name, title, value, min, max);
154 }
155 
156 void FitWithRooFit::initFsig(double value, double min, double max, const TString& name, const TString& title) {
157  fsig_ = std::make_unique<RooRealVar>(name, title, value, min, max);
159 }
160 
161 void FitWithRooFit::initA0(double value, double min, double max, const TString& name, const TString& title) {
162  a0_ = std::make_unique<RooRealVar>(name, title, value, min, max);
163  initVal_a0 = value;
164 }
165 
166 void FitWithRooFit::initA1(double value, double min, double max, const TString& name, const TString& title) {
167  a1_ = std::make_unique<RooRealVar>(name, title, value, min, max);
168  initVal_a1 = value;
169 }
170 
171 void FitWithRooFit::initA2(double value, double min, double max, const TString& name, const TString& title) {
172  a2_ = std::make_unique<RooRealVar>(name, title, value, min, max);
173  initVal_a2 = value;
174 }
175 
176 void FitWithRooFit::initA3(double value, double min, double max, const TString& name, const TString& title) {
177  a3_ = std::make_unique<RooRealVar>(name, title, value, min, max);
178  initVal_a3 = value;
179 }
180 
181 void FitWithRooFit::initA4(double value, double min, double max, const TString& name, const TString& title) {
182  a4_ = std::make_unique<RooRealVar>(name, title, value, min, max);
183  initVal_a4 = value;
184 }
185 
186 void FitWithRooFit::initA5(double value, double min, double max, const TString& name, const TString& title) {
187  a5_ = std::make_unique<RooRealVar>(name, title, value, min, max);
188  initVal_a5 = value;
189 }
190 
191 void FitWithRooFit::initA6(double value, double min, double max, const TString& name, const TString& title) {
192  a6_ = std::make_unique<RooRealVar>(name, title, value, min, max);
193  initVal_a6 = value;
194 }
195 
196 void FitWithRooFit::initAlpha(double value, double min, double max, const TString& name, const TString& title) {
197  alpha_ = std::make_unique<RooRealVar>(name, title, value, min, max);
199 }
200 
201 void FitWithRooFit::initN(double value, double min, double max, const TString& name, const TString& title) {
202  n_ = std::make_unique<RooRealVar>(name, title, value, min, max);
203  initVal_n = value;
204 }
205 
206 void FitWithRooFit::initFGCB(double value, double min, double max, const TString& name, const TString& title) {
207  fGCB_ = std::make_unique<RooRealVar>(name, title, value, min, max);
209 }
210 
212  auto initParam = [](std::unique_ptr<RooRealVar>& var, double val) {
213  if (var)
214  var->setVal(val);
215  };
216 
217  initParam(mean_, initVal_mean);
218  initParam(mean2_, initVal_mean2);
219  initParam(mean3_, initVal_mean3);
220  initParam(sigma_, initVal_sigma);
221  initParam(sigma2_, initVal_sigma2);
222  initParam(sigma3_, initVal_sigma3);
223  initParam(gamma_, initVal_gamma);
224  initParam(gaussFrac_, initVal_gaussFrac);
225  initParam(gaussFrac2_, initVal_gaussFrac2);
226  initParam(expCoeffa0_, initVal_expCoeffa0);
227  initParam(expCoeffa1_, initVal_expCoeffa1);
228  initParam(expCoeffa2_, initVal_expCoeffa2);
229  initParam(fsig_, initVal_fsig);
230  initParam(a0_, initVal_a0);
231  initParam(a1_, initVal_a1);
232  initParam(a2_, initVal_a2);
233  initParam(a3_, initVal_a3);
234  initParam(a4_, initVal_a4);
235  initParam(a5_, initVal_a5);
236  initParam(a6_, initVal_a6);
237  initParam(alpha_, initVal_alpha);
238  initParam(n_, initVal_n);
239  initParam(fGCB_, initVal_fGCB);
240 }
241 
243 std::unique_ptr<RooAbsPdf> FitWithRooFit::buildSignalModel(RooRealVar* x, const TString& signalType) {
244  if (signalType == "gaussian") {
245  // Fit a Gaussian p.d.f to the data
246  if ((mean_ == nullptr) || (sigma_ == nullptr)) {
247  edm::LogError("FitWithRooFit")
248  << "Error: one or more parameters are not initialized. Please be sure to initialize mean and sigma";
249  exit(1);
250  }
251  return std::make_unique<RooGaussian>("gauss", "gauss", *x, *mean_, *sigma_);
252  } else if (signalType == "doubleGaussian") {
253  // Fit with double gaussian
254  if ((mean_ == nullptr) || (sigma_ == nullptr) || (sigma2_ == nullptr)) {
255  edm::LogError("FitWithRooFit")
256  << "Error: one or more parameters are not initialized. Please be sure to initialize mean, sigma and sigma2";
257  exit(1);
258  }
259  RooGaussModel* gaussModel = new RooGaussModel("gaussModel", "gaussModel", *x, *mean_, *sigma_);
260  RooGaussModel* gaussModel2 = new RooGaussModel("gaussModel2", "gaussModel2", *x, *mean_, *sigma2_);
261  RooArgList components{*gaussModel, *gaussModel2};
262  auto out = std::make_unique<RooAddModel>("doubleGaussian", "double gaussian", components, *gaussFrac_);
263  out->addOwnedComponents(components);
264  return out;
265  } else if (signalType == "tripleGaussian") {
266  // Fit with triple gaussian
267  if ((mean_ == nullptr) || (mean2_ == nullptr) || (mean3_ == nullptr) || (sigma_ == nullptr) ||
268  (sigma2_ == nullptr) || (sigma3_ == nullptr)) {
269  edm::LogError("FitWithRooFit")
270  << "Error: one or more parameters are not initialized. Please be sure to initialize mean, mean2, "
271  "mean3, sigma, sigma2, sigma3";
272  exit(1);
273  }
274  RooGaussModel* gaussModel = new RooGaussModel("gaussModel", "gaussModel", *x, *mean_, *sigma_);
275  RooGaussModel* gaussModel2 = new RooGaussModel("gaussModel2", "gaussModel2", *x, *mean2_, *sigma2_);
276  RooGaussModel* gaussModel3 = new RooGaussModel("gaussModel3", "gaussModel3", *x, *mean3_, *sigma3_);
277  RooArgList components{*gaussModel, *gaussModel2, *gaussModel3};
278  auto out = std::make_unique<RooAddModel>(
279  "tripleGaussian", "triple gaussian", components, RooArgList{*gaussFrac_, *gaussFrac2_});
280  out->addOwnedComponents(components);
281  return out;
282  } else if (signalType == "breitWigner") {
283  // Fit a Breit-Wigner
284  if ((mean_ == nullptr) || (gamma_ == nullptr)) {
285  edm::LogError("FitWithRooFit")
286  << "Error: one or more parameters are not initialized. Please be sure to initialize mean and gamma";
287  exit(1);
288  }
289  return std::make_unique<RooBreitWigner>("breiWign", "breitWign", *x, *mean_, *gamma_);
290  } else if (signalType == "relBreitWigner") {
291  // Fit a relativistic Breit-Wigner
292  if ((mean_ == nullptr) || (gamma_ == nullptr)) {
293  edm::LogError("FitWithRooFit")
294  << "Error: one or more parameters are not initialized. Please be sure to initialize mean and gamma";
295  exit(1);
296  }
297  return std::make_unique<RooGenericPdf>("Relativistic Breit-Wigner",
298  "RBW",
299  "@0/(pow(@0*@0 - @1*@1,2) + @2*@2*@0*@0*@0*@0/(@1*@1))",
300  RooArgList{*x, *mean_, *gamma_});
301  } else if (signalType == "voigtian") {
302  // Fit a Voigtian
303  if ((mean_ == nullptr) || (sigma_ == nullptr) || (gamma_ == nullptr)) {
304  edm::LogError("FitWithRooFit")
305  << "Error: one or more parameters are not initialized. Please be sure to initialize mean, sigma and gamma";
306  exit(1);
307  }
308  return std::make_unique<RooVoigtian>("voigt", "voigt", *x, *mean_, *gamma_, *sigma_);
309  } else if (signalType == "crystalBall") {
310  // Fit a CrystalBall
311  if ((mean_ == nullptr) || (sigma_ == nullptr) || (alpha_ == nullptr) || (n_ == nullptr)) {
312  edm::LogError("FitWithRooFit")
313  << "Error: one or more parameters are not initialized. Please be sure to initialize mean, sigma, "
314  "alpha and n";
315  exit(1);
316  }
317  return std::make_unique<RooCBShape>("crystalBall", "crystalBall", *x, *mean_, *sigma_, *alpha_, *n_);
318  } else if (signalType == "breitWignerTimesCB") {
319  // Fit a Breit Wigner convoluted with a CrystalBall
320  if ((mean_ == nullptr) || (mean2_ == nullptr) || (sigma_ == nullptr) || (gamma_ == nullptr) ||
321  (alpha_ == nullptr) || (n_ == nullptr)) {
322  edm::LogError("FitWithRooFit")
323  << "Error: one or more parameters are not initialized. Please be sure to initialize mean, mean2, "
324  "sigma, gamma, alpha and n";
325  exit(1);
326  }
327  RooAbsPdf* bw = new RooBreitWigner("breiWigner", "breitWigner", *x, *mean_, *gamma_);
328  RooAbsPdf* cb = new RooCBShape("crystalBall", "crystalBall", *x, *mean2_, *sigma_, *alpha_, *n_);
329  auto out = std::make_unique<RooFFTConvPdf>("breitWignerTimesCB", "breitWignerTimesCB", *x, *bw, *cb);
330  out->addOwnedComponents({*bw, *cb});
331  return out;
332  } else if (signalType == "relBreitWignerTimesCB") {
333  // Fit a relativistic Breit Wigner convoluted with a CrystalBall
334  if ((mean_ == nullptr) || (mean2_ == nullptr) || (sigma_ == nullptr) || (gamma_ == nullptr) ||
335  (alpha_ == nullptr) || (n_ == nullptr)) {
336  edm::LogError("FitWithRooFit")
337  << "Error: one or more parameters are not initialized. Please be sure to initialize mean, mean2, "
338  "sigma, gamma, alpha and n";
339  exit(1);
340  }
341  RooGenericPdf* bw = new RooGenericPdf("Relativistic Breit-Wigner",
342  "RBW",
343  "@0/(pow(@0*@0 - @1*@1,2) + @2*@2*@0*@0*@0*@0/(@1*@1))",
344  {*x, *mean_, *gamma_});
345  RooAbsPdf* cb = new RooCBShape("crystalBall", "crystalBall", *x, *mean2_, *sigma_, *alpha_, *n_);
346  auto out = std::make_unique<RooFFTConvPdf>("relBreitWignerTimesCB", "relBreitWignerTimesCB", *x, *bw, *cb);
347  out->addOwnedComponents({*bw, *cb});
348  return out;
349  } else if (signalType == "gaussianPlusCrystalBall") {
350  // Fit a Gaussian + CrystalBall with the same mean
351  if ((mean_ == nullptr) || (sigma_ == nullptr) || (alpha_ == nullptr) || (n_ == nullptr) || (sigma2_ == nullptr) ||
352  (fGCB_ == nullptr)) {
353  edm::LogError("FitWithRooFit")
354  << "Error: one or more parameters are not initialized. Please be sure to initialize mean, sigma, "
355  "sigma2, alpha, n and fGCB";
356  exit(1);
357  }
358  RooAbsPdf* tempCB = new RooCBShape("crystalBall", "crystalBall", *x, *mean_, *sigma_, *alpha_, *n_);
359  RooAbsPdf* tempGaussian = new RooGaussian("gauss", "gauss", *x, *mean_, *sigma2_);
360  RooArgList components{*tempCB, *tempGaussian};
361 
362  auto out = std::make_unique<RooAddPdf>("gaussianPlusCrystalBall", "gaussianPlusCrystalBall", components, *fGCB_);
363  out->addOwnedComponents(components);
364  return out;
365  } else if (signalType == "voigtianPlusCrystalBall") {
366  // Fit a Voigtian + CrystalBall with the same mean
367  if ((mean_ == nullptr) || (sigma_ == nullptr) || (gamma_ == nullptr) || (alpha_ == nullptr) || (n_ == nullptr) ||
368  (sigma2_ == nullptr) || (fGCB_ == nullptr)) {
369  edm::LogError("FitWithRooFit")
370  << "Error: one or more parameters are not initialized. Please be sure to initialize mean, gamma, "
371  "sigma, sigma2, alpha, n and fGCB";
372  exit(1);
373  }
374  RooAbsPdf* tempVoigt = new RooVoigtian("voigt", "voigt", *x, *mean_, *gamma_, *sigma_);
375  RooAbsPdf* tempCB = new RooCBShape("crystalBall", "crystalBall", *x, *mean_, *sigma2_, *alpha_, *n_);
376  RooArgList components{*tempVoigt, *tempCB};
377 
378  auto out = std::make_unique<RooAddPdf>("voigtianPlusCrystalBall", "voigtianPlusCrystalBall", components, *fGCB_);
379  out->addOwnedComponents(components);
380  return out;
381  } else if (signalType == "breitWignerPlusCrystalBall") {
382  // Fit a Breit-Wigner + CrystalBall with the same mean
383  if ((mean_ == nullptr) || (gamma_ == nullptr) || (alpha_ == nullptr) || (n_ == nullptr) || (sigma2_ == nullptr) ||
384  (fGCB_ == nullptr)) {
385  edm::LogError("FitWithRooFit")
386  << "Error: one or more parameters are not initialized. Please be sure to initialize mean, gamma, "
387  "sigma, alpha, n and fGCB";
388  exit(1);
389  }
390  RooAbsPdf* tempBW = new RooBreitWigner("breitWign", "breitWign", *x, *mean_, *gamma_);
391  RooAbsPdf* tempCB = new RooCBShape("crystalBall", "crystalBall", *x, *mean_, *sigma2_, *alpha_, *n_);
392  RooArgList components{*tempCB, *tempBW};
393 
394  auto out =
395  std::make_unique<RooAddPdf>("breitWignerPlusCrystalBall", "breitWignerPlusCrystalBall", components, *fGCB_);
396  out->addOwnedComponents(components);
397  return out;
398  }
399 
400  else if (signalType != "") {
401  edm::LogError("FitWithRooFit") << "Unknown signal function: " << signalType << ". Signal will not be in the model";
402  }
403  return nullptr;
404 }
405 
407 std::unique_ptr<RooAbsPdf> FitWithRooFit::buildBackgroundModel(RooRealVar* x, const TString& backgroundType) {
408  if (backgroundType == "exponential") {
409  // Add an exponential for the background
410  if ((expCoeffa1_ == nullptr) || (fsig_ == nullptr)) {
411  edm::LogError("FitWithRooFit")
412  << "Error: one or more parameters are not initialized. Please be sure to initialize expCoeffa1 and fsig";
413  exit(1);
414  }
415  return std::make_unique<RooExponential>("exponential", "exponential", *x, *expCoeffa1_);
416  }
417 
418  if (backgroundType == "exponentialpol") {
419  // Add an exponential for the background
420  if ((expCoeffa0_ == nullptr) || (expCoeffa1_ == nullptr) || (expCoeffa2_ == nullptr) || (fsig_ == nullptr)) {
421  edm::LogError("FitWithRooFit")
422  << "Error: one or more parameters are not initialized. Please be sure to initialize expCoeff and fsig";
423  exit(1);
424  }
425  return std::make_unique<RooGenericPdf>("exponential",
426  "exponential",
427  "TMath::Exp(@1+@2*@0+@3*@0*@0)",
428  RooArgList{*x, *expCoeffa0_, *expCoeffa1_, *expCoeffa2_});
429  }
430 
431  else if (backgroundType == "chebychev0") {
432  // Add a linear background
433  if (a0_ == nullptr) {
434  edm::LogError("FitWithRooFit")
435  << "Error: one or more parameters are not initialized. Please be sure to initialize a0";
436  exit(1);
437  }
438  return std::make_unique<RooChebychev>("chebychev0", "chebychev0", *x, *a0_);
439  } else if (backgroundType == "chebychev1") {
440  // Add a 2nd order chebychev polynomial background
441  if ((a0_ == nullptr) || (a1_ == nullptr)) {
442  edm::LogError("FitWithRooFit")
443  << "Error: one or more parameters are not initialized. Please be sure to initialize a0 and a1";
444  exit(1);
445  }
446  return std::make_unique<RooChebychev>("chebychev1", "chebychev1", *x, RooArgList{*a0_, *a1_});
447  } else if (backgroundType == "chebychev3") {
448  // Add a 3rd order chebychev polynomial background
449  if ((a0_ == nullptr) || (a1_ == nullptr) || (a2_ == nullptr) || (a3_ == nullptr)) {
450  edm::LogError("FitWithRooFit")
451  << "Error: one or more parameters are not initialized. Please be sure to initialize a0, a1, a2 and a3";
452  exit(1);
453  }
454  return std::make_unique<RooChebychev>("3rdOrderPol", "3rdOrderPol", *x, RooArgList{*a0_, *a1_, *a2_, *a3_});
455  }
456 
457  else if (backgroundType == "chebychev6") {
458  // Add a 6th order chebychev polynomial background
459  if ((a0_ == nullptr) || (a1_ == nullptr) || (a2_ == nullptr) || (a3_ == nullptr) || (a4_ == nullptr) ||
460  (a5_ == nullptr) || (a6_ == nullptr)) {
461  edm::LogError("FitWithRooFit")
462  << "Error: one or more parameters are not initialized. Please be sure to initialize a0, a1, a2, a3, "
463  "a4, a5 and a6";
464  exit(1);
465  }
466  return std::make_unique<RooChebychev>(
467  "6thOrderPol", "6thOrderPol", *x, RooArgList{*a0_, *a1_, *a2_, *a3_, *a4_, *a5_, *a6_});
468  }
469  return nullptr;
470 }
471 
473 std::unique_ptr<RooAbsPdf> FitWithRooFit::buildModel(RooRealVar* x,
474  const TString& signalType,
475  const TString& backgroundType) {
476  std::unique_ptr<RooAbsPdf> model;
477 
478  std::unique_ptr<RooAbsPdf> signal = buildSignalModel(x, signalType);
479  std::unique_ptr<RooAbsPdf> background = buildBackgroundModel(x, backgroundType);
480 
481  if ((signal != nullptr) && (background != nullptr)) {
482  // Combine signal and background pdfs
483  edm::LogPrint("FitWithRooFit") << "Building model with signal and backgound";
484  RooArgList components{*signal.release(), *background.release()};
485  model = std::make_unique<RooAddPdf>("model", "model", components, *fsig_);
486  model->addOwnedComponents(components);
487  } else if (signal != nullptr) {
488  edm::LogPrint("FitWithRooFit") << "Building model with signal";
489  model = std::move(signal);
490  } else if (background != nullptr) {
491  edm::LogPrint("FitWithRooFit") << "Building model with backgound";
492  model = std::move(background);
493  } else {
494  edm::LogWarning("FitWithRooFit") << "Nothing to fit";
495  exit(0);
496  }
497  return model;
498 }
void initA0(double value, double min, double max, const TString &name="a0", const TString &title="a0")
double initVal_expCoeffa0
void initGamma(double value, double min, double max, const TString &name="gamma", const TString &title="gamma")
std::unique_ptr< RooRealVar > fsig_
void initA5(double value, double min, double max, const TString &name="a5", const TString &title="a5")
void initSigma(double value, double min, double max, const TString &name="sigma", const TString &title="sigma")
std::unique_ptr< RooRealVar > sigma3_
void initSigma2(double value, double min, double max, const TString &name="sigma2", const TString &title="sigma2")
std::unique_ptr< RooRealVar > mean_
std::unique_ptr< RooRealVar > mean2_
void initMean3(double value, double min, double max, const TString &name="mean3", const TString &title="mean3")
void initGaussFrac2(double value, double min, double max, const TString &name="GaussFrac2", const TString &title="GaussFrac2")
void initA1(double value, double min, double max, const TString &name="a1", const TString &title="a1")
double initVal_mean2
void initFGCB(double value, double min, double max, const TString &name="fGCB", const TString &title="fGCB")
void initA2(double value, double min, double max, const TString &name="a2", const TString &title="a2")
void initA6(double value, double min, double max, const TString &name="a6", const TString &title="a6")
Log< level::Error, false > LogError
std::unique_ptr< RooRealVar > a0_
void initMean2(double value, double min, double max, const TString &name="mean2", const TString &title="mean2")
std::unique_ptr< RooRealVar > a1_
std::unique_ptr< RooRealVar > mean3_
std::unique_ptr< RooRealVar > a5_
std::unique_ptr< RooRealVar > n_
std::unique_ptr< RooAbsPdf > buildSignalModel(RooRealVar *x, const TString &signalType)
Build the model for the specified signal type.
double initVal_expCoeffa1
double initVal_mean
double initVal_gamma
double initVal_alpha
double initVal_gaussFrac2
void reinitializeParameters()
std::unique_ptr< RooRealVar > a6_
double initVal_sigma3
void initN(double value, double min, double max, const TString &name="n", const TString &title="n")
std::unique_ptr< RooRealVar > sigma_
std::unique_ptr< RooAbsPdf > buildBackgroundModel(RooRealVar *x, const TString &backgroundType)
Build the model for the specified background type.
std::unique_ptr< RooRealVar > a4_
void initFsig(double value, double min, double max, const TString &name="fsig", const TString &title="signal fraction")
double initVal_sigma2
Definition: value.py:1
void initExpCoeffA1(double value, double min, double max, const TString &name="expCoeffa1", const TString &title="expCoeffa1")
Log< level::Warning, true > LogPrint
std::unique_ptr< RooRealVar > gaussFrac2_
std::unique_ptr< RooAbsPdf > buildModel(RooRealVar *x, const TString &signalType, const TString &backgroundType)
Build the model to fit.
double initVal_expCoeffa2
std::unique_ptr< RooRealVar > a2_
std::unique_ptr< RooRealVar > sigma2_
void initAlpha(double value, double min, double max, const TString &name="alpha", const TString &title="alpha")
void fit(TH1 *histo, const TString signalType, const TString backgroundType, double xMin=0., double xMax=0., bool sumW2Error=false)
std::unique_ptr< RooDataHist > importTH1(TH1 *histo, double xMin, double xMax)
Definition: FitWithRooFit.cc:6
std::unique_ptr< RooRealVar > gaussFrac_
std::unique_ptr< RooRealVar > gamma_
std::unique_ptr< RooRealVar > fGCB_
std::unique_ptr< RooRealVar > a3_
void initGaussFrac(double value, double min, double max, const TString &name="GaussFrac", const TString &title="GaussFrac")
std::unique_ptr< RooRealVar > expCoeffa0_
std::unique_ptr< RooRealVar > expCoeffa1_
double initVal_fGCB
void initExpCoeffA2(double value, double min, double max, const TString &name="expCoeffa2", const TString &title="expCoeffa2")
double initVal_fsig
float x
void initSigma3(double value, double min, double max, const TString &name="sigma3", const TString &title="sigma3")
void initA3(double value, double min, double max, const TString &name="a3", const TString &title="a3")
Log< level::Warning, false > LogWarning
void initExpCoeffA0(double value, double min, double max, const TString &name="expCoeffa0", const TString &title="expCoeffa0")
double initVal_gaussFrac
dh
Definition: cuy.py:354
std::unique_ptr< RooRealVar > expCoeffa2_
std::unique_ptr< RooRealVar > alpha_
double initVal_mean3
def move(src, dest)
Definition: eostools.py:511
void initMean(double value, double min, double max, const TString &name="mean", const TString &title="mean")
double initVal_sigma
def exit(msg="")
void initA4(double value, double min, double max, const TString &name="a4", const TString &title="a4")