File indexing completed on 2025-12-17 09:21:17
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 <fun4all/Fun4AllHistoManager.h> // for Fun4AllHistoManager
0010 #include <fun4all/Fun4AllReturnCodes.h> // for EVENT_OK, ABORTEVENT
0011 #include <fun4all/Fun4AllServer.h>
0012
0013 #include <phool/PHCompositeNode.h>
0014 #include <phool/PHPointerListIterator.h>
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::InitRun(PHCompositeNode *topNode)
0052 {
0053 createHistos();
0054
0055
0056
0057
0058 PHNodeIterator trkr_itr(topNode);
0059 PHCompositeNode *intt_node = dynamic_cast<PHCompositeNode *>(trkr_itr.findFirst("PHCompositeNode", "INTT"));
0060 if (!intt_node)
0061 {
0062 std::cout << PHWHERE << " No INTT node found, unregistering subsystem" << std::endl;
0063 Fun4AllServer::instance()->unregisterSubsystem(this);
0064 return Fun4AllReturnCodes::EVENT_OK;
0065 }
0066 PHNodeIterator intt_itr(intt_node);
0067 PHPointerListIterator<PHNode> iter(intt_itr.ls());
0068 PHNode *thisnode;
0069 while ((thisnode = iter()))
0070 {
0071 if (thisnode->getType() != "PHIODataNode")
0072 {
0073 continue;
0074 }
0075 PHIODataNode<InttRawHitContainer> *theNode = static_cast<PHIODataNode<InttRawHitContainer> *>(thisnode);
0076 if (theNode)
0077 {
0078 std::cout << PHWHERE << " Found INTT Raw hit container node " << theNode->getName() << std::endl;
0079
0080 if (theNode->getName() == "INTTEVENTHEADER" || thisnode->getName().find("G4") != std::string::npos)
0081 {
0082 continue;
0083 }
0084 auto *cont = (InttRawHitContainer *) theNode->getData();
0085 if (cont)
0086 {
0087 m_rawhit_containers.push_back(cont);
0088 }
0089 }
0090 }
0091 auto *hm = QAHistManagerDef::getHistoManager();
0092 assert(hm);
0093
0094 for (int felix = 0; felix < InttQa::kFelix_num; felix++)
0095 {
0096 std::string name = getHistoPrefix() + "intt" + std::to_string(felix);
0097 hist_fee_chip_chan_[felix] = dynamic_cast<TH3 *>(hm->getHisto(name));
0098
0099 std::string name_bco_event = name + "_ladder_bco_full_event_counter";
0100 hist_fee_bco_full_event_counter_[felix] = dynamic_cast<TH3 *>(hm->getHisto(name_bco_event));
0101
0102 std::string name_bco_event_diff = name + "_ladder_bco_full_event_counter_diff";
0103 hist_fee_bco_full_event_counter_diff_[felix] = dynamic_cast<TH3 *>(hm->getHisto(name_bco_event_diff));
0104
0105 std::string name_event_counter = name + "_event_counter";
0106 hist_event_counter_[felix] = dynamic_cast<TH1 *>(hm->getHisto(name_event_counter));
0107
0108 std::string name_event_counter_diff = name + "_event_counter_diff";
0109 hist_event_counter_diff_[felix] = dynamic_cast<TH1 *>(hm->getHisto(name_event_counter_diff));
0110 }
0111
0112 for (int felix = 0; felix < InttQa::kFelix_num; felix++)
0113 {
0114 for (int ladder = 0; ladder < InttQa::kFee_num; ladder++)
0115 {
0116 std::string name = getHistoPrefix() + "intt" + std::to_string(felix) + "_" + std::to_string(ladder);
0117 hist_hitmap_[felix][ladder] = dynamic_cast<TH2 *>(hm->getHisto(name));
0118 }
0119 }
0120
0121 hist_nhit_ = dynamic_cast<TH1 *>(hm->getHisto(std::string(getHistoPrefix() + "nhit")));
0122
0123 hist_nhit_south_ = dynamic_cast<TH1 *>(hm->getHisto(std::string(getHistoPrefix() + "nhit_south")));
0124 hist_nhit_north_ = dynamic_cast<TH1 *>(hm->getHisto(std::string(getHistoPrefix() + "nhit_north")));
0125 hist_pid_ = dynamic_cast<TH1 *>(hm->getHisto(std::string(getHistoPrefix() + "pid")));
0126 hist_adc_ = dynamic_cast<TH1 *>(hm->getHisto(std::string(getHistoPrefix() + "adc")));
0127 hist_bco_ = dynamic_cast<TH1 *>(hm->getHisto(std::string(getHistoPrefix() + "bco")));
0128 hist_bco_full_ = dynamic_cast<TH1 *>(hm->getHisto(std::string(getHistoPrefix() + "bco_full")));
0129
0130 return Fun4AllReturnCodes::EVENT_OK;
0131 }
0132
0133 void InttRawHitQA::createHistos()
0134 {
0135 auto *hm = QAHistManagerDef::getHistoManager();
0136 assert(hm);
0137
0138 for (int felix = 0; felix < InttQa::kFelix_num; felix++)
0139 {
0140 std::string name = getHistoPrefix() + "intt" + std::to_string(felix);
0141 std::string title = name + ";FELIX CH;Chip;Channel;Entries";
0142 {
0143 auto *h = new TH3D(name.c_str(), title.c_str(),
0144 InttQa::kFee_num, 0, InttQa::kFee_num,
0145 InttQa::kChip_num, 1, InttQa::kChip_num + 1,
0146 InttQa::kChan_num, 0, InttQa::kChan_num);
0147 hm->registerHisto(h);
0148 }
0149
0150 std::string name_bco_event = name + "_ladder_bco_full_event_counter";
0151 std::string title_bco_event = name + ";FELIX_CH;BCO full;Event Counter;Entries";
0152 {
0153 auto *h = new TH3D(name_bco_event.c_str(), title_bco_event.c_str(),
0154 InttQa::kFee_num, 0, InttQa::kFee_num,
0155 100, 0, std::pow(2, 40),
0156 1e4, 0, 1e7);
0157 hm->registerHisto(h);
0158 }
0159
0160 std::string name_bco_event_diff = name + "_ladder_bco_full_event_counter_diff";
0161 std::string title_bco_event_diff = name + ";FELIX_CH;#Delta BCO full;#Delta Event Counter;Entries";
0162 int max = 1000;
0163 {
0164 auto *h = new TH3D(name_bco_event_diff.c_str(), title_bco_event_diff.c_str(),
0165 InttQa::kFee_num, 0, InttQa::kFee_num,
0166 2 * max / 100, -max, max,
0167 2 * max / 100, -max, max);
0168 hm->registerHisto(h);
0169 }
0170
0171 std::string name_event_counter = name + "_event_counter";
0172 std::string title_event_counter = name_event_counter + ";Event Counter;Entries";
0173 {
0174 auto *h = new TH1D(name_event_counter.c_str(), title_event_counter.c_str(), 1e4, 0, 1e7);
0175 hm->registerHisto(h);
0176 }
0177
0178 std::string name_event_counter_diff = name + "_event_counter_diff";
0179 std::string title_event_counter_diff = name_event_counter_diff + ";Event Counter;Entries";
0180 {
0181 auto *h = new TH1D(name_event_counter_diff.c_str(), title_event_counter_diff.c_str(), 2 * max / 100, -max, max);
0182 hm->registerHisto(h);
0183 }
0184 }
0185
0186 for (int felix = 0; felix < InttQa::kFelix_num; felix++)
0187 {
0188 for (int ladder = 0; ladder < InttQa::kFee_num; ladder++)
0189 {
0190 std::string name = getHistoPrefix() + "intt" + std::to_string(felix) + "_" + std::to_string(ladder);
0191 std::string title = name + ";Chip;Channel;Entries";
0192 auto *h = new TH2I(name.c_str(), title.c_str(),
0193 InttQa::kChip_num, 1, InttQa::kChip_num + 1,
0194 InttQa::kChan_num, 0, InttQa::kChan_num);
0195 hm->registerHisto(h);
0196 }
0197 }
0198
0199 {
0200 auto *h = new TH1D(std::string(getHistoPrefix() + "nhit").c_str(), "#INTTRAWHIT per event;#hit;Entries", 1e4, 0, 1e4);
0201 InttQa::HistConfig(h);
0202 hm->registerHisto(h);
0203 }
0204
0205 {
0206 auto *h = new TH1D(std::string(getHistoPrefix() + "nhit_south").c_str(), "#INTTRAWHIT South;event;#hit", 1e4, 0, 1e7);
0207 InttQa::HistConfig(h);
0208 hm->registerHisto(h);
0209 }
0210
0211 {
0212 auto *h = new TH1D(std::string(getHistoPrefix() + "nhit_north").c_str(), "#INTTRAWHIT North;event;#hit", 1e4, 0, 1e7);
0213 InttQa::HistConfig(h);
0214 hm->registerHisto(h);
0215 }
0216
0217 {
0218 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);
0219 InttQa::HistConfig(h);
0220 hm->registerHisto(h);
0221 }
0222
0223 {
0224 auto *h = new TH1D(std::string(getHistoPrefix() + "adc").c_str(), "ADC distribution;ADC;Entries", 8, 0, 8);
0225 InttQa::HistConfig(h);
0226 hm->registerHisto(h);
0227 }
0228
0229 {
0230 auto *h = new TH1D(std::string(getHistoPrefix() + "bco").c_str(), "BCO distribution;BCO;Entries", InttQa::kBco_max + 10, -5, InttQa::kBco_max + 5);
0231 InttQa::HistConfig(h);
0232 hm->registerHisto(h);
0233 }
0234
0235 {
0236 auto *h = new TH1D(std::string(getHistoPrefix() + "bco_full").c_str(), "BCO full distribution;BCO full;Entries", 100, 0, std::pow(2, 40));
0237 InttQa::HistConfig(h);
0238 hm->registerHisto(h);
0239 }
0240 }
0241
0242 int InttRawHitQA::process_event(PHCompositeNode * )
0243 {
0244 for (auto &cont : m_rawhit_containers)
0245 {
0246 auto hits = GetHits(cont);
0247
0248 auto raw_hit_num = hits.size();
0249 hist_nhit_->Fill(raw_hit_num);
0250
0251
0252 if (raw_hit_num == 0 || cont == nullptr)
0253 {
0254 return Fun4AllReturnCodes::EVENT_OK;
0255 }
0256
0257 event_counter_by_myself_++;
0258
0259
0260
0261
0262 uint64_t bco_full = (cont->get_hit(0)->get_bco());
0263 hist_bco_full_->Fill(bco_full);
0264
0265
0266
0267
0268 uint32_t event_counter_ref = cont->get_hit(0)->get_event_counter();
0269
0270
0271
0272
0273
0274 for (auto *hit : hits)
0275 {
0276 int felix = hit->get_packetid() - InttQa::kFirst_pid;
0277 int felix_ch = hit->get_fee();
0278
0279
0280 int chip = hit->get_chip_id();
0281 if (chip > InttQa::kChip_num)
0282 {
0283 chip = chip - InttQa::kChip_num;
0284 }
0285
0286 int chan = hit->get_channel_id();
0287 auto adc = hit->get_adc();
0288 auto bco = hit->get_FPHX_BCO();
0289 int event_counter = hit->get_event_counter();
0290 if (is_first_event_ == true)
0291 {
0292 event_counter = 0;
0293 }
0294 else if (event_counter - previous_event_counter_ > 1000)
0295 {
0296 event_counter = -1;
0297 }
0298 else
0299 {
0300 last_event_counter_ = event_counter;
0301 }
0302
0303
0304
0305
0306
0307
0308
0309
0310
0311
0312
0313
0314
0315 hist_fee_chip_chan_[felix]->Fill(felix_ch, chip, chan);
0316
0317 hist_hitmap_[felix][felix_ch]->Fill(chip, chan);
0318
0319 hist_pid_->AddBinContent(felix + 1);
0320
0321 hist_adc_->Fill(adc);
0322 hist_bco_->Fill(bco);
0323
0324 hist_fee_bco_full_event_counter_[felix]
0325 ->Fill(felix_ch,
0326 hit->get_bco(),
0327 event_counter);
0328
0329 hist_fee_bco_full_event_counter_diff_[felix]
0330 ->Fill(felix_ch,
0331 hit->get_bco() - bco_full,
0332 event_counter - event_counter_ref);
0333
0334 hist_event_counter_[felix]
0335 ->Fill(event_counter);
0336
0337 hist_event_counter_diff_[felix]
0338 ->Fill(event_counter - event_counter_ref);
0339
0340 if (felix < 4)
0341 {
0342 hist_nhit_south_->Fill(event_counter);
0343 }
0344 else
0345 {
0346 hist_nhit_north_->Fill(event_counter);
0347 }
0348 }
0349
0350 is_first_event_ = false;
0351
0352 if (last_event_counter_ - previous_event_counter_ < 1000)
0353 {
0354 previous_event_counter_ = last_event_counter_;
0355 }
0356 else
0357 {
0358 previous_event_counter_ = -1;
0359 }
0360 }
0361
0362 return Fun4AllReturnCodes::EVENT_OK;
0363 }
0364
0365 int InttRawHitQA::End(PHCompositeNode * )
0366 {
0367 auto *hm = QAHistManagerDef::getHistoManager();
0368 assert(hm);
0369
0370 return Fun4AllReturnCodes::EVENT_OK;
0371 }
0372
0373 std::string InttRawHitQA::getHistoPrefix() const
0374 {
0375 return std::string("h_") + Name() + std::string("_");
0376 }