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