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