File indexing completed on 2025-06-30 07:55:46
0001 #include "TrackingEfficiency_processor.h"
0002
0003 #include <Acts/Definitions/TrackParametrization.hpp>
0004 #include <Acts/EventData/MultiTrajectoryHelpers.hpp>
0005 #include <ActsExamples/EventData/Trajectories.hpp>
0006 #include <JANA/JApplication.h>
0007 #include <JANA/JApplicationFwd.h>
0008 #include <JANA/JEvent.h>
0009 #include <JANA/Services/JGlobalRootLock.h>
0010 #include <Math/GenVector/Cartesian3D.h>
0011 #include <Math/GenVector/PxPyPzM4D.h>
0012 #include <Rtypes.h>
0013 #include <edm4eic/ReconstructedParticleCollection.h>
0014 #include <edm4hep/MCParticleCollection.h>
0015 #include <edm4hep/Vector3d.h>
0016 #include <edm4hep/Vector3f.h>
0017 #include <fmt/core.h>
0018 #include <spdlog/logger.h>
0019 #include <Eigen/Core>
0020 #include <cmath>
0021 #include <cstddef>
0022 #include <iterator>
0023 #include <map>
0024 #include <optional>
0025 #include <string>
0026 #include <vector>
0027
0028 #include "services/log/Log_service.h"
0029 #include "services/rootfile/RootFile_service.h"
0030
0031
0032
0033
0034 void TrackingEfficiency_processor::Init() {
0035 std::string plugin_name = ("tracking_efficiency");
0036
0037
0038 auto* app = GetApplication();
0039
0040
0041 auto root_file_service = app->GetService<RootFile_service>();
0042
0043
0044 auto globalRootLock = app->GetService<JGlobalRootLock>();
0045 globalRootLock->acquire_write_lock();
0046 auto* file = root_file_service->GetHistFile();
0047 globalRootLock->release_lock();
0048
0049
0050 m_dir_main = file->mkdir(plugin_name.c_str());
0051
0052
0053 m_log = app->GetService<Log_service>()->logger(plugin_name);
0054 }
0055
0056
0057
0058
0059 void TrackingEfficiency_processor::Process(const std::shared_ptr<const JEvent>& event) {
0060 using namespace ROOT;
0061
0062
0063
0064 const auto reco_particles =
0065 event->Get<edm4eic::ReconstructedParticle>("ReconstructedChargedParticles");
0066
0067 m_log->debug("Tracking reconstructed particles N={}: ", reco_particles.size());
0068 m_log->debug(" {:<5} {:>8} {:>8} {:>8} {:>8} {:>8}", "[i]", "[px]", "[py]", "[pz]", "[P]",
0069 "[P*3]");
0070
0071 for (std::size_t i = 0; i < reco_particles.size(); i++) {
0072 const auto& particle = *(reco_particles[i]);
0073
0074 double px = particle.getMomentum().x;
0075 double py = particle.getMomentum().y;
0076 double pz = particle.getMomentum().z;
0077
0078 ROOT::Math::Cartesian3D p(px, py, pz);
0079 m_log->debug(" {:<5} {:>8.2f} {:>8.2f} {:>8.2f} {:>8.2f} {:>8.2f}", i, px, py, pz, p.R(),
0080 p.R() * 3);
0081 }
0082
0083
0084
0085 auto acts_results = event->Get<ActsExamples::Trajectories>("CentralCKFActsTrajectories");
0086 m_log->debug("ACTS Trajectories( size: {} )", std::size(acts_results));
0087 m_log->debug("{:>10} {:>10} {:>10} {:>10} {:>10} {:>10} {:>12} {:>12} {:>12} {:>8}", "[loc 0]",
0088 "[loc 1]", "[phi]", "[theta]", "[q/p]", "[p]", "[err phi]", "[err th]", "[err q/p]",
0089 "[chi2]");
0090
0091
0092 for (const auto& traj : acts_results) {
0093
0094
0095
0096 const auto& mj = traj->multiTrajectory();
0097 const auto& trackTips = traj->tips();
0098 if (trackTips.empty()) {
0099 m_log->debug("Empty multiTrajectory.");
0100 continue;
0101 }
0102
0103 const auto& trackTip = trackTips.front();
0104
0105
0106 auto trajState = Acts::MultiTrajectoryHelpers::trajectoryState(mj, trackTip);
0107 if (traj->hasTrackParameters(trackTip)) {
0108 const auto& boundParam = traj->trackParameters(trackTip);
0109 const auto& parameter = boundParam.parameters();
0110 const auto& covariance = *boundParam.covariance();
0111 m_log->debug("{:>10.2f} {:>10.2f} {:>10.2f} {:>10.3f} {:>10.4f} {:>10.3f} {:>12.4e} "
0112 "{:>12.4e} {:>12.4e} {:>8.2f}",
0113 parameter[Acts::eBoundLoc0], parameter[Acts::eBoundLoc1],
0114 parameter[Acts::eBoundPhi], parameter[Acts::eBoundTheta],
0115 parameter[Acts::eBoundQOverP], 1.0 / parameter[Acts::eBoundQOverP],
0116 sqrt(covariance(Acts::eBoundPhi, Acts::eBoundPhi)),
0117 sqrt(covariance(Acts::eBoundTheta, Acts::eBoundTheta)),
0118 sqrt(covariance(Acts::eBoundQOverP, Acts::eBoundQOverP)), trajState.chi2Sum);
0119 }
0120 }
0121
0122
0123
0124 auto mc_particles = event->Get<edm4hep::MCParticle>("MCParticles");
0125 m_log->debug("MC particles N={}: ", mc_particles.size());
0126 m_log->debug(" {:<5} {:<6} {:<7} {:>8} {:>8} {:>8} {:>8}", "[i]", "status", "[PDG]", "[px]",
0127 "[py]", "[pz]", "[P]");
0128 for (std::size_t i = 0; i < mc_particles.size(); i++) {
0129 const auto* particle = mc_particles[i];
0130
0131
0132 if (particle->getGeneratorStatus() != 1) {
0133 continue;
0134 }
0135
0136 double px = particle->getMomentum().x;
0137 double py = particle->getMomentum().y;
0138 double pz = particle->getMomentum().z;
0139 ROOT::Math::PxPyPzM4D p4v(px, py, pz, particle->getMass());
0140 ROOT::Math::Cartesian3D p(px, py, pz);
0141 if (p.R() < 1) {
0142 continue;
0143 }
0144
0145 m_log->debug(" {:<5} {:<6} {:<7} {:>8.2f} {:>8.2f} {:>8.2f} {:>8.2f}", i,
0146 particle->getGeneratorStatus(), particle->getPDG(), px, py, pz, p.R());
0147 }
0148 }
0149
0150
0151
0152
0153 void TrackingEfficiency_processor::Finish() {
0154 fmt::print("OccupancyAnalysis::Finish() called\n");
0155
0156
0157
0158
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 }