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.

229 lines
6.2 KiB

  1. #include <cmath>
  2. #include <iostream>
  3. #include <string>
  4. #include <sstream>
  5. #include <fstream>
  6. #include <cstdlib>
  7. #include <vector>
  8. #include <map>
  9. #include <unistd.h>
  10. #include <TCanvas.h>
  11. #include <TChain.h>
  12. #include <TFile.h>
  13. #include <TGaxis.h>
  14. #include <TH1D.h>
  15. #include <TROOT.h>
  16. #include <TStyle.h>
  17. #include "../Paths.hpp"
  18. #include "../GlobalFunctions.hh"
  19. //Adds q2 and b mass branches into the evtGen tuple
  20. bool addQ2ToEvtGenMC(bool ReferenceChannel, bool PHSP){
  21. gStyle -> SetOptStat(0);
  22. LHCbStyle();
  23. gROOT->SetBatch(kTRUE);
  24. TH1::SetDefaultSumw2(kTRUE);
  25. TChain * ch = new TChain("MCDecayTreeTuple/MCDecayTree");
  26. string file = GetGenLevelFile(ReferenceChannel,PHSP);
  27. replace(file,".root","_orig.root");
  28. ch->Add(file.c_str());
  29. //Reading the initial tree
  30. unsigned int N = ch->GetEntries();
  31. if(N == 0){
  32. coutERROR("No entries found in TFiles");
  33. return 1;
  34. }
  35. else coutInfo("Processing " +to_string(N) + " events!");
  36. //Create new tree
  37. TTree * newTree = ch->CloneTree(0);
  38. newTree->SetName("DecayTree");
  39. double J_psi_M = 0;
  40. newTree->Branch("J_psi_M", &J_psi_M, "J_psi_M/D");
  41. double J_psi_PX = 0;
  42. newTree->Branch("J_psi_PX",&J_psi_PX,"J_psi_PX/D");
  43. double J_psi_PY = 0;
  44. newTree->Branch("J_psi_PY",&J_psi_PY,"J_psi_PY/D");
  45. double J_psi_PZ = 0;
  46. newTree->Branch("J_psi_PZ",&J_psi_PZ,"J_psi_PZ/D");
  47. double J_psi_PE = 0;
  48. newTree->Branch("J_psi_PE",&J_psi_PE,"J_psi_PE/D");
  49. double Q2 = 0;
  50. newTree->Branch("Q2", &Q2, "Q2/D");
  51. double B_plus_M = 0;
  52. newTree->Branch("B_plus_M", &B_plus_M, "B_plus_M/D");
  53. TFile *output = new TFile(GetGenLevelFile(ReferenceChannel,PHSP).c_str(),"RECREATE");
  54. double mm_pe, mm_px, mm_py, mm_pz;
  55. double mp_pe, mp_px, mp_py, mp_pz;
  56. double b_pe, b_px, b_py, b_pz;
  57. ch->SetBranchStatus("*",1);
  58. ch->SetBranchAddress("muplus_TRUEP_E", &mp_pe);
  59. ch->SetBranchAddress("muplus_TRUEP_X", &mp_px);
  60. ch->SetBranchAddress("muplus_TRUEP_Y", &mp_py);
  61. ch->SetBranchAddress("muplus_TRUEP_Z", &mp_pz);
  62. ch->SetBranchAddress("muminus_TRUEP_E", &mm_pe);
  63. ch->SetBranchAddress("muminus_TRUEP_X", &mm_px);
  64. ch->SetBranchAddress("muminus_TRUEP_Y", &mm_py);
  65. ch->SetBranchAddress("muminus_TRUEP_Z", &mm_pz);
  66. ch->SetBranchAddress("Bplus_TRUEP_E", &b_pe);
  67. ch->SetBranchAddress("Bplus_TRUEP_X", &b_px);
  68. ch->SetBranchAddress("Bplus_TRUEP_Y", &b_py);
  69. ch->SetBranchAddress("Bplus_TRUEP_Z", &b_pz);
  70. TLorentzVector tmp_plus;
  71. TLorentzVector tmp_minus;
  72. TLorentzVector tmp_J_psi;
  73. TLorentzVector tmp_Bplus;
  74. vector<double> double_E;
  75. for(unsigned int e = 0; e < N; e++){
  76. ch->GetEntry(e);
  77. tmp_plus.SetPxPyPzE(mp_px,mp_py,mp_pz,mp_pe);
  78. tmp_minus.SetPxPyPzE(mm_px,mm_py,mm_pz,mm_pe);
  79. tmp_Bplus.SetPxPyPzE(b_px,b_py,b_pz,b_pe);
  80. double_E.push_back(mm_pe);
  81. double_E.push_back(mp_pe);
  82. tmp_J_psi = tmp_plus + tmp_minus;
  83. J_psi_M = tmp_J_psi.M();
  84. J_psi_PX = tmp_J_psi.Px();
  85. J_psi_PY = tmp_J_psi.Py();
  86. J_psi_PZ = tmp_J_psi.Pz();
  87. J_psi_PE = tmp_J_psi.E();
  88. Q2 = tmp_J_psi.M2();
  89. B_plus_M = tmp_Bplus.M();
  90. newTree->Fill();
  91. }
  92. //Event number multiple candidates
  93. int counter = 0;
  94. for_indexed(auto firstNumber: double_E){
  95. for (vector<double>::iterator evt = double_E.begin()+i+1; evt < double_E.end(); evt++){
  96. if (double_E.at(i)==*evt)counter++;
  97. }
  98. }
  99. coutInfo("Total number of double candidates: " +to_string(counter));
  100. output->cd();
  101. newTree->Write("",TObject::kWriteDelete);
  102. output->Close();
  103. return 1;
  104. }
  105. void getQ2weights(){
  106. bool onlyCandOne = false;
  107. bool onlynTotOne = false;
  108. bool WriteToROOTfile = true;
  109. gStyle -> SetOptStat(0);
  110. LHCbStyle();
  111. gROOT->SetBatch(kTRUE);
  112. TH1::SetDefaultSumw2(kTRUE);
  113. double hmin = 0.1, hmax = 19.0;
  114. TChain * ch = new TChain("MCDecayTreeTuple/MCDecayTree");
  115. ch->Add("MCTruth300k_*.sim.root");
  116. unsigned int N = ch->GetEntries();
  117. if(N == 0){
  118. std::cout << "[ERROR]\tNo entries found in TFiles" << std::endl;
  119. return;
  120. }
  121. else
  122. std::cout << "[INFO]\t\tProcess " << N << " events!" << std::endl;
  123. double mm_pe, mm_px, mm_py, mm_pz;
  124. double mp_pe, mp_px, mp_py, mp_pz;
  125. ch->SetBranchAddress("muplus_TRUEP_E", &mp_pe);
  126. ch->SetBranchAddress("muplus_TRUEP_X", &mp_px);
  127. ch->SetBranchAddress("muplus_TRUEP_Y", &mp_py);
  128. ch->SetBranchAddress("muplus_TRUEP_Z", &mp_pz);
  129. ch->SetBranchAddress("muminus_TRUEP_E", &mm_pe);
  130. ch->SetBranchAddress("muminus_TRUEP_X", &mm_px);
  131. ch->SetBranchAddress("muminus_TRUEP_Y", &mm_py);
  132. ch->SetBranchAddress("muminus_TRUEP_Z", &mm_pz);
  133. unsigned int ncand, ntot;
  134. ch->SetBranchAddress("nCandidate", &ncand);
  135. ch->SetBranchAddress("totCandidates", &ntot);
  136. TH1D * h = new TH1D("q2_eff", "q2_eff; q^{2} [GeV^{2}]; events", 100, hmin, hmax);
  137. for(unsigned int e = 0; e < N; e++){
  138. if (e % (N/100) == 0)
  139. std::cout << std::setprecision(2) << std::fixed << "\r[" << round(double(e*100)/N) << "%]" << std::flush;
  140. ch->GetEntry(e);
  141. if(onlynTotOne){
  142. if(ntot != 1){
  143. continue;
  144. }
  145. }
  146. else if(onlyCandOne){
  147. if(ncand != 1){
  148. continue;
  149. }
  150. }
  151. h->Fill(((mm_pe+mp_pe)*(mm_pe+mp_pe)
  152. - (mm_px+mp_px)*(mm_px+mp_px)
  153. - (mm_py+mp_py)*(mm_py+mp_py)
  154. - (mm_pz+mp_pz)*(mm_pz+mp_pz))
  155. * 1.0e-6);
  156. }
  157. std::cout << std::endl;
  158. h->Scale(100.0 / h->Integral());
  159. TH1D * w = new TH1D("flatq2weights", "flatq2weights; q^{2} [GeV^{2}]; weights", 100, hmin, hmax);
  160. for(int i = 0; i < 100; i++)
  161. w->SetBinContent(i+1, 1.0);
  162. bool div_success = w->Divide(h);
  163. if(!div_success)
  164. std::cout << "Devition was not successfull!" << std::endl;
  165. TCanvas * cQ2 = new TCanvas("cQ2", "cQ2");
  166. cQ2->cd();
  167. h->GetYaxis()->SetRangeUser(0.5, 1.2);
  168. h->Draw();
  169. cQ2->SaveAs("q2efficiency.eps", "eps");
  170. w->GetYaxis()->SetRangeUser(0.0, 2.0);
  171. w->Draw();
  172. TLine * l = new TLine(hmin, 1.0, hmax, 1.0);
  173. l->SetLineColor(kGray+1);
  174. l->SetLineStyle(kDashed);
  175. l->Draw("SAME");
  176. cQ2->SaveAs("q2weights.eps", "eps");
  177. if(WriteToROOTfile){
  178. TFile * output = new TFile("/auto/data/dgerick/B2Kstmumu/FCNCfitter/config/phspweightq2.root", "UPDATE");
  179. output->cd();
  180. w->Write("",TObject::kWriteDelete);
  181. output->Close();
  182. }
  183. return;
  184. }