Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-08-05 08:11:23

0001 #include <iostream>
0002 #include "ReadEventBase.C"
0003 
0004 bool count_hist_remain(TH1F * hist_1, TH1F * hist_2)
0005 {
0006     int N_remain = 0;
0007     for (int i = 0; i < hist_1 -> GetNbinsX(); i++) { N_remain += abs(hist_1 -> GetBinContent(i+1) - hist_2 -> GetBinContent(i+1)); }
0008     return (N_remain == 0 && hist_1 -> GetEntries() > 1 && hist_1 -> GetStdDev() != 0 && hist_1 -> GetMean() != 0 && (hist_1 -> GetEntries() - hist_1 -> GetBinContent(1)) > 10 /*hist_1 -> GetBinContent(1)/double(hist_1 -> GetEntries()) < 0.8*/) ? true : false;
0009 }
0010 
0011 void CloneEvent_check(int study_run_index)
0012 {
0013     // int study_run_index = 4;
0014     vector<string> ZF_run_list = {"20864","20866","20867","20868","20869","20878","20880","20881","20883","20885"};
0015     string file_directory = "/sphenix/tg/tg01/commissioning/INTT/root_files";
0016     
0017     TFile * file_in;
0018     TTree * tree_in;
0019     ReadEventBase * data_in;
0020 
0021     int show_index = 0;
0022 
0023     int fNhits;
0024     int pid[100000];
0025     int module[100000];
0026     int chip_id[100000];
0027     int chan_id[100000];
0028     int adc[100000];
0029     int bco[100000];
0030     Long64_t bco_full[100000];
0031 
0032     TH1F * chan_dist_bco0 = new TH1F("","",128,0,128);
0033     TH1F * chan_dist_bco1 = new TH1F("","",128,0,128);
0034     TH1F * chan_dist_bco2 = new TH1F("","",128,0,128);
0035     TH1F * chan_dist_bco3 = new TH1F("","",128,0,128);
0036     TH1F * chan_dist_bco4 = new TH1F("","",128,0,128);
0037     TH1F * chan_dist_bco5 = new TH1F("","",128,0,128);
0038 
0039     for (int i = 0; i < 8; i++)
0040     {
0041         if ( gSystem->AccessPathName(Form("%s/beam_intt%i-000%s-0000_event_base.root", file_directory.c_str(), i, ZF_run_list[study_run_index].c_str())) == true  ) 
0042         {
0043             cout<<"file not exist: "<<Form("%s/beam_intt%i-000%s-0000_event_base.root", file_directory.c_str(), i, ZF_run_list[study_run_index].c_str())<<endl;
0044             continue;
0045         }
0046         file_in = TFile::Open(Form("%s/beam_intt%i-000%s-0000_event_base.root", file_directory.c_str(), i, ZF_run_list[study_run_index].c_str()));
0047         tree_in = (TTree*) file_in -> Get("tree");
0048         long long N_event = tree_in -> GetEntries();
0049         tree_in -> SetBranchStatus("*",0);
0050         tree_in -> SetBranchStatus("fNhits",1);
0051         // tree_in -> SetBranchStatus("fhitArray.pid",1);
0052         // tree_in -> SetBranchStatus("fhitArray.module",1);
0053         // tree_in -> SetBranchStatus("fhitArray.chip_id",1);
0054         tree_in -> SetBranchStatus("fhitArray.chan_id",1);
0055         // tree_in -> SetBranchStatus("fhitArray.adc",1);
0056         // tree_in -> SetBranchStatus("fhitArray.bco",1);
0057         // tree_in -> SetBranchStatus("fhitArray.bco_full",1);
0058 
0059         tree_in -> SetBranchAddress("fNhits",&fNhits);
0060         tree_in -> GetEntry(0); // note : actually I really don't know why this line is necessary.
0061         // tree_in -> SetBranchAddress("fhitArray.pid",&pid[0]);
0062         // tree_in -> SetBranchAddress("fhitArray.module",&module[0]);
0063         // tree_in -> SetBranchAddress("fhitArray.chip_id",&chip_id[0]);
0064         tree_in -> SetBranchAddress("fhitArray.chan_id",&chan_id[0]);
0065         // tree_in -> SetBranchAddress("fhitArray.adc",&adc[0]);
0066         // tree_in -> SetBranchAddress("fhitArray.bco",&bco[0]);
0067         // tree_in -> SetBranchAddress("fhitArray.bco_full",&bco_full[0]);
0068 
0069 
0070 
0071         // data_in = new ReadEventBase(tree_in);
0072         cout<<" "<<endl; 
0073         cout<<"check the run index: "<<ZF_run_list[study_run_index]<<endl;
0074         cout<<"checking the server ID: "<<i<<" N event : "<<N_event<<endl;
0075 
0076         double N_same_hit_evt_bco1 = 0;
0077         double N_same_hit_evt_bco2 = 0;
0078         double N_same_hit_evt_bco3 = 0;
0079         double N_same_hit_evt_bco4 = 0;
0080         double N_same_hit_evt_bco5 = 0;
0081 
0082         for (int event_i = 0; event_i < N_event; event_i++)
0083         {
0084             // if (event_i % 100 == 0) cout<<"event_i : "<<event_i<<endl;
0085 
0086             // data_in -> LoadTree(event_i);
0087             tree_in -> GetEntry(event_i);
0088             int fNhits_bco0 = fNhits;
0089             for (int hit_i = 0; hit_i < fNhits; hit_i++){
0090                 chan_dist_bco0 -> Fill(chan_id[hit_i]);
0091             }
0092             
0093             // note : ------------------ ------------------ ------------------ ------------------ ------------------ ------------------
0094             if (event_i == (tree_in -> GetEntries() - 1)) {
0095                 chan_dist_bco0 -> Reset("ICESM");
0096                 chan_dist_bco1 -> Reset("ICESM");
0097                 chan_dist_bco2 -> Reset("ICESM");
0098                 chan_dist_bco3 -> Reset("ICESM");
0099                 chan_dist_bco4 -> Reset("ICESM");
0100                 chan_dist_bco5 -> Reset("ICESM");
0101                 continue;
0102             }
0103             // data_in -> LoadTree(event_i+1);
0104             tree_in -> GetEntry(event_i+1);
0105             if (fNhits_bco0 == fNhits) { 
0106                 // N_same_hit_evt_bco1 += 1; 
0107                 for (int hit_i = 0; hit_i < fNhits; hit_i++) {chan_dist_bco1 -> Fill(chan_id[hit_i]);}
0108                 if (count_hist_remain(chan_dist_bco0, chan_dist_bco1)){
0109                     N_same_hit_evt_bco1 += 1;
0110                     cout<<"--- identical fNhits and chan distribution in bco gap 1, event_i : "<<event_i<<" with fNhits: "<<fNhits<<endl;
0111                 }
0112 
0113                 
0114             }
0115 
0116             // note : ------------------ ------------------ ------------------ ------------------ ------------------ ------------------
0117             if (event_i == (tree_in -> GetEntries() - 2)) {
0118                 chan_dist_bco0 -> Reset("ICESM");
0119                 chan_dist_bco1 -> Reset("ICESM");
0120                 chan_dist_bco2 -> Reset("ICESM");
0121                 chan_dist_bco3 -> Reset("ICESM");
0122                 chan_dist_bco4 -> Reset("ICESM");
0123                 chan_dist_bco5 -> Reset("ICESM");
0124                 continue;
0125             }
0126             // data_in -> LoadTree(event_i+2);
0127             tree_in -> GetEntry(event_i+2);
0128             if (fNhits_bco0 == fNhits) { 
0129                 // N_same_hit_evt_bco2 += 1; 
0130                 for (int hit_i = 0; hit_i < fNhits; hit_i++) {chan_dist_bco2 -> Fill(chan_id[hit_i]);}
0131                 if (count_hist_remain(chan_dist_bco0, chan_dist_bco2)){
0132                     N_same_hit_evt_bco2 += 1;
0133                     cout<<"--- identical fNhits and chan distribution in bco gap 2, event_i : "<<event_i<<" with fNhits: "<<fNhits<<endl;
0134                 }
0135 
0136 
0137             }
0138 
0139             // note : ------------------ ------------------ ------------------ ------------------ ------------------ ------------------
0140             if (event_i == (tree_in -> GetEntries() - 3)) {
0141                 chan_dist_bco0 -> Reset("ICESM");
0142                 chan_dist_bco1 -> Reset("ICESM");
0143                 chan_dist_bco2 -> Reset("ICESM");
0144                 chan_dist_bco3 -> Reset("ICESM");
0145                 chan_dist_bco4 -> Reset("ICESM");
0146                 chan_dist_bco5 -> Reset("ICESM");
0147                 continue;
0148             }
0149             // data_in -> LoadTree(event_i+3);
0150             tree_in -> GetEntry(event_i+3);
0151             if (fNhits_bco0 == fNhits) { 
0152                 // N_same_hit_evt_bco3 += 1; 
0153                 for (int hit_i = 0; hit_i < fNhits; hit_i++) {chan_dist_bco3 -> Fill(chan_id[hit_i]);}
0154                 if (count_hist_remain(chan_dist_bco0, chan_dist_bco3)){
0155                     N_same_hit_evt_bco3 += 1;
0156                     cout<<"--- identical fNhits and chan distribution in bco gap 3, event_i : "<<event_i<<" with fNhits: "<<fNhits<<endl;
0157                 }
0158 
0159 
0160             }
0161 
0162             // note : ------------------ ------------------ ------------------ ------------------ ------------------ ------------------
0163             if (event_i == (tree_in -> GetEntries() - 4)) {
0164                 chan_dist_bco0 -> Reset("ICESM");
0165                 chan_dist_bco1 -> Reset("ICESM");
0166                 chan_dist_bco2 -> Reset("ICESM");
0167                 chan_dist_bco3 -> Reset("ICESM");
0168                 chan_dist_bco4 -> Reset("ICESM");
0169                 chan_dist_bco5 -> Reset("ICESM");
0170                 continue;
0171             }
0172             // data_in -> LoadTree(event_i+4);
0173             tree_in -> GetEntry(event_i+4);
0174             if (fNhits_bco0 == fNhits) { 
0175                 // N_same_hit_evt_bco4 += 1; 
0176                 for (int hit_i = 0; hit_i < fNhits; hit_i++) {chan_dist_bco4 -> Fill(chan_id[hit_i]);}
0177                 if (count_hist_remain(chan_dist_bco0, chan_dist_bco4)){
0178                     N_same_hit_evt_bco4 += 1;
0179                     cout<<"--- identical fNhits and chan distribution in bco gap 4, event_i : "<<event_i<<" with fNhits: "<<fNhits<<endl;
0180                 }
0181 
0182 
0183             }
0184 
0185             // note : ------------------ ------------------ ------------------ ------------------ ------------------ ------------------
0186             if (event_i == (tree_in -> GetEntries() - 5)) {
0187                 chan_dist_bco0 -> Reset("ICESM");
0188                 chan_dist_bco1 -> Reset("ICESM");
0189                 chan_dist_bco2 -> Reset("ICESM");
0190                 chan_dist_bco3 -> Reset("ICESM");
0191                 chan_dist_bco4 -> Reset("ICESM");
0192                 chan_dist_bco5 -> Reset("ICESM");
0193                 continue;
0194             }
0195             // data_in -> LoadTree(event_i+5);
0196             tree_in -> GetEntry(event_i+5);
0197             if (fNhits_bco0 == fNhits) { 
0198                 // N_same_hit_evt_bco5 += 1; 
0199                 for (int hit_i = 0; hit_i < fNhits; hit_i++) {chan_dist_bco5 -> Fill(chan_id[hit_i]);}
0200                 if (count_hist_remain(chan_dist_bco0, chan_dist_bco5)){
0201                     N_same_hit_evt_bco5 += 1;
0202                     cout<<"--- identical fNhits and chan distribution in bco gap 5, event_i : "<<event_i<<" with fNhits: "<<fNhits<<endl;
0203                 }
0204 
0205 
0206             }
0207 
0208             // note : ------------------ ------------------ ------------------ ------------------ ------------------ ------------------
0209             chan_dist_bco0 -> Reset("ICESM");
0210             chan_dist_bco1 -> Reset("ICESM");
0211             chan_dist_bco2 -> Reset("ICESM");
0212             chan_dist_bco3 -> Reset("ICESM");
0213             chan_dist_bco4 -> Reset("ICESM");
0214             chan_dist_bco5 -> Reset("ICESM");
0215 
0216             // cout<<"event_i : "<<event_i<<endl;
0217             // cout<<"fNhits : "<<data_in -> fNhits<<endl;
0218         }
0219 
0220         cout<<"N_same_hit_evt_bco1 : "<<N_same_hit_evt_bco1<<endl;
0221         cout<<"N_same_hit_evt_bco2 : "<<N_same_hit_evt_bco2<<endl;
0222         cout<<"N_same_hit_evt_bco3 : "<<N_same_hit_evt_bco3<<endl;
0223         cout<<"N_same_hit_evt_bco4 : "<<N_same_hit_evt_bco4<<endl;
0224         cout<<"N_same_hit_evt_bco5 : "<<N_same_hit_evt_bco5<<endl;
0225 
0226         file_in -> Close();
0227     }
0228 
0229 }