Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-12-16 09:24:02

0001 // This file is part of the ACTS project.
0002 //
0003 // Copyright (C) 2016 CERN for the benefit of the ACTS project
0004 //
0005 // This Source Code Form is subject to the terms of the Mozilla Public
0006 // License, v. 2.0. If a copy of the MPL was not distributed with this
0007 // file, You can obtain one at https://mozilla.org/MPL/2.0/.
0008 
0009 #include "ActsExamples/Io/Root/RootTrackSummaryWriter.hpp"
0010 
0011 #include "Acts/Definitions/TrackParametrization.hpp"
0012 #include "Acts/EventData/VectorMultiTrajectory.hpp"
0013 #include "Acts/Surfaces/Surface.hpp"
0014 #include "Acts/TrackFitting/GsfOptions.hpp"
0015 #include "Acts/Utilities/Intersection.hpp"
0016 #include "Acts/Utilities/Result.hpp"
0017 #include "Acts/Utilities/detail/periodic.hpp"
0018 #include "ActsExamples/EventData/TruthMatching.hpp"
0019 #include "ActsExamples/Framework/AlgorithmContext.hpp"
0020 #include "ActsExamples/Framework/WriterT.hpp"
0021 #include "ActsExamples/Validation/TrackClassification.hpp"
0022 #include "ActsFatras/EventData/Barcode.hpp"
0023 
0024 #include <array>
0025 #include <cmath>
0026 #include <cstddef>
0027 #include <cstdint>
0028 #include <ios>
0029 #include <limits>
0030 #include <numbers>
0031 #include <optional>
0032 #include <ostream>
0033 #include <stdexcept>
0034 
0035 #include <TFile.h>
0036 #include <TTree.h>
0037 
0038 using Acts::VectorHelpers::eta;
0039 using Acts::VectorHelpers::perp;
0040 using Acts::VectorHelpers::phi;
0041 using Acts::VectorHelpers::theta;
0042 
0043 namespace ActsExamples {
0044 
0045 RootTrackSummaryWriter::RootTrackSummaryWriter(
0046     const RootTrackSummaryWriter::Config& config, Acts::Logging::Level level)
0047     : WriterT(config.inputTracks, "RootTrackSummaryWriter", level),
0048       m_cfg(config) {
0049   // tracks collection name is already checked by base ctor
0050   if (m_cfg.filePath.empty()) {
0051     throw std::invalid_argument("Missing output filename");
0052   }
0053   if (m_cfg.treeName.empty()) {
0054     throw std::invalid_argument("Missing tree name");
0055   }
0056 
0057   m_inputParticles.maybeInitialize(m_cfg.inputParticles);
0058   m_inputTrackParticleMatching.maybeInitialize(
0059       m_cfg.inputTrackParticleMatching);
0060 
0061   // Setup ROOT I/O
0062   auto path = m_cfg.filePath;
0063   m_outputFile = TFile::Open(path.c_str(), m_cfg.fileMode.c_str());
0064   if (m_outputFile == nullptr) {
0065     throw std::ios_base::failure("Could not open '" + path + "'");
0066   }
0067   m_outputFile->cd();
0068   m_outputTree = new TTree(m_cfg.treeName.c_str(), m_cfg.treeName.c_str());
0069   if (m_outputTree == nullptr) {
0070     throw std::bad_alloc();
0071   }
0072 
0073   // I/O parameters
0074   m_outputTree->Branch("event_nr", &m_eventNr);
0075   m_outputTree->Branch("track_nr", &m_trackNr);
0076 
0077   m_outputTree->Branch("nStates", &m_nStates);
0078   m_outputTree->Branch("nMeasurements", &m_nMeasurements);
0079   m_outputTree->Branch("nOutliers", &m_nOutliers);
0080   m_outputTree->Branch("nHoles", &m_nHoles);
0081   m_outputTree->Branch("nSharedHits", &m_nSharedHits);
0082   m_outputTree->Branch("chi2Sum", &m_chi2Sum);
0083   m_outputTree->Branch("NDF", &m_NDF);
0084   m_outputTree->Branch("measurementChi2", &m_measurementChi2);
0085   m_outputTree->Branch("outlierChi2", &m_outlierChi2);
0086   m_outputTree->Branch("measurementVolume", &m_measurementVolume);
0087   m_outputTree->Branch("measurementLayer", &m_measurementLayer);
0088   m_outputTree->Branch("outlierVolume", &m_outlierVolume);
0089   m_outputTree->Branch("outlierLayer", &m_outlierLayer);
0090 
0091   m_outputTree->Branch("nMajorityHits", &m_nMajorityHits);
0092   m_outputTree->Branch("majorityParticleId_vertex_primary",
0093                        &m_majorityParticleVertexPrimary);
0094   m_outputTree->Branch("majorityParticleId_vertex_secondary",
0095                        &m_majorityParticleVertexSecondary);
0096   m_outputTree->Branch("majorityParticleId_particle",
0097                        &m_majorityParticleParticle);
0098   m_outputTree->Branch("majorityParticleId_generation",
0099                        &m_majorityParticleGeneration);
0100   m_outputTree->Branch("majorityParticleId_sub_particle",
0101                        &m_majorityParticleSubParticle);
0102   m_outputTree->Branch("trackClassification", &m_trackClassification);
0103   m_outputTree->Branch("t_charge", &m_t_charge);
0104   m_outputTree->Branch("t_time", &m_t_time);
0105   m_outputTree->Branch("t_vx", &m_t_vx);
0106   m_outputTree->Branch("t_vy", &m_t_vy);
0107   m_outputTree->Branch("t_vz", &m_t_vz);
0108   m_outputTree->Branch("t_px", &m_t_px);
0109   m_outputTree->Branch("t_py", &m_t_py);
0110   m_outputTree->Branch("t_pz", &m_t_pz);
0111   m_outputTree->Branch("t_theta", &m_t_theta);
0112   m_outputTree->Branch("t_phi", &m_t_phi);
0113   m_outputTree->Branch("t_eta", &m_t_eta);
0114   m_outputTree->Branch("t_p", &m_t_p);
0115   m_outputTree->Branch("t_pT", &m_t_pT);
0116   m_outputTree->Branch("t_d0", &m_t_d0);
0117   m_outputTree->Branch("t_z0", &m_t_z0);
0118   m_outputTree->Branch("t_prodR", &m_t_prodR);
0119 
0120   m_outputTree->Branch("hasFittedParams", &m_hasFittedParams);
0121   m_outputTree->Branch("eLOC0_fit", &m_eLOC0_fit);
0122   m_outputTree->Branch("eLOC1_fit", &m_eLOC1_fit);
0123   m_outputTree->Branch("ePHI_fit", &m_ePHI_fit);
0124   m_outputTree->Branch("eTHETA_fit", &m_eTHETA_fit);
0125   m_outputTree->Branch("eQOP_fit", &m_eQOP_fit);
0126   m_outputTree->Branch("eT_fit", &m_eT_fit);
0127   m_outputTree->Branch("err_eLOC0_fit", &m_err_eLOC0_fit);
0128   m_outputTree->Branch("err_eLOC1_fit", &m_err_eLOC1_fit);
0129   m_outputTree->Branch("err_ePHI_fit", &m_err_ePHI_fit);
0130   m_outputTree->Branch("err_eTHETA_fit", &m_err_eTHETA_fit);
0131   m_outputTree->Branch("err_eQOP_fit", &m_err_eQOP_fit);
0132   m_outputTree->Branch("err_eT_fit", &m_err_eT_fit);
0133   m_outputTree->Branch("res_eLOC0_fit", &m_res_eLOC0_fit);
0134   m_outputTree->Branch("res_eLOC1_fit", &m_res_eLOC1_fit);
0135   m_outputTree->Branch("res_ePHI_fit", &m_res_ePHI_fit);
0136   m_outputTree->Branch("res_eTHETA_fit", &m_res_eTHETA_fit);
0137   m_outputTree->Branch("res_eQOP_fit", &m_res_eQOP_fit);
0138   m_outputTree->Branch("res_eT_fit", &m_res_eT_fit);
0139   m_outputTree->Branch("pull_eLOC0_fit", &m_pull_eLOC0_fit);
0140   m_outputTree->Branch("pull_eLOC1_fit", &m_pull_eLOC1_fit);
0141   m_outputTree->Branch("pull_ePHI_fit", &m_pull_ePHI_fit);
0142   m_outputTree->Branch("pull_eTHETA_fit", &m_pull_eTHETA_fit);
0143   m_outputTree->Branch("pull_eQOP_fit", &m_pull_eQOP_fit);
0144   m_outputTree->Branch("pull_eT_fit", &m_pull_eT_fit);
0145 
0146   if (m_cfg.writeGsfSpecific) {
0147     m_outputTree->Branch("max_material_fwd", &m_gsf_max_material_fwd);
0148     m_outputTree->Branch("sum_material_fwd", &m_gsf_sum_material_fwd);
0149   }
0150 
0151   if (m_cfg.writeCovMat) {
0152     // create one branch for every entry of covariance matrix
0153     // one block for every row of the matrix, every entry gets own branch
0154     m_outputTree->Branch("cov_eLOC0_eLOC0", &m_cov_eLOC0_eLOC0);
0155     m_outputTree->Branch("cov_eLOC0_eLOC1", &m_cov_eLOC0_eLOC1);
0156     m_outputTree->Branch("cov_eLOC0_ePHI", &m_cov_eLOC0_ePHI);
0157     m_outputTree->Branch("cov_eLOC0_eTHETA", &m_cov_eLOC0_eTHETA);
0158     m_outputTree->Branch("cov_eLOC0_eQOP", &m_cov_eLOC0_eQOP);
0159     m_outputTree->Branch("cov_eLOC0_eT", &m_cov_eLOC0_eT);
0160 
0161     m_outputTree->Branch("cov_eLOC1_eLOC0", &m_cov_eLOC1_eLOC0);
0162     m_outputTree->Branch("cov_eLOC1_eLOC1", &m_cov_eLOC1_eLOC1);
0163     m_outputTree->Branch("cov_eLOC1_ePHI", &m_cov_eLOC1_ePHI);
0164     m_outputTree->Branch("cov_eLOC1_eTHETA", &m_cov_eLOC1_eTHETA);
0165     m_outputTree->Branch("cov_eLOC1_eQOP", &m_cov_eLOC1_eQOP);
0166     m_outputTree->Branch("cov_eLOC1_eT", &m_cov_eLOC1_eT);
0167 
0168     m_outputTree->Branch("cov_ePHI_eLOC0", &m_cov_ePHI_eLOC0);
0169     m_outputTree->Branch("cov_ePHI_eLOC1", &m_cov_ePHI_eLOC1);
0170     m_outputTree->Branch("cov_ePHI_ePHI", &m_cov_ePHI_ePHI);
0171     m_outputTree->Branch("cov_ePHI_eTHETA", &m_cov_ePHI_eTHETA);
0172     m_outputTree->Branch("cov_ePHI_eQOP", &m_cov_ePHI_eQOP);
0173     m_outputTree->Branch("cov_ePHI_eT", &m_cov_ePHI_eT);
0174 
0175     m_outputTree->Branch("cov_eTHETA_eLOC0", &m_cov_eTHETA_eLOC0);
0176     m_outputTree->Branch("cov_eTHETA_eLOC1", &m_cov_eTHETA_eLOC1);
0177     m_outputTree->Branch("cov_eTHETA_ePHI", &m_cov_eTHETA_ePHI);
0178     m_outputTree->Branch("cov_eTHETA_eTHETA", &m_cov_eTHETA_eTHETA);
0179     m_outputTree->Branch("cov_eTHETA_eQOP", &m_cov_eTHETA_eQOP);
0180     m_outputTree->Branch("cov_eTHETA_eT", &m_cov_eTHETA_eT);
0181 
0182     m_outputTree->Branch("cov_eQOP_eLOC0", &m_cov_eQOP_eLOC0);
0183     m_outputTree->Branch("cov_eQOP_eLOC1", &m_cov_eQOP_eLOC1);
0184     m_outputTree->Branch("cov_eQOP_ePHI", &m_cov_eQOP_ePHI);
0185     m_outputTree->Branch("cov_eQOP_eTHETA", &m_cov_eQOP_eTHETA);
0186     m_outputTree->Branch("cov_eQOP_eQOP", &m_cov_eQOP_eQOP);
0187     m_outputTree->Branch("cov_eQOP_eT", &m_cov_eQOP_eT);
0188 
0189     m_outputTree->Branch("cov_eT_eLOC0", &m_cov_eT_eLOC0);
0190     m_outputTree->Branch("cov_eT_eLOC1", &m_cov_eT_eLOC1);
0191     m_outputTree->Branch("cov_eT_ePHI", &m_cov_eT_ePHI);
0192     m_outputTree->Branch("cov_eT_eTHETA", &m_cov_eT_eTHETA);
0193     m_outputTree->Branch("cov_eT_eQOP", &m_cov_eT_eQOP);
0194     m_outputTree->Branch("cov_eT_eT", &m_cov_eT_eT);
0195   }
0196 
0197   if (m_cfg.writeGx2fSpecific) {
0198     m_outputTree->Branch("nUpdatesGx2f", &m_nUpdatesGx2f);
0199   }
0200 }
0201 
0202 RootTrackSummaryWriter::~RootTrackSummaryWriter() {
0203   m_outputFile->Close();
0204 }
0205 
0206 ProcessCode RootTrackSummaryWriter::finalize() {
0207   m_outputFile->cd();
0208   m_outputTree->Write();
0209   m_outputFile->Close();
0210 
0211   if (m_cfg.writeCovMat) {
0212     ACTS_INFO("Wrote full covariance matrix to tree");
0213   }
0214   ACTS_INFO("Wrote parameters of tracks to tree '" << m_cfg.treeName << "' in '"
0215                                                    << m_cfg.filePath << "'");
0216 
0217   return ProcessCode::SUCCESS;
0218 }
0219 
0220 ProcessCode RootTrackSummaryWriter::writeT(const AlgorithmContext& ctx,
0221                                            const ConstTrackContainer& tracks) {
0222   // In case we do not have truth info, we bind to a empty collection
0223   const static SimParticleContainer emptyParticles;
0224   const static TrackParticleMatching emptyTrackParticleMatching;
0225 
0226   const auto& particles =
0227       m_inputParticles.isInitialized() ? m_inputParticles(ctx) : emptyParticles;
0228   const auto& trackParticleMatching =
0229       m_inputTrackParticleMatching.isInitialized()
0230           ? m_inputTrackParticleMatching(ctx)
0231           : emptyTrackParticleMatching;
0232 
0233   // For each particle within a track, how many hits did it contribute
0234   std::vector<ParticleHitCount> particleHitCounts;
0235 
0236   // Exclusive access to the tree while writing
0237   std::lock_guard<std::mutex> lock(m_writeMutex);
0238 
0239   // Get the event number
0240   m_eventNr = ctx.eventNumber;
0241 
0242   for (const auto& track : tracks) {
0243     m_trackNr.push_back(track.index());
0244 
0245     // Collect the trajectory summary info
0246     m_nStates.push_back(track.nTrackStates());
0247     m_nMeasurements.push_back(track.nMeasurements());
0248     m_nOutliers.push_back(track.nOutliers());
0249     m_nHoles.push_back(track.nHoles());
0250     m_nSharedHits.push_back(track.nSharedHits());
0251     m_chi2Sum.push_back(track.chi2());
0252     m_NDF.push_back(track.nDoF());
0253 
0254     {
0255       std::vector<double> measurementChi2;
0256       std::vector<std::uint32_t> measurementVolume;
0257       std::vector<std::uint32_t> measurementLayer;
0258       std::vector<double> outlierChi2;
0259       std::vector<std::uint32_t> outlierVolume;
0260       std::vector<std::uint32_t> outlierLayer;
0261       for (const auto& state : track.trackStatesReversed()) {
0262         const auto& geoID = state.referenceSurface().geometryId();
0263         const auto& volume = geoID.volume();
0264         const auto& layer = geoID.layer();
0265         if (state.typeFlags().test(Acts::TrackStateFlag::OutlierFlag)) {
0266           outlierChi2.push_back(state.chi2());
0267           outlierVolume.push_back(volume);
0268           outlierLayer.push_back(layer);
0269         } else if (state.typeFlags().test(
0270                        Acts::TrackStateFlag::MeasurementFlag)) {
0271           measurementChi2.push_back(state.chi2());
0272           measurementVolume.push_back(volume);
0273           measurementLayer.push_back(layer);
0274         }
0275       }
0276       m_measurementChi2.push_back(std::move(measurementChi2));
0277       m_measurementVolume.push_back(std::move(measurementVolume));
0278       m_measurementLayer.push_back(std::move(measurementLayer));
0279       m_outlierChi2.push_back(std::move(outlierChi2));
0280       m_outlierVolume.push_back(std::move(outlierVolume));
0281       m_outlierLayer.push_back(std::move(outlierLayer));
0282     }
0283 
0284     // Initialize the truth particle info
0285     ActsFatras::Barcode majorityParticleId{};
0286     TrackMatchClassification trackClassification =
0287         TrackMatchClassification::Unknown;
0288     unsigned int nMajorityHits = std::numeric_limits<unsigned int>::max();
0289     int t_charge = std::numeric_limits<int>::max();
0290     float t_time = NaNfloat;
0291     float t_vx = NaNfloat;
0292     float t_vy = NaNfloat;
0293     float t_vz = NaNfloat;
0294     float t_px = NaNfloat;
0295     float t_py = NaNfloat;
0296     float t_pz = NaNfloat;
0297     float t_theta = NaNfloat;
0298     float t_phi = NaNfloat;
0299     float t_eta = NaNfloat;
0300     float t_p = NaNfloat;
0301     float t_pT = NaNfloat;
0302     float t_d0 = NaNfloat;
0303     float t_z0 = NaNfloat;
0304     float t_qop = NaNfloat;
0305     float t_prodR = NaNfloat;
0306 
0307     // Get the perigee surface
0308     const Acts::Surface* pSurface =
0309         track.hasReferenceSurface() ? &track.referenceSurface() : nullptr;
0310 
0311     // Get the majority truth particle to this track
0312     auto match = trackParticleMatching.find(track.index());
0313     bool foundMajorityParticle = false;
0314     // Get the truth particle info
0315     if (match != trackParticleMatching.end() &&
0316         match->second.particle.has_value()) {
0317       // Get the barcode of the majority truth particle
0318       majorityParticleId = match->second.particle.value();
0319       trackClassification = match->second.classification;
0320       nMajorityHits = match->second.contributingParticles.front().hitCount;
0321 
0322       // Find the truth particle via the barcode
0323       auto ip = particles.find(majorityParticleId);
0324       if (ip != particles.end()) {
0325         foundMajorityParticle = true;
0326 
0327         const auto& particle = *ip;
0328         ACTS_VERBOSE("Find the truth particle with barcode "
0329                      << majorityParticleId << "=" << majorityParticleId.hash());
0330         // Get the truth particle info at vertex
0331         t_p = particle.absoluteMomentum();
0332         t_charge = static_cast<int>(particle.charge());
0333         t_time = particle.time();
0334         t_vx = particle.position().x();
0335         t_vy = particle.position().y();
0336         t_vz = particle.position().z();
0337         t_px = t_p * particle.direction().x();
0338         t_py = t_p * particle.direction().y();
0339         t_pz = t_p * particle.direction().z();
0340         t_theta = theta(particle.direction());
0341         t_phi = phi(particle.direction());
0342         t_eta = eta(particle.direction());
0343         t_pT = t_p * perp(particle.direction());
0344         t_qop = particle.qOverP();
0345         t_prodR = std::sqrt(t_vx * t_vx + t_vy * t_vy);
0346 
0347         if (pSurface != nullptr) {
0348           Acts::Intersection3D intersection =
0349               pSurface
0350                   ->intersect(ctx.geoContext, particle.position(),
0351                               particle.direction(),
0352                               Acts::BoundaryTolerance::Infinite())
0353                   .closest();
0354           auto position = intersection.position();
0355 
0356           // get the truth perigee parameter
0357           auto lpResult = pSurface->globalToLocal(ctx.geoContext, position,
0358                                                   particle.direction());
0359           if (lpResult.ok()) {
0360             t_d0 = lpResult.value()[Acts::BoundIndices::eBoundLoc0];
0361             t_z0 = lpResult.value()[Acts::BoundIndices::eBoundLoc1];
0362           } else {
0363             ACTS_ERROR("Global to local transformation did not succeed.");
0364           }
0365         }
0366       } else {
0367         ACTS_DEBUG("Truth particle with barcode "
0368                    << majorityParticleId << "=" << majorityParticleId.hash()
0369                    << " not found in the input collection!");
0370       }
0371     }
0372     if (!foundMajorityParticle) {
0373       ACTS_DEBUG("Truth particle for track " << track.tipIndex()
0374                                              << " not found!");
0375     }
0376 
0377     // Push the corresponding truth particle info for the track.
0378     // Always push back even if majority particle not found
0379     m_majorityParticleVertexPrimary.push_back(
0380         majorityParticleId.vertexPrimary());
0381     m_majorityParticleVertexSecondary.push_back(
0382         majorityParticleId.vertexSecondary());
0383     m_majorityParticleParticle.push_back(majorityParticleId.particle());
0384     m_majorityParticleGeneration.push_back(majorityParticleId.generation());
0385     m_majorityParticleSubParticle.push_back(majorityParticleId.subParticle());
0386     m_trackClassification.push_back(static_cast<int>(trackClassification));
0387     m_nMajorityHits.push_back(nMajorityHits);
0388     m_t_charge.push_back(t_charge);
0389     m_t_time.push_back(t_time);
0390     m_t_vx.push_back(t_vx);
0391     m_t_vy.push_back(t_vy);
0392     m_t_vz.push_back(t_vz);
0393     m_t_px.push_back(t_px);
0394     m_t_py.push_back(t_py);
0395     m_t_pz.push_back(t_pz);
0396     m_t_theta.push_back(t_theta);
0397     m_t_phi.push_back(t_phi);
0398     m_t_eta.push_back(t_eta);
0399     m_t_p.push_back(t_p);
0400     m_t_pT.push_back(t_pT);
0401     m_t_d0.push_back(t_d0);
0402     m_t_z0.push_back(t_z0);
0403     m_t_prodR.push_back(t_prodR);
0404 
0405     // Initialize the fitted track parameters info
0406     std::array<float, Acts::eBoundSize> param = {NaNfloat, NaNfloat, NaNfloat,
0407                                                  NaNfloat, NaNfloat, NaNfloat};
0408     std::array<float, Acts::eBoundSize> error = {NaNfloat, NaNfloat, NaNfloat,
0409                                                  NaNfloat, NaNfloat, NaNfloat};
0410 
0411     // get entries of covariance matrix. If no entry, return NaN
0412     auto getCov = [&](auto i, auto j) { return track.covariance()(i, j); };
0413 
0414     bool hasFittedParams = track.hasReferenceSurface();
0415     if (hasFittedParams) {
0416       const auto& parameter = track.parameters();
0417       for (unsigned int i = 0; i < Acts::eBoundSize; ++i) {
0418         param[i] = parameter[i];
0419       }
0420 
0421       for (unsigned int i = 0; i < Acts::eBoundSize; ++i) {
0422         double variance = getCov(i, i);
0423         error[i] = variance >= 0 ? std::sqrt(variance) : NaNfloat;
0424       }
0425     }
0426 
0427     std::array<float, Acts::eBoundSize> res = {NaNfloat, NaNfloat, NaNfloat,
0428                                                NaNfloat, NaNfloat, NaNfloat};
0429     std::array<float, Acts::eBoundSize> pull = {NaNfloat, NaNfloat, NaNfloat,
0430                                                 NaNfloat, NaNfloat, NaNfloat};
0431     if (foundMajorityParticle && hasFittedParams) {
0432       res = {param[Acts::eBoundLoc0] - t_d0,
0433              param[Acts::eBoundLoc1] - t_z0,
0434              Acts::detail::difference_periodic(
0435                  param[Acts::eBoundPhi], t_phi,
0436                  static_cast<float>(2 * std::numbers::pi)),
0437              param[Acts::eBoundTheta] - t_theta,
0438              param[Acts::eBoundQOverP] - t_qop,
0439              param[Acts::eBoundTime] - t_time};
0440 
0441       for (unsigned int i = 0; i < Acts::eBoundSize; ++i) {
0442         pull[i] = res[i] / error[i];
0443       }
0444     }
0445 
0446     // Push the fitted track parameters.
0447     // Always push back even if no fitted track parameters
0448     m_eLOC0_fit.push_back(param[Acts::eBoundLoc0]);
0449     m_eLOC1_fit.push_back(param[Acts::eBoundLoc1]);
0450     m_ePHI_fit.push_back(param[Acts::eBoundPhi]);
0451     m_eTHETA_fit.push_back(param[Acts::eBoundTheta]);
0452     m_eQOP_fit.push_back(param[Acts::eBoundQOverP]);
0453     m_eT_fit.push_back(param[Acts::eBoundTime]);
0454 
0455     m_res_eLOC0_fit.push_back(res[Acts::eBoundLoc0]);
0456     m_res_eLOC1_fit.push_back(res[Acts::eBoundLoc1]);
0457     m_res_ePHI_fit.push_back(res[Acts::eBoundPhi]);
0458     m_res_eTHETA_fit.push_back(res[Acts::eBoundTheta]);
0459     m_res_eQOP_fit.push_back(res[Acts::eBoundQOverP]);
0460     m_res_eT_fit.push_back(res[Acts::eBoundTime]);
0461 
0462     m_err_eLOC0_fit.push_back(error[Acts::eBoundLoc0]);
0463     m_err_eLOC1_fit.push_back(error[Acts::eBoundLoc1]);
0464     m_err_ePHI_fit.push_back(error[Acts::eBoundPhi]);
0465     m_err_eTHETA_fit.push_back(error[Acts::eBoundTheta]);
0466     m_err_eQOP_fit.push_back(error[Acts::eBoundQOverP]);
0467     m_err_eT_fit.push_back(error[Acts::eBoundTime]);
0468 
0469     m_pull_eLOC0_fit.push_back(pull[Acts::eBoundLoc0]);
0470     m_pull_eLOC1_fit.push_back(pull[Acts::eBoundLoc1]);
0471     m_pull_ePHI_fit.push_back(pull[Acts::eBoundPhi]);
0472     m_pull_eTHETA_fit.push_back(pull[Acts::eBoundTheta]);
0473     m_pull_eQOP_fit.push_back(pull[Acts::eBoundQOverP]);
0474     m_pull_eT_fit.push_back(pull[Acts::eBoundTime]);
0475 
0476     m_hasFittedParams.push_back(hasFittedParams);
0477 
0478     if (m_cfg.writeGsfSpecific) {
0479       using namespace Acts::GsfConstants;
0480       if (tracks.hasColumn(Acts::hashString(kFwdMaxMaterialXOverX0))) {
0481         m_gsf_max_material_fwd.push_back(
0482             track.template component<double>(kFwdMaxMaterialXOverX0));
0483       } else {
0484         m_gsf_max_material_fwd.push_back(NaNfloat);
0485       }
0486 
0487       if (tracks.hasColumn(Acts::hashString(kFwdSumMaterialXOverX0))) {
0488         m_gsf_sum_material_fwd.push_back(
0489             track.template component<double>(kFwdSumMaterialXOverX0));
0490       } else {
0491         m_gsf_sum_material_fwd.push_back(NaNfloat);
0492       }
0493     }
0494 
0495     if (m_cfg.writeCovMat) {
0496       // write all entries of covariance matrix to output file
0497       // one branch for every entry of the matrix.
0498       m_cov_eLOC0_eLOC0.push_back(getCov(0, 0));
0499       m_cov_eLOC0_eLOC1.push_back(getCov(0, 1));
0500       m_cov_eLOC0_ePHI.push_back(getCov(0, 2));
0501       m_cov_eLOC0_eTHETA.push_back(getCov(0, 3));
0502       m_cov_eLOC0_eQOP.push_back(getCov(0, 4));
0503       m_cov_eLOC0_eT.push_back(getCov(0, 5));
0504 
0505       m_cov_eLOC1_eLOC0.push_back(getCov(1, 0));
0506       m_cov_eLOC1_eLOC1.push_back(getCov(1, 1));
0507       m_cov_eLOC1_ePHI.push_back(getCov(1, 2));
0508       m_cov_eLOC1_eTHETA.push_back(getCov(1, 3));
0509       m_cov_eLOC1_eQOP.push_back(getCov(1, 4));
0510       m_cov_eLOC1_eT.push_back(getCov(1, 5));
0511 
0512       m_cov_ePHI_eLOC0.push_back(getCov(2, 0));
0513       m_cov_ePHI_eLOC1.push_back(getCov(2, 1));
0514       m_cov_ePHI_ePHI.push_back(getCov(2, 2));
0515       m_cov_ePHI_eTHETA.push_back(getCov(2, 3));
0516       m_cov_ePHI_eQOP.push_back(getCov(2, 4));
0517       m_cov_ePHI_eT.push_back(getCov(2, 5));
0518 
0519       m_cov_eTHETA_eLOC0.push_back(getCov(3, 0));
0520       m_cov_eTHETA_eLOC1.push_back(getCov(3, 1));
0521       m_cov_eTHETA_ePHI.push_back(getCov(3, 2));
0522       m_cov_eTHETA_eTHETA.push_back(getCov(3, 3));
0523       m_cov_eTHETA_eQOP.push_back(getCov(3, 4));
0524       m_cov_eTHETA_eT.push_back(getCov(3, 5));
0525 
0526       m_cov_eQOP_eLOC0.push_back(getCov(4, 0));
0527       m_cov_eQOP_eLOC1.push_back(getCov(4, 1));
0528       m_cov_eQOP_ePHI.push_back(getCov(4, 2));
0529       m_cov_eQOP_eTHETA.push_back(getCov(4, 3));
0530       m_cov_eQOP_eQOP.push_back(getCov(4, 4));
0531       m_cov_eQOP_eT.push_back(getCov(4, 5));
0532 
0533       m_cov_eT_eLOC0.push_back(getCov(5, 0));
0534       m_cov_eT_eLOC1.push_back(getCov(5, 1));
0535       m_cov_eT_ePHI.push_back(getCov(5, 2));
0536       m_cov_eT_eTHETA.push_back(getCov(5, 3));
0537       m_cov_eT_eQOP.push_back(getCov(5, 4));
0538       m_cov_eT_eT.push_back(getCov(5, 5));
0539     }
0540 
0541     if (m_cfg.writeGx2fSpecific) {
0542       if (tracks.hasColumn(Acts::hashString("Gx2fnUpdateColumn"))) {
0543         int nUpdate = static_cast<int>(
0544             track.template component<std::uint32_t,
0545                                      Acts::hashString("Gx2fnUpdateColumn")>());
0546         m_nUpdatesGx2f.push_back(nUpdate);
0547       } else {
0548         m_nUpdatesGx2f.push_back(-1);
0549       }
0550     }
0551   }
0552 
0553   // fill the variables
0554   m_outputTree->Fill();
0555 
0556   m_trackNr.clear();
0557   m_nStates.clear();
0558   m_nMeasurements.clear();
0559   m_nOutliers.clear();
0560   m_nHoles.clear();
0561   m_nSharedHits.clear();
0562   m_chi2Sum.clear();
0563   m_NDF.clear();
0564   m_measurementChi2.clear();
0565   m_outlierChi2.clear();
0566   m_measurementVolume.clear();
0567   m_measurementLayer.clear();
0568   m_outlierVolume.clear();
0569   m_outlierLayer.clear();
0570 
0571   m_nMajorityHits.clear();
0572   m_majorityParticleVertexPrimary.clear();
0573   m_majorityParticleVertexSecondary.clear();
0574   m_majorityParticleParticle.clear();
0575   m_majorityParticleGeneration.clear();
0576   m_majorityParticleSubParticle.clear();
0577   m_trackClassification.clear();
0578   m_t_charge.clear();
0579   m_t_time.clear();
0580   m_t_vx.clear();
0581   m_t_vy.clear();
0582   m_t_vz.clear();
0583   m_t_px.clear();
0584   m_t_py.clear();
0585   m_t_pz.clear();
0586   m_t_theta.clear();
0587   m_t_phi.clear();
0588   m_t_p.clear();
0589   m_t_pT.clear();
0590   m_t_eta.clear();
0591   m_t_d0.clear();
0592   m_t_z0.clear();
0593   m_t_prodR.clear();
0594 
0595   m_hasFittedParams.clear();
0596   m_eLOC0_fit.clear();
0597   m_eLOC1_fit.clear();
0598   m_ePHI_fit.clear();
0599   m_eTHETA_fit.clear();
0600   m_eQOP_fit.clear();
0601   m_eT_fit.clear();
0602   m_err_eLOC0_fit.clear();
0603   m_err_eLOC1_fit.clear();
0604   m_err_ePHI_fit.clear();
0605   m_err_eTHETA_fit.clear();
0606   m_err_eQOP_fit.clear();
0607   m_err_eT_fit.clear();
0608   m_res_eLOC0_fit.clear();
0609   m_res_eLOC1_fit.clear();
0610   m_res_ePHI_fit.clear();
0611   m_res_eTHETA_fit.clear();
0612   m_res_eQOP_fit.clear();
0613   m_res_eT_fit.clear();
0614   m_pull_eLOC0_fit.clear();
0615   m_pull_eLOC1_fit.clear();
0616   m_pull_ePHI_fit.clear();
0617   m_pull_eTHETA_fit.clear();
0618   m_pull_eQOP_fit.clear();
0619   m_pull_eT_fit.clear();
0620 
0621   m_gsf_max_material_fwd.clear();
0622   m_gsf_sum_material_fwd.clear();
0623 
0624   if (m_cfg.writeCovMat) {
0625     m_cov_eLOC0_eLOC0.clear();
0626     m_cov_eLOC0_eLOC1.clear();
0627     m_cov_eLOC0_ePHI.clear();
0628     m_cov_eLOC0_eTHETA.clear();
0629     m_cov_eLOC0_eQOP.clear();
0630     m_cov_eLOC0_eT.clear();
0631 
0632     m_cov_eLOC1_eLOC0.clear();
0633     m_cov_eLOC1_eLOC1.clear();
0634     m_cov_eLOC1_ePHI.clear();
0635     m_cov_eLOC1_eTHETA.clear();
0636     m_cov_eLOC1_eQOP.clear();
0637     m_cov_eLOC1_eT.clear();
0638 
0639     m_cov_ePHI_eLOC0.clear();
0640     m_cov_ePHI_eLOC1.clear();
0641     m_cov_ePHI_ePHI.clear();
0642     m_cov_ePHI_eTHETA.clear();
0643     m_cov_ePHI_eQOP.clear();
0644     m_cov_ePHI_eT.clear();
0645 
0646     m_cov_eTHETA_eLOC0.clear();
0647     m_cov_eTHETA_eLOC1.clear();
0648     m_cov_eTHETA_ePHI.clear();
0649     m_cov_eTHETA_eTHETA.clear();
0650     m_cov_eTHETA_eQOP.clear();
0651     m_cov_eTHETA_eT.clear();
0652 
0653     m_cov_eQOP_eLOC0.clear();
0654     m_cov_eQOP_eLOC1.clear();
0655     m_cov_eQOP_ePHI.clear();
0656     m_cov_eQOP_eTHETA.clear();
0657     m_cov_eQOP_eQOP.clear();
0658     m_cov_eQOP_eT.clear();
0659 
0660     m_cov_eT_eLOC0.clear();
0661     m_cov_eT_eLOC1.clear();
0662     m_cov_eT_ePHI.clear();
0663     m_cov_eT_eTHETA.clear();
0664     m_cov_eT_eQOP.clear();
0665     m_cov_eT_eT.clear();
0666   }
0667 
0668   m_nUpdatesGx2f.clear();
0669 
0670   return ProcessCode::SUCCESS;
0671 }
0672 
0673 }  // namespace ActsExamples