Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-12-19 09:18:43

0001 #include "MvtxHitMapWriter.h"
0002 
0003 #include <mvtx/MvtxHitMap.h>
0004 #include <mvtx/MvtxPixelDefs.h>
0005 
0006 #include <ffarawobjects/MvtxRawEvtHeader.h>
0007 #include <ffarawobjects/MvtxRawHit.h>
0008 #include <ffarawobjects/MvtxRawHitContainer.h>
0009 
0010 #include <fun4all/Fun4AllReturnCodes.h>
0011 #include <fun4all/PHTFileServer.h>
0012 
0013 #include <phool/PHCompositeNode.h>
0014 #include <phool/getClass.h>
0015 #include <phool/phool.h>
0016 
0017 #include <TTree.h>
0018 
0019 #include <algorithm>
0020 #include <cstdint>
0021 #include <cstdlib>
0022 #include <iostream>
0023 #include <string>
0024 #include <utility>
0025 #include <vector>
0026 
0027 // MvtxHitMapWriter class
0028 //==============================================================================
0029 int MvtxHitMapWriter::InitRun(PHCompositeNode* /*topNode*/)
0030 {
0031   // Create the hit map
0032   delete m_hit_map;  // make cppcheck happy
0033   m_hit_map = new MvtxHitMap();
0034   m_hit_map->clear();
0035 
0036   std::cout << "MvtxHitMapWriter::InitRun - Writing output to " << m_outputfile << std::endl;
0037 
0038   // Create the output file and trees
0039   PHTFileServer::open(m_outputfile, "RECREATE");
0040   // main tree
0041   m_tree_info = new TTree("info", "MVTX Hit Map Info");
0042   m_tree_info->Branch("num_strobes", &m_num_strobes, "num_strobes/i");
0043   m_tree_info->Branch("nhits_total", &m_nhits_total, "nhits_total/i");
0044   m_tree_info->Branch("num_fired_pixels", &m_num_fired_pixels, "num_fired_pixels/i");
0045 
0046   // hit map tree
0047   m_tree = new TTree("hit_map", "MVTX Hit Map");
0048   m_tree->Branch("pixels", &m_pixels);
0049   m_tree->Branch("nhits", &m_nhits);
0050 
0051   m_num_strobes = 0;
0052   m_nhits_total = 0;
0053   m_num_fired_pixels = 0;
0054 
0055   return Fun4AllReturnCodes::EVENT_OK;
0056 }
0057 
0058 int MvtxHitMapWriter::get_nodes(PHCompositeNode* topNode)
0059 {
0060   // get dst nodes
0061   m_mvtx_raw_event_header = findNode::getClass<MvtxRawEvtHeader>(topNode, "MVTXRAWEVTHEADER");
0062   if (!m_mvtx_raw_event_header)
0063   {
0064     std::cout << PHWHERE << "::" << __func__ << ": Could not get MVTXRAWEVTHEADER from Node Tree" << std::endl;
0065     exit(1);
0066   }
0067   if (Verbosity() > 2)
0068   {
0069     m_mvtx_raw_event_header->identify();
0070   }
0071 
0072   m_mvtx_raw_hit_container = findNode::getClass<MvtxRawHitContainer>(topNode, "MVTXRAWHIT");
0073   if (!m_mvtx_raw_hit_container)
0074   {
0075     std::cout << PHWHERE << "::" << __func__ << ": Could not get MVTXRAWHIT from Node Tree" << std::endl;
0076     exit(1);
0077   }
0078   if (Verbosity() > 2)
0079   {
0080     m_mvtx_raw_hit_container->identify();
0081   }
0082 
0083   return Fun4AllReturnCodes::EVENT_OK;
0084 }
0085 
0086 int MvtxHitMapWriter::process_event(PHCompositeNode* topNode)
0087 {
0088   // get the nodes
0089   if (m_num_strobes % 100000 == 0)
0090   {
0091     std::cout << "MvtxHitMapWriter::process_event - processing strobe number: " << m_num_strobes << std::endl;
0092   }
0093 
0094   get_nodes(topNode);
0095 
0096   // get lls from the MVTX raw event header
0097   for (unsigned int ihit = 0; ihit < m_mvtx_raw_hit_container->get_nhits(); ihit++)
0098   {
0099     // get this hit
0100     auto *mvtx_hit = m_mvtx_raw_hit_container->get_hit(ihit);
0101     if (!mvtx_hit)
0102     {
0103       std::cout << PHWHERE << "::" << __func__ << ": Could not get MVTX hit from container. Hit index: " << ihit << std::endl;
0104       continue;
0105     }
0106 
0107     // get the hit info
0108     uint64_t strobe = mvtx_hit->get_bco();
0109     if (strobe > m_last_strobe)
0110     {
0111       m_last_strobe = strobe;
0112       m_num_strobes++;
0113     }
0114 
0115     // generate pixel key
0116     MvtxPixelDefs::pixelkey this_pixelkey = MvtxPixelDefs::gen_pixelkey(mvtx_hit);
0117 
0118     // add the hit to the hit map
0119     m_hit_map->add_hit(this_pixelkey, 1);
0120   }
0121 
0122   return Fun4AllReturnCodes::EVENT_OK;
0123 }
0124 
0125 int MvtxHitMapWriter::FillTree()
0126 {
0127   // get the hit map
0128   MvtxHitMap::pixel_hit_vector_t pixel_hit_vector = m_hit_map->get_pixel_hit_vector();
0129   std::sort(pixel_hit_vector.begin(), pixel_hit_vector.end(), [](const MvtxHitMap::pixel_hits_pair_t& a, const MvtxHitMap::pixel_hits_pair_t& b)
0130             { return a.second > b.second; });
0131 
0132   m_num_fired_pixels = pixel_hit_vector.size();
0133   if (m_num_fired_pixels == 0)
0134   {
0135     std::cout << "MvtxHitMapWriter::FillTree - No pixels with hits" << std::endl;
0136     return Fun4AllReturnCodes::EVENT_OK;
0137   }
0138 
0139   m_nhits_total = 0;
0140   m_pixels.clear();
0141   m_nhits.clear();
0142 
0143   for (auto & it : pixel_hit_vector)
0144   {
0145     m_pixels.push_back(it.first);
0146     m_nhits.push_back(it.second);
0147     m_nhits_total += it.second;
0148   }
0149 
0150   m_tree->Fill();
0151   m_tree_info->Fill();
0152 
0153   std::cout << "MvtxHitMapWriter::FillTree - Number of pixels with hits: " << m_num_fired_pixels << std::endl;
0154   std::cout << "MvtxHitMapWriter::FillTree - Total number of hits: " << m_nhits_total << std::endl;
0155   std::cout << "MvtxHitMapWriter::FillTree - Number of strobes: " << m_num_strobes << std::endl;
0156 
0157   return Fun4AllReturnCodes::EVENT_OK;
0158 }
0159 
0160 int MvtxHitMapWriter::End(PHCompositeNode* /*topNode*/)
0161 {
0162   // Fill the current mask tree
0163   FillTree();
0164   std::cout << "MvtxHitMapWriter::End - Writing output to " << m_outputfile << std::endl;
0165   PHTFileServer::cd(m_outputfile);
0166   PHTFileServer::write(m_outputfile);
0167   return Fun4AllReturnCodes::EVENT_OK;
0168 }