File indexing completed on 2026-05-07 08:04:14
0001 #include "TrackingEfficiency_processor.h"
0002
0003 #include <Acts/Definitions/TrackParametrization.hpp>
0004 #include <Acts/EventData/TrackContainer.hpp>
0005 #include <Acts/EventData/TrackProxy.hpp>
0006 #include <Acts/EventData/VectorMultiTrajectory.hpp>
0007 #include <Acts/EventData/VectorTrackContainer.hpp>
0008 #include <ActsExamples/EventData/Track.hpp>
0009 #include <JANA/JApplication.h>
0010 #include <JANA/JApplicationFwd.h>
0011 #include <JANA/JEvent.h>
0012 #include <JANA/Services/JGlobalRootLock.h>
0013 #include <Math/GenVector/Cartesian3D.h>
0014 #include <Math/GenVector/PxPyPzM4D.h>
0015 #include <edm4eic/ReconstructedParticleCollection.h>
0016 #include <edm4hep/MCParticleCollection.h>
0017 #include <edm4hep/Vector3d.h>
0018 #include <edm4hep/Vector3f.h>
0019 #include <fmt/format.h>
0020 #include <spdlog/logger.h>
0021 #include <Eigen/Core>
0022 #include <cassert>
0023 #include <cmath>
0024 #include <cstddef>
0025 #include <map>
0026 #include <string>
0027 #include <vector>
0028
0029 #include "services/log/Log_service.h"
0030 #include "services/rootfile/RootFile_service.h"
0031
0032
0033
0034
0035 void TrackingEfficiency_processor::Init() {
0036 std::string plugin_name = ("tracking_efficiency");
0037
0038
0039 auto* app = GetApplication();
0040
0041
0042 auto root_file_service = app->GetService<RootFile_service>();
0043
0044
0045 auto globalRootLock = app->GetService<JGlobalRootLock>();
0046 globalRootLock->acquire_write_lock();
0047 auto* file = root_file_service->GetHistFile();
0048 globalRootLock->release_lock();
0049
0050
0051 m_dir_main = file->mkdir(plugin_name.c_str());
0052
0053
0054 m_log = app->GetService<Log_service>()->logger(plugin_name);
0055 }
0056
0057
0058
0059
0060 void TrackingEfficiency_processor::Process(const std::shared_ptr<const JEvent>& event) {
0061 using namespace ROOT;
0062
0063
0064
0065 const auto* reco_particles =
0066 event->GetCollection<edm4eic::ReconstructedParticle>("ReconstructedChargedParticles");
0067
0068 m_log->debug("Tracking reconstructed particles N={}: ", reco_particles->size());
0069 m_log->debug(" {:<5} {:>8} {:>8} {:>8} {:>8} {:>8}", "[i]", "[px]", "[py]", "[pz]", "[P]",
0070 "[P*3]");
0071
0072 for (std::size_t i = 0; i < reco_particles->size(); i++) {
0073 const auto& particle = (*reco_particles)[i];
0074
0075 double px = particle.getMomentum().x;
0076 double py = particle.getMomentum().y;
0077 double pz = particle.getMomentum().z;
0078
0079 ROOT::Math::Cartesian3D p(px, py, pz);
0080 m_log->debug(" {:<5} {:>8.2f} {:>8.2f} {:>8.2f} {:>8.2f} {:>8.2f}", i, px, py, pz, p.R(),
0081 p.R() * 3);
0082 }
0083
0084
0085
0086 auto acts_track_states =
0087 event->Get<Acts::ConstVectorMultiTrajectory>("CentralCKFActsTrackStates");
0088 auto acts_tracks = event->Get<Acts::ConstVectorTrackContainer>("CentralCKFActsTracks");
0089 m_log->debug("ACTS Tracks( track states size: {}, tracks size: {} )", acts_track_states.size(),
0090 acts_tracks.size());
0091 m_log->debug("{:>10} {:>10} {:>10} {:>10} {:>10} {:>10} {:>12} {:>12} {:>12} {:>8} {:>8}",
0092 "[loc 0]", "[loc 1]", "[phi]", "[theta]", "[q/p]", "[p]", "[err phi]", "[err th]",
0093 "[err q/p]", "[chi2]", "[ndf]");
0094
0095
0096 if (!acts_track_states.empty() && !acts_tracks.empty()) {
0097 assert(acts_track_states.front() != nullptr &&
0098 "ConstVectorMultiTrajectory pointer should not be null");
0099 assert(acts_tracks.front() != nullptr &&
0100 "ConstVectorTrackContainer pointer should not be null");
0101
0102
0103 auto trackStateContainer =
0104 std::make_shared<Acts::ConstVectorMultiTrajectory>(*acts_track_states.front());
0105 auto trackContainer = std::make_shared<Acts::ConstVectorTrackContainer>(*acts_tracks.front());
0106 ActsExamples::ConstTrackContainer track_container(trackContainer, trackStateContainer);
0107
0108 for (const auto& track : track_container) {
0109
0110 const auto& parameter = track.parameters();
0111 const auto& covariance = track.covariance();
0112 auto chi2 = track.chi2();
0113 auto ndf = track.nDoF();
0114
0115 m_log->debug("{:>10.2f} {:>10.2f} {:>10.2f} {:>10.3f} {:>10.4f} {:>10.3f} {:>12.4e} "
0116 "{:>12.4e} {:>12.4e} {:>8.2f} {:<6}",
0117 parameter[Acts::eBoundLoc0], parameter[Acts::eBoundLoc1],
0118 parameter[Acts::eBoundPhi], parameter[Acts::eBoundTheta],
0119 parameter[Acts::eBoundQOverP], 1.0 / parameter[Acts::eBoundQOverP],
0120 sqrt(covariance(Acts::eBoundPhi, Acts::eBoundPhi)),
0121 sqrt(covariance(Acts::eBoundTheta, Acts::eBoundTheta)),
0122 sqrt(covariance(Acts::eBoundQOverP, Acts::eBoundQOverP)), chi2, ndf);
0123 }
0124 }
0125
0126
0127
0128 const auto* mc_particles = event->GetCollection<edm4hep::MCParticle>("MCParticles");
0129 m_log->debug("MC particles N={}: ", mc_particles->size());
0130 m_log->debug(" {:<5} {:<6} {:<7} {:>8} {:>8} {:>8} {:>8}", "[i]", "status", "[PDG]", "[px]",
0131 "[py]", "[pz]", "[P]");
0132 for (std::size_t i = 0; i < mc_particles->size(); i++) {
0133 const auto& particle = (*mc_particles)[i];
0134
0135
0136 if (particle.getGeneratorStatus() != 1) {
0137 continue;
0138 }
0139
0140 double px = particle.getMomentum().x;
0141 double py = particle.getMomentum().y;
0142 double pz = particle.getMomentum().z;
0143 ROOT::Math::PxPyPzM4D p4v(px, py, pz, particle.getMass());
0144 ROOT::Math::Cartesian3D p(px, py, pz);
0145 if (p.R() < 1) {
0146 continue;
0147 }
0148
0149 m_log->debug(" {:<5} {:<6} {:<7} {:>8.2f} {:>8.2f} {:>8.2f} {:>8.2f}", i,
0150 particle.getGeneratorStatus(), particle.getPDG(), px, py, pz, p.R());
0151 }
0152 }
0153
0154
0155
0156
0157 void TrackingEfficiency_processor::Finish() {
0158 fmt::print("OccupancyAnalysis::Finish() called\n");
0159
0160
0161
0162
0163
0164
0165
0166
0167
0168
0169
0170
0171
0172
0173
0174
0175
0176
0177
0178
0179
0180
0181
0182
0183
0184
0185
0186
0187
0188
0189
0190 }