ROOT Analysis for the Inclusive Detachted Dilepton Trigger Lines
You can not select more than 25 topics Topics must start with a letter or number, can include dashes ('-') and can be up to 35 characters long.

557 lines
19 KiB

7 months ago
7 months ago
7 months ago
7 months ago
7 months ago
7 months ago
7 months ago
7 months ago
7 months ago
7 months ago
7 months ago
7 months ago
7 months ago
7 months ago
7 months ago
7 months ago
7 months ago
7 months ago
  1. #ifndef BASIC_ANALYSIS
  2. #define BASIC_ANALYSIS
  3. #include <string>
  4. #include <iostream>
  5. #include <cmath>
  6. #include <algorithm>
  7. #include <filesystem>
  8. #include <string_view>
  9. #include "TTree.h"
  10. #include "RooHist.h"
  11. #include "constants.h"
  12. std::pair<double, double> DivWithErr(double x, double dx, double y, double dy) {
  13. double err_x = (1 / y) * dx;
  14. double err_y = - (x / (y*y)) * dy;
  15. return std::make_pair(x / y, TMath::Sqrt((TMath::Sq(err_x) + TMath::Sq(err_y))));
  16. }
  17. class FourVect
  18. {
  19. private:
  20. Float_t px, py, pz, energy;
  21. Float_t *PtrPX()
  22. {
  23. return &px;
  24. }
  25. Float_t *PtrPY()
  26. {
  27. return &py;
  28. }
  29. Float_t *PtrPZ()
  30. {
  31. return &pz;
  32. }
  33. Float_t *PtrEnergy()
  34. {
  35. return &energy;
  36. }
  37. public:
  38. static FourVect *Init(TTree *tree, const char *particle)
  39. {
  40. FourVect *v = new FourVect{};
  41. tree->SetBranchAddress(TString::Format("%s_PX", particle).Data(), v->PtrPX());
  42. tree->SetBranchAddress(TString::Format("%s_PY", particle).Data(), v->PtrPY());
  43. tree->SetBranchAddress(TString::Format("%s_PZ", particle).Data(), v->PtrPZ());
  44. tree->SetBranchAddress(TString::Format("%s_ENERGY", particle).Data(), v->PtrEnergy());
  45. return v;
  46. }
  47. TLorentzVector LorentzVector()
  48. {
  49. return TLorentzVector(px, py, pz, energy);
  50. }
  51. TLorentzVector LorentzVector(double sub_mass)
  52. {
  53. TVector3 momentum(px, py, pz);
  54. float energy = TMath::Sqrt(TMath::Sq(sub_mass) + momentum.Mag2());
  55. return TLorentzVector(momentum, energy);
  56. }
  57. std::string ToString()
  58. {
  59. return TString::Format("(%f, %f, %f, %f)", px, py, pz, energy).Data();
  60. }
  61. };
  62. struct FittedParam
  63. {
  64. std::string title;
  65. std::string name;
  66. double value;
  67. double error;
  68. bool at_limit;
  69. bool constant;
  70. int decimals;
  71. FittedParam(RooRealVar var, int decimals)
  72. {
  73. this->title = var.GetTitle();
  74. this->name = var.GetName();
  75. double val = var.getVal();
  76. this->value = val;
  77. this->constant = var.isConstant();
  78. if (!this->constant)
  79. {
  80. this->error = var.getError();
  81. this->at_limit = ((val == var.getMin()) || (val == var.getMax()));
  82. }
  83. else
  84. {
  85. this->error = 0.;
  86. }
  87. this->decimals = decimals;
  88. }
  89. FittedParam(std::string name, std::string title, double value, double err, int decimals)
  90. {
  91. this->title = title;
  92. this->name = name;
  93. this->value = value;
  94. this->error = err;
  95. this->decimals = decimals;
  96. this->constant = false;
  97. this->at_limit = false;
  98. }
  99. std::string ToString() const
  100. {
  101. if (this->constant)
  102. {
  103. return TString::Format("%s = %.*f", this->title.c_str(), this->decimals, this->value).Data();
  104. }
  105. else
  106. {
  107. return TString::Format("%s = %.*f #pm %.*f", this->title.c_str(), this->decimals, this->value, this->decimals, this->error).Data();
  108. }
  109. }
  110. };
  111. struct ShapeParamters
  112. {
  113. double alpha_left;
  114. double n_left;
  115. double alpha_right;
  116. double n_right;
  117. double sigma_lr;
  118. };
  119. struct RooFitSummary
  120. {
  121. RooPlot *fit_histogram;
  122. RooPlot *pull_histogram;
  123. std::map<std::string, std::string> pdf_names;
  124. std::pair<double, double> signal_yield;
  125. std::pair<double, double> background_yield;
  126. std::vector<FittedParam> fitted_params;
  127. ShapeParamters shape_parameters;
  128. };
  129. void DrawInDefaultCanvas(TH1 *histogram, const char *folder, double margin_left = 0.1, Option_t *option = "")
  130. {
  131. std::filesystem::create_directory(TString::Format("output_files/analysis/%s", folder).Data());
  132. TString name = TString::Format("%s_canvas", histogram->GetName());
  133. TCanvas *c = new TCanvas(name, histogram->GetName(), 0, 0, 800, 600);
  134. c->SetLeftMargin(margin_left);
  135. histogram->SetStats(0);
  136. histogram->Draw(option);
  137. c->Draw();
  138. c->SaveAs(TString::Format("output_files/analysis/%s/%s.pdf", folder, name.Data()).Data());
  139. }
  140. void DrawInDefaultCanvasStacked(std::vector<TH1D *> histograms, std::vector<Color_t> colors, std::vector<Style_t> fill_style, const char *folder, double margin_left = 0.1, Option_t *option = "")
  141. {
  142. std::filesystem::create_directory(TString::Format("output_files/analysis/%s", folder).Data());
  143. TString name = TString::Format("%s_stack_canvas", histograms[0]->GetName());
  144. TCanvas *c = new TCanvas(name, histograms[0]->GetName(), 0, 0, 800, 600);
  145. c->SetLeftMargin(margin_left);
  146. std::string drwopt_2 = std::string(option).empty() ? "SAME HIST" : TString::Format("%s SAME HIST", option).Data();
  147. for (size_t i = 0; i < histograms.size(); i++)
  148. {
  149. const Double_t scaling_factor = 1.;
  150. auto hist_clone = (TH1 *)histograms[i]->Clone(TString::Format("%s_clone", histograms[i]->GetName()));
  151. hist_clone->Scale(scaling_factor / histograms[i]->GetMaximum());
  152. hist_clone->SetLineColor(colors[i]);
  153. hist_clone->SetMaximum(scaling_factor + (scaling_factor * 0.05));
  154. hist_clone->SetMinimum(0.);
  155. hist_clone->SetStats(0);
  156. if (fill_style[i] != 0)
  157. {
  158. hist_clone->SetFillStyle(fill_style[i]);
  159. hist_clone->SetFillColor(colors[i]);
  160. }
  161. if (i > 0)
  162. {
  163. hist_clone->Draw(drwopt_2.c_str());
  164. }
  165. else
  166. {
  167. hist_clone->Draw(drwopt_2.c_str());
  168. }
  169. }
  170. c->BuildLegend(0.55, 0.79, 0.87, 0.89);
  171. c->Draw();
  172. c->SaveAs(TString::Format("output_files/analysis/%s/%s.pdf", folder, name.Data()).Data());
  173. }
  174. void DrawInDefaultCanvas(RooPlot *rooPlot, const char *folder, double margin_left = 0, Option_t *option = "")
  175. {
  176. std::filesystem::create_directory(TString::Format("output_files/analysis/%s", folder).Data());
  177. TString name = TString::Format("%s_canvas", rooPlot->GetName());
  178. TCanvas *c = new TCanvas(name, rooPlot->GetName(), 0, 0, 800, 600);
  179. c->SetLeftMargin(margin_left);
  180. rooPlot->Draw(option);
  181. c->Draw();
  182. c->SaveAs(TString::Format("output_files/analysis/%s/%s.pdf", folder, name.Data()).Data());
  183. }
  184. void PrintStyles(RooPlot *plt)
  185. {
  186. auto Xaxis = plt->GetXaxis();
  187. auto Yaxis = plt->GetYaxis();
  188. std::cout << "## Styles for plot " << plt->GetName() << std::endl;
  189. std::cout << "## X axis\n"
  190. << "Label Size: " << Xaxis->GetLabelSize() << ",\n"
  191. << "Label Offset: " << Xaxis->GetLabelOffset() << ",\n"
  192. << "Title Size: " << Xaxis->GetTitleSize() << ",\n"
  193. << "Title Offset: " << Xaxis->GetTitleOffset() << "." << std::endl;
  194. std::cout << "## Y axis\n"
  195. << "Label Size: " << Yaxis->GetLabelSize() << ",\n"
  196. << "Label Offset: " << Yaxis->GetLabelOffset() << ",\n"
  197. << "Title Size: " << Yaxis->GetTitleSize() << ",\n"
  198. << "Title Offset: " << Yaxis->GetTitleOffset() << "." << std::endl;
  199. }
  200. void DrawInDefaultCanvas(RooFitSummary fitSummary, const char *folder)
  201. {
  202. std::filesystem::create_directory(TString::Format("output_files/analysis/%s", folder).Data());
  203. TString name = TString::Format("%s_canvas", fitSummary.fit_histogram->GetName());
  204. TCanvas *c = new TCanvas(name, fitSummary.fit_histogram->GetTitle(), 0, 0, 800, 600);
  205. Double_t xlow, ylow, xup, yup;
  206. c->GetPad(0)->GetPadPar(xlow, ylow, xup, yup);
  207. c->Divide(1, 2);
  208. Double_t upPad_ylow = ylow + 0.25 * (yup - ylow);
  209. Double_t dwPad_yup = ylow + 0.25 * (yup - ylow);
  210. TVirtualPad *upPad = c->GetPad(1);
  211. upPad->SetPad(xlow, upPad_ylow, xup, yup);
  212. TVirtualPad *dwPad = c->GetPad(2);
  213. dwPad->SetPad(xlow, ylow, xup, dwPad_yup);
  214. dwPad->SetTopMargin(0);
  215. Double_t upPad_area = (yup - upPad_ylow) * (xup - xlow);
  216. Double_t dwPad_area = (dwPad_yup - ylow) * (xup - xlow);
  217. // std::cout << "### UP AREA: " << upPad_area << " LOW AREA: " << dwPad_area << std::endl;
  218. const Double_t pull_title_size = 0.09;
  219. const Double_t pull_label_size = 0.09;
  220. Double_t fit_title_size = (pull_title_size * dwPad_area) / upPad_area;
  221. Double_t fit_label_size = (pull_label_size * dwPad_area) / upPad_area;
  222. // canvas->Update();
  223. c->cd(1);
  224. // TPad *pull_pad = new TPad(TString::Format("%s_canv_pull", fitSummary.fit_histogram->GetName()), "Pull Pad", 0., 0., 1., 0.3);
  225. // pull_pad->Draw();
  226. // pull_pad->SetTopMargin(0.001);
  227. // pull_pad->SetBottomMargin(0.3);
  228. // pull_pad->SetGrid();
  229. // TPad *fit_pad = new TPad(TString::Format("%s_canv_fit", fitSummary.fit_histogram->GetName()), "Fit Pad", 0., 0.3, 1., 1.);
  230. // fit_pad->Draw();
  231. // fit_pad->SetBottomMargin(0.001);
  232. // fit_pad->cd();
  233. auto fit_Xaxis = fitSummary.fit_histogram->GetXaxis();
  234. auto fit_Yaxis = fitSummary.fit_histogram->GetYaxis();
  235. fit_Xaxis->SetLabelOffset(0.02);
  236. fit_Xaxis->SetLabelSize(fit_label_size);
  237. fit_Xaxis->SetTitleOffset(1.4);
  238. fit_Xaxis->SetTitleSize(fit_title_size);
  239. fit_Yaxis->SetLabelOffset(0.01);
  240. fit_Yaxis->SetLabelSize(fit_label_size);
  241. fit_Yaxis->SetTitleSize(fit_title_size);
  242. fit_Yaxis->SetTitleOffset(0);
  243. fitSummary.fit_histogram->Draw();
  244. TLegend *leg1 = new TLegend(0.55, 0.45, 0.87, 0.89);
  245. leg1->SetFillColor(kWhite);
  246. leg1->SetLineColor(kWhite);
  247. for (const auto &[name, title] : fitSummary.pdf_names)
  248. {
  249. leg1->AddEntry(fitSummary.fit_histogram->findObject(name.c_str()), title.c_str(), "LP");
  250. }
  251. for (const auto &par : fitSummary.fitted_params)
  252. {
  253. if (!par.constant)
  254. {
  255. leg1->AddEntry((TObject *)0, par.ToString().c_str(), "");
  256. }
  257. }
  258. leg1->Draw();
  259. c->cd(2);
  260. Double_t pull_min = fitSummary.pull_histogram->GetMinimum();
  261. Double_t pull_max = fitSummary.pull_histogram->GetMaximum();
  262. // std::cout << "### (" << fitSummary.pull_histogram->GetName() << ") PULL MIN: " << pull_min << " PULL MAX: " << pull_max << std::endl;
  263. double bound = 0;
  264. if (TMath::Abs(pull_min) > TMath::Abs(pull_max))
  265. {
  266. bound = TMath::Abs(pull_min);
  267. }
  268. else
  269. {
  270. bound = TMath::Abs(pull_max);
  271. }
  272. auto pull_Xaxis = fitSummary.pull_histogram->GetXaxis();
  273. auto pull_Yaxis = fitSummary.pull_histogram->GetYaxis();
  274. pull_Xaxis->SetLabelOffset(0.02);
  275. pull_Xaxis->SetLabelSize(pull_label_size);
  276. pull_Xaxis->SetTitleSize(pull_title_size);
  277. pull_Xaxis->SetTitle("");
  278. pull_Yaxis->SetLabelOffset(0.01);
  279. pull_Yaxis->SetLabelSize(pull_label_size);
  280. pull_Yaxis->SetTitleSize(pull_title_size);
  281. pull_Yaxis->SetTitleOffset(0.45);
  282. pull_Yaxis->SetTitle("Pull Distribution");
  283. fitSummary.pull_histogram->SetTitle("");
  284. fitSummary.pull_histogram->SetMaximum(bound);
  285. fitSummary.pull_histogram->SetMinimum(-bound);
  286. fitSummary.pull_histogram->Draw();
  287. auto line_mid = new TLine(B_MASS_VAR_MIN, 0, B_MASS_VAR_MAX, 0);
  288. line_mid->SetLineColor(kOrange + 4);
  289. line_mid->Draw();
  290. auto line_up = new TLine(B_MASS_VAR_MIN, bound / 2, B_MASS_VAR_MAX, bound / 2);
  291. line_up->SetLineStyle(kDashed);
  292. line_up->SetLineColor(kOrange + 4);
  293. line_up->Draw();
  294. auto line_low = new TLine(B_MASS_VAR_MIN, -bound / 2, B_MASS_VAR_MAX, -bound / 2);
  295. line_low->SetLineStyle(kDashed);
  296. line_low->SetLineColor(kOrange + 4);
  297. line_low->Draw();
  298. c->Draw();
  299. c->SaveAs(TString::Format("output_files/analysis/%s/%s.pdf", folder, name.Data()).Data());
  300. }
  301. void PrintProgress(const char *title, unsigned int total, unsigned int every, unsigned int current)
  302. {
  303. if ((current + 1) % every == 0 || current + 1 == total)
  304. {
  305. std::cout << "["
  306. << title
  307. << "] Processed event: " << current + 1 << " (" << TString::Format("%.2f", ((double)(current + 1) / (double)total) * 100.) << "%)" << std::endl;
  308. }
  309. }
  310. RooPlot *CreateRooFitHistogram(TH1D *hist)
  311. {
  312. RooRealVar roo_var_mass("var_mass", "B Mass Variable", B_MASS_VAR_MIN, B_MASS_VAR_MAX);
  313. roo_var_mass.setRange("fitting_range", B_MASS_VAR_MIN, B_MASS_VAR_MAX);
  314. RooDataHist roohist_B_M("roohist_B_M", "B Mass Histogram", roo_var_mass, RooFit::Import(*hist));
  315. RooPlot *roo_frame_mass = roo_var_mass.frame(RooFit::Title(hist->GetTitle()), RooFit::Name(TString::Format("%s_rplt", hist->GetName())));
  316. roohist_B_M.plotOn(roo_frame_mass, RooFit::Binning(B_MASS_HIST_BINS), RooFit::Name("B Mass Distribution"));
  317. roo_frame_mass->GetXaxis()->SetTitle(hist->GetXaxis()->GetTitle());
  318. return roo_frame_mass;
  319. }
  320. RooFitSummary CreateRooDataSetAndFitCB(TTree *dataSet, TString var_name, TString xAxis, bool hasExpBkg, bool useExtShape, ShapeParamters extShape, bool const_sigma = false, Double_t fit_low = 0, Double_t fit_up = 0)
  321. {
  322. auto suffix_name = [name = dataSet->GetName()](const char *text)
  323. {
  324. return TString::Format("%s_%s", text, name);
  325. };
  326. Double_t fitRangeUp = fit_low == 0 ? B_MASS_VAR_MIN : fit_low;
  327. Double_t fitRangeLow = fit_up == 0 ? B_MASS_VAR_MAX : fit_up;
  328. RooRealVar roo_var_mass(var_name, "B Mass Variable", B_MASS_VAR_MIN, B_MASS_VAR_MAX);
  329. const char *fitting_range_name = "fitting_range";
  330. roo_var_mass.setRange(fitting_range_name, fitRangeUp, fitRangeLow);
  331. TString dataset_name = suffix_name("roodataset_B_M");
  332. // RooDataHist roodataset_B_M(hist_name, "B Mass Histogram", roo_var_mass, RooFit::Import(*hist));
  333. RooDataSet roodataset_B_M(dataset_name, "B Mass Data Set", roo_var_mass, RooFit::Import(*dataSet));
  334. RooPlot *roo_frame_mass = roo_var_mass.frame(RooFit::Title(dataSet->GetTitle()), RooFit::Name(TString::Format("%s_rplt", dataSet->GetName())));
  335. roodataset_B_M.plotOn(roo_frame_mass, RooFit::Name(dataset_name));
  336. roo_frame_mass->GetXaxis()->SetTitle(xAxis);
  337. // Crystal Ball for Signal
  338. RooRealVar var_mass_x0(suffix_name("var_mass_x0"), "#mu", 5278., 5170., 5500.);
  339. RooRealVar var_mass_sigmaLR(suffix_name("var_mass_sigmaLR"), "#sigma_{LR}", 16., 5., 40.);
  340. RooRealVar var_mass_alphaL(suffix_name("var_mass_alphaL"), "#alpha_{L}", 2., 0., 4.);
  341. RooRealVar var_mass_nL(suffix_name("var_mass_nL"), "n_{L}", 5., 0., 15.);
  342. RooRealVar var_mass_alphaR(suffix_name("var_mass_alphaR"), "#alpha_{R}", 2., 0., 4.);
  343. RooRealVar var_mass_nR(suffix_name("var_mass_nR"), "n_{R}", 5., 0., 15.);
  344. if (useExtShape)
  345. {
  346. if (extShape.alpha_left != 0.)
  347. {
  348. var_mass_alphaL.setConstant(true);
  349. var_mass_alphaL.setVal(extShape.alpha_left);
  350. }
  351. if (extShape.n_left != 0.)
  352. {
  353. var_mass_nL.setConstant(true);
  354. var_mass_nL.setVal(extShape.n_left);
  355. }
  356. if (extShape.alpha_right != 0.)
  357. {
  358. var_mass_alphaR.setConstant(true);
  359. var_mass_alphaR.setVal(extShape.alpha_right);
  360. }
  361. if (extShape.n_right != 0.)
  362. {
  363. var_mass_nR.setConstant(true);
  364. var_mass_nR.setVal(extShape.n_right);
  365. }
  366. if (extShape.sigma_lr != 0. && const_sigma)
  367. {
  368. var_mass_sigmaLR.setConstant(true);
  369. var_mass_sigmaLR.setVal(extShape.sigma_lr);
  370. }
  371. }
  372. TString signal_name = suffix_name("sig_cb");
  373. RooCrystalBall sig_cb(signal_name, "CB Signal", roo_var_mass, var_mass_x0, var_mass_sigmaLR, var_mass_alphaL, var_mass_nL, var_mass_alphaR, var_mass_nR);
  374. std::vector<FittedParam> fitted_params{};
  375. std::map<std::string, std::string> pdf_names{};
  376. pdf_names[signal_name.Data()] = sig_cb.getTitle().Data();
  377. TString pull_compare_name{};
  378. double sig_val = 0.;
  379. double sig_err = 0.;
  380. double bkg_val = 0.;
  381. double bkg_err = 0.;
  382. if (hasExpBkg)
  383. {
  384. // Exponential for Background
  385. RooRealVar var_mass_bkg_c(suffix_name("var_mass_bkg_c"), "#lambda", -0.0014, -0.004, -0.000);
  386. TString background_name = suffix_name("bgk_exp");
  387. RooExponential bkg_exp(background_name, "Exp Background", roo_var_mass, var_mass_bkg_c);
  388. pdf_names[background_name.Data()] = bkg_exp.getTitle().Data();
  389. RooRealVar var_mass_nsig(suffix_name("nsig"), "N_{Sig}", 0., dataSet->GetEntries());
  390. RooRealVar var_mass_nbkg(suffix_name("nbkg"), "N_{Bkg}", 0., dataSet->GetEntries());
  391. TString fitted_pdf_name = suffix_name("sigplusbkg");
  392. RooAddPdf fitted_pdf(fitted_pdf_name, "Sig + Bkg", RooArgList(sig_cb, bkg_exp), RooArgList(var_mass_nsig, var_mass_nbkg));
  393. pdf_names[fitted_pdf_name.Data()] = fitted_pdf.getTitle().Data();
  394. RooFitResult *fitres = fitted_pdf.fitTo(roodataset_B_M, RooFit::Save(), RooFit::PrintLevel(1), RooFit::Range(fitting_range_name));
  395. // sigplusbkg.plotOn(roo_frame_mass, RooFit::VisualizeError(*fitres, 1), RooFit::FillColor(kOrange + 1), RooFit::FillStyle(3144));
  396. fitted_pdf.plotOn(roo_frame_mass, RooFit::LineColor(kRed), RooFit::LineStyle(kSolid), RooFit::Range(fitting_range_name), RooFit::Name(fitted_pdf_name));
  397. fitted_pdf.plotOn(roo_frame_mass, RooFit::Components(RooArgSet(bkg_exp)), RooFit::LineColor(kBlue - 7), RooFit::LineStyle(kDashed), RooFit::Range(fitting_range_name), RooFit::Name(background_name));
  398. fitted_pdf.plotOn(roo_frame_mass, RooFit::Components(RooArgSet(sig_cb)), RooFit::LineColor(kRed - 7), RooFit::LineStyle(kDashed), RooFit::Range(fitting_range_name), RooFit::Name(signal_name));
  399. fitted_params.push_back(FittedParam(var_mass_nsig, 2));
  400. fitted_params.push_back(FittedParam(var_mass_nbkg, 2));
  401. sig_val = var_mass_nsig.getVal();
  402. sig_err = var_mass_nsig.getError();
  403. bkg_val = var_mass_nbkg.getVal();
  404. bkg_err = var_mass_nbkg.getError();
  405. double significance_val = sig_val / TMath::Sqrt(sig_val + bkg_val);
  406. double err_prop_sig = (sig_val + 2 * bkg_val) / (2 * TMath::Power((sig_val + bkg_val), (3 / 2)));
  407. double err_prop_bkg = -sig_val / (2 * TMath::Power((sig_val + bkg_val), (3 / 2)));
  408. double significance_err = TMath::Sqrt(TMath::Sq(err_prop_sig * sig_err) + TMath::Sq(err_prop_bkg * bkg_err));
  409. auto sig_over_bkg = DivWithErr(sig_val, sig_err, bkg_val, bkg_err);
  410. fitted_params.push_back(FittedParam("significance", "N_{Sig}/#sqrt{N_{Sig} + N_{Bkg}}", significance_val, significance_err, 2));
  411. fitted_params.push_back(FittedParam("sig_over_bkg", "N_{Sig}/N_{Bkg}", sig_over_bkg.first, sig_over_bkg.second, 2));
  412. fitted_params.push_back(FittedParam(var_mass_bkg_c, 5));
  413. pull_compare_name = fitted_pdf_name;
  414. }
  415. else
  416. {
  417. RooFitResult *fitres = sig_cb.fitTo(roodataset_B_M, RooFit::Save(), RooFit::PrintLevel(1), RooFit::Range(fitting_range_name));
  418. // sigplusbkg.plotOn(roo_frame_mass, RooFit::VisualizeError(*fitres, 1), RooFit::FillColor(kOrange + 1), RooFit::FillStyle(3144));
  419. sig_cb.plotOn(roo_frame_mass, RooFit::LineColor(kRed), RooFit::LineStyle(kSolid), RooFit::Range(fitting_range_name), RooFit::Name(signal_name));
  420. pull_compare_name = signal_name;
  421. }
  422. fitted_params.push_back(FittedParam(var_mass_x0, 2));
  423. fitted_params.push_back(FittedParam(var_mass_sigmaLR, 2));
  424. fitted_params.push_back(FittedParam(var_mass_alphaL, 2));
  425. fitted_params.push_back(FittedParam(var_mass_nL, 2));
  426. fitted_params.push_back(FittedParam(var_mass_alphaR, 2));
  427. fitted_params.push_back(FittedParam(var_mass_nR, 2));
  428. RooPlot *roo_frame_pull = roo_var_mass.frame(RooFit::Title("Pull Distribution"), RooFit::Name(TString::Format("%s_pull_rplt", dataSet->GetName())));
  429. RooHist *pull_hist = roo_frame_mass->pullHist(dataset_name, pull_compare_name);
  430. roo_frame_pull->addPlotable(pull_hist, "P");
  431. // Double_t pull_min = pull_hist->GetMinimum();
  432. // Double_t pull_max = pull_hist->GetMaximum();
  433. // std::cout << "##### (" << pull_hist->GetName() << ") PULL MIN: " << pull_min << " PULL MAX: " << pull_max << std::endl;
  434. return RooFitSummary{
  435. roo_frame_mass,
  436. roo_frame_pull,
  437. pdf_names,
  438. std::make_pair(sig_val, sig_err),
  439. std::make_pair(bkg_val, bkg_err),
  440. fitted_params,
  441. ShapeParamters{
  442. var_mass_alphaL.getVal(),
  443. var_mass_nL.getVal(),
  444. var_mass_alphaR.getVal(),
  445. var_mass_nR.getVal(),
  446. var_mass_sigmaLR.getVal()}};
  447. }
  448. bool InRange(double value, double center, double low_intvl, double up_intvl)
  449. {
  450. return center - low_intvl < value && value < center + up_intvl;
  451. }
  452. bool InRange(double value, double center, double intvl)
  453. {
  454. return InRange(value, center, intvl, intvl);
  455. }
  456. #endif