Angular analysis of B+->K*+(K+pi0)mumu
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.

230 lines
9.0 KiB

  1. //Renata Kopecna
  2. #include <TCanvas.h>
  3. #include <TPaveText.h>
  4. #include <TTree.h>
  5. #include <TBranch.h>
  6. #include <TF1.h>
  7. #include <TH1D.h>
  8. #include <TStyle.h>
  9. #include <string>
  10. #include <vector>
  11. #include <event.hh>
  12. #include <parameters.hh>
  13. #include <funcs.hh>
  14. #include <pdf.hh>
  15. #include <options.hh>
  16. #include <multifit.hh>
  17. using namespace std;
  18. using namespace fcnc;
  19. multifit::multifit(options* o):
  20. opts(o)
  21. {
  22. output = new TFile(("result" + opts->name + ".root").c_str(), "RECREATE");
  23. };
  24. multifit::~multifit(){
  25. output->Write();
  26. output->Close();
  27. delete output;
  28. };
  29. void multifit::fit_gaussian(vector<double>& values, double& gauss_mean, double& sigma_gauss_mean, double& gauss_width, double& sigma_gauss_width, double& chi_squared) const
  30. {
  31. double max = *max_element(std::begin(values), std::end(values));
  32. double min = *min_element(std::begin(values), std::end(values));
  33. double mean=0.0;
  34. double rms=0.0;
  35. TTree* t = new TTree("tree", "pull values");
  36. double value;
  37. TBranch *branch = t->Branch("value",&value,"value/D");
  38. for (unsigned int i = 0; i < values.size(); i++) {
  39. value = values.at(i);
  40. branch->Fill();
  41. t->Fill();
  42. mean += value/values.size();
  43. }
  44. for (unsigned int i = 0; i < values.size(); i++){
  45. value = values.at(i);
  46. rms += (value - mean)*(value - mean)/values.size();
  47. }
  48. rms = sqrt(rms);
  49. TF1* f1 = new TF1("f1", "gaus(0)/sqrt(2.0*3.1415926536)/abs([2])", min-5.0, max+5.0);
  50. f1->SetParameters(1.0,mean,rms);
  51. //f1->SetParLimits(0, 1.0, 1.0);
  52. f1->FixParameter(0, 1.0);
  53. t->UnbinnedFit("f1", "value");
  54. gauss_mean = f1->GetParameter(1);
  55. sigma_gauss_mean = f1->GetParError(1);
  56. gauss_width = f1->GetParameter(2);
  57. sigma_gauss_width = f1->GetParError(2);
  58. delete f1;
  59. delete t;
  60. };
  61. void multifit::update_pull(parameter* p, vector<double>& values, vector<double>& errors, double& gauss_mean, double& sigma_gauss_mean, double& gauss_width, double& sigma_gauss_width, double& chi_squared, string appendix) const
  62. {
  63. unsigned int runs=values.size();
  64. string parname(p->get_name());
  65. string descr(p->get_description());
  66. gStyle->SetOptStat(0);
  67. gStyle->SetOptFit(0);
  68. //finally the pull histo
  69. TH1D* pull_histo = new TH1D((parname + "_pull_").c_str(),
  70. (descr + " pull distribution;" + "(" + descr
  71. + "^{fitted}-" + descr
  72. + "^{generated})/#sigma").c_str(),
  73. 50, -5.0, 5.0);//sensible for pulls
  74. //fill histos
  75. //this is done every time because the time is negligible
  76. double start_value = p->get_start_value();
  77. vector<double> pull_values;
  78. for (unsigned int k=0; k<runs; k++){
  79. if (errors.at(k) != 0.0){
  80. double pull = (values.at(k)-start_value)/errors.at(k);
  81. pull_histo->Fill(pull);
  82. pull_values.push_back(pull);
  83. }
  84. }
  85. fit_gaussian(pull_values, gauss_mean, sigma_gauss_mean, gauss_width, sigma_gauss_width, chi_squared);
  86. TCanvas *pull_canvas = new TCanvas((parname + "_pull_" + appendix).c_str(), "Bs Likelihood Analysis: Pull distribution",1600,1200);
  87. pull_canvas->cd();
  88. pull_histo->Draw();
  89. TF1 *g = new TF1("g","[0]*exp(-(x-[1])*(x-[1])/2.0/[2]/[2])",-5.0,5.0);
  90. //for normalization: multiply with dx*histint = width/Nbins*histint
  91. g->SetParameters(1.0/sqrt(2.0*TMath::Pi())/gauss_width*10.0/50.0*pull_histo->Integral(), gauss_mean, gauss_width);
  92. g->Draw("same");
  93. TPaveText* text = new TPaveText(0.6,0.7,0.88,0.88,"NDC");
  94. text->SetFillColor(0);
  95. std::ostringstream stream_mean;
  96. stream_mean << "Mean: " << fixed << setprecision(3) << gauss_mean << "#pm" << sigma_gauss_mean;
  97. std::ostringstream stream_width;
  98. stream_width << "Width: " << fixed << setprecision(3) << gauss_width << "#pm" << sigma_gauss_width;
  99. text->AddText(stream_mean.str().c_str());
  100. text->AddText(stream_width.str().c_str());
  101. text->Draw("same");
  102. string afilename; //TODO: move to paths
  103. if (appendix != "") afilename = "plots/pull_" + parname + "_" + appendix + ".eps";
  104. else afilename = "plots/pull_" + parname + ".eps";
  105. if (opts->write_eps)pull_canvas->Print(afilename.c_str(), "eps");
  106. output->cd();
  107. pull_canvas->Write();
  108. delete text;
  109. delete g;
  110. delete pull_histo;
  111. delete pull_canvas;
  112. }
  113. void multifit::update_value(parameter * p, vector<double>& values, vector<double>& errors, double& gauss_mean, double& sigma_gauss_mean, double& gauss_width, double& sigma_gauss_width, double& chi_squared, string appendix) const
  114. {
  115. unsigned int runs=values.size();
  116. string parname(p->get_name());
  117. string descr(p->get_description());
  118. //first find min and max values
  119. double max = *max_element(std::begin(values), std::end(values));
  120. double min = *min_element(std::begin(values), std::end(values));
  121. double width = max-min;
  122. double hist_min = min - 0.5*width;
  123. double hist_max = max + 0.5*width;
  124. gStyle->SetOptStat(0);
  125. gStyle->SetOptFit(0);
  126. TH1D* value_histo = new TH1D((parname+"_values").c_str(),
  127. (descr + ";" + descr).c_str(),
  128. 50, hist_min, hist_max);
  129. for (unsigned int k=0; k<runs; k++) {
  130. value_histo->Fill(values.at(k));
  131. }
  132. fit_gaussian(values, gauss_mean, sigma_gauss_mean, gauss_width, sigma_gauss_width, chi_squared);
  133. TCanvas *value_canvas = new TCanvas((parname + "_value_" + appendix).c_str(), "Bs Likelihood Analysis: Value distribution",1600,1200);
  134. value_canvas->cd();
  135. value_histo->Draw();
  136. TF1 *g = new TF1("g","[0]*exp(-(x-[1])*(x-[1])/2.0/[2]/[2])",-5.0,5.0);
  137. g->SetParameters(1.0/sqrt(2.0*TMath::Pi())/gauss_width*(hist_max-hist_min)/50.0*value_histo->Integral(), gauss_mean, gauss_width);
  138. g->Draw("same");
  139. TPaveText* text = new TPaveText(0.6,0.7,0.88,0.88,"NDC");
  140. text->SetFillColor(0);
  141. std::ostringstream stream_mean;
  142. stream_mean << "Mean: " << fixed << setprecision(3) << gauss_mean << "#pm" << sigma_gauss_mean;
  143. std::ostringstream stream_width;
  144. stream_width << "Width: " << fixed << setprecision(3) << gauss_width << "#pm" << sigma_gauss_width;
  145. text->AddText(stream_mean.str().c_str());
  146. text->AddText(stream_width.str().c_str());
  147. text->Draw("same");
  148. //TODO: move
  149. string afilename;
  150. if (appendix != "") afilename = "plots/value_" + parname + "_" + appendix + ".eps";
  151. else afilename = "plots/value_" + parname + ".eps";
  152. if (opts->write_eps) value_canvas->Print(afilename.c_str(), "eps");
  153. output->cd();
  154. value_canvas->Write();
  155. delete text;
  156. delete value_histo;
  157. delete value_canvas;
  158. delete g;
  159. }
  160. void multifit::update_error(parameter* p, vector<double>& values, vector<double>& errors, double& gauss_mean, double& sigma_gauss_mean, double& gauss_width, double& sigma_gauss_width, double& chi_squared, string appendix) const
  161. {
  162. unsigned int runs=values.size();
  163. string parname(p->get_name());
  164. string descr(p->get_description());
  165. //then find min and max errors
  166. double max = *max_element(std::begin(errors), std::end(errors));
  167. double min = *min_element(std::begin(errors), std::end(errors));
  168. double width = max - min;
  169. double hist_min = min - 0.5*width;
  170. double hist_max = max + 0.5*width;
  171. gStyle->SetOptStat(0);
  172. gStyle->SetOptFit(0);
  173. TH1D* error_histo = new TH1D((parname+"_errors").c_str(),
  174. ("#sigma(" + descr + ");#sigma(" + descr + ")").c_str(),
  175. 100, hist_min, hist_max);
  176. for (unsigned int k=0; k<runs; k++){
  177. error_histo->Fill(errors.at(k));
  178. }
  179. fit_gaussian(errors, gauss_mean, sigma_gauss_mean, gauss_width, sigma_gauss_width, chi_squared);
  180. TCanvas *error_canvas = new TCanvas((parname + "_error_" + appendix).c_str(), "Bs Likelihood Analysis: Error distribution",1600,1200);
  181. error_canvas->cd();
  182. error_histo->Draw();
  183. TF1 *g = new TF1("g","[0]*exp(-(x-[1])*(x-[1])/2.0/[2]/[2])",-5.0,5.0);
  184. g->SetParameters(1.0/sqrt(2.0*TMath::Pi())/gauss_width*(hist_max-hist_min)/100.0*error_histo->Integral(), gauss_mean, gauss_width);
  185. g->Draw("same");
  186. TPaveText* text = new TPaveText(0.6,0.7,0.88,0.88,"NDC");
  187. text->SetFillColor(0);
  188. std::ostringstream stream_mean;
  189. stream_mean << "Mean: " << fixed << setprecision(3) << gauss_mean << "#pm" << sigma_gauss_mean;
  190. std::ostringstream stream_width;
  191. stream_width << "Width: " << fixed << setprecision(3) << gauss_width << "#pm" << sigma_gauss_width;
  192. text->AddText(stream_mean.str().c_str());
  193. text->AddText(stream_width.str().c_str());
  194. text->Draw("same");
  195. //TODO: move
  196. string afilename;
  197. if (appendix != "") afilename = "plots/error_" + parname + "_" + appendix + ".eps";
  198. else afilename = "plots/error_" + parname + ".eps";
  199. if (opts->write_eps) error_canvas->Print(afilename.c_str(), "eps");
  200. output->cd();
  201. error_canvas->Write();
  202. delete text;
  203. delete error_histo;
  204. delete error_canvas;
  205. delete g;
  206. }