Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-12-16 09:17:31

0001 #include <TFile.h>
0002 #include <TH1D.h>
0003 #include <TH2D.h>
0004 
0005 void combine_1Dhist(string histname, TFile* f_MB, TFile* f_Jet5GeV, TFile* f_Jet10GeV, TFile* f_Jet15GeV, TFile* f_Jet20GeV, TFile* f_Jet30GeV, TFile* f_Jet50GeV, TFile* f_Jet70GeV, TFile* f_out);
0006 void combine_2Dhist(string histname, TFile* f_MB, TFile* f_Jet5GeV, TFile* f_Jet10GeV, TFile* f_Jet15GeV, TFile* f_Jet20GeV, TFile* f_Jet30GeV, TFile* f_Jet50GeV, TFile* f_Jet70GeV, TFile* f_out);
0007 
0008 const double mb_cross_section = 4.197e+10;
0009 const double jet5_cross_section = 1.369e+08;
0010 const double jet10_cross_section = 3.997e+06;
0011 const double jet15_cross_section = 4.073e+05;
0012 const double jet20_cross_section = 6.218e+04;
0013 const double jet30_cross_section = 2.502e+03;
0014 const double jet50_cross_section = 7.269;
0015 const double jet70_cross_section = 1.034e-02;
0016 
0017 void get_combinedoutput(int radius_index = 4) {
0018 
0019   TFile* f_MB = new TFile(Form("output_sim_hadd/output_sim_r0%d_MB.root", radius_index), "READ");
0020   TFile* f_Jet5GeV = new TFile(Form("output_sim_hadd/output_sim_r0%d_Jet5GeV.root", radius_index), "READ");
0021   TFile* f_Jet10GeV = new TFile(Form("output_sim_hadd/output_sim_r0%d_Jet10GeV.root", radius_index), "READ");
0022   TFile* f_Jet15GeV = new TFile(Form("output_sim_hadd/output_sim_r0%d_Jet15GeV.root", radius_index), "READ");
0023   TFile* f_Jet20GeV = new TFile(Form("output_sim_hadd/output_sim_r0%d_Jet20GeV.root", radius_index), "READ");
0024   TFile* f_Jet30GeV = new TFile(Form("output_sim_hadd/output_sim_r0%d_Jet30GeV.root", radius_index), "READ");
0025   TFile* f_Jet50GeV = new TFile(Form("output_sim_hadd/output_sim_r0%d_Jet50GeV.root", radius_index), "READ");
0026   TFile* f_Jet70GeV = new TFile(Form("output_sim_hadd/output_sim_r0%d_Jet70GeV.root", radius_index), "READ");
0027   if (!f_MB || !f_Jet5GeV || !f_Jet10GeV || !f_Jet15GeV || !f_Jet20GeV || !f_Jet30GeV || !f_Jet50GeV || !f_Jet70GeV) {
0028     std::cout << "Error: cannot open one or more input files." << std::endl;
0029     return;
0030   }
0031   TFile* f_combined = new TFile(Form("output_sim_r0%d.root", radius_index), "RECREATE");
0032 
0033   combine_1Dhist("h_recojet_pt_record", f_MB, f_Jet5GeV, f_Jet10GeV, f_Jet15GeV, f_Jet20GeV, f_Jet30GeV, f_Jet50GeV, f_Jet70GeV, f_combined);
0034   combine_1Dhist("h_recojet_pt_passcut_record", f_MB, f_Jet5GeV, f_Jet10GeV, f_Jet15GeV, f_Jet20GeV, f_Jet30GeV, f_Jet50GeV, f_Jet70GeV, f_combined);
0035   combine_1Dhist("h_truthjet_pt_record", f_MB, f_Jet5GeV, f_Jet10GeV, f_Jet15GeV, f_Jet20GeV, f_Jet30GeV, f_Jet50GeV, f_Jet70GeV, f_combined);
0036   combine_1Dhist("h_truthjet_pt_passcut_record", f_MB, f_Jet5GeV, f_Jet10GeV, f_Jet15GeV, f_Jet20GeV, f_Jet30GeV, f_Jet50GeV, f_Jet70GeV, f_combined);
0037 
0038   combine_1Dhist("h_truth_nominal", f_MB, f_Jet5GeV, f_Jet10GeV, f_Jet15GeV, f_Jet20GeV, f_Jet30GeV, f_Jet50GeV, f_Jet70GeV, f_combined);
0039   combine_1Dhist("h_measure_nominal", f_MB, f_Jet5GeV, f_Jet10GeV, f_Jet15GeV, f_Jet20GeV, f_Jet30GeV, f_Jet50GeV, f_Jet70GeV, f_combined);
0040   combine_2Dhist("h_respmatrix_nominal", f_MB, f_Jet5GeV, f_Jet10GeV, f_Jet15GeV, f_Jet20GeV, f_Jet30GeV, f_Jet50GeV, f_Jet70GeV, f_combined);
0041   combine_1Dhist("h_fake_nominal", f_MB, f_Jet5GeV, f_Jet10GeV, f_Jet15GeV, f_Jet20GeV, f_Jet30GeV, f_Jet50GeV, f_Jet70GeV, f_combined);
0042   combine_1Dhist("h_miss_nominal", f_MB, f_Jet5GeV, f_Jet10GeV, f_Jet15GeV, f_Jet20GeV, f_Jet30GeV, f_Jet50GeV, f_Jet70GeV, f_combined);
0043 
0044   combine_1Dhist("h_truth_jesup", f_MB, f_Jet5GeV, f_Jet10GeV, f_Jet15GeV, f_Jet20GeV, f_Jet30GeV, f_Jet50GeV, f_Jet70GeV, f_combined);
0045   combine_1Dhist("h_measure_jesup", f_MB, f_Jet5GeV, f_Jet10GeV, f_Jet15GeV, f_Jet20GeV, f_Jet30GeV, f_Jet50GeV, f_Jet70GeV, f_combined);
0046   combine_2Dhist("h_respmatrix_jesup", f_MB, f_Jet5GeV, f_Jet10GeV, f_Jet15GeV, f_Jet20GeV, f_Jet30GeV, f_Jet50GeV, f_Jet70GeV, f_combined);
0047   combine_1Dhist("h_fake_jesup", f_MB, f_Jet5GeV, f_Jet10GeV, f_Jet15GeV, f_Jet20GeV, f_Jet30GeV, f_Jet50GeV, f_Jet70GeV, f_combined);
0048   combine_1Dhist("h_miss_jesup", f_MB, f_Jet5GeV, f_Jet10GeV, f_Jet15GeV, f_Jet20GeV, f_Jet30GeV, f_Jet50GeV, f_Jet70GeV, f_combined);
0049 
0050   combine_1Dhist("h_truth_jesdown", f_MB, f_Jet5GeV, f_Jet10GeV, f_Jet15GeV, f_Jet20GeV, f_Jet30GeV, f_Jet50GeV, f_Jet70GeV, f_combined);
0051   combine_1Dhist("h_measure_jesdown", f_MB, f_Jet5GeV, f_Jet10GeV, f_Jet15GeV, f_Jet20GeV, f_Jet30GeV, f_Jet50GeV, f_Jet70GeV, f_combined);
0052   combine_2Dhist("h_respmatrix_jesdown", f_MB, f_Jet5GeV, f_Jet10GeV, f_Jet15GeV, f_Jet20GeV, f_Jet30GeV, f_Jet50GeV, f_Jet70GeV, f_combined);
0053   combine_1Dhist("h_fake_jesdown", f_MB, f_Jet5GeV, f_Jet10GeV, f_Jet15GeV, f_Jet20GeV, f_Jet30GeV, f_Jet50GeV, f_Jet70GeV, f_combined);
0054   combine_1Dhist("h_miss_jesdown", f_MB, f_Jet5GeV, f_Jet10GeV, f_Jet15GeV, f_Jet20GeV, f_Jet30GeV, f_Jet50GeV, f_Jet70GeV, f_combined);
0055 
0056   combine_1Dhist("h_truth_jerup", f_MB, f_Jet5GeV, f_Jet10GeV, f_Jet15GeV, f_Jet20GeV, f_Jet30GeV, f_Jet50GeV, f_Jet70GeV, f_combined);
0057   combine_1Dhist("h_measure_jerup", f_MB, f_Jet5GeV, f_Jet10GeV, f_Jet15GeV, f_Jet20GeV, f_Jet30GeV, f_Jet50GeV, f_Jet70GeV, f_combined);
0058   combine_2Dhist("h_respmatrix_jerup", f_MB, f_Jet5GeV, f_Jet10GeV, f_Jet15GeV, f_Jet20GeV, f_Jet30GeV, f_Jet50GeV, f_Jet70GeV, f_combined);
0059   combine_1Dhist("h_fake_jerup", f_MB, f_Jet5GeV, f_Jet10GeV, f_Jet15GeV, f_Jet20GeV, f_Jet30GeV, f_Jet50GeV, f_Jet70GeV, f_combined);
0060   combine_1Dhist("h_miss_jerup", f_MB, f_Jet5GeV, f_Jet10GeV, f_Jet15GeV, f_Jet20GeV, f_Jet30GeV, f_Jet50GeV, f_Jet70GeV, f_combined);
0061 
0062   combine_1Dhist("h_truth_jerdown", f_MB, f_Jet5GeV, f_Jet10GeV, f_Jet15GeV, f_Jet20GeV, f_Jet30GeV, f_Jet50GeV, f_Jet70GeV, f_combined);
0063   combine_1Dhist("h_measure_jerdown", f_MB, f_Jet5GeV, f_Jet10GeV, f_Jet15GeV, f_Jet20GeV, f_Jet30GeV, f_Jet50GeV, f_Jet70GeV, f_combined);
0064   combine_2Dhist("h_respmatrix_jerdown", f_MB, f_Jet5GeV, f_Jet10GeV, f_Jet15GeV, f_Jet20GeV, f_Jet30GeV, f_Jet50GeV, f_Jet70GeV, f_combined);
0065   combine_1Dhist("h_fake_jerdown", f_MB, f_Jet5GeV, f_Jet10GeV, f_Jet15GeV, f_Jet20GeV, f_Jet30GeV, f_Jet50GeV, f_Jet70GeV, f_combined);
0066   combine_1Dhist("h_miss_jerdown", f_MB, f_Jet5GeV, f_Jet10GeV, f_Jet15GeV, f_Jet20GeV, f_Jet30GeV, f_Jet50GeV, f_Jet70GeV, f_combined);
0067 
0068   combine_1Dhist("h_fullclosure_measure", f_MB, f_Jet5GeV, f_Jet10GeV, f_Jet15GeV, f_Jet20GeV, f_Jet30GeV, f_Jet50GeV, f_Jet70GeV, f_combined);
0069   combine_1Dhist("h_fullclosure_truth", f_MB, f_Jet5GeV, f_Jet10GeV, f_Jet15GeV, f_Jet20GeV, f_Jet30GeV, f_Jet50GeV, f_Jet70GeV, f_combined);
0070   combine_2Dhist("h_fullclosure_respmatrix", f_MB, f_Jet5GeV, f_Jet10GeV, f_Jet15GeV, f_Jet20GeV, f_Jet30GeV, f_Jet50GeV, f_Jet70GeV, f_combined);
0071   combine_1Dhist("h_fullclosure_fake", f_MB, f_Jet5GeV, f_Jet10GeV, f_Jet15GeV, f_Jet20GeV, f_Jet30GeV, f_Jet50GeV, f_Jet70GeV, f_combined);
0072   combine_1Dhist("h_fullclosure_miss", f_MB, f_Jet5GeV, f_Jet10GeV, f_Jet15GeV, f_Jet20GeV, f_Jet30GeV, f_Jet50GeV, f_Jet70GeV, f_combined);
0073 
0074   combine_1Dhist("h_halfclosure_measure", f_MB, f_Jet5GeV, f_Jet10GeV, f_Jet15GeV, f_Jet20GeV, f_Jet30GeV, f_Jet50GeV, f_Jet70GeV, f_combined);
0075   combine_1Dhist("h_halfclosure_truth", f_MB, f_Jet5GeV, f_Jet10GeV, f_Jet15GeV, f_Jet20GeV, f_Jet30GeV, f_Jet50GeV, f_Jet70GeV, f_combined);
0076   combine_2Dhist("h_halfclosure_respmatrix", f_MB, f_Jet5GeV, f_Jet10GeV, f_Jet15GeV, f_Jet20GeV, f_Jet30GeV, f_Jet50GeV, f_Jet70GeV, f_combined);
0077   combine_1Dhist("h_halfclosure_fake", f_MB, f_Jet5GeV, f_Jet10GeV, f_Jet15GeV, f_Jet20GeV, f_Jet30GeV, f_Jet50GeV, f_Jet70GeV, f_combined);
0078   combine_1Dhist("h_halfclosure_miss", f_MB, f_Jet5GeV, f_Jet10GeV, f_Jet15GeV, f_Jet20GeV, f_Jet30GeV, f_Jet50GeV, f_Jet70GeV, f_combined);
0079 
0080   combine_1Dhist("h_halfclosure_measure_check", f_MB, f_Jet5GeV, f_Jet10GeV, f_Jet15GeV, f_Jet20GeV, f_Jet30GeV, f_Jet50GeV, f_Jet70GeV, f_combined);
0081   combine_1Dhist("h_halfclosure_truth_check", f_MB, f_Jet5GeV, f_Jet10GeV, f_Jet15GeV, f_Jet20GeV, f_Jet30GeV, f_Jet50GeV, f_Jet70GeV, f_combined);
0082   combine_2Dhist("h_halfclosure_respmatrix_check", f_MB, f_Jet5GeV, f_Jet10GeV, f_Jet15GeV, f_Jet20GeV, f_Jet30GeV, f_Jet50GeV, f_Jet70GeV, f_combined);
0083   combine_1Dhist("h_halfclosure_fake_check", f_MB, f_Jet5GeV, f_Jet10GeV, f_Jet15GeV, f_Jet20GeV, f_Jet30GeV, f_Jet50GeV, f_Jet70GeV, f_combined);
0084   combine_1Dhist("h_halfclosure_miss_check", f_MB, f_Jet5GeV, f_Jet10GeV, f_Jet15GeV, f_Jet20GeV, f_Jet30GeV, f_Jet50GeV, f_Jet70GeV, f_combined);
0085 
0086   f_combined->Close();
0087   std::cout << "Combined output file created: " << Form("output_sim_r0%d.root", radius_index) << std::endl;
0088 }
0089 
0090 void combine_1Dhist(string histname, TFile* f_MB, TFile* f_Jet5GeV, TFile* f_Jet10GeV, TFile* f_Jet15GeV, TFile* f_Jet20GeV, TFile* f_Jet30GeV, TFile* f_Jet50GeV, TFile* f_Jet70GeV, TFile* f_out) {
0091   TH1D *h_MB_event = (TH1D*)f_MB->Get("h_event"); int mb_event = h_MB_event->GetBinContent(1); double mb_scale = mb_cross_section / (double)mb_event;
0092   TH1D *h_Jet5GeV_event = (TH1D*)f_Jet5GeV->Get("h_event"); int jet5_event = h_Jet5GeV_event->GetBinContent(1); double jet5_scale = jet5_cross_section / (double)jet5_event;
0093   TH1D *h_Jet10GeV_event = (TH1D*)f_Jet10GeV->Get("h_event"); int jet10_event = h_Jet10GeV_event->GetBinContent(1); double jet10_scale = jet10_cross_section / (double)jet10_event;
0094   TH1D *h_Jet15GeV_event = (TH1D*)f_Jet15GeV->Get("h_event"); int jet15_event = h_Jet15GeV_event->GetBinContent(1); double jet15_scale = jet15_cross_section / (double)jet15_event;
0095   TH1D *h_Jet20GeV_event = (TH1D*)f_Jet20GeV->Get("h_event"); int jet20_event = h_Jet20GeV_event->GetBinContent(1); double jet20_scale = jet20_cross_section / (double)jet20_event;
0096   TH1D *h_Jet30GeV_event = (TH1D*)f_Jet30GeV->Get("h_event"); int jet30_event = h_Jet30GeV_event->GetBinContent(1); double jet30_scale = jet30_cross_section / (double)jet30_event;
0097   TH1D *h_Jet50GeV_event = (TH1D*)f_Jet50GeV->Get("h_event"); int jet50_event = h_Jet50GeV_event->GetBinContent(1); double jet50_scale = jet50_cross_section / (double)jet50_event;
0098   TH1D *h_Jet70GeV_event = (TH1D*)f_Jet70GeV->Get("h_event"); int jet70_event = h_Jet70GeV_event->GetBinContent(1); double jet70_scale = jet70_cross_section / (double)jet70_event;
0099 
0100   TH1D* h_MB_forcombine = (TH1D*)f_MB->Get(histname.c_str());
0101   TH1D* h_Jet5GeV_forcombine = (TH1D*)f_Jet5GeV->Get(histname.c_str());
0102   TH1D* h_Jet10GeV_forcombine = (TH1D*)f_Jet10GeV->Get(histname.c_str());
0103   TH1D* h_Jet15GeV_forcombine = (TH1D*)f_Jet15GeV->Get(histname.c_str());
0104   TH1D* h_Jet20GeV_forcombine = (TH1D*)f_Jet20GeV->Get(histname.c_str());
0105   TH1D* h_Jet30GeV_forcombine = (TH1D*)f_Jet30GeV->Get(histname.c_str());
0106   TH1D* h_Jet50GeV_forcombine = (TH1D*)f_Jet50GeV->Get(histname.c_str());
0107   TH1D* h_Jet70GeV_forcombine = (TH1D*)f_Jet70GeV->Get(histname.c_str());
0108 
0109   TH1D* h_combined = (TH1D*)h_MB_forcombine->Clone(histname.c_str());
0110   h_combined->Scale(mb_scale);
0111   h_combined->Add(h_Jet5GeV_forcombine, jet5_scale);
0112   h_combined->Add(h_Jet10GeV_forcombine, jet10_scale);
0113   h_combined->Add(h_Jet15GeV_forcombine, jet15_scale);
0114   h_combined->Add(h_Jet20GeV_forcombine, jet20_scale);
0115   h_combined->Add(h_Jet30GeV_forcombine, jet30_scale);
0116   h_combined->Add(h_Jet50GeV_forcombine, jet50_scale);
0117   h_combined->Add(h_Jet70GeV_forcombine, jet70_scale);
0118 
0119   f_out->cd();
0120   h_combined->Write();
0121 }
0122 
0123 void combine_2Dhist(string histname, TFile* f_MB, TFile* f_Jet5GeV, TFile* f_Jet10GeV, TFile* f_Jet15GeV, TFile* f_Jet20GeV, TFile* f_Jet30GeV, TFile* f_Jet50GeV, TFile* f_Jet70GeV, TFile* f_out) {
0124   TH1D *h_MB_event = (TH1D*)f_MB->Get("h_event"); int mb_event = h_MB_event->GetBinContent(1); double mb_scale = mb_cross_section / (double)mb_event;
0125   TH1D *h_Jet5GeV_event = (TH1D*)f_Jet5GeV->Get("h_event"); int jet5_event = h_Jet5GeV_event->GetBinContent(1); double jet5_scale = jet5_cross_section / (double)jet5_event;
0126   TH1D *h_Jet10GeV_event = (TH1D*)f_Jet10GeV->Get("h_event"); int jet10_event = h_Jet10GeV_event->GetBinContent(1); double jet10_scale = jet10_cross_section / (double)jet10_event;
0127   TH1D *h_Jet15GeV_event = (TH1D*)f_Jet15GeV->Get("h_event"); int jet15_event = h_Jet15GeV_event->GetBinContent(1); double jet15_scale = jet15_cross_section / (double)jet15_event;
0128   TH1D *h_Jet20GeV_event = (TH1D*)f_Jet20GeV->Get("h_event"); int jet20_event = h_Jet20GeV_event->GetBinContent(1); double jet20_scale = jet20_cross_section / (double)jet20_event;
0129   TH1D *h_Jet30GeV_event = (TH1D*)f_Jet30GeV->Get("h_event"); int jet30_event = h_Jet30GeV_event->GetBinContent(1); double jet30_scale = jet30_cross_section / (double)jet30_event;
0130   TH1D *h_Jet50GeV_event = (TH1D*)f_Jet50GeV->Get("h_event"); int jet50_event = h_Jet50GeV_event->GetBinContent(1); double jet50_scale = jet50_cross_section / (double)jet50_event;
0131   TH1D *h_Jet70GeV_event = (TH1D*)f_Jet70GeV->Get("h_event"); int jet70_event = h_Jet70GeV_event->GetBinContent(1); double jet70_scale = jet70_cross_section / (double)jet70_event;
0132 
0133   TH2D* h_MB_forcombine = (TH2D*)f_MB->Get(histname.c_str());
0134   TH2D* h_Jet5GeV_forcombine = (TH2D*)f_Jet5GeV->Get(histname.c_str());
0135   TH2D* h_Jet10GeV_forcombine = (TH2D*)f_Jet10GeV->Get(histname.c_str());
0136   TH2D* h_Jet15GeV_forcombine = (TH2D*)f_Jet15GeV->Get(histname.c_str());
0137   TH2D* h_Jet20GeV_forcombine = (TH2D*)f_Jet20GeV->Get(histname.c_str());
0138   TH2D* h_Jet30GeV_forcombine = (TH2D*)f_Jet30GeV->Get(histname.c_str());
0139   TH2D* h_Jet50GeV_forcombine = (TH2D*)f_Jet50GeV->Get(histname.c_str());
0140   TH2D* h_Jet70GeV_forcombine = (TH2D*)f_Jet70GeV->Get(histname.c_str());
0141 
0142   TH2D* h_combined = (TH2D*)h_MB_forcombine->Clone(histname.c_str());
0143   h_combined->Scale(mb_scale);
0144   h_combined->Add(h_Jet5GeV_forcombine, jet5_scale);
0145   h_combined->Add(h_Jet10GeV_forcombine, jet10_scale);
0146   h_combined->Add(h_Jet15GeV_forcombine, jet15_scale);
0147   h_combined->Add(h_Jet20GeV_forcombine, jet20_scale);
0148   h_combined->Add(h_Jet30GeV_forcombine, jet30_scale);
0149   h_combined->Add(h_Jet50GeV_forcombine, jet50_scale);
0150   h_combined->Add(h_Jet70GeV_forcombine, jet70_scale);
0151 
0152   f_out->cd();
0153   h_combined->Write();
0154 }