File indexing completed on 2025-04-03 08:04:50
0001
0002
0003
0004
0005
0006
0007
0008
0009 #include "DDG4/Geant4EventAction.h"
0010 #include "DD4hep/Printout.h"
0011 #include "DDG4/Geant4Kernel.h"
0012
0013
0014 #include "G4Event.hh"
0015 #include "G4TrajectoryContainer.hh"
0016 #include "G4VTrajectory.hh"
0017 #include "G4VTrajectoryPoint.hh"
0018 #include "G4SystemOfUnits.hh"
0019 #include "G4UnitsTable.hh"
0020 #include "G4AttValue.hh"
0021 #include "G4AttDef.hh"
0022 #include "G4RichTrajectory.hh"
0023 #include "G4RichTrajectoryPoint.hh"
0024 #include "G4SmoothTrajectory.hh"
0025 #include "G4SmoothTrajectoryPoint.hh"
0026 #include "G4Trajectory.hh"
0027 #include "G4TrajectoryPoint.hh"
0028
0029
0030 #include <fstream>
0031 #include <iostream>
0032 #include <vector>
0033 #include <sstream>
0034 #include <string>
0035 #include <map>
0036 #include <cmath>
0037 #include <limits>
0038 #include <fmt/core.h>
0039 #include <fmt/format.h>
0040
0041
0042 namespace dd4hep {
0043
0044
0045 namespace sim {
0046
0047
0048
0049
0050
0051
0052
0053
0054
0055
0056
0057
0058
0059
0060
0061
0062
0063
0064
0065
0066 class FirebirdTrajectoryWriterEventAction : public Geant4EventAction {
0067 protected:
0068
0069 std::string m_outputFile {"trajectories.firebird.json"};
0070
0071
0072 std::string m_componentName {"Geant4Trajectories"};
0073
0074
0075 bool m_saveOptical {false};
0076
0077
0078 bool m_onlyPrimary {false};
0079
0080
0081 bool m_vertexCut {false};
0082
0083
0084 double m_vertexZMin {-5000};
0085
0086
0087 double m_vertexZMax {5000};
0088
0089
0090 bool m_stepCut {false};
0091
0092
0093 double m_stepZMin {-5000};
0094
0095
0096 double m_stepZMax {5000};
0097
0098
0099 double m_stepRMax {5000};
0100
0101
0102 double m_minMomentum {150};
0103
0104
0105 double m_maxMomentum {1e6};
0106
0107
0108 double m_minTrackLength {0};
0109
0110
0111 std::vector<int> m_saveParticles {};
0112
0113
0114 bool m_requireRichTrajectory {true};
0115
0116
0117 bool m_verboseTimeExtraction {false};
0118
0119
0120 long m_totalTrajectories {0};
0121 long m_filteredTrajectories {0};
0122 long m_savedTrajectories {0};
0123 long m_trajectoryWithoutTime {0};
0124 long m_stepsFiltered {0};
0125
0126
0127 std::vector<std::string> m_entries;
0128
0129
0130 bool isValidForJSON(G4double value) const {
0131 return !std::isinf(value) && !std::isnan(value);
0132 }
0133
0134
0135 G4double getSafeValue(G4double value, G4double defaultValue = 0.0) const {
0136 return isValidForJSON(value) ? value : defaultValue;
0137 }
0138
0139
0140 G4double calculateR(G4ThreeVector position) const {
0141 return std::sqrt(position.x() * position.x() + position.y() * position.y());
0142 }
0143
0144
0145 G4double extractTimeFromPoint(G4VTrajectoryPoint* point, int pointIndex) {
0146
0147 G4RichTrajectoryPoint* richPoint = dynamic_cast<G4RichTrajectoryPoint*>(point);
0148 if (!richPoint) {
0149 if (m_requireRichTrajectory) {
0150 if (m_verboseTimeExtraction) {
0151 warning("[firebird-writer] Point %d is not a rich trajectory point, cannot extract time", pointIndex);
0152 }
0153 m_trajectoryWithoutTime++;
0154 return -1.0;
0155 }
0156
0157
0158 return pointIndex * 0.1 * CLHEP::ns;
0159 }
0160
0161 std::vector<G4AttValue>* attValues = richPoint->CreateAttValues();
0162 if (!attValues) {
0163 if (m_verboseTimeExtraction) {
0164 warning("[firebird-writer] Point %d has no attribute values", pointIndex);
0165 }
0166 return -1.0;
0167 }
0168
0169
0170 std::string timeAttName = (pointIndex == 0) ? "PreT" : "PostT";
0171 G4double extractedTime = -1.0;
0172
0173 for (const auto& attValue : *attValues) {
0174 if (attValue.GetName() == timeAttName) {
0175 std::string valueStr = attValue.GetValue();
0176
0177
0178 std::istringstream iss(valueStr);
0179 G4double timeValue;
0180 std::string unit;
0181
0182
0183 iss >> timeValue >> unit;
0184
0185
0186 if (!unit.empty()) {
0187 if (unit == "ns") {
0188 extractedTime = timeValue * CLHEP::ns;
0189 } else if (unit == "s") {
0190 extractedTime = timeValue * CLHEP::s;
0191 } else if (unit == "ms") {
0192 extractedTime = timeValue * CLHEP::ms;
0193 } else if (unit == "us" || unit == "µs") {
0194 extractedTime = timeValue * CLHEP::microsecond;
0195 } else {
0196
0197 extractedTime = timeValue * CLHEP::ns;
0198 }
0199 } else {
0200
0201 extractedTime = timeValue * CLHEP::ns;
0202 }
0203
0204 if (m_verboseTimeExtraction) {
0205 info("[firebird-writer] Extracted time %s = %f ns from point %d",
0206 timeAttName.c_str(), extractedTime / CLHEP::ns, pointIndex);
0207 }
0208
0209
0210 break;
0211 }
0212 }
0213
0214
0215 delete attValues;
0216
0217
0218
0219 if (extractedTime < 0 && m_requireRichTrajectory) {
0220 if (m_verboseTimeExtraction) {
0221 warning("[firebird-writer] Could not find %s in point %d", timeAttName.c_str(), pointIndex);
0222 }
0223 m_trajectoryWithoutTime++;
0224 return -1.0;
0225 }
0226
0227
0228 if (extractedTime < 0) {
0229 return pointIndex * 0.1 * CLHEP::ns;
0230 }
0231
0232 return extractedTime;
0233 }
0234
0235 public:
0236
0237 FirebirdTrajectoryWriterEventAction(Geant4Context* context, const std::string& name = "FirebirdTrajectoryWriterEventAction")
0238 : Geant4EventAction(context, name)
0239 {
0240 declareProperty("OutputFile", m_outputFile);
0241 declareProperty("ComponentName", m_componentName);
0242 declareProperty("SaveOptical", m_saveOptical);
0243 declareProperty("OnlyPrimary", m_onlyPrimary);
0244 declareProperty("VertexCut", m_vertexCut);
0245 declareProperty("VertexZMin", m_vertexZMin);
0246 declareProperty("VertexZMax", m_vertexZMax);
0247 declareProperty("StepCut", m_stepCut);
0248 declareProperty("StepZMin", m_stepZMin);
0249 declareProperty("StepZMax", m_stepZMax);
0250 declareProperty("StepRMax", m_stepRMax);
0251 declareProperty("MomentumMin", m_minMomentum);
0252 declareProperty("MomentumMax", m_maxMomentum);
0253 declareProperty("TrackLengthMin", m_minTrackLength);
0254 declareProperty("SaveParticles", m_saveParticles);
0255 declareProperty("RequireRichTrajectory", m_requireRichTrajectory);
0256 declareProperty("VerboseTimeExtraction", m_verboseTimeExtraction);
0257
0258
0259 info("[firebird-writer] Trajectory filtering configuration:");
0260 info("[firebird-writer] OutputFile: %s", m_outputFile.c_str());
0261 info("[firebird-writer] ComponentName: %s", m_componentName.c_str());
0262 info("[firebird-writer] SaveOptical: %s", m_saveOptical ? "true" : "false");
0263 info("[firebird-writer] OnlyPrimary: %s", m_onlyPrimary ? "true" : "false");
0264 info("[firebird-writer] VertexCut: %s", m_vertexCut ? "true" : "false");
0265 info("[firebird-writer] VertexZMin: %.2f mm", m_vertexZMin);
0266 info("[firebird-writer] VertexZMax: %.2f mm", m_vertexZMax);
0267 info("[firebird-writer] StepCut: %s", m_stepCut ? "true" : "false");
0268 info("[firebird-writer] StepZMin: %.2f mm", m_stepZMin);
0269 info("[firebird-writer] StepZMax: %.2f mm", m_stepZMax);
0270 info("[firebird-writer] StepRMax: %.2f mm", m_stepRMax);
0271 info("[firebird-writer] MinMomentum: %.3f MeV/c", m_minMomentum);
0272 info("[firebird-writer] MaxMomentum: %.3f MeV/c", m_maxMomentum);
0273 info("[firebird-writer] MinTrackLength: %.2f mm", m_minTrackLength);
0274 info("[firebird-writer] RequireRichTrajectory: %s", m_requireRichTrajectory ? "true" : "false");
0275 info("[firebird-writer] VerboseTimeExtraction: %s", m_verboseTimeExtraction ? "true" : "false");
0276
0277 if (!m_saveParticles.empty()) {
0278 std::stringstream ss;
0279 for (size_t i = 0; i < m_saveParticles.size(); ++i) {
0280 if (i > 0) ss << ", ";
0281 ss << m_saveParticles[i];
0282 }
0283 info("[firebird-writer] SaveParticles: %s", ss.str().c_str());
0284 } else {
0285 info("[firebird-writer] SaveParticles: [all]");
0286 }
0287 }
0288
0289
0290 virtual ~FirebirdTrajectoryWriterEventAction() {
0291
0292 if (!m_entries.empty()) {
0293 try {
0294
0295 std::ofstream output(m_outputFile);
0296 if (output.is_open()) {
0297
0298 output << fmt::format(R"({{"type":"firebird-dex-json","version":"0.03","origin":{{"file":"{}","entries_count":{}}},)",
0299 m_outputFile, m_entries.size());
0300
0301
0302 output << "\"entries\":[";
0303
0304
0305 for (size_t i = 0; i < m_entries.size(); ++i) {
0306 output << m_entries[i];
0307 if (i < m_entries.size() - 1) {
0308 output << ",";
0309 }
0310 }
0311
0312
0313 output << "]}";
0314 output.close();
0315
0316 info("[firebird-writer] Successfully wrote JSON trajectories to: %s", m_outputFile.c_str());
0317 } else {
0318 error("[firebird-writer] Failed to open output file: %s", m_outputFile.c_str());
0319 }
0320 } catch (const std::exception& e) {
0321 error("[firebird-writer] Error writing JSON file: %s", e.what());
0322 }
0323 } else {
0324 warning("[firebird-writer] No events were processed. Output file not created.");
0325 }
0326
0327
0328 info("[firebird-writer] Trajectory filtering statistics:");
0329 info("[firebird-writer] Total trajectories processed: %ld", m_totalTrajectories);
0330 info("[firebird-writer] Filtered (skipped) trajectories: %ld (%0.1f%%)",
0331 m_filteredTrajectories,
0332 m_totalTrajectories > 0 ? (m_filteredTrajectories * 100.0 / m_totalTrajectories) : 0.0);
0333 info("[firebird-writer] Saved trajectories: %ld (%0.1f%%)",
0334 m_savedTrajectories,
0335 m_totalTrajectories > 0 ? (m_savedTrajectories * 100.0 / m_totalTrajectories) : 0.0);
0336 if (m_stepCut) {
0337 info("[firebird-writer] Steps filtered due to position limits: %ld", m_stepsFiltered);
0338 }
0339 if (m_requireRichTrajectory) {
0340 info("[firebird-writer] Trajectories without proper time information: %ld", m_trajectoryWithoutTime);
0341 }
0342 }
0343
0344
0345 bool passesFilters(G4VTrajectory* trajectory) {
0346
0347 int pdgCode = trajectory->GetPDGEncoding();
0348 std::string particleName = trajectory->GetParticleName();
0349 int parentID = trajectory->GetParentID();
0350 G4ThreeVector momentum = trajectory->GetInitialMomentum();
0351 double p = momentum.mag() / CLHEP::MeV;
0352
0353
0354 if (particleName == "opticalphoton" && m_saveOptical) {
0355 return true;
0356 }
0357
0358
0359 if (m_onlyPrimary && parentID != 0) {
0360 return false;
0361 }
0362
0363
0364 if (p < m_minMomentum || p > m_maxMomentum) {
0365 return false;
0366 }
0367
0368
0369 if (!m_saveParticles.empty()) {
0370 bool found = false;
0371 for (int code : m_saveParticles) {
0372 if (pdgCode == code) {
0373 found = true;
0374 break;
0375 }
0376 }
0377 if (!found) {
0378 return false;
0379 }
0380 }
0381
0382
0383 if (m_minTrackLength > 0) {
0384
0385 G4int n_points = trajectory->GetPointEntries();
0386
0387 if (n_points <= 1) {
0388
0389 return false;
0390 }
0391
0392
0393 double trackLength = 0;
0394 G4ThreeVector prevPos = trajectory->GetPoint(0)->GetPosition();
0395
0396 for (G4int i = 1; i < n_points; i++) {
0397 G4ThreeVector pos = trajectory->GetPoint(i)->GetPosition();
0398 trackLength += (pos - prevPos).mag();
0399 prevPos = pos;
0400 }
0401
0402
0403 if (trackLength / CLHEP::mm < m_minTrackLength) {
0404 return false;
0405 }
0406 }
0407
0408
0409 if (m_vertexCut) {
0410
0411 G4int n_points = trajectory->GetPointEntries();
0412
0413 if (n_points > 0) {
0414 G4ThreeVector vertex = trajectory->GetPoint(0)->GetPosition();
0415 double vz = vertex.z() / CLHEP::mm;
0416
0417 if (vz < m_vertexZMin || vz > m_vertexZMax) {
0418 return false;
0419 }
0420 }
0421 }
0422
0423
0424 if (m_requireRichTrajectory) {
0425 G4RichTrajectory* richTrajectory = dynamic_cast<G4RichTrajectory*>(trajectory);
0426 if (!richTrajectory) {
0427 if (m_verboseTimeExtraction) {
0428 warning("[firebird-writer] Trajectory is not a rich trajectory, skipping");
0429 }
0430 return false;
0431 }
0432
0433
0434 G4int n_points = trajectory->GetPointEntries();
0435 if (n_points > 0) {
0436 G4double time = extractTimeFromPoint(trajectory->GetPoint(0), 0);
0437 if (time < 0) {
0438 if (m_verboseTimeExtraction) {
0439 warning("[firebird-writer] First point of trajectory has no time information, skipping");
0440 }
0441 return false;
0442 }
0443 }
0444 }
0445
0446
0447 return true;
0448 }
0449
0450
0451 std::string generateTrackParams(G4VTrajectory* trajectory) {
0452
0453 G4ThreeVector momentum = trajectory->GetInitialMomentum();
0454 G4double p = momentum.mag();
0455
0456
0457 if (p < 1e-10) {
0458 p = 1e-10;
0459 }
0460
0461
0462 int pdgCode = trajectory->GetPDGEncoding();
0463 std::string particleName = trajectory->GetParticleName();
0464
0465
0466 G4double charge = trajectory->GetCharge();
0467
0468
0469 G4double theta = momentum.theta();
0470 G4double phi = momentum.phi();
0471
0472
0473 G4double qOverP = charge / (p / CLHEP::GeV);
0474
0475
0476 G4double px = momentum.x() / CLHEP::MeV;
0477 G4double py = momentum.y() / CLHEP::MeV;
0478 G4double pz = momentum.z() / CLHEP::MeV;
0479
0480
0481 G4ThreeVector vertex(0, 0, 0);
0482 G4double time = 0.0;
0483
0484 if (trajectory->GetPointEntries() > 0) {
0485 G4VTrajectoryPoint* point = trajectory->GetPoint(0);
0486 vertex = point->GetPosition();
0487
0488
0489 time = extractTimeFromPoint(point, 0);
0490 if (time < 0) {
0491 time = 0.0;
0492 }
0493 time /= CLHEP::ns;
0494 }
0495
0496
0497 G4double vx = vertex.x() / CLHEP::mm;
0498 G4double vy = vertex.y() / CLHEP::mm;
0499 G4double vz = vertex.z() / CLHEP::mm;
0500
0501
0502 px = getSafeValue(px);
0503 py = getSafeValue(py);
0504 pz = getSafeValue(pz);
0505 vx = getSafeValue(vx);
0506 vy = getSafeValue(vy);
0507 vz = getSafeValue(vz);
0508 theta = getSafeValue(theta);
0509 phi = getSafeValue(phi);
0510 qOverP = getSafeValue(qOverP);
0511 time = getSafeValue(time);
0512
0513
0514 G4double locA = 0.0;
0515 G4double locB = 0.0;
0516
0517
0518
0519 return fmt::format("[{},\"{}\",{},{},{},{},{},{},{},{},{},{},{},{},{}]",
0520 pdgCode, particleName, charge,
0521 px, py, pz,
0522 vx, vy, vz,
0523 theta, phi, qOverP,
0524 locA, locB, time);
0525 }
0526
0527
0528 std::string processTrajectoryPoints(G4VTrajectory* trajectory) {
0529 G4int n_points = trajectory->GetPointEntries();
0530
0531 if (n_points == 0) {
0532 return "[]";
0533 }
0534
0535 std::string pointsStr = "[";
0536 bool firstPoint = true;
0537
0538
0539 for (G4int i = 0; i < n_points; i++) {
0540 G4VTrajectoryPoint* point = trajectory->GetPoint(i);
0541 G4ThreeVector position = point->GetPosition();
0542
0543
0544 if (m_stepCut) {
0545 double z = position.z() / CLHEP::mm;
0546 double r = calculateR(position) / CLHEP::mm;
0547
0548
0549 if (z < m_stepZMin || z > m_stepZMax || r > m_stepRMax) {
0550 m_stepsFiltered++;
0551 continue;
0552 }
0553 }
0554
0555
0556 G4double time = extractTimeFromPoint(point, i);
0557 if (time < 0) {
0558
0559 time = i * 0.1 * CLHEP::ns;
0560 }
0561 time /= CLHEP::ns;
0562
0563
0564 double x = position.x() / CLHEP::mm;
0565 double y = position.y() / CLHEP::mm;
0566 double z = position.z() / CLHEP::mm;
0567
0568
0569 x = getSafeValue(x);
0570 y = getSafeValue(y);
0571 z = getSafeValue(z);
0572 time = getSafeValue(time);
0573
0574
0575 if (!firstPoint) {
0576 pointsStr += ",";
0577 } else {
0578 firstPoint = false;
0579 }
0580
0581 pointsStr += fmt::format("[{},{},{},{},{}]", x, y, z, time, 0);
0582 }
0583
0584 pointsStr += "]";
0585 return pointsStr;
0586 }
0587
0588
0589 virtual void begin(const G4Event* ) override {
0590
0591 }
0592
0593
0594 virtual void end(const G4Event* event) override {
0595 G4TrajectoryContainer* trajectoryContainer = event->GetTrajectoryContainer();
0596 if (!trajectoryContainer) {
0597 warning("[firebird-writer] No trajectory container found for event %d", event->GetEventID());
0598 return;
0599 }
0600
0601 G4int n_trajectories = trajectoryContainer->entries();
0602 if (n_trajectories == 0) {
0603 warning("[firebird-writer] No trajectories found for event %d", event->GetEventID());
0604 return;
0605 }
0606
0607
0608 m_totalTrajectories += n_trajectories;
0609
0610
0611 int filtered_event = 0;
0612 int saved_event = 0;
0613
0614
0615 std::string eventEntry = fmt::format(R"({{"id":{},"components":[)", event->GetEventID());
0616
0617
0618 eventEntry += fmt::format(R"({{"name":"{}","type":"TrackerLinePointTrajectory",)", m_componentName);
0619
0620
0621 eventEntry += R"("originType":["G4VTrajectory","G4VTrajectoryPoint"],)";
0622
0623
0624 eventEntry += R"("paramColumns":["pdg","type","charge","px","py","pz","vx","vy","vz","theta","phi","q_over_p","loc_a","loc_b","time"],)";
0625
0626
0627 eventEntry += R"("pointColumns":["x","y","z","t","aux"],)";
0628
0629
0630 eventEntry += R"("lines":[)";
0631
0632
0633 bool firstLine = true;
0634 for (G4int i = 0; i < n_trajectories; i++) {
0635 G4VTrajectory* trajectory = (*trajectoryContainer)[i];
0636
0637
0638 if (!passesFilters(trajectory)) {
0639 filtered_event++;
0640 continue;
0641 }
0642
0643
0644 std::string pointsStr = processTrajectoryPoints(trajectory);
0645
0646
0647 if (pointsStr == "[]") {
0648 filtered_event++;
0649 continue;
0650 }
0651
0652
0653 saved_event++;
0654
0655
0656 if (!firstLine) {
0657 eventEntry += ",";
0658 } else {
0659 firstLine = false;
0660 }
0661
0662
0663 eventEntry += "{\"points\":";
0664
0665
0666 eventEntry += pointsStr;
0667
0668
0669 eventEntry += ",\"params\":";
0670 eventEntry += generateTrackParams(trajectory);
0671
0672
0673 eventEntry += "}";
0674 }
0675
0676
0677 eventEntry += "]}]}";
0678
0679
0680 if (saved_event > 0) {
0681 m_entries.push_back(eventEntry);
0682 }
0683
0684
0685 m_filteredTrajectories += filtered_event;
0686 m_savedTrajectories += saved_event;
0687
0688
0689 info("[firebird-writer] Event %d: processed %d trajectories, filtered %d, saved %d",
0690 event->GetEventID(), n_trajectories, filtered_event, saved_event);
0691 }
0692 };
0693 }
0694 }
0695
0696
0697 #include "DDG4/Factories.h"
0698 DECLARE_GEANT4ACTION(FirebirdTrajectoryWriterEventAction)