Back to home page

sPhenix code displayed by LXR

 
 

    


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

0001 #include "InttRawHitQA.h"
0002 #include "InttQaCommon.h"  // for HistConfig, kFee_num
0003 
0004 #include <qautils/QAHistManagerDef.h>
0005 
0006 #include <ffarawobjects/InttRawHit.h>           // for InttRawHit
0007 #include <ffarawobjects/InttRawHitContainer.h>  // for InttRawHitContainer
0008 
0009 #include <phool/PHCompositeNode.h>
0010 #include <phool/PHPointerListIterator.h>
0011 #include <fun4all/Fun4AllHistoManager.h>  // for Fun4AllHistoManager
0012 #include <fun4all/Fun4AllReturnCodes.h>   // for EVENT_OK, ABORTEVENT
0013 #include <fun4all/Fun4AllServer.h>
0014 
0015 #include <phool/getClass.h>  // for getClass
0016 #include <phool/phool.h>     // for PHWHERE
0017 
0018 #include <TH1.h>
0019 #include <TH2.h>
0020 #include <TH3.h>
0021 
0022 #include <cassert>
0023 #include <cmath>
0024 #include <cstdint>
0025 #include <iostream>
0026 #include <string>
0027 
0028 //____________________________________________________________________________..
0029 InttRawHitQA::InttRawHitQA(const std::string &name)
0030   : SubsysReco(name)
0031 {
0032 }
0033 
0034 std::vector<InttRawHit *> InttRawHitQA::GetHits(InttRawHitContainer* container)
0035 {
0036   std::vector<InttRawHit *> hits;
0037   if(container == nullptr)
0038   {
0039     return hits;
0040   }
0041   auto raw_hit_num = container->get_nhits();
0042   for (unsigned int i = 0; i < raw_hit_num; i++)
0043   {
0044     auto hit = container->get_hit(i);
0045     hits.push_back(hit);
0046   }
0047 
0048   return hits;
0049 }
0050 
0051 int InttRawHitQA::Init(PHCompositeNode * /*unused*/)
0052 {
0053   return Fun4AllReturnCodes::EVENT_OK;
0054 }
0055 
0056 int InttRawHitQA::InitRun(PHCompositeNode *topNode)
0057 {
0058   createHistos();
0059 
0060   /////////////////////////////////////////////////////////////////////////
0061   // INTT raw hit
0062   /////////////////////////////////////////////////////////////////////////
0063   PHNodeIterator trkr_itr(topNode);
0064   PHCompositeNode *intt_node = dynamic_cast<PHCompositeNode *>(
0065       trkr_itr.findFirst("PHCompositeNode", "INTT"));  
0066   if(!intt_node)
0067   {
0068     std::cout << PHWHERE << " No INTT node found, unregistering subsystem" << std::endl;
0069     Fun4AllServer::instance()->unregisterSubsystem(this);
0070     return Fun4AllReturnCodes::EVENT_OK;
0071   }
0072   PHNodeIterator intt_itr(intt_node);
0073   PHPointerListIterator<PHNode> iter(intt_itr.ls());
0074   PHNode *thisnode;
0075   while((thisnode = iter()))
0076   {
0077     if(thisnode->getType() !="PHIODataNode")
0078     {
0079       continue;
0080     }
0081     PHIODataNode<InttRawHitContainer> *theNode = static_cast<PHIODataNode<InttRawHitContainer> *>(thisnode);
0082     if(theNode)
0083     {
0084       std::cout << PHWHERE << " Found INTT Raw hit container node " << theNode->getName() << std::endl;
0085       // skip if theNode->getName() is INTTEVENTHEADER
0086       if(theNode->getName() == "INTTEVENTHEADER")
0087       {
0088         continue;
0089       }
0090       auto cont = (InttRawHitContainer*)theNode->getData();
0091       if(cont)
0092       {
0093         m_rawhit_containers.push_back(cont);
0094       }
0095     }
0096   }
0097   auto hm = QAHistManagerDef::getHistoManager();
0098   assert(hm);
0099 
0100   for (int felix = 0; felix < InttQa::kFelix_num; felix++)
0101   {
0102     std::string name = getHistoPrefix() + "intt" + std::to_string(felix);
0103     hist_fee_chip_chan_[felix] = dynamic_cast<TH3D *>(hm->getHisto(name.c_str()));
0104 
0105     std::string name_bco_event = name + "_ladder_bco_full_event_counter";
0106     hist_fee_bco_full_event_counter_[felix] = dynamic_cast<TH3D *>(hm->getHisto(name_bco_event.c_str()));
0107 
0108     std::string name_bco_event_diff = name + "_ladder_bco_full_event_counter_diff";
0109     hist_fee_bco_full_event_counter_diff_[felix] = dynamic_cast<TH3D *>(hm->getHisto(name_bco_event_diff.c_str()));
0110 
0111     std::string name_event_counter = name + "_event_counter";
0112     hist_event_counter_[felix] = dynamic_cast<TH1D *>(hm->getHisto(name_event_counter.c_str()));
0113 
0114     std::string name_event_counter_diff = name + "_event_counter_diff";
0115     hist_event_counter_diff_[felix] = dynamic_cast<TH1D *>(hm->getHisto(name_event_counter_diff.c_str()));
0116   }
0117 
0118   for (int felix = 0; felix < InttQa::kFelix_num; felix++)
0119   {
0120     for (int ladder = 0; ladder < InttQa::kFee_num; ladder++)
0121     {
0122       std::string name = getHistoPrefix() + "intt" + std::to_string(felix) + "_" + std::to_string(ladder);
0123       hist_hitmap_[felix][ladder] = dynamic_cast<TH2I *>(hm->getHisto(name.c_str()));
0124     }
0125   }
0126 
0127   hist_nhit_ = dynamic_cast<TH1D *>(hm->getHisto(std::string(getHistoPrefix() + "nhit").c_str()));
0128 
0129   hist_nhit_south_ = dynamic_cast<TH1D *>(hm->getHisto(std::string(getHistoPrefix() + "nhit_south").c_str()));
0130   hist_nhit_north_ = dynamic_cast<TH1D *>(hm->getHisto(std::string(getHistoPrefix() + "nhit_north").c_str()));
0131   hist_pid_ = dynamic_cast<TH1D *>(hm->getHisto(std::string(getHistoPrefix() + "pid").c_str()));
0132   hist_adc_ = dynamic_cast<TH1D *>(hm->getHisto(std::string(getHistoPrefix() + "adc").c_str()));
0133   hist_bco_ = dynamic_cast<TH1D *>(hm->getHisto(std::string(getHistoPrefix() + "bco").c_str()));
0134   hist_bco_full_ = dynamic_cast<TH1D *>(hm->getHisto(std::string(getHistoPrefix() + "bco_full").c_str()));
0135 
0136   return Fun4AllReturnCodes::EVENT_OK;
0137 }
0138 
0139 void InttRawHitQA::createHistos()
0140 {
0141   auto hm = QAHistManagerDef::getHistoManager();
0142   assert(hm);
0143 
0144   for (int felix = 0; felix < InttQa::kFelix_num; felix++)
0145   {
0146     std::string name = getHistoPrefix() + "intt" + std::to_string(felix);
0147     std::string title = name + ";FELIX CH;Chip;Channel;Entries";
0148     {
0149       auto h = new TH3D(name.c_str(), title.c_str(),
0150                         InttQa::kFee_num, 0, InttQa::kFee_num,
0151                         InttQa::kChip_num, 1, InttQa::kChip_num + 1,
0152                         InttQa::kChan_num, 0, InttQa::kChan_num);
0153       hm->registerHisto(h);
0154     }
0155 
0156     std::string name_bco_event = name + "_ladder_bco_full_event_counter";
0157     std::string title_bco_event = name + ";FELIX_CH;BCO full;Event Counter;Entries";
0158     {
0159       auto h = new TH3D(name_bco_event.c_str(), title_bco_event.c_str(),
0160                         InttQa::kFee_num, 0, InttQa::kFee_num,
0161                         100, 0, std::pow(2, 40),
0162                         1e4, 0, 1e7);
0163       hm->registerHisto(h);
0164     }
0165 
0166     std::string name_bco_event_diff = name + "_ladder_bco_full_event_counter_diff";
0167     std::string title_bco_event_diff = name + ";FELIX_CH;#Delta BCO full;#Delta Event Counter;Entries";
0168     int max = 1000;
0169     {
0170       auto h = new TH3D(name_bco_event_diff.c_str(), title_bco_event_diff.c_str(),
0171                         InttQa::kFee_num, 0, InttQa::kFee_num,
0172                         2 * max / 100, -max, max,
0173                         2 * max / 100, -max, max);
0174       hm->registerHisto(h);
0175     }
0176 
0177     std::string name_event_counter = name + "_event_counter";
0178     std::string title_event_counter = name_event_counter + ";Event Counter;Entries";
0179     {
0180       auto h = new TH1D(name_event_counter.c_str(), title_event_counter.c_str(), 1e4, 0, 1e7);
0181       hm->registerHisto(h);
0182     }
0183 
0184     std::string name_event_counter_diff = name + "_event_counter_diff";
0185     std::string title_event_counter_diff = name_event_counter_diff + ";Event Counter;Entries";
0186     {
0187       auto h = new TH1D(name_event_counter_diff.c_str(), title_event_counter_diff.c_str(), 2 * max / 100, -max, max);
0188       hm->registerHisto(h);
0189     }
0190   }
0191 
0192   for (int felix = 0; felix < InttQa::kFelix_num; felix++)
0193   {
0194     for (int ladder = 0; ladder < InttQa::kFee_num; ladder++)
0195     {
0196       std::string name = getHistoPrefix() + "intt" + std::to_string(felix) + "_" + std::to_string(ladder);
0197       std::string title = name + ";Chip;Channel;Entries";
0198       auto h = new TH2I(name.c_str(), title.c_str(),
0199                         InttQa::kChip_num, 1, InttQa::kChip_num + 1,
0200                         InttQa::kChan_num, 0, InttQa::kChan_num);
0201       hm->registerHisto(h);
0202     }
0203   }
0204 
0205   {
0206     auto h = new TH1D(std::string(getHistoPrefix() + "nhit").c_str(), "#INTTRAWHIT per event;#hit;Entries", 1e4, 0, 1e4);
0207     InttQa::HistConfig(h);
0208     hm->registerHisto(h);
0209   }
0210 
0211   {
0212     auto h = new TH1D(std::string(getHistoPrefix() + "nhit_south").c_str(), "#INTTRAWHIT South;event;#hit", 1e4, 0, 1e7);
0213     InttQa::HistConfig(h);
0214     hm->registerHisto(h);
0215   }
0216 
0217   {
0218     auto h = new TH1D(std::string(getHistoPrefix() + "nhit_north").c_str(), "#INTTRAWHIT North;event;#hit", 1e4, 0, 1e7);
0219     InttQa::HistConfig(h);
0220     hm->registerHisto(h);
0221   }
0222 
0223   {
0224     auto h = new TH1D(std::string(getHistoPrefix() + "pid").c_str(), "Packet ID distribution;pid;Entries", InttQa::kFelix_num, InttQa::kFirst_pid, InttQa::kFirst_pid + InttQa::kFelix_num);
0225     InttQa::HistConfig(h);
0226     hm->registerHisto(h);
0227   }
0228 
0229   {
0230     auto h = new TH1D(std::string(getHistoPrefix() + "adc").c_str(), "ADC distribution;ADC;Entries", 8, 0, 8);
0231     InttQa::HistConfig(h);
0232     hm->registerHisto(h);
0233   }
0234 
0235   {
0236     auto h = new TH1D(std::string(getHistoPrefix() + "bco").c_str(), "BCO distribution;BCO;Entries", InttQa::kBco_max + 10, -5, InttQa::kBco_max + 5);
0237     InttQa::HistConfig(h);
0238     hm->registerHisto(h);
0239   }
0240 
0241   {
0242     auto h = new TH1D(std::string(getHistoPrefix() + "bco_full").c_str(), "BCO full distribution;BCO full;Entries", 100, 0, std::pow(2, 40));
0243     InttQa::HistConfig(h);
0244     hm->registerHisto(h);
0245   }
0246 }
0247 
0248 int InttRawHitQA::process_event(PHCompositeNode * /*unused*/)
0249 {
0250   for(auto& cont : m_rawhit_containers)
0251   {
0252    
0253   auto hits = GetHits(cont);
0254 
0255   auto raw_hit_num = hits.size();
0256   hist_nhit_->Fill(raw_hit_num);
0257 
0258   // if no raw hit is found, skip this event
0259   if (raw_hit_num == 0 || cont == nullptr)
0260   {
0261     return Fun4AllReturnCodes::EVENT_OK;
0262   }
0263 
0264   event_counter_by_myself_++;
0265 
0266   //////////////////////////////////////////////////////////////////
0267   // processes for each event                                     //
0268   //////////////////////////////////////////////////////////////////
0269   uint64_t bco_full = (cont->get_hit(0)->get_bco());
0270   hist_bco_full_->Fill(bco_full);
0271 
0272   //////////////////////////////////////////////////////////////////
0273   // primary raw hit sweep to get some reference values           //
0274   //////////////////////////////////////////////////////////////////
0275   uint32_t event_counter_ref = cont->get_hit(0)->get_event_counter();
0276 
0277   //////////////////////////////////////////////////////////////////
0278   // processes for each raw hit                                   //
0279   //////////////////////////////////////////////////////////////////
0280   // loop over all raw hits
0281   for (auto hit : hits)
0282   {
0283     int felix = hit->get_packetid() - InttQa::kFirst_pid;
0284     int felix_ch = hit->get_fee();
0285 
0286     // uint16_t InttRawHit::get_chip_id
0287     int chip = hit->get_chip_id();
0288     if (chip > InttQa::kChip_num)
0289     {
0290       chip = chip - InttQa::kChip_num;
0291     }
0292 
0293     int chan = hit->get_channel_id();
0294     auto adc = hit->get_adc();
0295     auto bco = hit->get_FPHX_BCO();
0296     int event_counter = hit->get_event_counter();  // uint32_t IttRawHit::get_event_counter()
0297     if (is_first_event_ == true)
0298     {
0299       event_counter = 0;
0300     }
0301     else if (event_counter - previous_event_counter_ > 1000)
0302     {
0303       event_counter = -1;  // it means bad
0304     }
0305     else
0306     {
0307       last_event_counter_ = event_counter;
0308     }
0309 
0310     /*
0311         int bco_diff = 0;
0312         if( (bco_full & 0x7f )  > bco )
0313           bco_diff = int(bco_full & 0x7f ) - bco;
0314         else
0315           bco_diff = int(bco_full & 0x7f ) + (128 - bco);
0316     */
0317 
0318     //////////////////////////////////////////////////////////////////
0319     // Filling hists                                                //
0320     //////////////////////////////////////////////////////////////////
0321 
0322     hist_fee_chip_chan_[felix]->Fill(felix_ch, chip, chan);
0323 
0324     hist_hitmap_[felix][felix_ch]->Fill(chip, chan);
0325 
0326     hist_pid_->AddBinContent(felix + 1);
0327 
0328     hist_adc_->Fill(adc);
0329     hist_bco_->Fill(bco);
0330 
0331     hist_fee_bco_full_event_counter_[felix]
0332         ->Fill(felix_ch,  // chip, chan );
0333                hit->get_bco(),
0334                event_counter);
0335 
0336     hist_fee_bco_full_event_counter_diff_[felix]
0337         ->Fill(felix_ch,
0338                hit->get_bco() - bco_full,
0339                event_counter - event_counter_ref);
0340 
0341     hist_event_counter_[felix]
0342         ->Fill(event_counter);
0343 
0344     hist_event_counter_diff_[felix]
0345         ->Fill(event_counter - event_counter_ref);
0346 
0347     if (felix < 4)
0348     {
0349       hist_nhit_south_->Fill(event_counter);
0350     }
0351     else
0352     {
0353       hist_nhit_north_->Fill(event_counter);
0354     }
0355   }
0356 
0357   is_first_event_ = false;
0358 
0359   if (last_event_counter_ - previous_event_counter_ < 1000)
0360   {
0361     previous_event_counter_ = last_event_counter_;  // in the case of reasonable event counter
0362   }
0363   else
0364   {
0365     previous_event_counter_ = -1;  // in the case of a crazy event counter
0366   }
0367   }
0368   // cout << "-------------------------------------------------" << std::endl;
0369   return Fun4AllReturnCodes::EVENT_OK;
0370 }
0371 
0372 int InttRawHitQA::ResetEvent(PHCompositeNode * /*unused*/)
0373 {
0374   // Intitialize for Clone hit counter
0375   return Fun4AllReturnCodes::EVENT_OK;
0376 }
0377 
0378 int InttRawHitQA::EndRun(const int /*unused*/)
0379 {
0380   auto hm = QAHistManagerDef::getHistoManager();
0381   assert(hm);
0382 
0383   return Fun4AllReturnCodes::EVENT_OK;
0384 }
0385 
0386 int InttRawHitQA::End(PHCompositeNode * /*unused*/)
0387 {
0388   return Fun4AllReturnCodes::EVENT_OK;
0389 }
0390 
0391 std::string InttRawHitQA::getHistoPrefix() const { return std::string("h_") + Name() + std::string("_"); }
0392 
0393 int InttRawHitQA::Reset(PHCompositeNode * /*unused*/)
0394 {
0395   return Fun4AllReturnCodes::EVENT_OK;
0396 }