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.

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