File indexing completed on 2025-12-16 09:16:13
0001 #include "TFile.h"
0002 #include "TH1F.h"
0003 #include "TH2F.h"
0004
0005 #include "../util/binning.h"
0006
0007 void get_phiyield(TH1F* h_yield, TH1D* h_mass, float xval)
0008 {
0009 TF1* phi_massfit = new TF1("phi_fit","[0]*TMath::Voigt(x-[1],[2],[3])+pol1(4)",1.,1.1);
0010 phi_massfit->SetParameter(0,1000.);
0011 phi_massfit->SetParameter(1,1.02);
0012 phi_massfit->SetParameter(2,0.01);
0013 phi_massfit->SetParameter(3,0.01);
0014 phi_massfit->SetParameter(4,1000.);
0015 phi_massfit->SetParameter(5,-0.);
0016
0017 h_mass->Fit(phi_massfit,"R");
0018
0019 TF1* phi_signal = new TF1("phi_signal","[0]*TMath::Voigt(x-[1],[2],[3])",1.,1.1);
0020 phi_signal->SetParameter(0,phi_massfit->GetParameter(0));
0021 phi_signal->SetParError(0,phi_massfit->GetParError(0));
0022 phi_signal->SetParameter(1,phi_massfit->GetParameter(1));
0023 phi_signal->SetParError(1,phi_massfit->GetParError(1));
0024 phi_signal->SetParameter(2,phi_massfit->GetParameter(2));
0025 phi_signal->SetParError(2,phi_massfit->GetParError(2));
0026 phi_signal->SetParameter(3,phi_massfit->GetParameter(3));
0027 phi_signal->SetParError(3,phi_massfit->GetParError(3));
0028
0029 float phi_yield = phi_signal->Integral(0.95,1.1)/h_mass->GetBinWidth(2);
0030 float phi_yielderror = phi_signal->IntegralError(0.95,1.1)/h_mass->GetBinWidth(2);
0031
0032 h_yield->SetBinContent(h_yield->FindBin(xval),phi_yield);
0033 h_yield->SetBinError(h_yield->FindBin(xval),phi_yielderror);
0034 delete phi_massfit;
0035 delete phi_signal;
0036 }
0037
0038 void get_Ksyield(TH1F* h_yield, TH1D* h_mass, float xval)
0039 {
0040 TF1* ks_massfit = new TF1("ks_fit","gaus(0)+pol1(3)",0.4,0.6);
0041 ks_massfit->SetParameter(0,1000.);
0042 ks_massfit->SetParameter(1,0.5);
0043 ks_massfit->SetParameter(2,0.1);
0044 ks_massfit->SetParameter(3,100.);
0045 ks_massfit->SetParameter(4,-0.1);
0046
0047 h_mass->Fit(ks_massfit,"R");
0048
0049 TF1* ks_signal = new TF1("ks_signal","gaus(0)",0.4,0.6);
0050 ks_signal->SetParameter(0,ks_massfit->GetParameter(0));
0051 ks_signal->SetParError(0,ks_massfit->GetParError(0));
0052 ks_signal->SetParameter(1,ks_massfit->GetParameter(1));
0053 ks_signal->SetParError(1,ks_massfit->GetParError(1));
0054 ks_signal->SetParameter(2,ks_massfit->GetParameter(2));
0055 ks_signal->SetParError(2,ks_massfit->GetParError(2));
0056
0057 float ks_yield = ks_signal->Integral(0.4,0.6)/h_mass->GetBinWidth(2);
0058 float ks_yielderror = ks_signal->IntegralError(0.4,0.6)/h_mass->GetBinWidth(2);
0059
0060 h_yield->SetBinContent(h_yield->FindBin(xval),ks_yield);
0061 h_yield->SetBinError(h_yield->FindBin(xval),ks_yielderror);
0062 delete ks_massfit;
0063 delete ks_signal;
0064 }
0065
0066 void get_lambdayield(TH1F* h_yield, TH1D* h_mass, float xval)
0067 {
0068 TF1* lambda_massfit = new TF1("lambda_fit","gaus(0)+pol2(3)",1.08,1.2);
0069 lambda_massfit->SetParameter(0,1000.);
0070 lambda_massfit->SetParameter(1,1.11);
0071 lambda_massfit->SetParameter(2,0.1);
0072 lambda_massfit->SetParameter(3,100.);
0073 lambda_massfit->SetParameter(4,0.1);
0074 lambda_massfit->SetParameter(5,-0.1);
0075
0076 h_mass->Fit(lambda_massfit,"R");
0077
0078 TF1* lambda_signal = new TF1("lambda_signal","gaus(0)",0.4,0.6);
0079 lambda_signal->SetParameter(0,lambda_massfit->GetParameter(0));
0080 lambda_signal->SetParError(0,lambda_massfit->GetParError(0));
0081 lambda_signal->SetParameter(1,lambda_massfit->GetParameter(1));
0082 lambda_signal->SetParError(1,lambda_massfit->GetParError(1));
0083 lambda_signal->SetParameter(2,lambda_massfit->GetParameter(2));
0084 lambda_signal->SetParError(2,lambda_massfit->GetParError(2));
0085
0086 float lambda_yield = lambda_signal->Integral(1.05,1.2)/h_mass->GetBinWidth(2);
0087 float lambda_yielderror = lambda_signal->IntegralError(1.05,1.2)/h_mass->GetBinWidth(2);
0088
0089 h_yield->SetBinContent(h_yield->FindBin(xval),lambda_yield);
0090 h_yield->SetBinError(h_yield->FindBin(xval),lambda_yielderror);
0091 delete lambda_massfit;
0092 delete lambda_signal;
0093 }
0094
0095 void calculate_ratios()
0096 {
0097 gStyle->SetOptStat(1);
0098
0099 TFile* Ks_file = TFile::Open("../mass_histograms/data_Kshort_mass.root");
0100 TFile* phi_file = TFile::Open("../mass_histograms/data_phi_mass.root");
0101 TFile* lambda_file = TFile::Open("../mass_histograms/data_lambda_mass.root");
0102
0103 TH1D* Ks_mass = (TH1D*)Ks_file->Get("Kshort_mass");
0104 TH2F* Ks_massvspt = (TH2F*)Ks_file->Get("Kshort_mass_vspt");
0105 TH2F* Ks_massvsntrk = (TH2F*)Ks_file->Get("Kshort_mass_vsntrk");
0106 TH2F* Ks_massvsy = (TH2F*)Ks_file->Get("Kshort_mass_vsy");
0107 TH2F* Ks_massvsphi = (TH2F*)Ks_file->Get("Kshort_mass_vsphi");
0108
0109 TH1D* phi_mass = (TH1D*)phi_file->Get("phi_mass");
0110 TH2F* phi_massvspt = (TH2F*)phi_file->Get("phi_mass_vspt");
0111 TH2F* phi_massvsntrk = (TH2F*)phi_file->Get("phi_mass_vsntrk");
0112 TH2F* phi_massvsy = (TH2F*)phi_file->Get("phi_mass_vsy");
0113 TH2F* phi_massvsphi = (TH2F*)phi_file->Get("phi_mass_vsphi");
0114
0115 TH1D* lambda_mass = (TH1D*)lambda_file->Get("lambda_mass");
0116 TH2F* lambda_massvspt = (TH2F*)lambda_file->Get("lambda_mass_vspt");
0117 TH2F* lambda_massvsntrk = (TH2F*)lambda_file->Get("lambda_mass_vsntrk");
0118 TH2F* lambda_massvsy = (TH2F*)lambda_file->Get("lambda_mass_vsy");
0119 TH2F* lambda_massvsphi = (TH2F*)lambda_file->Get("lambda_mass_vsphi");
0120
0121 TH1F* integrated_phi_yield = new TH1F("integrated_phi_yield","Integrated #{phi} yield",1,0.,1.);
0122 TH1F* integrated_Ks_yield = new TH1F("integrated_Ks_yield","Integrated K_{s}^{0} yield",1,0.,1.);
0123 TH1F* integrated_lambda_yield = new TH1F("integrated_lambda_yield","Integrated #{Lambda} yield",1,0.,1.);
0124
0125 TH1F* phi_yield_vspt = makeHistogram("phi_yield","#{phi} yield",BinInfo::final_pt_bins);
0126 TH1F* Ks_yield_vspt = makeHistogram("Ks_yield","K_{s}^{0} yield",BinInfo::final_pt_bins);
0127 TH1F* lambda_yield_vspt = makeHistogram("lambda_yield","#{Lambda} yield",BinInfo::final_pt_bins);
0128
0129 TH1F* phi_yield_vsy = makeHistogram("phi_yield","#{phi} yield",BinInfo::final_rapidity_bins);
0130 TH1F* Ks_yield_vsy = makeHistogram("Ks_yield","K_{s}^{0} yield",BinInfo::final_rapidity_bins);
0131 TH1F* lambda_yield_vsy = makeHistogram("lambda_yield","#{Lambda} yield",BinInfo::final_rapidity_bins);
0132
0133 TH1F* phi_yield_vsphi = makeHistogram("phi_yield","#{phi} yield",BinInfo::final_phi_bins);
0134 TH1F* Ks_yield_vsphi = makeHistogram("Ks_yield","K_{s}^{0} yield",BinInfo::final_phi_bins);
0135 TH1F* lambda_yield_vsphi = makeHistogram("lambda_yield","#{Lambda} yield",BinInfo::final_phi_bins);
0136
0137 TH1F* phi_yield_vsntrack = makeHistogram("phi_yield","#{phi} yield",BinInfo::final_ntrack_bins);
0138 TH1F* Ks_yield_vsntrack = makeHistogram("Ks_yield","K_{s}^{0} yield",BinInfo::final_ntrack_bins);
0139 TH1F* lambda_yield_vsntrack = makeHistogram("lambda_yield","#{Lambda} yield",BinInfo::final_ntrack_bins);
0140
0141 get_phiyield(integrated_phi_yield,phi_mass,0.5);
0142 get_Ksyield(integrated_Ks_yield,Ks_mass,0.5);
0143 get_lambdayield(integrated_lambda_yield,lambda_mass,0.5);
0144
0145 for(int i=1;i<=BinInfo::pt_bins.nBins;i++)
0146 {
0147 TH1D* phi_mass_slice = phi_massvspt->ProjectionX("phi_slice",i,i+1,"e");
0148 TH1D* Ks_mass_slice = Ks_massvspt->ProjectionX("Ks_slice",i,i+1,"e");
0149 TH1D* lambda_mass_slice = lambda_massvspt->ProjectionX("lambda_slice",i,i+1,"e");
0150
0151 get_phiyield(phi_yield_vspt,phi_mass_slice,phi_yield_vspt->GetBinCenter(i));
0152 get_Ksyield(Ks_yield_vspt,Ks_mass_slice,Ks_yield_vspt->GetBinCenter(i));
0153 get_lambdayield(lambda_yield_vspt,lambda_mass_slice,lambda_yield_vspt->GetBinCenter(i));
0154 delete phi_mass_slice;
0155 delete Ks_mass_slice;
0156 delete lambda_mass_slice;
0157 }
0158
0159 for(int i=1;i<=BinInfo::rapidity_bins.nBins;i++)
0160 {
0161 TH1D* phi_mass_slice = phi_massvsy->ProjectionX("phi_slice",i,i+1,"e");
0162 TH1D* Ks_mass_slice = Ks_massvsy->ProjectionX("Ks_slice",i,i+1,"e");
0163 TH1D* lambda_mass_slice = lambda_massvsy->ProjectionX("lambda_slice",i,i+1,"e");
0164
0165 get_phiyield(phi_yield_vsy,phi_mass_slice,phi_yield_vsy->GetBinCenter(i));
0166 get_Ksyield(Ks_yield_vsy,Ks_mass_slice,Ks_yield_vsy->GetBinCenter(i));
0167 get_lambdayield(lambda_yield_vsy,lambda_mass_slice,lambda_yield_vsy->GetBinCenter(i));
0168 delete phi_mass_slice;
0169 delete Ks_mass_slice;
0170 delete lambda_mass_slice;
0171 }
0172
0173 for(int i=1;i<=BinInfo::phi_bins.nBins;i++)
0174 {
0175 TH1D* phi_mass_slice = phi_massvsphi->ProjectionX("phi_slice",i,i+1,"e");
0176 TH1D* Ks_mass_slice = Ks_massvsphi->ProjectionX("Ks_slice",i,i+1,"e");
0177 TH1D* lambda_mass_slice = lambda_massvsphi->ProjectionX("lambda_slice",i,i+1,"e");
0178
0179 get_phiyield(phi_yield_vsphi,phi_mass_slice,phi_yield_vsphi->GetBinCenter(i));
0180 get_Ksyield(Ks_yield_vsphi,Ks_mass_slice,Ks_yield_vsphi->GetBinCenter(i));
0181 get_lambdayield(lambda_yield_vsphi,lambda_mass_slice,lambda_yield_vsphi->GetBinCenter(i));
0182 delete phi_mass_slice;
0183 delete Ks_mass_slice;
0184 delete lambda_mass_slice;
0185 }
0186
0187 for(int i=1;i<=BinInfo::ntrack_bins.nBins;i++)
0188 {
0189 TH1D* phi_mass_slice = phi_massvsntrk->ProjectionX("phi_slice",i,i+1,"e");
0190 TH1D* Ks_mass_slice = Ks_massvsntrk->ProjectionX("Ks_slice",i,i+1,"e");
0191 TH1D* lambda_mass_slice = lambda_massvsntrk->ProjectionX("lambda_slice",i,i+1,"e");
0192
0193 get_phiyield(phi_yield_vsntrack,phi_mass_slice,phi_yield_vsntrack->GetBinCenter(i));
0194 get_Ksyield(Ks_yield_vsntrack,Ks_mass_slice,Ks_yield_vsntrack->GetBinCenter(i));
0195 get_lambdayield(lambda_yield_vsntrack,lambda_mass_slice,lambda_yield_vsntrack->GetBinCenter(i));
0196 delete phi_mass_slice;
0197 delete Ks_mass_slice;
0198 delete lambda_mass_slice;
0199 }
0200
0201 TH1F* lambdaKsratio = new TH1F("integrated_lambdaKs_ratio","Integrated #{Lambda}/K_{s}^{0} Ratio",1,0.,1.);
0202 TH1F* phiKsratio = new TH1F("integrated_phiKs_ratio","Integrated #{phi}/K_{s}^{0} Ratio",1,0.,1.);
0203
0204 lambdaKsratio->Divide(integrated_lambda_yield,integrated_Ks_yield);
0205 phiKsratio->Divide(integrated_phi_yield,integrated_Ks_yield);
0206
0207 TH1F* lambdaKsratio_vspt = makeHistogram("lambdaKsratio","#{Lambda}/K_{s}^{0} Ratio",BinInfo::final_pt_bins);
0208 TH1F* phiKsratio_vspt = makeHistogram("phiKsratio","#{phi}/K_{s}^{0} Ratio",BinInfo::final_pt_bins);
0209
0210 lambdaKsratio_vspt->Divide(lambda_yield_vspt,Ks_yield_vspt);
0211 phiKsratio_vspt->Divide(phi_yield_vspt,Ks_yield_vspt);
0212
0213 TH1F* lambdaKsratio_vsy = makeHistogram("lambdaKsratio","#{Lambda}/K_{s}^{0} Ratio",BinInfo::final_rapidity_bins);
0214 TH1F* phiKsratio_vsy = makeHistogram("phiKsratio","#{phi}/K_{s}^{0} Ratio",BinInfo::final_rapidity_bins);
0215
0216 lambdaKsratio_vsy->Divide(lambda_yield_vsy,Ks_yield_vsy);
0217 phiKsratio_vsy->Divide(phi_yield_vsy,Ks_yield_vsy);
0218
0219 TH1F* lambdaKsratio_vsphi = makeHistogram("lambdaKsratio","#{Lambda}/K_{s}^{0} Ratio",BinInfo::final_phi_bins);
0220 TH1F* phiKsratio_vsphi = makeHistogram("phiKsratio","#{phi}/K_{s}^{0} Ratio",BinInfo::final_phi_bins);
0221
0222 lambdaKsratio_vsphi->Divide(lambda_yield_vsphi,Ks_yield_vsphi);
0223 phiKsratio_vsphi->Divide(phi_yield_vsphi,Ks_yield_vsphi);
0224
0225 TH1F* lambdaKsratio_vsntrack = makeHistogram("lambdaKsratio","#{Lambda}/K_{s}^{0} Ratio",BinInfo::final_ntrack_bins);
0226 TH1F* phiKsratio_vsntrack = makeHistogram("phiKsratio","#{phi}/K_{s}^{0} Ratio",BinInfo::final_ntrack_bins);
0227
0228 lambdaKsratio_vsntrack->Divide(lambda_yield_vsntrack,Ks_yield_vsntrack);
0229 phiKsratio_vsntrack->Divide(phi_yield_vsntrack,Ks_yield_vsntrack);
0230
0231 TFile* fout = new TFile("fits.root","RECREATE");
0232 Ks_mass->Write();
0233 phi_mass->Write();
0234 lambda_mass->Write();
0235
0236 integrated_phi_yield->Write();
0237 integrated_Ks_yield->Write();
0238 integrated_lambda_yield->Write();
0239
0240 phi_yield_vspt->Write();
0241 Ks_yield_vspt->Write();
0242 lambda_yield_vspt->Write();
0243
0244 phi_yield_vsy->Write();
0245 Ks_yield_vsy->Write();
0246 lambda_yield_vsy->Write();
0247
0248 phi_yield_vsphi->Write();
0249 Ks_yield_vsphi->Write();
0250 lambda_yield_vsphi->Write();
0251
0252 phi_yield_vsntrack->Write();
0253 Ks_yield_vsntrack->Write();
0254 lambda_yield_vsntrack->Write();
0255
0256 lambdaKsratio->Write();
0257 phiKsratio->Write();
0258
0259 lambdaKsratio_vspt->Write();
0260 phiKsratio_vspt->Write();
0261
0262 lambdaKsratio_vsy->Write();
0263 phiKsratio_vsy->Write();
0264
0265 lambdaKsratio_vsphi->Write();
0266 phiKsratio_vsphi->Write();
0267
0268 lambdaKsratio_vsntrack->Write();
0269 phiKsratio_vsntrack->Write();
0270 }