Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-12-17 09:21:21

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