Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-08-06 08:18:42

0001 #include "InttStreamQA.h"
0002 
0003 /// Fun4All includes
0004 #include <fun4all/Fun4AllReturnCodes.h>
0005 #include <fun4all/Fun4AllServer.h>
0006 #include <phool/PHCompositeNode.h>
0007 #include <phool/getClass.h>
0008 
0009 #include <ffarawobjects/Gl1Packet.h>
0010 #include <ffarawobjects/InttRawHit.h>
0011 #include <ffarawobjects/InttRawHitContainer.h>
0012 #include <phool/PHPointerListIterator.h>
0013 #include <qautils/QAHistManagerDef.h>
0014 
0015 #include <TFile.h>
0016 #include <TH1.h>
0017 #include <TH2.h>
0018 
0019 #include <array>
0020 #include <map>
0021 #include <set>
0022 #include <sstream>
0023 
0024 /**
0025  * Constructor of module
0026  */
0027 InttStreamQA::InttStreamQA(const std::string& name)
0028   : SubsysReco(name)
0029 {
0030 }
0031 
0032 int InttStreamQA::InitRun(PHCompositeNode* topNode)
0033 {
0034   if (Verbosity() > 5)
0035   {
0036     std::cout << "Beginning InitRun in InttStreamQA" << std::endl;
0037   }
0038 
0039   PHNodeIterator trkr_itr(topNode);
0040   PHCompositeNode *intt_node = dynamic_cast<PHCompositeNode *>(
0041       trkr_itr.findFirst("PHCompositeNode", "INTT"));  
0042   if(!intt_node)
0043   {
0044     std::cout << PHWHERE << " No INTT node found, exit" << std::endl;
0045     return Fun4AllReturnCodes::ABORTRUN;
0046   }
0047   PHNodeIterator intt_itr(intt_node);
0048   PHPointerListIterator<PHNode> iter(intt_itr.ls());
0049   PHNode *thisnode;
0050   while((thisnode = iter()))
0051   {
0052     if(thisnode->getType() !="PHIODataNode")
0053     {
0054       continue;
0055     }
0056     PHIODataNode<InttRawHitContainer> *theNode = static_cast<PHIODataNode<InttRawHitContainer> *>(thisnode);
0057     if(theNode)
0058     {
0059       std::cout << PHWHERE << " Found INTT Raw hit container node " << theNode->getName() << std::endl;
0060       auto cont = (InttRawHitContainer*)theNode->getData();
0061       if(cont)
0062       {
0063         m_rawhit_containers.push_back(cont);
0064       }
0065     }
0066   }
0067 
0068   createHistos();
0069 
0070   return Fun4AllReturnCodes::EVENT_OK;
0071 }
0072 
0073 /**
0074  * Main workhorse function where each event is looped over and
0075  * data from each event is collected from the node tree for analysis
0076  */
0077 int InttStreamQA::process_event(PHCompositeNode* topNode)
0078 {
0079   static int event = 0;
0080   static int nskip = 0;
0081 
0082   Gl1Packet* gl1 = findNode::getClass<Gl1Packet>(topNode, "GL1RAWHIT");
0083 
0084   uint64_t bco_gl1 = (gl1 != nullptr) ? gl1->getBCO() : std::numeric_limits<uint64_t>::max();
0085   // int      evt_gl1  = (gl1 !=nullptr) ? gl1->getEvtSequence() : -1;
0086   int bunch_gl1 = (gl1 != nullptr) ? gl1->getBunchNumber() : -1;
0087 
0088   // uint64_t trig =  (gl1!=nullptr) ? gl1->getLiveVector() : 0;
0089   // uint64_t trig =  (gl1!=nullptr) ? gl1->getScaledVector() : 0;
0090 
0091 for(auto& rawhitmap : m_rawhit_containers)
0092 {
0093   if (rawhitmap == nullptr)
0094   {
0095     std::cout << PHWHERE << "rawhit  is null" << std::endl;
0096     return Fun4AllReturnCodes::EVENT_OK;
0097   }
0098 
0099   uint64_t bcointt = (rawhitmap->get_nhits() > 0)
0100                          ? rawhitmap->get_hit(0)->get_bco()
0101                          : std::numeric_limits<uint64_t>::max();
0102 
0103   static uint64_t prebcointt = 0;
0104 
0105   if (bco_gl1 == std::numeric_limits<uint64_t>::max())
0106   {
0107     std::cout << "StreamQA bco is max. not valid" << std::endl;
0108   }
0109   if (bcointt == std::numeric_limits<uint64_t>::max())
0110   {
0111     std::cout << "StreamQA inttbco is max. no intt data valid. skip nth: " << nskip << std::endl;
0112     nskip++;
0113   }
0114 
0115   // to 40bit
0116   bco_gl1 &= 0xFFFFFFFFFFULL;
0117   bcointt &= 0xFFFFFFFFFFULL;
0118 
0119   uint64_t bcodiff = bcointt - prebcointt;
0120 
0121   int evtcnt = (rawhitmap->get_nhits() > 0)
0122                    ? rawhitmap->get_hit(0)->get_event_counter()
0123                    : std::numeric_limits<int>::max();
0124 
0125   //////////////////////////////
0126 
0127   int64_t diff_inttgl1 = bcointt - bco_gl1;
0128   h_bcointtgl1_diff->Fill(diff_inttgl1);
0129 
0130   if (Verbosity() > 5)
0131   {
0132     std::cout << event << "   bco : intt:0x" << std::hex << bcointt << ", diff 0x" << bcodiff
0133               << " : " << std::dec << " " << evtcnt << " " << diff_inttgl1 << " gl1:0x" << std::hex << bco_gl1 << std::dec << std::endl;
0134   }
0135 
0136   event++;
0137 
0138   h_bunch_gl1->Fill(bunch_gl1);
0139 
0140   std::array<std::set<uint>,8> vUnique;
0141   std::array<std::map<uint, int>,8> vchipbco;  // key, ihit
0142 
0143   // loop rawhits to remove copy hit
0144   uint nhit = rawhitmap->get_nhits();
0145   for (uint ihit = 0; ihit < nhit; ihit++)
0146   {
0147     InttRawHit* hit = rawhitmap->get_hit(ihit);
0148 
0149     int ifelix = hit->get_packetid() - 3001;
0150     uint bco = hit->get_FPHX_BCO();                         // 7bit
0151     uint64_t bcofull = (hit->get_bco() & 0xFFFFFFFFFFULL);  // 7bit
0152 
0153     uint ladder = hit->get_fee();               // 0-13 (4bit)
0154     uint chip = (hit->get_chip_id() - 1) % 26;  // 0-26 (5bit)
0155     uint chan = hit->get_channel_id();          // 0-127 (7bit)
0156     uint adc = hit->get_adc();                  // 3bit
0157 
0158     // check the difference between strobeBCO(bcofull) and gl1BCO
0159     int64_t bcogl1diff = bcofull - bco_gl1;
0160     // std::cout<<"bco-gl1diff " <<bcogl1diff<<std::endl;
0161     h_bcogl1diff_felix[ifelix]->Fill(bcogl1diff);
0162 
0163     // lad[25-22]+chip[21-17]+chan[16-10]+adc[9-7]+bco[6-0]
0164     uint key = ((ladder & 0xFU) << 22U) | ((chip & 0x1FU) << 17U) | ((chan & 0x7FU) << 10U) | ((adc & 0x7U) << 7U) | (bco & 0x7FU);
0165 
0166     auto ret = vUnique[ifelix].insert(key);
0167     if (ret.second)  // insertion successfull --> entry did not exist
0168     {
0169       uint chipbcokey = ((ladder & 0xFU) << 22U) | ((chip & 0x1FU) << 17U) | (bco & 0x7FU);
0170       vchipbco[ifelix].insert(std::make_pair(chipbcokey, ihit));  // no ADC info
0171 
0172       h_bco[ifelix]->Fill(ladder * 26 + chip + 0.5, bco + 0.5);
0173       h_hit[ifelix]->Fill(ladder * 26 + chip + 0.5, chan + 0.5);
0174     }
0175 
0176     // std::cout<<"    hit : "<<ihit<<" "<<ifelix<<" 0x"<<std::hex<<bco<<std::dec<<std::endl;
0177   }
0178 
0179   ////////////////////////////
0180   // felix by felix analysis
0181   std::array<std::map<int, int>,8> vbcodiff_felix;
0182   for (int ifelix = 0; ifelix < 8; ifelix++)
0183   {
0184     for (auto val : vchipbco[ifelix])
0185     {
0186       uint bco = (val.first) & 0x7FU;
0187       h_bco_felix[ifelix]->Fill(bco);
0188 
0189       InttRawHit* hit = rawhitmap->get_hit(val.second);
0190 
0191       uint64_t bcofull = (hit->get_bco() & 0xFFFFFFFFFFULL);  // 7bit
0192 
0193       // stream mode
0194       // uint64_t bcofull_reco = bco + bcointt;
0195       int bcofull_reco = bco + bcofull;
0196 
0197       int bcointtgl1_diff = bcofull_reco - bco_gl1;
0198       int bcointthit_diff = bcofull_reco - bcointt;
0199 
0200       h_bcoreco_diff[ifelix]->Fill(bcointtgl1_diff);
0201       h_bcorecointt_diff[ifelix]->Fill(bcointthit_diff);
0202 
0203       h_bunch_bco[ifelix]->Fill(bco, bunch_gl1);
0204       h_bunch_evt_bcodiff[ifelix]->Fill(bcointtgl1_diff, bunch_gl1);
0205 
0206       auto bco_itr = vbcodiff_felix[ifelix].find(bcointtgl1_diff);
0207       if (bco_itr == vbcodiff_felix[ifelix].end())
0208       {
0209         vbcodiff_felix[ifelix].insert(std::make_pair(bcointtgl1_diff, 1));
0210         h_bcoreco_evt_diff[ifelix]->Fill(bcointtgl1_diff);  // event by event for each felix
0211       }
0212       else
0213       {
0214         bco_itr->second += 1;
0215       }
0216     }
0217   }
0218 
0219   ////////////////////////////
0220   // all felix combined analysis
0221   std::map<int, int> vbcodiff_all;
0222 
0223   for (auto& ifelix : vbcodiff_felix)
0224   {
0225     for (auto& val : ifelix)
0226     {
0227       int bcointtgl1_diff = val.first;
0228       // int count           = val.second;
0229 
0230       // std::cout<<"             recobco diff : "<<bcointtgl1_diff<<" "<<count<<" "<<ifelix<<std::endl;
0231 
0232       auto bco_all_itr = vbcodiff_all.find(bcointtgl1_diff);
0233       if (bco_all_itr == vbcodiff_all.end())
0234       {
0235         vbcodiff_all.insert(std::make_pair(bcointtgl1_diff, 1));
0236         h_bcoreco_evt_diff_all->Fill(bcointtgl1_diff);
0237 
0238         if (bcointtgl1_diff == 23)  // should use variable
0239         {
0240           h_bunch_all->Fill(bunch_gl1);
0241         }
0242       }
0243       else
0244       {
0245         bco_all_itr->second += 1;
0246       }
0247     }
0248   }
0249 
0250   prebcointt = bcointt;
0251 }
0252   return Fun4AllReturnCodes::EVENT_OK;
0253 }
0254 
0255 /**
0256  * End the module and finish any data collection. Clean up any remaining
0257  * loose ends
0258  */
0259 int InttStreamQA::End(PHCompositeNode* /*topNode*/)
0260 {
0261   if (Verbosity() > 1)
0262   {
0263     std::cout << "Ending InttStreamQA analysis package" << std::endl;
0264   }
0265 
0266   /*
0267     TFile *froot = TFile::Open("streamrecof4a.root","recreate");
0268     for(int i=0; i<8; i++){
0269       h_bco[i]->Write();
0270       h_hit[i]->Write();
0271       h_bco_felix[i]->Write();
0272       //h_bunch[i]->Write();
0273       h_bcoreco_diff[i]->Write();
0274       h_bcoreco_evt_diff[i]->Write();
0275       h_bcorecointt_diff[i]->Write();
0276       h_bcogl1diff_felix[i]->Write();
0277 
0278       //h_bunch_strb[i]->Write();
0279       h_bunch_evt_bcodiff[i]->Write();
0280       h_bunch_bco[i]->Write();
0281       //h_bcoprediff[i]->Write();
0282     }
0283     h_bcoreco_evt_diff_all->Write();
0284     h_bunch_all->Write();
0285     h_bunch_gl1->Write();
0286     h_bcointtgl1_diff->Write();
0287     froot->Close();
0288   */
0289 
0290   return 0;
0291 }
0292 
0293 std::string InttStreamQA::getHistoPrefix() const
0294 {
0295   return std::string("h_") + Name() + std::string("_");
0296 }
0297 
0298 void InttStreamQA::createHistos()
0299 {
0300   auto hm = QAHistManagerDef::getHistoManager();
0301   assert(hm);
0302 
0303   std::string sname, stitle;
0304   for (int i = 0; i < 8; i++)
0305   {
0306     sname = (getHistoPrefix() + "bco_" + std::to_string(i));
0307     stitle = ("bco_" + std::to_string(i) + ";bco;ladder*chip");
0308     h_bco[i] = new TH2F(sname.c_str(), stitle.c_str(), 26 * 14, 0, 26 * 14, 140, -7, 133);
0309     hm->registerHisto(h_bco[i]);
0310 
0311     sname = (getHistoPrefix() + "hit_" + std::to_string(i));
0312     stitle = ("hit_" + std::to_string(i) + ";bco;ladder*chip");
0313     h_hit[i] = new TH2F(sname.c_str(), stitle.c_str(), 26 * 14, 0, 26 * 14, 128, 0, 128);
0314     hm->registerHisto(h_hit[i]);
0315 
0316     // RecoBCO - GL1BCO
0317     sname = (getHistoPrefix() + "bcoreco_diff_" + std::to_string(i));
0318     stitle = ("bcoreco diff_" + std::to_string(i));
0319     h_bcoreco_diff[i] = new TH1F(sname.c_str(), stitle.c_str(), 540, -270, 270);
0320     hm->registerHisto(h_bcoreco_diff[i]);
0321 
0322     sname = (getHistoPrefix() + "bcoreco_evt_diff_" + std::to_string(i));
0323     stitle = ("bcoreco evt diff_" + std::to_string(i));
0324     h_bcoreco_evt_diff[i] = new TH1F(sname.c_str(), stitle.c_str(), 540, -270, 270);
0325     hm->registerHisto(h_bcoreco_evt_diff[i]);
0326 
0327     // RecoBCO - StrobeBCO, same as FPHX BCO
0328     sname = (getHistoPrefix() + "bcorecointt_diff_" + std::to_string(i));
0329     stitle = ("bcoreco intt diff_" + std::to_string(i));
0330     h_bcorecointt_diff[i] = new TH1F(sname.c_str(), stitle.c_str(), 540, -270, 270);
0331     hm->registerHisto(h_bcorecointt_diff[i]);
0332 
0333     // ChipBCO
0334     sname = (getHistoPrefix() + "bco_felix_" + std::to_string(i));
0335     stitle = ("bco_felix_" + std::to_string(i));
0336     h_bco_felix[i] = new TH1F(sname.c_str(), stitle.c_str(), 128, 0, 128);
0337     hm->registerHisto(h_bco_felix[i]);
0338 
0339     // StrobeBCO - GL1BCO
0340     sname = (getHistoPrefix() + "bcogl1diff_felix_" + std::to_string(i));
0341     stitle = ("bcogl1diff_felix_" + std::to_string(i));
0342     h_bcogl1diff_felix[i] = new TH1F(sname.c_str(), stitle.c_str(), 1024, -512, 512);
0343     hm->registerHisto(h_bcogl1diff_felix[i]);
0344 
0345     // h_bunch[i] = new TH1F((getHistoPrefix()+Form("bunch_%d", i)).c_str(), Form("bunch @ trigger_%d", i), 150, -15, 135);
0346     // hm->registerHisto(h_bunch[i]);
0347 
0348     // h_bunch_strb[i] = new TH1F((getHistoPrefix()+Form("bunch_strb_%d", i)).c_str(), Form("bunch @ strobe_%d", i), 150, -15, 135);
0349     // hm->registerHisto(h_bunch_strb[i]);
0350 
0351     // RecoBCO - Gl1BCO, vs  bunch, to see the peak
0352     sname = (getHistoPrefix() + "bunch_evt_bcodiff_" + std::to_string(i));
0353     stitle = ("bunch @ strobe_" + std::to_string(i));
0354     h_bunch_evt_bcodiff[i] = new TH2F(sname.c_str(), stitle.c_str(), 750, -250, 500, 150, -15, 135);
0355     hm->registerHisto(h_bunch_evt_bcodiff[i]);
0356 
0357     // ChipBCO vs Bunch to check the linear correlation
0358     sname = (getHistoPrefix() + "bunch_bco_" + std::to_string(i));
0359     stitle = ("bunch vs BCO " + std::to_string(i));
0360     h_bunch_bco[i] = new TH2F(sname.c_str(), stitle.c_str(), 150, 0, 150, 150, -15, 135);
0361     hm->registerHisto(h_bunch_bco[i]);
0362 
0363     // StrobeBCO - Prev StrobeBCO
0364     // h_bcoprediff[i] = new TH1F((getHistoPrefix()+Form("bcoprediff_%d", i)).c_str(), Form("BCO - PreBCO %d", i), 1000, 0, 1000);
0365   }
0366 
0367   h_bunch_all = new TH1F((getHistoPrefix() + "bunch_all").c_str(), "bunch @ evt all felix", 150, -15, 135);
0368   h_bunch_gl1 = new TH1F((getHistoPrefix() + "bunch_gl1").c_str(), "bunch @ gl1", 150, -15, 135);
0369   hm->registerHisto(h_bunch_all);
0370   hm->registerHisto(h_bunch_gl1);
0371 
0372   h_bcoreco_evt_diff_all = new TH1F((getHistoPrefix() + "h_bcoreco_evt_diff_all").c_str(), "bcoreco evt diff_all", 540, -270, 270);
0373   h_bcointtgl1_diff = new TH1F((getHistoPrefix() + "h_bcointtgl1_diff").c_str(), "bco intt gl1 diff_", 540, -270, 270);
0374   hm->registerHisto(h_bcoreco_evt_diff_all);
0375   hm->registerHisto(h_bcointtgl1_diff);
0376 }