Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-10-17 08:20:45

0001 #include "InttQa.h"
0002 
0003 // #include <ffarawobjects/Gl1Packet.h>
0004 // #include <trackbase/InttEventHeader.h> // why in trackbase and not ffarawobjects??????????
0005 #include <ffarawobjects/InttRawHit.h>
0006 #include <ffarawobjects/InttRawHitContainer.h>
0007 
0008 #include <fun4all/Fun4AllHistoManager.h>  // required by QAHistManagerDef
0009 #include <fun4all/Fun4AllReturnCodes.h>
0010 #include <fun4all/Fun4AllServer.h>
0011 
0012 #include <phool/PHCompositeNode.h>
0013 #include <phool/PHPointerListIterator.h>
0014 #include <phool/getClass.h>
0015 #include <phool/phool.h>  // PHWHERE
0016 
0017 #include <qautils/QAHistManagerDef.h>
0018 
0019 #include <TH1D.h>
0020 #include <TH2D.h>
0021 
0022 #include <format> // -std=c++20
0023 #include <boost/format.hpp> // std::format,vformat are a pain to use with non-const strings
0024 #include <exception>
0025 #include <vector>
0026 
0027 InttQa::InttQa (
0028     std::string const& name
0029 ) : SubsysReco(name) {
0030 
0031     auto hm = QAHistManagerDef::getHistoManager();
0032 
0033     // The numeric values are important for naming the histograms,
0034     // which is why I do not use range-based for loops here
0035     for (int felix_server = 0; felix_server < n_felix_servers; ++felix_server) {
0036 
0037         m_felix_server_hit_distribution[felix_server] = new TH2D (
0038             std::format("{}_hit_distribution_server{:01d}", m_prefix, felix_server).c_str(),
0039             std::format("Hit Distribution intt{:01d}", felix_server).c_str(),
0040             n_chips, -0.5, n_chips-0.5,
0041             n_felix_channels, -0.5, n_felix_channels-0.5
0042         );
0043         hm->registerHisto(m_felix_server_hit_distribution[felix_server]);
0044 
0045         //...
0046 
0047         for (int felix_channel = 0; felix_channel < n_felix_channels; ++felix_channel) {
0048 
0049             m_felix_channel_bco_distribution[felix_server][felix_channel] = new TH1D (
0050                 std::format("{}_bco_distribution_server{:01d}_channel{:02d}", m_prefix, felix_server, felix_channel).c_str(),
0051                 std::format("BCO Distribution intt{:01d}, Felix Channel {:02d}", felix_server, felix_channel).c_str(),
0052                 n_bcos, -0.5, n_bcos-0.5
0053             );
0054             hm->registerHisto(m_felix_channel_bco_distribution[felix_server][felix_channel]);
0055 
0056             m_felix_channel_hit_distribution[felix_server][felix_channel] = new TH2D (
0057                 std::format("{}_hit_distribution_server{:01d}_channel{:02d}", m_prefix, felix_server, felix_channel).c_str(),
0058                 std::format("Hit Distribution intt{:01d}, Felix Channel {:02d}", felix_server, felix_channel).c_str(),
0059                 n_channels, -0.5, n_channels-0.5,
0060                 n_chips, -0.5, n_chips-0.5
0061             );
0062             hm->registerHisto(m_felix_channel_hit_distribution[felix_server][felix_channel]);
0063 
0064             // ...
0065 
0066             for (int chip = 0; chip < n_chips; ++chip) {
0067 
0068                 m_chip_hit_distribution[felix_server][felix_channel][chip] = new TH1D (
0069                     std::format("{}_hit_distribution_server{:01d}_channel{:02d}_chip{:02d}", m_prefix, felix_server, felix_channel, chip).c_str(),
0070                     std::format("Hit Distribution intt{:01d}, Felix Channel {:02d}, Chip{:02d}", felix_server, felix_channel, chip).c_str(),
0071                     n_bcos, -0.5, n_bcos-0.5
0072                 );
0073                 hm->registerHisto(m_chip_hit_distribution[felix_server][felix_channel][chip]);
0074 
0075                 m_chip_adc_distribution[felix_server][felix_channel][chip] = new TH1D (
0076                     std::format("{}_adc_distribution_server{:01d}_channel{:02d}_chip{:02d}", m_prefix, felix_server, felix_channel, chip).c_str(),
0077                     std::format("ADC Distribution intt{:01d}, Felix Channel {:02d}, Chip{:02d}", felix_server, felix_channel, chip).c_str(),
0078                     n_adcs, -0.5, n_adcs-0.5
0079                 );
0080                 hm->registerHisto(m_chip_adc_distribution[felix_server][felix_channel][chip]);
0081 
0082                 //...
0083             }
0084         }
0085     }
0086 }
0087 
0088 int
0089 InttQa::InitRun (
0090     PHCompositeNode* top_node
0091 ) {
0092     // In case this is called multiple times per instance lifetime
0093     m_intt_raw_hit_containers.clear();
0094 
0095     PHNodeIterator dst_itr(top_node);
0096 
0097     PHCompositeNode* intt_node = dynamic_cast<PHCompositeNode*>(dst_itr.findFirst("INTT"));
0098     if (!intt_node) {
0099         if (Verbosity()) {
0100             std::cout
0101                 << PHWHERE
0102                 << "\tNo 'INTT' node found on the NodeTree, Unregistering\n"
0103                 << std::flush;
0104         }
0105         Fun4AllServer::instance()->unregisterSubsystem(this);
0106         return Fun4AllReturnCodes::EVENT_OK;
0107     }
0108     PHNodeIterator intt_itr(intt_node);
0109 
0110     // For the record, if you're reading this, I do not like this loop structure either
0111     // I think this is moot b/c in practice most DST will only have one node, 'INTTRAWHIT'
0112     PHPointerListIterator<PHNode> next_intt_node(intt_itr.ls());
0113     for (PHNode* itr_node; (itr_node = next_intt_node());) {
0114         auto intt_raw_hit_container_node = static_cast<PHIODataNode<InttRawHitContainer>*>(itr_node);
0115         if (!intt_raw_hit_container_node) continue;
0116         if (Verbosity()) {
0117             std::cout
0118                 << PHWHERE
0119                 << "Got InttRawHitContainerNode: '" << intt_raw_hit_container_node->getName() << "'"
0120                 << std::endl;
0121         }
0122 
0123         auto intt_raw_hit_container = dynamic_cast<InttRawHitContainer*>(intt_raw_hit_container_node->getData());
0124         if (!intt_raw_hit_container) continue; 
0125         m_intt_raw_hit_containers.push_back(intt_raw_hit_container);
0126     }
0127 
0128     if (m_intt_raw_hit_containers.empty()) {
0129         if (Verbosity()) {
0130             std::cout
0131                 << PHWHERE
0132                 << "No InttRawHitContainers found on the NodeTree, Unregistering"
0133                 << std::endl;
0134         }
0135         Fun4AllServer::instance()->unregisterSubsystem(this);
0136         return Fun4AllReturnCodes::EVENT_OK;
0137     }
0138 
0139     return Fun4AllReturnCodes::EVENT_OK;
0140 }
0141 
0142 int
0143 InttQa::process_event (
0144     PHCompositeNode* // top_node
0145 ) {
0146     for (auto const& intt_raw_hit_container : m_intt_raw_hit_containers) {
0147         for (unsigned int hit_index{0}; hit_index < intt_raw_hit_container->get_nhits(); ++hit_index) {
0148 
0149             auto hit = intt_raw_hit_container->get_hit(hit_index);
0150 
0151             int felix_server = hit->get_packetid() - 3001; // Only place this literal is used
0152             int felix_channel = hit->get_fee();
0153             int chip = (hit->get_chip_id() + n_chips - 1) % n_chips; // Hardware is base 1 index, Offline is base 0 index
0154             int channel = hit->get_channel_id();
0155 
0156             // Fine for triggered case, for streaming we will need the GL1
0157             // (and a usage of InttOdbcQuery to check if the data is streaming)
0158             // int bco_diff = hit->get_FPHX_BCO() + hit->get_bco() - gl1_bco 
0159 
0160             int bco_diff = (hit->get_FPHX_BCO() - int(hit->get_bco() & 0x7FU) + n_bcos) % n_bcos;
0161             int adc = hit->get_adc();
0162 
0163             m_felix_server_hit_distribution[felix_server]->Fill(chip, felix_channel);
0164             m_felix_channel_bco_distribution[felix_server][felix_channel]->Fill(bco_diff);
0165             m_felix_channel_hit_distribution[felix_server][felix_channel]->Fill(channel, chip);
0166             m_chip_hit_distribution[felix_server][felix_channel][chip]->Fill(channel);
0167             m_chip_adc_distribution[felix_server][felix_channel][chip]->Fill(adc);
0168         }
0169     }
0170 
0171     return Fun4AllReturnCodes::EVENT_OK;
0172 }
0173