Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-08-06 08:12:25

0001 #include "sPhenixStyle.C"
0002 #include "../ana_map_folder/ana_map_v1.h"
0003 namespace ana_map_version = ANA_MAP_V3;
0004 
0005 #include "../ReadINTTZ/ReadINTTZ.C"
0006 #include "../ReadINTTZ/ReadINTTZCombine.C"
0007 
0008 double gaus_func(double *x, double *par)
0009 {
0010     // note : par[0] : size
0011     // note : par[1] : mean
0012     // note : par[2] : width
0013     // note : par[3] : offset 
0014     return par[0] * TMath::Gaus(x[0],par[1],par[2]) + par[3];
0015 }
0016 
0017 void zvtx_comp()
0018 {
0019     string private_gen_zvtx_file = "/sphenix/user/ChengWei/INTT/INTT_commissioning/ZeroField/20869/condor_zvtx_folder/file_merged.root";
0020     string F4A_gen_zvtx_file = "/sphenix/user/ChengWei/INTT/INTT_commissioning/ZeroField/F4A_20869/2024_05_07/folder_Data_CombinedNtuple_Run20869_HotDead_BCO_ADC_Survey_test/evt_vtxZ/complete_file/merged_file.root";    
0021     string output_folder = "/sphenix/user/ChengWei/INTT/INTT_commissioning/ZeroField/F4A_20869/2024_05_07/folder_Data_CombinedNtuple_Run20869_HotDead_BCO_ADC_Survey_test/evt_vtxZ/complete_file";
0022     int Nclus_cut = 1000;
0023 
0024     TFile * file_private_gen_zvtx = TFile::Open(private_gen_zvtx_file.c_str());
0025     TTree * tree_private_gen_zvtx = (TTree *)file_private_gen_zvtx -> Get("tree_Z");
0026     ReadINTTZ * private_gen_zvtx = new ReadINTTZ(tree_private_gen_zvtx);
0027 
0028     TFile * file_F4A_gen_zvtx = TFile::Open(F4A_gen_zvtx_file.c_str());
0029     TTree * tree_F4A_gen_zvtx = (TTree *)file_F4A_gen_zvtx -> Get("tree_Z");
0030     ReadINTTZCombine * F4A_gen_zvtx = new ReadINTTZCombine(tree_F4A_gen_zvtx);
0031 
0032     TH1F * zvtx_diff = new TH1F("","zvtx_diff;Private - F4A [cm];Entry",100,-2,2);
0033     TH2F * zvtx_correlation = new TH2F("","zvtx_correlation;Private [cm];F4A [cm]",100,-50,10,100,-50,10);
0034     TH2F * zvtx_diff_nclus = new TH2F("","zvtx_diff_nclus;NClus;Private - F4A [cm]",100,0,9000,100,-2,2);
0035     TH2F * zvtx_diff_privateZ = new TH2F("","zvtx_diff_privateZ;Private reco. Z [cm];Private - F4A [cm]",200,-50,10,200,-2,2);
0036     TH2F * F4AMBDz_privateINTTz_diff_nclus = new TH2F("","F4AMBDz_privateINTTz_diff_nclus;NClus;Private INTT recoZ - F4A MBD recoZ [cm]",200,0,9000,200,-10,10);
0037     TH2F * F4AMBDz_privateINTTz_diff_F4AMBDz = new TH2F("","F4AMBDz_privateINTTz_diff_F4AMBDz;F4A MBD recoZ [cm];Private INTT recoZ - F4A MBD recoZ [cm]",200,-50,10,200,-10,10);
0038 
0039     string unit_text = "cm";
0040     double unit_correction = 0.1;
0041 
0042     int private_gen_eID_offset = 0;
0043 
0044     for (int event_i = 0; event_i < 400000; event_i ++)
0045     {
0046         private_gen_zvtx -> LoadTree(event_i + 1 + private_gen_eID_offset );
0047         private_gen_zvtx -> GetEntry(event_i + 1 + private_gen_eID_offset );
0048 
0049         F4A_gen_zvtx -> LoadTree(event_i);
0050         F4A_gen_zvtx -> GetEntry(event_i);
0051 
0052         // note : MC, but the event is not the minimum bias event
0053         if (F4A_gen_zvtx -> is_min_bias_wozdc == 0) {continue;}
0054 
0055         // note : to get rid of the nan value
0056         if (F4A_gen_zvtx -> Centrality_float != F4A_gen_zvtx -> Centrality_float) {continue;}
0057         if (F4A_gen_zvtx -> MBD_reco_z       != F4A_gen_zvtx -> MBD_reco_z)       {continue;}
0058 
0059         if (F4A_gen_zvtx -> nclu_inner <= 0)   {continue;}
0060         if (F4A_gen_zvtx -> nclu_outer <= 0)   {continue;}
0061         if ((F4A_gen_zvtx -> nclu_inner + F4A_gen_zvtx -> nclu_outer) < Nclus_cut) {continue;}
0062 
0063         double MBD_charge_assy = (F4A_gen_zvtx -> MBD_north_charge_sum - F4A_gen_zvtx -> MBD_south_charge_sum) / (F4A_gen_zvtx -> MBD_north_charge_sum + F4A_gen_zvtx -> MBD_south_charge_sum);
0064         if ( MBD_charge_assy < -1 * ana_map_version::MBD_assy_ratio_cut || MBD_charge_assy >  ana_map_version::MBD_assy_ratio_cut) {continue;}
0065 
0066         if (F4A_gen_zvtx -> LB_Gaus_Width_width < ana_map_version::INTT_zvtx_recohist_gaus_fit_width_cut_l || F4A_gen_zvtx -> LB_Gaus_Width_width > ana_map_version::INTT_zvtx_recohist_gaus_fit_width_cut_r) {continue;}
0067         if (F4A_gen_zvtx -> LB_cut_peak_width <   ana_map_version::INTT_zvtx_recohist_cutgroup_width_cut_l || F4A_gen_zvtx -> LB_cut_peak_width >   ana_map_version::INTT_zvtx_recohist_cutgroup_width_cut_r) {continue;} 
0068         
0069         if (private_gen_zvtx -> LB_Gaus_Width_width < ana_map_version::INTT_zvtx_recohist_gaus_fit_width_cut_l || private_gen_zvtx -> LB_Gaus_Width_width > ana_map_version::INTT_zvtx_recohist_gaus_fit_width_cut_r) {continue;}
0070         if (private_gen_zvtx -> LB_cut_peak_width <   ana_map_version::INTT_zvtx_recohist_cutgroup_width_cut_l || private_gen_zvtx -> LB_cut_peak_width >   ana_map_version::INTT_zvtx_recohist_cutgroup_width_cut_r) {continue;} 
0071 
0072         // note : the INTTF4A_gen_zvtx, MBDF4A_gen_zvtx diff cut
0073         // if ( (F4A_gen_zvtx -> LB_Gaus_Mean_mean - F4A_gen_zvtx -> MBD_reco_z) < ana_map_version::INTTz_MBDz_diff_cut_l || (F4A_gen_zvtx -> LB_Gaus_Mean_mean - F4A_gen_zvtx -> MBD_reco_z) > ana_map_version::INTTz_MBDz_diff_cut_r ) {continue;}
0074 
0075 
0076         if (private_gen_zvtx -> bco_full  != F4A_gen_zvtx -> bco_full) {
0077             cout<<endl;
0078             cout<<"private_gen_eID_offset: "<<private_gen_eID_offset<<endl;
0079             cout<<"at the event : "<<event_i<<", the bco_full is not the same"<<endl;
0080             cout<<"private_gen_zvtx -> bco_full : "<<private_gen_zvtx -> bco_full<<endl;
0081             cout<<"F4A_gen_zvtx -> bco_full : "<<F4A_gen_zvtx -> bco_full<<endl;
0082             private_gen_eID_offset += 1;
0083             continue;
0084         }
0085 
0086         double private_z = private_gen_zvtx -> LB_Gaus_Mean_mean * unit_correction; 
0087         double F4A_z     = F4A_gen_zvtx -> LB_Gaus_Mean_mean * unit_correction;
0088         double F4A_MBD_z = F4A_gen_zvtx -> MBD_reco_z * unit_correction;
0089         double zvtx_diff_value = private_z - F4A_z;
0090         double nclus = F4A_gen_zvtx -> nclu_inner + F4A_gen_zvtx -> nclu_outer;
0091 
0092         zvtx_diff -> Fill(zvtx_diff_value);
0093         zvtx_correlation -> Fill(private_z, F4A_z);
0094         zvtx_diff_nclus -> Fill(nclus, zvtx_diff_value);
0095         zvtx_diff_privateZ -> Fill(private_z, zvtx_diff_value);
0096         F4AMBDz_privateINTTz_diff_nclus -> Fill(nclus, private_z - F4A_MBD_z);
0097         F4AMBDz_privateINTTz_diff_F4AMBDz -> Fill(F4A_MBD_z, private_z - F4A_MBD_z);
0098     }
0099 
0100     SetsPhenixStyle();
0101 
0102     TLatex * ltx = new TLatex();
0103     ltx->SetNDC();
0104     ltx->SetTextSize(0.045);
0105     ltx->SetTextAlign(31);
0106 
0107     TLatex * draw_text = new TLatex();
0108     draw_text -> SetNDC();
0109     draw_text -> SetTextSize(0.03);
0110 
0111     // double legend_text_size = 0.03;
0112     // double legend_upper_y = 0.9;
0113     // TLegend * legend = new TLegend(
0114     //     0.6, // note : x1
0115     //     legend_upper_y - (input_file_full_info.size()) * (legend_text_size + 0.005), // note : y1
0116     //     0.8, // note : x2
0117     //     legend_upper_y // note : y2
0118     // );
0119     // // legend -> SetMargin(0);
0120     // legend->SetTextSize(legend_text_size);
0121 
0122     TCanvas * c1 = new TCanvas("c1","",950,800);
0123 
0124     TLine * coord_line = new TLine();
0125     coord_line -> SetLineWidth(1);
0126     coord_line -> SetLineColor(16);
0127     coord_line -> SetLineStyle(2);
0128 
0129     TF1 * gaus_fit = new TF1("gaus_fit", gaus_func, zvtx_diff->GetXaxis()->GetXmin(), zvtx_diff->GetXaxis()->GetXmax(),4);
0130     gaus_fit -> SetLineColor(2);
0131     gaus_fit -> SetLineWidth(2);
0132     gaus_fit -> SetNpx(1000);
0133 
0134     TF1 * linear_fit = new TF1("linear_fit","pol1", -500, 500);
0135     linear_fit -> SetLineColor(2);
0136     linear_fit -> SetLineWidth(1);
0137     linear_fit -> SetLineStyle(9);
0138 
0139     c1 -> cd();
0140     gaus_fit -> SetParameters(zvtx_diff -> GetBinContent( zvtx_diff -> GetMaximumBin() ), zvtx_diff -> GetBinCenter( zvtx_diff -> GetMaximumBin() ), zvtx_diff -> GetStdDev() * 0.6, 0);
0141     gaus_fit -> SetParLimits(0,0,100000);  // note : size 
0142     gaus_fit -> SetParLimits(2,0,10000);   // note : Width
0143     gaus_fit -> SetParLimits(3,0,10000);   // note : offset
0144     zvtx_diff -> Fit(gaus_fit, "N", "", zvtx_diff -> GetBinCenter( zvtx_diff -> GetMaximumBin() ) - (0.4 * zvtx_diff -> GetStdDev() ), zvtx_diff -> GetBinCenter( zvtx_diff -> GetMaximumBin() ) + (0.4 * zvtx_diff -> GetStdDev() ) );
0145     gaus_fit -> SetRange( gaus_fit->GetParameter(1) - gaus_fit->GetParameter(2) * 0.5, gaus_fit->GetParameter(1) + gaus_fit->GetParameter(2) * 0.5 ); 
0146     
0147     zvtx_diff -> Draw("hist");
0148     gaus_fit -> Draw("lsame");
0149     draw_text -> DrawLatex(0.21, 0.83+0.04, Form("Entries: %.0f", zvtx_diff -> GetEntries()));
0150     draw_text -> DrawLatex(0.21, 0.79+0.04, Form("Mean: %.3f %s", zvtx_diff -> GetMean(), unit_text.c_str() ));
0151     draw_text -> DrawLatex(0.21, 0.75+0.04, Form("Width: %.3f %s", zvtx_diff -> GetStdDev(), unit_text.c_str() ));
0152     draw_text -> DrawLatex(0.21, 0.75, Form("Gaus mean: %.3f %s",gaus_fit -> GetParameter(1), unit_text.c_str() ));
0153     draw_text -> DrawLatex(0.21, 0.69, Form("NClus > %i", Nclus_cut));
0154     ltx->DrawLatex(1 - gPad->GetRightMargin(), 1 - gPad->GetTopMargin() + 0.01, Form("#it{#bf{sPHENIX}} Work-in-progress"));
0155 
0156     c1 -> Print(Form("%s/comp_w_private_zvtx_diff.pdf", output_folder.c_str()));
0157     c1 -> Clear();
0158 
0159 
0160     c1 -> cd();
0161     zvtx_correlation -> Draw("colz0"); 
0162     ltx->DrawLatex(1 - gPad->GetRightMargin(), 1 - gPad->GetTopMargin() + 0.01, Form("#it{#bf{sPHENIX}} Work-in-progress"));
0163     zvtx_correlation -> Fit(linear_fit, "NQ");
0164     linear_fit -> Draw("l same");
0165     draw_text -> DrawLatex(0.21, 0.83, Form("NClus > %i", Nclus_cut));
0166     draw_text -> DrawLatex(0.21, 0.87, Form("Y = %.2fx + %.2f", linear_fit -> GetParameter(1), linear_fit -> GetParameter(0)));
0167     c1 -> Print(Form("%s/comp_w_private_zvtx_correlation.pdf",output_folder.c_str()));
0168     c1 -> Clear();
0169 
0170     c1 -> cd();
0171     zvtx_diff_nclus -> Draw("colz0"); 
0172     draw_text -> DrawLatex(0.21, 0.87, Form("NClus > %i", Nclus_cut));
0173     ltx->DrawLatex(1 - gPad->GetRightMargin(), 1 - gPad->GetTopMargin() + 0.01, Form("#it{#bf{sPHENIX}} Work-in-progress"));
0174     c1 -> Print(Form("%s/comp_w_private_zvtx_diff_nclus.pdf",output_folder.c_str()));
0175     c1 -> Clear();
0176 
0177     c1 -> cd();
0178     zvtx_diff_privateZ -> Draw("colz0"); 
0179     draw_text -> DrawLatex(0.21, 0.87, Form("NClus > %i", Nclus_cut));
0180     ltx->DrawLatex(1 - gPad->GetRightMargin(), 1 - gPad->GetTopMargin() + 0.01, Form("#it{#bf{sPHENIX}} Work-in-progress"));
0181     c1 -> Print(Form("%s/comp_w_private_zvtx_diff_privateZ.pdf",output_folder.c_str()));
0182     c1 -> Clear();
0183 
0184     c1 -> cd();
0185     F4AMBDz_privateINTTz_diff_nclus -> Draw("colz0"); 
0186     draw_text -> DrawLatex(0.21, 0.87, Form("NClus > %i", Nclus_cut));
0187     ltx->DrawLatex(1 - gPad->GetRightMargin(), 1 - gPad->GetTopMargin() + 0.01, Form("#it{#bf{sPHENIX}} Work-in-progress"));
0188     c1 -> Print(Form("%s/comp_w_private_F4AMBDz_privateINTTz_diff_nclus.pdf",output_folder.c_str()));
0189     c1 -> Clear();
0190 
0191     c1 -> cd();
0192     F4AMBDz_privateINTTz_diff_F4AMBDz -> Draw("colz0"); 
0193     draw_text -> DrawLatex(0.21, 0.87, Form("NClus > %i", Nclus_cut));
0194     ltx->DrawLatex(1 - gPad->GetRightMargin(), 1 - gPad->GetTopMargin() + 0.01, Form("#it{#bf{sPHENIX}} Work-in-progress"));
0195     c1 -> Print(Form("%s/comp_w_private_F4AMBDz_privateINTTz_diff_F4AMBDz.pdf",output_folder.c_str()));
0196     c1 -> Clear();
0197 
0198 
0199 }