Back to home page

sPhenix code displayed by LXR

 
 

    


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

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 <intt/InttMapping.h>
0009 #include <intt/InttOdbcQuery.h>
0010 
0011 #include <fun4all/Fun4AllHistoManager.h>  // required by QAHistManagerDef
0012 #include <fun4all/Fun4AllReturnCodes.h>
0013 #include <fun4all/Fun4AllServer.h>
0014 
0015 #include <phool/PHCompositeNode.h>
0016 #include <phool/PHPointerListIterator.h>
0017 #include <phool/getClass.h>
0018 #include <phool/phool.h>  // PHWHERE
0019 
0020 #include <qautils/QAHistManagerDef.h>
0021 
0022 #include <TH1D.h>
0023 #include <TH2D.h>
0024 
0025 #include <format> // -std=c++20
0026 #include <boost/format.hpp> // std::format,vformat are a pain to use with non-const strings
0027 #include <exception>
0028 #include <vector>
0029 
0030 InttQa::InttQa (
0031     std::string const& name
0032 ) : SubsysReco(name) {
0033 
0034     auto* hm = QAHistManagerDef::getHistoManager();
0035 
0036     m_is_streaming_hist = new TH1I (
0037         std::format("{}_is_streaming", m_prefix).c_str(),
0038         "GetBinContent(1) == 1 for streaming, GetBinContent(1) == 0 for triggered",
0039         1, 0.5, 1.5
0040     );
0041     hm->registerHisto(m_is_streaming_hist);
0042 
0043     // The numeric values are important for naming the histograms,
0044     // which is why I do not use range-based for loops here
0045     for (int felix_server = 0; felix_server < n_felix_servers; ++felix_server) {
0046 
0047         m_felix_server_hit_distribution[felix_server] = new TH2D (
0048             std::format("{}_hit_distribution_server{:01d}", m_prefix, felix_server).c_str(),
0049             std::format("Hit Distribution intt{:01d}", felix_server).c_str(),
0050             n_chips, -0.5, n_chips-0.5,
0051             n_felix_channels, -0.5, n_felix_channels-0.5
0052         );
0053         hm->registerHisto(m_felix_server_hit_distribution[felix_server]);
0054 
0055         //...
0056 
0057         for (int felix_channel = 0; felix_channel < n_felix_channels; ++felix_channel) {
0058 
0059             m_felix_channel_bco_distribution[felix_server][felix_channel] = new TH1D (
0060                 std::format("{}_bco_distribution_server{:01d}_channel{:02d}", m_prefix, felix_server, felix_channel).c_str(),
0061                 std::format("BCO Distribution intt{:01d}, Felix Channel {:02d}", felix_server, felix_channel).c_str(),
0062                 n_bcos, -0.5, n_bcos-0.5
0063             );
0064             hm->registerHisto(m_felix_channel_bco_distribution[felix_server][felix_channel]);
0065 
0066             m_felix_channel_hit_distribution[felix_server][felix_channel] = new TH2D (
0067                 std::format("{}_hit_distribution_server{:01d}_channel{:02d}", m_prefix, felix_server, felix_channel).c_str(),
0068                 std::format("Hit Distribution intt{:01d}, Felix Channel {:02d}", felix_server, felix_channel).c_str(),
0069                 n_channels, -0.5, n_channels-0.5,
0070                 n_chips, -0.5, n_chips-0.5
0071             );
0072             hm->registerHisto(m_felix_channel_hit_distribution[felix_server][felix_channel]);
0073 
0074             // ...
0075 
0076             for (int chip = 0; chip < n_chips; ++chip) {
0077 
0078                 m_chip_hit_distribution[felix_server][felix_channel][chip] = new TH1D (
0079                     std::format("{}_hit_distribution_server{:01d}_channel{:02d}_chip{:02d}", m_prefix, felix_server, felix_channel, chip).c_str(),
0080                     std::format("Hit Distribution intt{:01d}, Felix Channel {:02d}, Chip{:02d}", felix_server, felix_channel, chip).c_str(),
0081                     n_bcos, -0.5, n_bcos-0.5
0082                 );
0083                 hm->registerHisto(m_chip_hit_distribution[felix_server][felix_channel][chip]);
0084 
0085                 m_chip_adc_distribution[felix_server][felix_channel][chip] = new TH1D (
0086                     std::format("{}_adc_distribution_server{:01d}_channel{:02d}_chip{:02d}", m_prefix, felix_server, felix_channel, chip).c_str(),
0087                     std::format("ADC Distribution intt{:01d}, Felix Channel {:02d}, Chip{:02d}", felix_server, felix_channel, chip).c_str(),
0088                     n_adcs, -0.5, n_adcs-0.5
0089                 );
0090                 hm->registerHisto(m_chip_adc_distribution[felix_server][felix_channel][chip]);
0091 
0092                 //...
0093             }
0094         }
0095     }
0096 
0097     for (int barrel = 0; barrel < n_barrels; ++barrel) {
0098 
0099         m_barrel_hit_distribution[barrel] = new TH2D (
0100             std::format("{}_hit_distribution_barrel_{:01d}", m_prefix, barrel).c_str(),
0101             std::format("Hit Distribution {} Barrel", (barrel ? "Outer" : "Inner")).c_str(),
0102             n_chips, -0.5, n_chips-0.5,
0103             n_ladders[barrel], -3.1416 * (1.0 + 1.0 / n_ladders[barrel]), +3.1416 * (1.0 - 1.0 / n_ladders[barrel])
0104         );
0105 
0106         hm->registerHisto(m_barrel_hit_distribution[barrel]);
0107     }
0108 }
0109 
0110 void
0111 InttQa::query (
0112 ) {
0113     int runnumber = Fun4AllServer::instance()->RunNumber();
0114     if (runnumber == 0) { // default value in Fun4AllServer.h
0115         std::cout
0116             << PHWHERE
0117             << "Runnumber is default value (0) and was probably initialized, Unregistering"
0118             << std::endl;
0119         Fun4AllServer::instance()->unregisterSubsystem(this);
0120         return;
0121     }
0122 
0123     InttOdbcQuery intt_query;
0124     if (intt_query.Query(runnumber) != 0) {
0125         std::cout
0126             << PHWHERE
0127             << "Database query unsuccessful, Unregistering"
0128             << std::endl;
0129         Fun4AllServer::instance()->unregisterSubsystem(this);
0130         return;
0131     }
0132 
0133     m_is_streaming = intt_query.IsStreaming();
0134     m_is_streaming_hist->SetBinContent(1, m_is_streaming ? 1 : 0);
0135 
0136     if (Verbosity()) {
0137         std::cout
0138             << PHWHERE
0139             << " m_is_streaming: " << m_is_streaming
0140             << " bin contents: " << m_is_streaming_hist->GetBinContent(1)
0141             << std::endl;
0142     }
0143 
0144 }
0145 
0146 int
0147 InttQa::get_nodes (
0148     PHCompositeNode* top_node
0149 ) {
0150     // In case this is called multiple times per instance lifetime
0151     m_intt_raw_hit_containers.clear();
0152 
0153     PHNodeIterator dst_itr(top_node);
0154 
0155     PHCompositeNode* intt_node = dynamic_cast<PHCompositeNode*>(dst_itr.findFirst("INTT"));
0156     if (!intt_node) {
0157         std::cout
0158             << PHWHERE
0159             << "\tNo 'INTT' node found on the NodeTree, Unregistering"
0160             << std::endl;
0161         Fun4AllServer::instance()->unregisterSubsystem(this);
0162         return Fun4AllReturnCodes::EVENT_OK;
0163     }
0164     PHNodeIterator intt_itr(intt_node);
0165 
0166     // For the record, if you're reading this, I do not like this loop structure either
0167     // I think this is moot b/c in practice most DST will only have one node, 'INTTRAWHIT'
0168     PHPointerListIterator<PHNode> next_intt_node(intt_itr.ls());
0169     for (PHNode* itr_node; (itr_node = next_intt_node());) {
0170         auto* intt_raw_hit_container_node = static_cast<PHIODataNode<InttRawHitContainer>*>(itr_node);
0171         if (!intt_raw_hit_container_node) { continue; }
0172         if (Verbosity()) {
0173             std::cout
0174                 << PHWHERE
0175                 << "Got InttRawHitContainerNode: '" << intt_raw_hit_container_node->getName() << "'"
0176                 << std::endl;
0177         }
0178 
0179         auto* intt_raw_hit_container = dynamic_cast<InttRawHitContainer*>(intt_raw_hit_container_node->getData());
0180         if (!intt_raw_hit_container) { continue; }
0181         m_intt_raw_hit_containers.push_back(intt_raw_hit_container);
0182     }
0183 
0184     if (m_intt_raw_hit_containers.empty()) {
0185         std::cout
0186             << PHWHERE
0187             << "No InttRawHitContainers found on the NodeTree, Unregistering"
0188             << std::endl;
0189         Fun4AllServer::instance()->unregisterSubsystem(this);
0190         return Fun4AllReturnCodes::EVENT_OK;
0191     }
0192 
0193     if (m_is_streaming) { // need GL1
0194         m_gl1_packet = findNode::getClass<Gl1Packet>(top_node, "GL1RAWHIT");
0195         if (!m_gl1_packet) {
0196             std::cout
0197                 << PHWHERE
0198                 << "INTT is streaming, but no GL1 found on node tree, Unregistering"
0199                 << std::endl;
0200             Fun4AllServer::instance()->unregisterSubsystem(this);
0201         }
0202         return Fun4AllReturnCodes::EVENT_OK;
0203     }
0204 
0205     return Fun4AllReturnCodes::EVENT_OK;
0206 }
0207 
0208 int
0209 InttQa::InitRun (
0210     PHCompositeNode* top_node
0211 ) {
0212     query();
0213     return get_nodes(top_node);
0214 }
0215 
0216 int
0217 InttQa::process_event (
0218     PHCompositeNode* /*unused*/
0219 ) {
0220     for (auto const& intt_raw_hit_container : m_intt_raw_hit_containers) {
0221         for (unsigned int hit_index{0}; hit_index < intt_raw_hit_container->get_nhits(); ++hit_index) {
0222 
0223             auto* hit = intt_raw_hit_container->get_hit(hit_index);
0224             InttNameSpace::RawData_s raw = InttNameSpace::RawFromHit(hit);
0225 
0226             int bco_diff{0};
0227             int adc = hit->get_adc();
0228 
0229             if (m_is_streaming) {
0230                 uint64_t gl1_bco = m_gl1_packet->getBCO() & 0xFFFFFFFFFFU;
0231                 uint64_t intt_bco = hit->get_FPHX_BCO() & 0xFFFFFFFFFFU;
0232                 bco_diff = int(intt_bco - gl1_bco) + hit->get_bco();
0233             } else {
0234                 bco_diff = (hit->get_FPHX_BCO() - int(hit->get_bco() & 0x7FU) + n_bcos) % n_bcos;
0235             }
0236 
0237             m_felix_server_hit_distribution[raw.felix_server]->Fill(raw.chip, raw.felix_channel);
0238             m_felix_channel_bco_distribution[raw.felix_server][raw.felix_channel]->Fill(bco_diff);
0239             m_felix_channel_hit_distribution[raw.felix_server][raw.felix_channel]->Fill(raw.channel, raw.chip);
0240             m_chip_hit_distribution[raw.felix_server][raw.felix_channel][raw.chip]->Fill(raw.channel);
0241             m_chip_adc_distribution[raw.felix_server][raw.felix_channel][raw.chip]->Fill(adc);
0242 
0243             InttNameSpace::Offline_s offline = InttNameSpace::ToOffline(raw);
0244             int layer = offline.layer - 3;
0245             int barrel = layer / 2;
0246 
0247             // Every other layer is staggered
0248             double phi = 6.2832 * (2.0 * offline.ladder_phi - (layer % 2)) / n_ladders[barrel];
0249             while (3.1416 * (1.0 - 1.0 / n_ladders[barrel]) < phi) { phi -= 6.2832; }
0250 
0251             double z_index = offline.strip_y;
0252             switch (offline.ladder_z) {
0253             case 0:
0254                 z_index += 5;
0255                 break;
0256             case 2:
0257                 z_index += 13;
0258                 break;
0259             case 3:
0260                 z_index += 21;
0261                 break;
0262             default:
0263                 break;
0264             }
0265             m_barrel_hit_distribution[barrel]->Fill(z_index, phi);
0266         }
0267     }
0268 
0269     return Fun4AllReturnCodes::EVENT_OK;
0270 }
0271