File indexing completed on 2025-08-06 08:18:45
0001 #include "MvtxRawHitQA.h"
0002
0003 #include <qautils/QAHistManagerDef.h>
0004 #include <qautils/QAUtil.h>
0005
0006 #include <fun4all/Fun4AllHistoManager.h>
0007 #include <fun4all/Fun4AllReturnCodes.h>
0008 #include <fun4all/Fun4AllServer.h>
0009
0010 #include <phool/PHPointerListIterator.h>
0011 #include <phool/PHCompositeNode.h>
0012 #include <phool/getClass.h>
0013
0014 #include <TMath.h>
0015 #include <TH2.h>
0016 #include <TProfile.h>
0017 #include <TProfile2D.h>
0018 #include <cassert>
0019
0020 #include "TGaxis.h"
0021
0022
0023 MvtxRawHitQA::MvtxRawHitQA(const std::string &name)
0024 : SubsysReco(name)
0025 {
0026 }
0027
0028
0029 int MvtxRawHitQA::InitRun(PHCompositeNode *topNode)
0030 {
0031 createHistos();
0032
0033 TGaxis::SetMaxDigits(3);
0034
0035 PHNodeIterator trkr_itr(topNode);
0036 PHCompositeNode *mvtx_node = dynamic_cast<PHCompositeNode *>(
0037 trkr_itr.findFirst("PHCompositeNode", "MVTX"));
0038 if(!mvtx_node)
0039 {
0040 std::cout << PHWHERE << " No MVTX node found, exit" << std::endl;
0041 Fun4AllServer::instance()->unregisterSubsystem(this);
0042 return Fun4AllReturnCodes::EVENT_OK;
0043 }
0044 PHNodeIterator mvtx_itr(mvtx_node);
0045 PHPointerListIterator<PHNode> iter(mvtx_itr.ls());
0046 PHNode *thisnode;
0047 while((thisnode = iter()))
0048 {
0049 if(thisnode->getType() !="PHIODataNode")
0050 {
0051 continue;
0052 }
0053
0054 if((thisnode->getName()).find("HEADER") != std::string::npos)
0055 {
0056 continue;
0057 }
0058 PHIODataNode<MvtxRawHitContainer> *theNode = static_cast<PHIODataNode<MvtxRawHitContainer> *>(thisnode);
0059 if(theNode)
0060 {
0061 std::cout << PHWHERE << " Found Mvtx Raw hit container node " << theNode->getName() << std::endl;
0062 auto cont = (MvtxRawHitContainer*)theNode->getData();
0063 if(cont)
0064 {
0065 m_rawhit_containers.push_back(cont);
0066 }
0067 }
0068 }
0069
0070 auto hm = QAHistManagerDef::getHistoManager();
0071 assert(hm);
0072
0073 h_nhits_layer0 = dynamic_cast<TH1 *>(hm->getHisto(std::string(getHistoPrefix() + "nhits_layer0").c_str()));
0074 h_nhits_layer1 = dynamic_cast<TH1 *>(hm->getHisto(std::string(getHistoPrefix() + "nhits_layer1").c_str()));
0075 h_nhits_layer2 = dynamic_cast<TH1 *>(hm->getHisto(std::string(getHistoPrefix() + "nhits_layer2").c_str()));
0076
0077 h_bco = dynamic_cast<TH1 *>(hm->getHisto(std::string(getHistoPrefix() + "bco").c_str()));
0078 h_strobe_bc = dynamic_cast<TH1 *>(hm->getHisto(std::string(getHistoPrefix() + "strobe_bc").c_str()));
0079 h_chip_bc = dynamic_cast<TH1 *>(hm->getHisto(std::string(getHistoPrefix() + "chip_bc").c_str()));
0080
0081 h_nhits_stave_chip_layer0 = dynamic_cast<TH2 *>(hm->getHisto(std::string(getHistoPrefix() + "nhits_stave_chip_layer0").c_str()));
0082 h_nhits_stave_chip_layer1 = dynamic_cast<TH2 *>(hm->getHisto(std::string(getHistoPrefix() + "nhits_stave_chip_layer1").c_str()));
0083 h_nhits_stave_chip_layer2 = dynamic_cast<TH2 *>(hm->getHisto(std::string(getHistoPrefix() + "nhits_stave_chip_layer2").c_str()));
0084
0085 return Fun4AllReturnCodes::EVENT_OK;
0086 }
0087
0088
0089 int MvtxRawHitQA::process_event(PHCompositeNode * )
0090 {
0091
0092 std::vector < MvtxRawHit* > hits;
0093 std::vector < uint64_t > bcos;
0094 std::vector < uint32_t > strobe_bcs;
0095 std::vector < uint32_t > chip_bcs;
0096 std::vector < uint8_t > layers;
0097 std::vector < uint8_t > staves;
0098 std::vector < uint8_t > chips;
0099 std::vector < uint16_t > rows;
0100 std::vector < uint16_t > cols;
0101
0102 hits.clear();
0103 bcos.clear();
0104 strobe_bcs.clear();
0105 chip_bcs.clear();
0106 layers.clear();
0107 staves.clear();
0108 chips.clear();
0109 rows.clear();
0110 cols.clear();
0111
0112 unsigned int raw_hit_num = 0;
0113 unsigned int total_raw_hit_num = 0;
0114 int nhit_layer0=0;
0115 int nhit_layer1=0;
0116 int nhit_layer2=0;
0117
0118 for(auto& rawhitcont : m_rawhit_containers)
0119 {
0120 if (rawhitcont)
0121 {
0122 raw_hit_num = rawhitcont->get_nhits();
0123 total_raw_hit_num += raw_hit_num;
0124
0125 for (unsigned int i = 0; i < raw_hit_num; i++)
0126 {
0127 auto hit = rawhitcont->get_hit(i);
0128 auto bco = hit->get_bco();
0129 auto strobe_bc = hit->get_strobe_bc();
0130 auto chip_bc = hit->get_chip_bc();
0131 auto layer = hit->get_layer_id();
0132 auto stave = hit->get_stave_id();
0133 auto chip = hit->get_chip_id();
0134 auto row = hit->get_row();
0135 auto col = hit->get_col();
0136 hits.push_back( hit );
0137 bcos.push_back( bco );
0138 strobe_bcs.push_back( strobe_bc );
0139 chip_bcs.push_back( chip_bc );
0140 layers.push_back( layer );
0141 staves.push_back( stave );
0142 chips.push_back( chip );
0143 rows.push_back( row );
0144 cols.push_back( col );
0145 }
0146 }
0147
0148 }
0149
0150
0151 if( total_raw_hit_num == 0 )
0152 {
0153 return Fun4AllReturnCodes::EVENT_OK;
0154 }
0155
0156 for (int i=0; i<(int)total_raw_hit_num; i++)
0157 {
0158 if (layers[i]==0)
0159 {
0160 nhit_layer0++;
0161 h_nhits_stave_chip_layer0->Fill(chips[i],staves[i]);
0162 }
0163 if (layers[i]==1)
0164 {
0165 nhit_layer1++;
0166 h_nhits_stave_chip_layer1->Fill(chips[i],staves[i]);
0167 }
0168 if (layers[i]==2)
0169 {
0170 nhit_layer2++;
0171 h_nhits_stave_chip_layer2->Fill(chips[i],staves[i]);
0172 }
0173 }
0174
0175 h_nhits_layer0->Fill(nhit_layer0);
0176 h_nhits_layer1->Fill(nhit_layer1);
0177 h_nhits_layer2->Fill(nhit_layer2);
0178
0179 h_bco->Fill(bcos[0]);
0180 h_strobe_bc->Fill(strobe_bcs[0]);
0181 h_chip_bc->Fill(chip_bcs[0]);
0182
0183 return Fun4AllReturnCodes::EVENT_OK;
0184 }
0185
0186
0187 int MvtxRawHitQA::EndRun(const int )
0188 {
0189 auto hm = QAHistManagerDef::getHistoManager();
0190 assert(hm);
0191
0192 return Fun4AllReturnCodes::EVENT_OK;
0193 }
0194
0195
0196 int MvtxRawHitQA::End(PHCompositeNode * ) { return Fun4AllReturnCodes::EVENT_OK; }
0197
0198 std::string MvtxRawHitQA::getHistoPrefix() const { return std::string("h_") + Name() + std::string("_"); }
0199
0200 void MvtxRawHitQA::createHistos()
0201 {
0202 auto hm = QAHistManagerDef::getHistoManager();
0203 assert(hm);
0204
0205 {
0206 auto h = new TH1F(std::string(getHistoPrefix() + "nhits_layer0").c_str(), "Number of hits in layer 0;Number of hits;Entries",100,0,25000);
0207 hm->registerHisto(h);
0208 }
0209
0210 {
0211 auto h = new TH1F(std::string(getHistoPrefix() + "nhits_layer1").c_str(), "Number of hits in layer 1;Number of hits;Entries",100,0,25000);
0212 hm->registerHisto(h);
0213 }
0214
0215 {
0216 auto h = new TH1F(std::string(getHistoPrefix() + "nhits_layer2").c_str(), "Number of hits in layer 2;Number of hits;Entries",100,0,25000);
0217 hm->registerHisto(h);
0218 }
0219
0220 {
0221 auto h = new TH1F(std::string(getHistoPrefix() + "bco").c_str(), "BCO distribution;BCO;Entries",100,0,TMath::Power( 2, 40 ));
0222 hm->registerHisto(h);
0223 }
0224
0225 {
0226 auto h = new TH1F(std::string(getHistoPrefix() + "strobe_bc").c_str(), "Strobe BC distribution;Strobe BC;Entries",100,0,4000);
0227 hm->registerHisto(h);
0228 }
0229
0230 {
0231 auto h = new TH1F(std::string(getHistoPrefix() + "chip_bc").c_str(), "Chip BC distribution;Chip BC;Entries",100,0,500);
0232 hm->registerHisto(h);
0233 }
0234
0235 {
0236 auto h = new TH2F(std::string(getHistoPrefix() + "nhits_stave_chip_layer0").c_str(), "Hitmap in layer 0;ChipID;StaveID;Number of hits",9,0,9,12,0,12);
0237 hm->registerHisto(h);
0238 }
0239
0240 {
0241 auto h = new TH2F(std::string(getHistoPrefix() + "nhits_stave_chip_layer1").c_str(), "Hitmap in layer 1;ChipID;StaveID;Number of hits",9,0,9,16,0,16);
0242 hm->registerHisto(h);
0243 }
0244
0245 {
0246 auto h = new TH2F(std::string(getHistoPrefix() + "nhits_stave_chip_layer2").c_str(), "Hitmap in layer 2;ChipID;StaveID;Number of hits",9,0,9,20,0,20);
0247 hm->registerHisto(h);
0248 }
0249 }