|
| 1 | +//////////////////////////////////////////////////////////////////////// |
| 2 | +// Class: CRTRateAnalysis |
| 3 | +// Plugin Type: analyzer (Unknown Unknown) |
| 4 | +// File: CRTRateAnalysis_module.cc |
| 5 | +// |
| 6 | +// Generated at Tue Jan 7 05:28:06 2025 by Henry Lay using cetskelgen |
| 7 | +// from cetlib version 3.18.02. |
| 8 | +//////////////////////////////////////////////////////////////////////// |
| 9 | + |
| 10 | +#include "art/Framework/Core/EDAnalyzer.h" |
| 11 | +#include "art/Framework/Core/ModuleMacros.h" |
| 12 | +#include "art/Framework/Principal/Event.h" |
| 13 | +#include "art/Framework/Principal/Handle.h" |
| 14 | +#include "art/Framework/Principal/Run.h" |
| 15 | +#include "art/Framework/Principal/SubRun.h" |
| 16 | +#include "canvas/Utilities/InputTag.h" |
| 17 | +#include "fhiclcpp/ParameterSet.h" |
| 18 | +#include "messagefacility/MessageLogger/MessageLogger.h" |
| 19 | +#include "art_root_io/TFileService.h" |
| 20 | + |
| 21 | +#include "TTree.h" |
| 22 | + |
| 23 | +#include "canvas/Persistency/Common/FindOneP.h" |
| 24 | +#include "canvas/Persistency/Common/FindManyP.h" |
| 25 | + |
| 26 | +#include "artdaq-core/Data/RawEvent.hh" |
| 27 | + |
| 28 | +#include "sbnobj/SBND/CRT/FEBData.hh" |
| 29 | +#include "sbnobj/SBND/CRT/CRTSpacePoint.hh" |
| 30 | +#include "sbnobj/SBND/CRT/CRTCluster.hh" |
| 31 | +#include "sbnobj/SBND/CRT/CRTBlob.hh" |
| 32 | + |
| 33 | +#include "sbndcode/Geometry/GeometryWrappers/CRTGeoService.h" |
| 34 | +#include "sbndcode/ChannelMaps/CRT/CRTChannelMapService.h" |
| 35 | + |
| 36 | +namespace sbnd { |
| 37 | + namespace crt { |
| 38 | + class CRTRateAnalysis; |
| 39 | + } |
| 40 | +} |
| 41 | + |
| 42 | + |
| 43 | +class sbnd::crt::CRTRateAnalysis : public art::EDAnalyzer { |
| 44 | +public: |
| 45 | + explicit CRTRateAnalysis(fhicl::ParameterSet const& p); |
| 46 | + // The compiler-generated destructor is fine for non-base |
| 47 | + // classes without bare pointers or other resource use. |
| 48 | + |
| 49 | + // Plugins should not be copied or assigned. |
| 50 | + CRTRateAnalysis(CRTRateAnalysis const&) = delete; |
| 51 | + CRTRateAnalysis(CRTRateAnalysis&&) = delete; |
| 52 | + CRTRateAnalysis& operator=(CRTRateAnalysis const&) = delete; |
| 53 | + CRTRateAnalysis& operator=(CRTRateAnalysis&&) = delete; |
| 54 | + |
| 55 | + // Required functions. |
| 56 | + void analyze(art::Event const& e) override; |
| 57 | + |
| 58 | + void ResetEventVars(); |
| 59 | + void ResetRawVars(); |
| 60 | + void ResetSpacePointVars(); |
| 61 | + void ResetBlobVars(); |
| 62 | + |
| 63 | +private: |
| 64 | + |
| 65 | + art::ServiceHandle<CRTGeoService> fCRTGeoService; |
| 66 | + art::ServiceHandle<SBND::CRTChannelMapService> fCRTChannelMapService; |
| 67 | + |
| 68 | + std::string fFEBDataModuleLabel, fCRTSpacePointModuleLabel, |
| 69 | + fCRTBlobModuleLabel, fDAQHeaderModuleLabel, fDAQHeaderInstanceLabel; |
| 70 | + |
| 71 | + TTree *fEventTree, *fRawTree, *fSpacePointTree, *fBlobTree; |
| 72 | + |
| 73 | + int32_t _run, _subrun, _event; |
| 74 | + uint32_t _event_header_ts; |
| 75 | + |
| 76 | + int16_t _tagger, _module; |
| 77 | + uint16_t _max_channel; |
| 78 | + |
| 79 | + int16_t _n_hits; |
| 80 | + double _x, _y, _z, _pe; |
| 81 | + |
| 82 | + int16_t _n_total_sps, _n_bottom_sps, _n_bottom_sps_twohits, _n_south_sps, |
| 83 | + _n_north_sps, _n_west_sps, _n_east_sps, _n_toplow_sps, _n_tophigh_sps; |
| 84 | +}; |
| 85 | + |
| 86 | + |
| 87 | +sbnd::crt::CRTRateAnalysis::CRTRateAnalysis(fhicl::ParameterSet const& p) |
| 88 | + : EDAnalyzer{p} |
| 89 | + , fFEBDataModuleLabel(p.get<std::string>("FEBDataModuleLabel")) |
| 90 | + , fCRTSpacePointModuleLabel(p.get<std::string>("CRTSpacePointModuleLabel")) |
| 91 | + , fCRTBlobModuleLabel(p.get<std::string>("CRTBlobModuleLabel")) |
| 92 | + , fDAQHeaderModuleLabel(p.get<std::string>("DAQHeaderModuleLabel")) |
| 93 | + , fDAQHeaderInstanceLabel(p.get<std::string>("DAQHeaderInstanceLabel")) |
| 94 | +{ |
| 95 | + art::ServiceHandle<art::TFileService> fs; |
| 96 | + |
| 97 | + fEventTree = fs->make<TTree>("events", ""); |
| 98 | + fEventTree->Branch("run", &_run); |
| 99 | + fEventTree->Branch("subrun", &_subrun); |
| 100 | + fEventTree->Branch("event", &_event); |
| 101 | + fEventTree->Branch("event_header_ts", &_event_header_ts); |
| 102 | + |
| 103 | + fRawTree = fs->make<TTree>("readouts", ""); |
| 104 | + fRawTree->Branch("run", &_run); |
| 105 | + fRawTree->Branch("subrun", &_subrun); |
| 106 | + fRawTree->Branch("event", &_event); |
| 107 | + fRawTree->Branch("event_header_ts", &_event_header_ts); |
| 108 | + fRawTree->Branch("tagger", &_tagger); |
| 109 | + fRawTree->Branch("module", &_module); |
| 110 | + fRawTree->Branch("max_channel", &_max_channel); |
| 111 | + |
| 112 | + fSpacePointTree = fs->make<TTree>("spacepoints", ""); |
| 113 | + fSpacePointTree->Branch("run", &_run); |
| 114 | + fSpacePointTree->Branch("subrun", &_subrun); |
| 115 | + fSpacePointTree->Branch("event", &_event); |
| 116 | + fSpacePointTree->Branch("event_header_ts", &_event_header_ts); |
| 117 | + fSpacePointTree->Branch("tagger", &_tagger); |
| 118 | + fSpacePointTree->Branch("n_hits", &_n_hits); |
| 119 | + fSpacePointTree->Branch("x", &_x); |
| 120 | + fSpacePointTree->Branch("y", &_y); |
| 121 | + fSpacePointTree->Branch("z", &_z); |
| 122 | + fSpacePointTree->Branch("pe", &_pe); |
| 123 | + |
| 124 | + fBlobTree = fs->make<TTree>("blobs", ""); |
| 125 | + fBlobTree->Branch("run", &_run); |
| 126 | + fBlobTree->Branch("subrun", &_subrun); |
| 127 | + fBlobTree->Branch("event", &_event); |
| 128 | + fBlobTree->Branch("event_header_ts", &_event_header_ts); |
| 129 | + fBlobTree->Branch("n_total_sps", &_n_total_sps); |
| 130 | + fBlobTree->Branch("n_bottom_sps", &_n_bottom_sps); |
| 131 | + fBlobTree->Branch("n_bottom_sps_twohits", &_n_bottom_sps_twohits); |
| 132 | + fBlobTree->Branch("n_south_sps", &_n_south_sps); |
| 133 | + fBlobTree->Branch("n_north_sps", &_n_north_sps); |
| 134 | + fBlobTree->Branch("n_west_sps", &_n_west_sps); |
| 135 | + fBlobTree->Branch("n_east_sps", &_n_east_sps); |
| 136 | + fBlobTree->Branch("n_toplow_sps", &_n_toplow_sps); |
| 137 | + fBlobTree->Branch("n_tophigh_sps", &_n_tophigh_sps); |
| 138 | + fBlobTree->Branch("pe", &_pe); |
| 139 | +} |
| 140 | + |
| 141 | +void sbnd::crt::CRTRateAnalysis::analyze(art::Event const& e) |
| 142 | +{ |
| 143 | + ResetEventVars(); |
| 144 | + |
| 145 | + _run = e.id().run(); |
| 146 | + _subrun = e.id().subRun(); |
| 147 | + _event = e.id().event(); |
| 148 | + |
| 149 | + art::Handle<artdaq::detail::RawEventHeader> DAQHeaderHandle; |
| 150 | + e.getByLabel(fDAQHeaderModuleLabel, fDAQHeaderInstanceLabel, DAQHeaderHandle); |
| 151 | + |
| 152 | + if(DAQHeaderHandle.isValid()) |
| 153 | + { |
| 154 | + artdaq::RawEvent rawHeaderEvent = artdaq::RawEvent(*DAQHeaderHandle); |
| 155 | + uint64_t raw_ts = rawHeaderEvent.timestamp(); |
| 156 | + _event_header_ts = raw_ts / static_cast<uint32_t>(1e9); |
| 157 | + } |
| 158 | + |
| 159 | + fEventTree->Fill(); |
| 160 | + |
| 161 | + art::Handle<std::vector<FEBData>> FEBDataHandle; |
| 162 | + e.getByLabel(fFEBDataModuleLabel, FEBDataHandle); |
| 163 | + |
| 164 | + std::vector<art::Ptr<FEBData>> FEBDataVec; |
| 165 | + art::fill_ptr_vector(FEBDataVec, FEBDataHandle); |
| 166 | + |
| 167 | + for(auto const& data : FEBDataVec) |
| 168 | + { |
| 169 | + ResetRawVars(); |
| 170 | + _tagger = fCRTGeoService->AuxDetIndexToTaggerEnum(data->Mac5()); |
| 171 | + _module = data->Mac5(); |
| 172 | + |
| 173 | + int max_adc = -1, max_ch = -1; |
| 174 | + |
| 175 | + for(int ch = 0; ch < 32; ++ch) |
| 176 | + { |
| 177 | + int adc = data->ADC(ch); |
| 178 | + |
| 179 | + if(adc > max_adc) |
| 180 | + { |
| 181 | + max_adc = adc; |
| 182 | + max_ch = ch; |
| 183 | + } |
| 184 | + } |
| 185 | + |
| 186 | + _max_channel = fCRTChannelMapService->ConstructOfflineChannelIDFromOfflineModuleIDAndOfflineLocalChannel(_module, max_ch); |
| 187 | + |
| 188 | + fRawTree->Fill(); |
| 189 | + } |
| 190 | + |
| 191 | + art::Handle<std::vector<CRTSpacePoint>> CRTSpacePointHandle; |
| 192 | + e.getByLabel(fCRTSpacePointModuleLabel, CRTSpacePointHandle); |
| 193 | + |
| 194 | + std::vector<art::Ptr<CRTSpacePoint>> CRTSpacePointVec; |
| 195 | + art::fill_ptr_vector(CRTSpacePointVec, CRTSpacePointHandle); |
| 196 | + |
| 197 | + art::FindOneP<CRTCluster> spacePointsToClusters(CRTSpacePointHandle, e, fCRTSpacePointModuleLabel); |
| 198 | + |
| 199 | + for(auto const& spacePoint : CRTSpacePointVec) |
| 200 | + { |
| 201 | + ResetSpacePointVars(); |
| 202 | + |
| 203 | + const art::Ptr<CRTCluster> cluster = spacePointsToClusters.at(spacePoint.key()); |
| 204 | + |
| 205 | + _tagger = cluster->Tagger(); |
| 206 | + _n_hits = cluster->NHits(); |
| 207 | + _x = spacePoint->X(); |
| 208 | + _y = spacePoint->Y(); |
| 209 | + _z = spacePoint->Z(); |
| 210 | + _pe = spacePoint->PE(); |
| 211 | + |
| 212 | + fSpacePointTree->Fill(); |
| 213 | + } |
| 214 | + |
| 215 | + art::Handle<std::vector<CRTBlob>> CRTBlobHandle; |
| 216 | + e.getByLabel(fCRTBlobModuleLabel, CRTBlobHandle); |
| 217 | + |
| 218 | + std::vector<art::Ptr<CRTBlob>> CRTBlobVec; |
| 219 | + art::fill_ptr_vector(CRTBlobVec, CRTBlobHandle); |
| 220 | + |
| 221 | + art::FindManyP<CRTSpacePoint> blobsToSpacePoints(CRTBlobHandle, e, fCRTBlobModuleLabel); |
| 222 | + |
| 223 | + for(auto const& blob : CRTBlobVec) |
| 224 | + { |
| 225 | + ResetBlobVars(); |
| 226 | + |
| 227 | + _n_total_sps = blob->TotalSpacePoints(); |
| 228 | + _n_bottom_sps = blob->SpacePointsInTagger(kBottomTagger); |
| 229 | + _n_south_sps = blob->SpacePointsInTagger(kSouthTagger); |
| 230 | + _n_north_sps = blob->SpacePointsInTagger(kNorthTagger); |
| 231 | + _n_west_sps = blob->SpacePointsInTagger(kWestTagger); |
| 232 | + _n_east_sps = blob->SpacePointsInTagger(kEastTagger); |
| 233 | + _n_toplow_sps = blob->SpacePointsInTagger(kTopLowTagger); |
| 234 | + _n_tophigh_sps = blob->SpacePointsInTagger(kTopHighTagger); |
| 235 | + _pe = blob->PE(); |
| 236 | + |
| 237 | + _n_bottom_sps_twohits = 0; |
| 238 | + |
| 239 | + const std::vector<art::Ptr<CRTSpacePoint>> spacePoints = blobsToSpacePoints.at(blob.key()); |
| 240 | + |
| 241 | + for(auto const& spacePoint : spacePoints) |
| 242 | + { |
| 243 | + const art::Ptr<CRTCluster> cluster = spacePointsToClusters.at(spacePoint.key()); |
| 244 | + |
| 245 | + if(cluster->Tagger() == kBottomTagger && cluster->NHits() > 1) |
| 246 | + ++_n_bottom_sps_twohits; |
| 247 | + } |
| 248 | + |
| 249 | + fBlobTree->Fill(); |
| 250 | + } |
| 251 | +} |
| 252 | + |
| 253 | +void sbnd::crt::CRTRateAnalysis::ResetEventVars() |
| 254 | +{ |
| 255 | + _run = -1; _subrun = -1; _event = -1; |
| 256 | + _event_header_ts = std::numeric_limits<uint32_t>::max(); |
| 257 | +} |
| 258 | + |
| 259 | +void sbnd::crt::CRTRateAnalysis::ResetRawVars() |
| 260 | +{ |
| 261 | + _tagger = -1; _module = -1; |
| 262 | + _max_channel = std::numeric_limits<uint16_t>::max(); |
| 263 | +} |
| 264 | + |
| 265 | +void sbnd::crt::CRTRateAnalysis::ResetSpacePointVars() |
| 266 | +{ |
| 267 | + _tagger = -1; |
| 268 | + |
| 269 | + _n_hits = -1; |
| 270 | + |
| 271 | + _x = std::numeric_limits<double>::lowest(); |
| 272 | + _y = std::numeric_limits<double>::lowest(); |
| 273 | + _z = std::numeric_limits<double>::lowest(); |
| 274 | + _pe = std::numeric_limits<double>::lowest(); |
| 275 | +} |
| 276 | + |
| 277 | +void sbnd::crt::CRTRateAnalysis::ResetBlobVars() |
| 278 | +{ |
| 279 | + _n_total_sps = -1; |
| 280 | + _n_bottom_sps = -1; |
| 281 | + _n_bottom_sps_twohits = -1; |
| 282 | + _n_south_sps = -1; |
| 283 | + _n_north_sps = -1; |
| 284 | + _n_west_sps = -1; |
| 285 | + _n_east_sps = -1; |
| 286 | + _n_toplow_sps = -1; |
| 287 | + _n_tophigh_sps = -1; |
| 288 | + |
| 289 | + _pe = std::numeric_limits<double>::lowest(); |
| 290 | +} |
| 291 | + |
| 292 | +DEFINE_ART_MODULE(sbnd::crt::CRTRateAnalysis) |
0 commit comments