Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-08-03 08:13:31

0001 #include "TFile.h"
0002 #include "TH1F.h"
0003 #include "TH2F.h"
0004 #include "TObject.h"
0005 #include "TTree.h"
0006 #include <iostream>
0007 
0008 void merge_trees(TFile *outputfile, TTree *tree1, TTree *tree2)
0009 {
0010     outputfile->cd();
0011     TObjArray *branches = tree2->GetListOfBranches();
0012     TIter next(branches);
0013     while (TObject *obj = next())
0014     {
0015         TBranch *branch = dynamic_cast<TBranch *>(obj);
0016         branch->SetTree(tree1);
0017         tree1->GetListOfBranches()->Add(branch);
0018         tree1->GetListOfLeaves()->Add(branch->GetLeaf(branch->GetName()));
0019         tree2->GetListOfBranches()->Remove(branch);
0020     }
0021     tree1->Write("", TObject::kWriteDelete);
0022 }
0023 
0024 void sync_mbd_intt_DST(const char *outputfname, const char *f_mbd_name, const char *t_mbd_name, const char *f_intt_name, const char *t_intt_name, int Nevent)
0025 {
0026     TFile *f_mbd = TFile::Open(f_mbd_name);
0027     TTree *t_mbd = dynamic_cast<TTree *>(f_mbd->Get(t_mbd_name));
0028     int t_mbd_Nevt = t_mbd->GetEntries();
0029     UShort_t femclk;
0030     t_mbd->SetBranchAddress("femclk", &femclk);
0031     std::cout << "--> MBD tree entries:" << t_mbd_Nevt << std::endl;
0032 
0033     TFile *f_intt = TFile::Open(f_intt_name);
0034     TTree *t_intt = dynamic_cast<TTree *>(f_intt->Get(t_intt_name));
0035     int t_intt_Nevt = t_intt->GetEntries();
0036     uint64_t INTT_BCO;
0037     t_intt->SetBranchAddress("INTT_BCO", &INTT_BCO);
0038     std::cout << "--> INTT tree entries:" << t_intt_Nevt << std::endl;
0039 
0040     TFile *outfile = TFile::Open(outputfname, "recreate");
0041     TTree *t_intt_clone = t_intt->CloneTree(0);
0042     TTree *t_mbd_clone = t_mbd->CloneTree(0);
0043 
0044     int run_event = (Nevent > 0 && Nevent < t_mbd_Nevt) ? Nevent : t_mbd_Nevt;
0045 
0046     TH2F *intt_mbd_bco_presync = new TH2F("intt_mbd_bco_presync", "INTT - MBD", 100, 0, run_event * 1.5, 100, -10, 100000);
0047     intt_mbd_bco_presync->GetXaxis()->SetTitle("evt");
0048     intt_mbd_bco_presync->GetYaxis()->SetTitle("clock_diff");
0049 
0050     TH1F *intt_mbd_bco_presync_1D = new TH1F("intt_mbd_bco_presync_1D", "INTT - MBD", 100, -10, 100000);
0051     intt_mbd_bco_presync_1D->GetXaxis()->SetTitle("clock_diff");
0052     intt_mbd_bco_presync_1D->GetYaxis()->SetTitle("entry");
0053 
0054     TH2F *intt_mbd_bco_postsync = new TH2F("intt_mbd_bco_postsync", "INTT - MBD", 100, 0, run_event * 1.5, 100, -10, 100000);
0055     intt_mbd_bco_postsync->GetXaxis()->SetTitle("evt");
0056     intt_mbd_bco_postsync->GetYaxis()->SetTitle("clock_diff");
0057 
0058     TH1F *intt_mbd_bco_postsync_1D = new TH1F("intt_mbd_bco_postsync_1D", "INTT - MBD", 100, -10, 100000);
0059     intt_mbd_bco_postsync_1D->GetXaxis()->SetTitle("clock_diff");
0060     intt_mbd_bco_postsync_1D->GetYaxis()->SetTitle("entry");
0061 
0062     int prev_mbdclk = 0;
0063     int prev_inttbcofull = 0;
0064 
0065     int mbd_evt_offset = 0;
0066     int intt_evt_offset = 0;
0067 
0068     std::cout << "--> start to synchronize the events" << std::endl;
0069     for (int entry = 0; entry < run_event; entry++)
0070     {
0071         t_mbd->GetEntry(entry + mbd_evt_offset);
0072         t_intt->GetEntry(entry + intt_evt_offset);
0073 
0074         int mbdclk = femclk;
0075         int inttbcofull = INTT_BCO;
0076         int inttbco16bit = inttbcofull & 0xFFFF;
0077         int mbd_prvdiff = (mbdclk - prev_mbdclk) & 0xFFFF;
0078         int intt_prvdiff = inttbcofull - prev_inttbcofull;
0079 
0080         if (entry % 1000 == 0)
0081         {
0082             std::cout << std::dec << entry << " mbdclk=" << mbdclk << " inttbco16bit=" << inttbco16bit << " (mbd-intt) " << ((mbdclk - inttbco16bit) & 0xFFFF) << " (MBD, femclk-prev) " << mbd_prvdiff
0083                       << " (INTT, bco-prev) " << intt_prvdiff << " intt_evt_offset " << intt_evt_offset << std::endl;
0084         }
0085 
0086         intt_mbd_bco_presync->Fill(entry, (mbdclk - inttbco16bit) & 0xFFFF);
0087         intt_mbd_bco_presync_1D->Fill((mbdclk - inttbco16bit) & 0xFFFF);
0088 
0089         prev_mbdclk = mbdclk;
0090         prev_inttbcofull = inttbcofull;
0091 
0092         bool is_synced = false;
0093 
0094         while (entry + 1 + intt_evt_offset < t_intt_Nevt)
0095         {
0096             t_intt->GetEntry(entry + 1 + intt_evt_offset);
0097             int next_inttbco16bit = INTT_BCO & 0xFFFF;
0098 
0099             t_mbd->GetEntry(entry + 1);
0100             int next_mbdclk = femclk;
0101 
0102             if (((next_mbdclk - next_inttbco16bit) & 0xFFFF) != ((mbdclk - inttbco16bit) & 0xFFFF))
0103             {
0104                 std::cout << "Not synced: intt entry(entry+1+intt_evt_offset)=" << entry + 1 + intt_evt_offset << ", mbd entry(entry+1)=" << entry + 1 << ", intt_evt_offset=" << intt_evt_offset
0105                           << ", ((next_mbdclk - next_inttbco16bit) & 0xFFFF)=" << ((next_mbdclk - next_inttbco16bit) & 0xFFFF)
0106                           << ", ((mbdclk - inttbco16bit) & 0xFFFF)=" << ((mbdclk - inttbco16bit) & 0xFFFF)
0107                           << ", ((next_mbdclk - next_inttbco16bit) & 0xFFFF)-((mbdclk - inttbco16bit) & 0xFFFF)=" << ((next_mbdclk - next_inttbco16bit) & 0xFFFF) - ((mbdclk - inttbco16bit) & 0xFFFF)
0108                           << std::endl;
0109 
0110                 intt_evt_offset += 1;
0111             }
0112             else
0113             {
0114                 is_synced = true;
0115                 break;
0116             }
0117         }
0118 
0119         if (is_synced)
0120         {
0121             intt_mbd_bco_postsync->Fill(entry, (mbdclk - inttbco16bit) & 0xFFFF);
0122             intt_mbd_bco_postsync_1D->Fill((mbdclk - inttbco16bit) & 0xFFFF);
0123             t_intt_clone->Fill();
0124             t_mbd_clone->Fill();
0125         }
0126     }
0127 
0128     std::cout << "--> INTT Nevent post-sync: " << t_intt_clone->GetEntries() << std::endl << "--> MBD Nevent post-sync: " << t_mbd_clone->GetEntries() << std::endl;
0129 
0130     merge_trees(outfile, t_intt_clone, t_mbd_clone);
0131 
0132     intt_mbd_bco_presync->Write("", TObject::kWriteDelete);
0133     intt_mbd_bco_presync_1D->Write("", TObject::kWriteDelete);
0134     intt_mbd_bco_postsync->Write("", TObject::kWriteDelete);
0135     intt_mbd_bco_postsync_1D->Write("", TObject::kWriteDelete);
0136     f_intt->Close();
0137     f_mbd->Close();
0138     outfile->Close();
0139 
0140     std::cout << "--> Synchronization done" << std::endl
0141               << "--> Merged file:" << outputfname << ", Tree name: " << t_intt_name << std::endl
0142               << "--> TH2F plot: intt_mbd_bco_pre/postsync, can be checked in the merged file" << std::endl
0143               << "--> TH1F plot: intt_mbd_bco_pre/postsync_1D, can be checked in the merged file" << std::endl;
0144 }
0145 
0146 void intt_mbd_evt_combiner(int Nevent = -1,                         //
0147                            const char *outputfname = "output.root", //
0148                            const char *f_intt_name = "intt.root",   //
0149                            const char *t_intt_name = "EventTree",   //
0150                            const char *f_mbd_name = "mbd.root",     //
0151                            const char *t_mbd_name = "EventTree")
0152 {
0153     gBenchmark->Start("InttMbdEvtCombiner");
0154 
0155     sync_mbd_intt_DST(outputfname, f_mbd_name, t_mbd_name, f_intt_name, t_intt_name, Nevent);
0156 
0157     gBenchmark->Show("InttMbdEvtCombiner");
0158 }