File indexing completed on 2025-12-17 09:21:53
0001
0002
0003
0004
0005
0006 #include "TrackEvaluation.h"
0007
0008 #include <g4detectors/PHG4CylinderGeomContainer.h>
0009 #include <g4detectors/PHG4TpcGeom.h>
0010 #include <g4detectors/PHG4TpcGeomContainer.h>
0011
0012 #include <g4main/PHG4Hit.h>
0013 #include <g4main/PHG4HitContainer.h>
0014 #include <g4main/PHG4Particle.h>
0015 #include <g4main/PHG4TruthInfoContainer.h>
0016 #include <g4main/PHG4VtxPoint.h>
0017
0018 #include <micromegas/CylinderGeomMicromegas.h>
0019 #include <micromegas/MicromegasDefs.h>
0020
0021 #include <trackbase/ActsGeometry.h>
0022 #include <trackbase/ClusterErrorPara.h>
0023 #include <trackbase/TrkrCluster.h>
0024 #include <trackbase/TrkrClusterContainer.h>
0025 #include <trackbase/TrkrClusterHitAssoc.h>
0026 #include <trackbase/TrkrDefs.h>
0027 #include <trackbase/TrkrHit.h>
0028 #include <trackbase/TrkrHitSet.h>
0029 #include <trackbase/TrkrHitSetContainer.h>
0030 #include <trackbase/TrkrHitTruthAssoc.h>
0031
0032 #include <trackbase_historic/SvtxTrack.h>
0033 #include <trackbase_historic/SvtxTrackMap.h>
0034 #include <trackbase_historic/SvtxTrackState.h>
0035 #include <trackbase_historic/TrackSeed.h>
0036
0037 #include <fun4all/Fun4AllReturnCodes.h>
0038
0039 #include <phool/PHCompositeNode.h>
0040 #include <phool/PHIODataNode.h>
0041 #include <phool/PHNode.h>
0042 #include <phool/PHNodeIterator.h>
0043 #include <phool/PHObject.h>
0044 #include <phool/getClass.h>
0045
0046 #include <TVector3.h>
0047
0048 #include <algorithm>
0049 #include <cassert>
0050 #include <cmath>
0051 #include <cstdint>
0052 #include <cstdlib>
0053 #include <iostream>
0054 #include <limits>
0055 #include <numeric>
0056 #include <vector>
0057
0058
0059 namespace
0060 {
0061
0062
0063 template <class T>
0064 class range_adaptor
0065 {
0066 public:
0067 explicit range_adaptor(const T& range)
0068 : m_range(range)
0069 {
0070 }
0071 const typename T::first_type& begin() { return m_range.first; }
0072 const typename T::second_type& end() { return m_range.second; }
0073
0074 private:
0075 T m_range;
0076 };
0077
0078
0079 template <class T>
0080 constexpr T square(const T& x)
0081 {
0082 return x * x;
0083 }
0084
0085
0086 template <class T>
0087 constexpr T get_r(T x, T y)
0088 {
0089 return std::sqrt(square(x) + square(y));
0090 }
0091
0092
0093 template <class T>
0094 T get_pt(T px, T py)
0095 {
0096 return std::sqrt(square(px) + square(py));
0097 }
0098
0099
0100 template <class T>
0101 T get_p(T px, T py, T pz)
0102 {
0103 return std::sqrt(square(px) + square(py) + square(pz));
0104 }
0105
0106
0107 template <class T>
0108 T get_eta(T p, T pz)
0109 {
0110 return std::log((p + pz) / (p - pz)) / 2;
0111 }
0112
0113
0114 struct interpolation_data_t
0115 {
0116 using list = std::vector<interpolation_data_t>;
0117 double x() const { return position.x(); }
0118 double y() const { return position.y(); }
0119 double z() const { return position.z(); }
0120
0121 double px() const { return momentum.x(); }
0122 double py() const { return momentum.y(); }
0123 double pz() const { return momentum.z(); }
0124
0125 TVector3 position;
0126 TVector3 momentum;
0127 double weight = 1;
0128 };
0129
0130
0131 template <double (interpolation_data_t::*accessor)() const>
0132 double average(const interpolation_data_t::list& hits)
0133 {
0134
0135
0136 double sw = 0;
0137 double swx = 0;
0138
0139 bool valid(false);
0140 for (const auto& hit : hits)
0141 {
0142 const double x = (hit.*accessor)();
0143 if (std::isnan(x))
0144 {
0145 continue;
0146 }
0147
0148 const double w = hit.weight;
0149 if (w <= 0)
0150 {
0151 continue;
0152 }
0153
0154 valid = true;
0155 sw += w;
0156 swx += w * x;
0157 }
0158
0159 if (!valid)
0160 {
0161 return std::numeric_limits<double>::quiet_NaN();
0162 }
0163 return swx / sw;
0164 }
0165
0166
0167 std::vector<TrkrDefs::cluskey> get_cluster_keys(SvtxTrack* track)
0168 {
0169 std::vector<TrkrDefs::cluskey> out;
0170 for (const auto& seed : {track->get_silicon_seed(), track->get_tpc_seed()})
0171 {
0172 if (seed)
0173 {
0174 std::copy(seed->begin_cluster_keys(), seed->end_cluster_keys(), std::back_inserter(out));
0175 }
0176 }
0177
0178 return out;
0179 }
0180
0181
0182 inline int is_primary(PHG4Particle* particle)
0183 {
0184 return particle->get_parent_id() == 0;
0185 }
0186
0187
0188 int64_t get_mask(SvtxTrack* track)
0189 {
0190 const auto cluster_keys = get_cluster_keys(track);
0191 return std::accumulate(cluster_keys.begin(), cluster_keys.end(), int64_t(0),
0192 [](int64_t value, const TrkrDefs::cluskey& key)
0193 {
0194 return TrkrDefs::getLayer(key) < 64 ? (uint64_t) value | (1ULL << TrkrDefs::getLayer(key)) : value;
0195 });
0196 }
0197
0198
0199 template <int type>
0200 int get_clusters(SvtxTrack* track)
0201 {
0202 const auto cluster_keys = get_cluster_keys(track);
0203 return std::count_if(cluster_keys.begin(), cluster_keys.end(),
0204 [](const TrkrDefs::cluskey& key)
0205 { return TrkrDefs::getTrkrId(key) == type; });
0206 }
0207
0208
0209 TrackEvaluationContainerv1::TrackStruct create_track(SvtxTrack* track)
0210 {
0211 TrackEvaluationContainerv1::TrackStruct trackStruct;
0212
0213 trackStruct.charge = track->get_charge();
0214 trackStruct.nclusters = track->size_cluster_keys();
0215 trackStruct.mask = get_mask(track);
0216 trackStruct.nclusters_mvtx = get_clusters<TrkrDefs::mvtxId>(track);
0217 trackStruct.nclusters_intt = get_clusters<TrkrDefs::inttId>(track);
0218 trackStruct.nclusters_tpc = get_clusters<TrkrDefs::tpcId>(track);
0219 trackStruct.nclusters_micromegas = get_clusters<TrkrDefs::micromegasId>(track);
0220
0221 trackStruct.chisquare = track->get_chisq();
0222 trackStruct.ndf = track->get_ndf();
0223
0224 trackStruct.x = track->get_x();
0225 trackStruct.y = track->get_y();
0226 trackStruct.z = track->get_z();
0227 trackStruct.r = get_r(trackStruct.x, trackStruct.y);
0228 trackStruct.phi = std::atan2(trackStruct.y, trackStruct.x);
0229
0230 trackStruct.px = track->get_px();
0231 trackStruct.py = track->get_py();
0232 trackStruct.pz = track->get_pz();
0233 trackStruct.pt = get_pt(trackStruct.px, trackStruct.py);
0234 trackStruct.p = get_p(trackStruct.px, trackStruct.py, trackStruct.pz);
0235 trackStruct.eta = get_eta(trackStruct.p, trackStruct.pz);
0236
0237 return trackStruct;
0238 }
0239
0240
0241 void add_cluster_size(TrackEvaluationContainerv1::ClusterStruct& cluster, TrkrCluster* trk_clus)
0242 {
0243 cluster.size = (unsigned char) trk_clus->getSize();
0244 cluster.phi_size = trk_clus->getPhiSize();
0245 cluster.z_size = trk_clus->getZSize();
0246 cluster.ovlp = (unsigned char) trk_clus->getOverlap();
0247 cluster.edge = (unsigned char) trk_clus->getEdge();
0248 cluster.adc = trk_clus->getAdc();
0249 cluster.max_adc = trk_clus->getMaxAdc();
0250 }
0251
0252
0253 void add_cluster_energy(TrackEvaluationContainerv1::ClusterStruct& cluster, TrkrDefs::cluskey clus_key,
0254 TrkrClusterHitAssoc* cluster_hit_map,
0255 TrkrHitSetContainer* hitsetcontainer)
0256 {
0257
0258 if (!(cluster_hit_map && hitsetcontainer))
0259 {
0260 return;
0261 }
0262
0263
0264 const auto detId = TrkrDefs::getTrkrId(clus_key);
0265 if (detId != TrkrDefs::micromegasId)
0266 {
0267 return;
0268 }
0269
0270 const auto hitset_key = TrkrDefs::getHitSetKeyFromClusKey(clus_key);
0271 auto *const hitset = hitsetcontainer->findHitSet(hitset_key);
0272 if (!hitset)
0273 {
0274 return;
0275 }
0276
0277 const auto range = cluster_hit_map->getHits(clus_key);
0278 cluster.energy_max = 0;
0279 cluster.energy_sum = 0;
0280
0281 for (const auto& pair : range_adaptor(range))
0282 {
0283 auto *const hit = hitset->getHit(pair.second);
0284 if (hit)
0285 {
0286 const auto energy = hit->getEnergy();
0287 cluster.energy_sum += energy;
0288 cluster.energy_max = std::max<double>(energy, cluster.energy_max);
0289 }
0290 }
0291 }
0292
0293
0294 void add_truth_information(TrackEvaluationContainerv1::TrackStruct& track, PHG4Particle* particle, PHG4TruthInfoContainer* truthinfo)
0295 {
0296 if (particle)
0297 {
0298 PHG4VtxPoint* vtx = truthinfo->GetVtx(particle->get_vtx_id());
0299 track.is_primary = is_primary(particle);
0300 track.pid = particle->get_pid();
0301 track.gtrackID = particle->get_track_id();
0302 track.truth_t = vtx->get_t();
0303 track.truth_px = particle->get_px();
0304 track.truth_py = particle->get_py();
0305 track.truth_pz = particle->get_pz();
0306 track.truth_pt = get_pt(track.truth_px, track.truth_py);
0307 track.truth_p = get_p(track.truth_px, track.truth_py, track.truth_pz);
0308 track.truth_eta = get_eta(track.truth_p, track.truth_pz);
0309 }
0310 }
0311
0312
0313 double line_circle_intersection(const TVector3& p0, const TVector3& p1, double radius)
0314 {
0315 const double A = square(p1.x() - p0.x()) + square(p1.y() - p0.y());
0316 const double B = 2 * p0.x() * (p1.x() - p0.x()) + 2 * p0.y() * (p1.y() - p0.y());
0317 const double C = square(p0.x()) + square(p0.y()) - square(radius);
0318 const double delta = square(B) - 4 * A * C;
0319 if (delta < 0)
0320 {
0321 return -1;
0322 }
0323
0324
0325 const double tup = (-B + std::sqrt(delta)) / (2 * A);
0326 if (tup >= 0 && tup < 1)
0327 {
0328 return tup;
0329 }
0330
0331
0332 const double tdn = (-B - sqrt(delta)) / (2 * A);
0333 if (tdn >= 0 && tdn < 1)
0334 {
0335 return tdn;
0336 }
0337
0338
0339 return -1;
0340 }
0341
0342
0343 [[maybe_unused]] inline std::ostream& operator<<(std::ostream& out, const TrackEvaluationContainerv1::ClusterStruct& cluster)
0344 {
0345 out << "ClusterStruct" << std::endl;
0346 out << " cluster: (" << cluster.x << "," << cluster.y << "," << cluster.z << ")" << std::endl;
0347 out << " track: (" << cluster.trk_x << "," << cluster.trk_y << "," << cluster.trk_z << ")" << std::endl;
0348 out << " truth: (" << cluster.truth_x << "," << cluster.truth_y << "," << cluster.truth_z << ")" << std::endl;
0349 return out;
0350 }
0351
0352 [[maybe_unused]] inline std::ostream& operator<<(std::ostream& out, const TVector3& position)
0353 {
0354 out << "(" << position.x() << ", " << position.y() << ", " << position.z() << ")";
0355 return out;
0356 }
0357
0358 }
0359
0360
0361 TrackEvaluation::TrackEvaluation(const std::string& name)
0362 : SubsysReco(name)
0363 {
0364 }
0365
0366
0367 int TrackEvaluation::Init(PHCompositeNode* topNode)
0368 {
0369
0370 PHNodeIterator iter(topNode);
0371 auto *dstNode = dynamic_cast<PHCompositeNode*>(iter.findFirst("PHCompositeNode", "DST"));
0372 if (!dstNode)
0373 {
0374 std::cout << "TrackEvaluation::Init - DST Node missing" << std::endl;
0375 return Fun4AllReturnCodes::ABORTEVENT;
0376 }
0377
0378
0379 iter = PHNodeIterator(dstNode);
0380 auto *evalNode = dynamic_cast<PHCompositeNode*>(iter.findFirst("PHCompositeNode", "EVAL"));
0381 if (!evalNode)
0382 {
0383
0384 std::cout << "TrackEvaluation::Init - EVAL node missing - creating" << std::endl;
0385 evalNode = new PHCompositeNode("EVAL");
0386 dstNode->addNode(evalNode);
0387 }
0388
0389 auto *newNode = new PHIODataNode<PHObject>(new TrackEvaluationContainerv1, "TrackEvaluationContainer", "PHObject");
0390 evalNode->addNode(newNode);
0391
0392 return Fun4AllReturnCodes::EVENT_OK;
0393 }
0394
0395
0396 int TrackEvaluation::InitRun(PHCompositeNode* )
0397 {
0398 return Fun4AllReturnCodes::EVENT_OK;
0399 }
0400
0401
0402 int TrackEvaluation::process_event(PHCompositeNode* topNode)
0403 {
0404
0405 auto res = load_nodes(topNode);
0406 if (res != Fun4AllReturnCodes::EVENT_OK)
0407 {
0408 return res;
0409 }
0410
0411
0412 std::cout << "start..." << std::endl;
0413 if (m_container)
0414 {
0415 m_container->Reset();
0416 }
0417 std::cout << "event..." << std::endl;
0418 if (m_flags & EvalEvent)
0419 {
0420 evaluate_event();
0421 }
0422 std::cout << "clusters..." << std::endl;
0423 if (m_flags & EvalClusters)
0424 {
0425 evaluate_clusters();
0426 }
0427 std::cout << "tracks..." << std::endl;
0428 if (m_flags & EvalTracks)
0429 {
0430 evaluate_tracks();
0431 }
0432 std::cout << "end..." << std::endl;
0433
0434 m_g4hit_map.clear();
0435 return Fun4AllReturnCodes::EVENT_OK;
0436 }
0437
0438
0439 int TrackEvaluation::End(PHCompositeNode* )
0440 {
0441 return Fun4AllReturnCodes::EVENT_OK;
0442 }
0443
0444
0445 int TrackEvaluation::load_nodes(PHCompositeNode* topNode)
0446 {
0447
0448 m_tGeometry = findNode::getClass<ActsGeometry>(topNode, "ActsGeometry");
0449 assert(m_tGeometry);
0450
0451
0452 m_track_map = findNode::getClass<SvtxTrackMap>(topNode, "SvtxTrackMap");
0453
0454
0455 m_cluster_map = findNode::getClass<TrkrClusterContainer>(topNode, "CORRECTED_TRKR_CLUSTER");
0456 if (!m_cluster_map)
0457 {
0458 m_cluster_map = findNode::getClass<TrkrClusterContainer>(topNode, "TRKR_CLUSTER");
0459 }
0460
0461
0462 m_cluster_hit_map = findNode::getClass<TrkrClusterHitAssoc>(topNode, "TRKR_CLUSTERHITASSOC");
0463
0464
0465 m_hit_truth_map = findNode::getClass<TrkrHitTruthAssoc>(topNode, "TRKR_HITTRUTHASSOC");
0466
0467
0468 m_container = findNode::getClass<TrackEvaluationContainerv1>(topNode, "TrackEvaluationContainer");
0469
0470
0471 m_hitsetcontainer = findNode::getClass<TrkrHitSetContainer>(topNode, "TRKR_HITSET");
0472
0473
0474 m_g4hits_tpc = findNode::getClass<PHG4HitContainer>(topNode, "G4HIT_TPC");
0475 m_g4hits_intt = findNode::getClass<PHG4HitContainer>(topNode, "G4HIT_INTT");
0476 m_g4hits_mvtx = findNode::getClass<PHG4HitContainer>(topNode, "G4HIT_MVTX");
0477 m_g4hits_micromegas = findNode::getClass<PHG4HitContainer>(topNode, "G4HIT_MICROMEGAS");
0478
0479
0480 m_g4truthinfo = findNode::getClass<PHG4TruthInfoContainer>(topNode, "G4TruthInfo");
0481
0482
0483 m_tpc_geom_container = findNode::getClass<PHG4TpcGeomContainer>(topNode, "TPCGEOMCONTAINER");
0484 assert(m_tpc_geom_container);
0485
0486
0487 m_micromegas_geom_container = findNode::getClass<PHG4CylinderGeomContainer>(topNode, "CYLINDERGEOM_MICROMEGAS_FULL");
0488
0489 return Fun4AllReturnCodes::EVENT_OK;
0490 }
0491
0492
0493 void TrackEvaluation::evaluate_event()
0494 {
0495 if (!m_container)
0496 {
0497 return;
0498 }
0499
0500
0501 TrackEvaluationContainerv1::EventStruct event;
0502 if (m_cluster_map)
0503 {
0504 for (const auto& hitsetkey : m_cluster_map->getHitSetKeys())
0505 {
0506 const auto trkrId = TrkrDefs::getTrkrId(hitsetkey);
0507 const auto layer = TrkrDefs::getLayer(hitsetkey);
0508 assert(layer < TrackEvaluationContainerv1::EventStruct::max_layer);
0509
0510
0511 const auto clusters = m_cluster_map->getClusters(hitsetkey);
0512 const int nclusters = std::distance(clusters.first, clusters.second);
0513
0514 switch (trkrId)
0515 {
0516 case TrkrDefs::mvtxId:
0517 event.nclusters_mvtx += nclusters;
0518 break;
0519 case TrkrDefs::inttId:
0520 event.nclusters_intt += nclusters;
0521 break;
0522 case TrkrDefs::tpcId:
0523 event.nclusters_tpc += nclusters;
0524 break;
0525 case TrkrDefs::micromegasId:
0526 event.nclusters_micromegas += nclusters;
0527 break;
0528 default:
0529 std::cout << "unknown trkrId: " << trkrId << std::endl;
0530 exit(1);
0531 }
0532
0533 event.nclusters[layer] += nclusters;
0534 }
0535 }
0536
0537
0538 m_container->addEvent(event);
0539 }
0540
0541
0542 void TrackEvaluation::evaluate_clusters()
0543 {
0544 if (!(m_cluster_map && m_hitsetcontainer && m_container))
0545 {
0546 return;
0547 }
0548
0549
0550 m_container->clearClusters();
0551 SvtxTrack* track = nullptr;
0552
0553 for (const auto& hitsetkey : m_cluster_map->getHitSetKeys())
0554 {
0555 for (const auto& [key, cluster] : range_adaptor(m_cluster_map->getClusters(hitsetkey)))
0556 {
0557
0558 auto cluster_struct = create_cluster(key, cluster, track);
0559 add_cluster_size(cluster_struct, cluster);
0560 add_cluster_energy(cluster_struct, key, m_cluster_hit_map, m_hitsetcontainer);
0561
0562
0563 const auto g4hits = find_g4hits(key);
0564 const bool is_micromegas(TrkrDefs::getTrkrId(key) == TrkrDefs::micromegasId);
0565 if (is_micromegas)
0566 {
0567 const int tileid = MicromegasDefs::getTileId(key);
0568 add_truth_information_micromegas(cluster_struct, tileid, g4hits);
0569 }
0570 else
0571 {
0572 add_truth_information(cluster_struct, g4hits);
0573 }
0574
0575
0576 m_container->addCluster(cluster_struct);
0577 }
0578 }
0579 }
0580
0581
0582 void TrackEvaluation::evaluate_tracks()
0583 {
0584 if (!(m_track_map && m_cluster_map && m_container))
0585 {
0586 return;
0587 }
0588
0589
0590 m_container->clearTracks();
0591 for (const auto& [track_id, track] : *m_track_map)
0592 {
0593 auto track_struct = create_track(track);
0594
0595
0596 const auto [id, contributors] = get_max_contributor(track);
0597 track_struct.contributors = contributors;
0598
0599
0600 auto *particle = m_g4truthinfo->GetParticle(id);
0601 track_struct.embed = get_embed(particle);
0602 ::add_truth_information(track_struct, particle, m_g4truthinfo);
0603
0604
0605 auto state_iter = track->begin_states();
0606
0607 for (const auto& cluster_key : get_cluster_keys(track))
0608 {
0609 auto *cluster = m_cluster_map->findCluster(cluster_key);
0610 if (!cluster)
0611 {
0612 std::cout << "TrackEvaluation::evaluate_tracks - unable to find cluster for key " << cluster_key << std::endl;
0613 continue;
0614 }
0615
0616 auto cluster_struct = create_cluster(cluster_key, cluster, track);
0617 add_cluster_size(cluster_struct, cluster);
0618 add_cluster_energy(cluster_struct, cluster_key, m_cluster_hit_map, m_hitsetcontainer);
0619
0620 const auto g4hits = find_g4hits(cluster_key);
0621 const bool is_micromegas(TrkrDefs::getTrkrId(cluster_key) == TrkrDefs::micromegasId);
0622 if (is_micromegas)
0623 {
0624 const int tileid = MicromegasDefs::getTileId(cluster_key);
0625 add_truth_information_micromegas(cluster_struct, tileid, g4hits);
0626 }
0627 else
0628 {
0629 add_truth_information(cluster_struct, g4hits);
0630 }
0631
0632
0633
0634 const auto radius(cluster_struct.r);
0635 float dr_min = -1;
0636 for (auto iter = state_iter; iter != track->end_states(); ++iter)
0637 {
0638 const auto dr = std::abs(radius - get_r(iter->second->get_x(), iter->second->get_y()));
0639 if (dr_min < 0 || dr < dr_min)
0640 {
0641 state_iter = iter;
0642 dr_min = dr;
0643 }
0644 else
0645 {
0646 break;
0647 }
0648 }
0649
0650
0651 if (is_micromegas)
0652 {
0653 const int tileid = MicromegasDefs::getTileId(cluster_key);
0654 add_trk_information_micromegas(cluster_struct, tileid, state_iter->second);
0655 }
0656 else
0657 {
0658 add_trk_information(cluster_struct, state_iter->second);
0659 }
0660
0661
0662 track_struct.clusters.push_back(cluster_struct);
0663 }
0664 m_container->addTrack(track_struct);
0665 }
0666 }
0667
0668
0669 TrackEvaluation::G4HitSet TrackEvaluation::find_g4hits(TrkrDefs::cluskey cluster_key) const
0670 {
0671
0672 if (!(m_cluster_hit_map && m_hit_truth_map))
0673 {
0674 return G4HitSet();
0675 }
0676
0677
0678 auto map_iter = m_g4hit_map.lower_bound(cluster_key);
0679 if (map_iter != m_g4hit_map.end() && cluster_key == map_iter->first)
0680 {
0681 return map_iter->second;
0682 }
0683
0684
0685 G4HitSet out;
0686 const auto hitset_key = TrkrDefs::getHitSetKeyFromClusKey(cluster_key);
0687 for (const auto& [first, hit_key] : range_adaptor(m_cluster_hit_map->getHits(cluster_key)))
0688 {
0689
0690 TrkrHitTruthAssoc::MMap g4hit_map;
0691 m_hit_truth_map->getG4Hits(hitset_key, hit_key, g4hit_map);
0692
0693
0694 for (auto& truth_iter : g4hit_map)
0695 {
0696 const auto g4hit_key = truth_iter.second.second;
0697 PHG4Hit* g4hit = nullptr;
0698
0699 switch (TrkrDefs::getTrkrId(hitset_key))
0700 {
0701 case TrkrDefs::mvtxId:
0702 if (m_g4hits_mvtx)
0703 {
0704 g4hit = m_g4hits_mvtx->findHit(g4hit_key);
0705 }
0706 break;
0707
0708 case TrkrDefs::inttId:
0709 if (m_g4hits_intt)
0710 {
0711 g4hit = m_g4hits_intt->findHit(g4hit_key);
0712 }
0713 break;
0714
0715 case TrkrDefs::tpcId:
0716 if (m_g4hits_tpc)
0717 {
0718 g4hit = m_g4hits_tpc->findHit(g4hit_key);
0719 }
0720 break;
0721
0722 case TrkrDefs::micromegasId:
0723 if (m_g4hits_micromegas)
0724 {
0725 g4hit = m_g4hits_micromegas->findHit(g4hit_key);
0726 }
0727 break;
0728
0729 default:
0730 break;
0731 }
0732
0733 if (g4hit)
0734 {
0735 out.insert(g4hit);
0736 }
0737 else
0738 {
0739 std::cout << "TrackEvaluation::find_g4hits - g4hit not found " << g4hit_key << std::endl;
0740 }
0741 }
0742 }
0743
0744
0745 return m_g4hit_map.insert(map_iter, std::make_pair(cluster_key, std::move(out)))->second;
0746 }
0747
0748
0749 std::pair<int, int> TrackEvaluation::get_max_contributor(SvtxTrack* track) const
0750 {
0751 if (!(m_track_map && m_cluster_map && m_g4truthinfo))
0752 {
0753 return {0, 0};
0754 }
0755
0756
0757 using IdMap = std::map<int, int>;
0758 IdMap contributor_map;
0759
0760
0761 for (const auto& cluster_key : get_cluster_keys(track))
0762 {
0763 for (const auto& hit : find_g4hits(cluster_key))
0764 {
0765 const int trkid = hit->get_trkid();
0766 auto iter = contributor_map.lower_bound(trkid);
0767 if (iter == contributor_map.end() || iter->first != trkid)
0768 {
0769 contributor_map.insert(iter, std::make_pair(trkid, 1));
0770 }
0771 else
0772 {
0773 ++iter->second;
0774 }
0775 }
0776 }
0777
0778 if (contributor_map.empty())
0779 {
0780 return {0, 0};
0781 }
0782
0783 return *std::max_element(
0784 contributor_map.cbegin(), contributor_map.cend(),
0785 [](const IdMap::value_type& first, const IdMap::value_type& second)
0786 { return first.second < second.second; });
0787 }
0788
0789
0790 int TrackEvaluation::get_embed(PHG4Particle* particle) const
0791 {
0792 return (m_g4truthinfo && particle) ? m_g4truthinfo->isEmbeded(particle->get_primary_id()) : 0;
0793 }
0794
0795
0796 TrackEvaluationContainerv1::ClusterStruct TrackEvaluation::create_cluster(TrkrDefs::cluskey key, TrkrCluster* cluster, SvtxTrack* track) const
0797 {
0798 TrackSeed* si_seed = nullptr;
0799 TrackSeed* tpc_seed = nullptr;
0800 if (track != nullptr)
0801 {
0802 si_seed = track->get_silicon_seed();
0803 tpc_seed = track->get_tpc_seed();
0804 }
0805
0806 Acts::Vector3 global;
0807 global = m_tGeometry->getGlobalPosition(key, cluster);
0808 TrackEvaluationContainerv1::ClusterStruct cluster_struct;
0809 cluster_struct.layer = TrkrDefs::getLayer(key);
0810 cluster_struct.x = global.x();
0811 cluster_struct.y = global.y();
0812 cluster_struct.z = global.z();
0813 cluster_struct.r = get_r(cluster_struct.x, cluster_struct.y);
0814 cluster_struct.phi = std::atan2(cluster_struct.y, cluster_struct.x);
0815 cluster_struct.phi_error = 0.0;
0816 cluster_struct.z_error = 0.0;
0817 cluster_struct.trk_alpha = 0.0;
0818 cluster_struct.trk_beta = 0.0;
0819
0820 ClusterErrorPara ClusErrPara;
0821
0822 if (track != nullptr)
0823 {
0824 float r = cluster_struct.r;
0825
0826 if (cluster_struct.layer >= 7)
0827 {
0828 auto para_errors_mm = ClusterErrorPara::get_clusterv5_modified_error(cluster, r, key);
0829
0830 cluster_struct.phi_error = cluster->getRPhiError() / cluster_struct.r;
0831 cluster_struct.z_error = cluster->getZError();
0832 cluster_struct.para_phi_error = sqrt(para_errors_mm.first) / cluster_struct.r;
0833 cluster_struct.para_z_error = sqrt(para_errors_mm.second);
0834
0835 cluster_struct.trk_radius = 1.0 / tpc_seed->get_qOverR();
0836 cluster_struct.trk_alpha = (r * r) / (2 * r * std::abs(1.0 / tpc_seed->get_qOverR()));
0837 cluster_struct.trk_beta = std::atan(tpc_seed->get_slope());
0838 }
0839 else
0840 {
0841 auto para_errors_mvtx = ClusErrPara.get_cluster_error(cluster, r, key, si_seed->get_qOverR(), si_seed->get_slope());
0842 cluster_struct.phi_error = sqrt(para_errors_mvtx.first) / cluster_struct.r;
0843 cluster_struct.z_error = sqrt(para_errors_mvtx.second);
0844
0845 cluster_struct.trk_radius = 1.0 / tpc_seed->get_qOverR();
0846 cluster_struct.trk_alpha = (r * r) / (2 * r * std::abs(1.0 / tpc_seed->get_qOverR()));
0847 cluster_struct.trk_beta = std::atan(si_seed->get_slope());
0848 }
0849 }
0850
0851 return cluster_struct;
0852 }
0853
0854
0855 void TrackEvaluation::add_trk_information(TrackEvaluationContainerv1::ClusterStruct& cluster, SvtxTrackState* state) const
0856 {
0857
0858 const auto trk_r = get_r(state->get_x(), state->get_y());
0859 const auto dr = cluster.r - trk_r;
0860 const auto trk_drdt = (state->get_x() * state->get_px() + state->get_y() * state->get_py()) / trk_r;
0861 const auto trk_dxdr = state->get_px() / trk_drdt;
0862 const auto trk_dydr = state->get_py() / trk_drdt;
0863 const auto trk_dzdr = state->get_pz() / trk_drdt;
0864
0865
0866 cluster.trk_x = state->get_x() + dr * trk_dxdr;
0867 cluster.trk_y = state->get_y() + dr * trk_dydr;
0868 cluster.trk_z = state->get_z() + dr * trk_dzdr;
0869 cluster.trk_r = get_r(cluster.trk_x, cluster.trk_y);
0870 cluster.trk_phi = std::atan2(cluster.trk_y, cluster.trk_x);
0871
0872
0873 cluster.trk_px = state->get_px();
0874 cluster.trk_py = state->get_py();
0875 cluster.trk_pz = state->get_pz();
0876
0877
0878
0879
0880
0881 const auto cosphi(std::cos(cluster.trk_phi));
0882 const auto sinphi(std::sin(cluster.trk_phi));
0883 const auto trk_pphi = -state->get_px() * sinphi + state->get_py() * cosphi;
0884 const auto trk_pr = state->get_px() * cosphi + state->get_py() * sinphi;
0885 const auto trk_pz = state->get_pz();
0886 cluster.trk_alpha = std::atan2(trk_pphi, trk_pr);
0887 cluster.trk_beta = std::atan2(trk_pz, trk_pr);
0888 cluster.trk_phi_error = state->get_phi_error();
0889 cluster.trk_z_error = state->get_z_error();
0890 }
0891
0892
0893 void TrackEvaluation::add_trk_information_micromegas(TrackEvaluationContainerv1::ClusterStruct& cluster, int tileid, SvtxTrackState* state) const
0894 {
0895
0896 const auto layer = cluster.layer;
0897 auto *const layergeom = dynamic_cast<CylinderGeomMicromegas*>(m_micromegas_geom_container->GetLayerGeom(layer));
0898 assert(layergeom);
0899
0900
0901 const TVector3 cluster_world(cluster.x, cluster.y, cluster.z);
0902 const auto cluster_local = layergeom->get_local_from_world_coords(tileid, m_tGeometry, cluster_world);
0903
0904
0905 TVector3 track_world(state->get_x(), state->get_y(), state->get_z());
0906 TVector3 track_local = layergeom->get_local_from_world_coords(tileid, m_tGeometry, track_world);
0907
0908
0909 const TVector3 direction_world(state->get_px(), state->get_py(), state->get_pz());
0910 const TVector3 direction_local = layergeom->get_local_from_world_vect(tileid, m_tGeometry, direction_world);
0911
0912
0913 const auto delta_z = cluster_local.z() - track_local.z();
0914 track_local += TVector3(
0915 delta_z * direction_local.x() / direction_local.z(),
0916 delta_z * direction_local.y() / direction_local.z(),
0917 delta_z);
0918
0919
0920 track_world = layergeom->get_world_from_local_coords(tileid, m_tGeometry, track_local);
0921
0922
0923 cluster.trk_x = track_world.x();
0924 cluster.trk_y = track_world.y();
0925 cluster.trk_z = track_world.z();
0926 cluster.trk_r = get_r(cluster.trk_x, cluster.trk_y);
0927 cluster.trk_phi = std::atan2(cluster.trk_y, cluster.trk_x);
0928
0929
0930 cluster.trk_px = state->get_px();
0931 cluster.trk_py = state->get_py();
0932 cluster.trk_pz = state->get_pz();
0933
0934
0935
0936
0937
0938 const auto cosphi(std::cos(cluster.trk_phi));
0939 const auto sinphi(std::sin(cluster.trk_phi));
0940 const auto trk_pphi = -state->get_px() * sinphi + state->get_py() * cosphi;
0941 const auto trk_pr = state->get_px() * cosphi + state->get_py() * sinphi;
0942 const auto trk_pz = state->get_pz();
0943 cluster.trk_alpha = std::atan2(trk_pphi, trk_pr);
0944 cluster.trk_beta = std::atan2(trk_pz, trk_pr);
0945 cluster.trk_phi_error = state->get_phi_error();
0946 cluster.trk_z_error = state->get_z_error();
0947 }
0948
0949
0950 void TrackEvaluation::add_truth_information(TrackEvaluationContainerv1::ClusterStruct& cluster, const std::set<PHG4Hit*>& g4hits) const
0951 {
0952
0953 cluster.truth_size = g4hits.size();
0954
0955
0956 const auto layer = cluster.layer;
0957 const bool is_tpc(layer >= 7 && layer < 55);
0958 const PHG4TpcGeom* layergeom = is_tpc ? m_tpc_geom_container->GetLayerCellGeom(layer) : nullptr;
0959 const auto rin = layergeom ? layergeom->get_radius() - layergeom->get_thickness() / 2 : 0;
0960 const auto rout = layergeom ? layergeom->get_radius() + layergeom->get_thickness() / 2 : 0;
0961
0962
0963 interpolation_data_t::list hits;
0964 for (const auto& g4hit : g4hits)
0965 {
0966 interpolation_data_t::list tmp_hits;
0967 const auto weight = g4hit->get_edep();
0968 for (int i = 0; i < 2; ++i)
0969 {
0970 const TVector3 g4hit_world(g4hit->get_x(i), g4hit->get_y(i), g4hit->get_z(i));
0971 const TVector3 momentum_world(g4hit->get_px(i), g4hit->get_py(i), g4hit->get_pz(i));
0972 tmp_hits.push_back({.position = g4hit_world, .momentum = momentum_world, .weight = weight});
0973 }
0974
0975 if (is_tpc)
0976 {
0977
0978
0979 auto r0 = get_r(tmp_hits[0].x(), tmp_hits[0].y());
0980 auto r1 = get_r(tmp_hits[1].x(), tmp_hits[1].y());
0981 if (r0 > r1)
0982 {
0983 std::swap(tmp_hits[0], tmp_hits[1]);
0984 std::swap(r0, r1);
0985 }
0986
0987
0988 if (r1 <= rin || r0 >= rout)
0989 {
0990 continue;
0991 }
0992
0993
0994 const auto dr_old = r1 - r0;
0995
0996
0997 if (r0 < rin)
0998 {
0999 const auto t = line_circle_intersection(tmp_hits[0].position, tmp_hits[1].position, rin);
1000 if (t < 0)
1001 {
1002 continue;
1003 }
1004
1005 tmp_hits[0].position = tmp_hits[0].position * (1. - t) + tmp_hits[1].position * t;
1006 tmp_hits[0].momentum = tmp_hits[0].momentum * (1. - t) + tmp_hits[1].momentum * t;
1007 r0 = rin;
1008 }
1009
1010 if (r1 > rout)
1011 {
1012 const auto t = line_circle_intersection(tmp_hits[0].position, tmp_hits[1].position, rout);
1013 if (t < 0)
1014 {
1015 continue;
1016 }
1017
1018 tmp_hits[1].position = tmp_hits[0].position * (1. - t) + tmp_hits[1].position * t;
1019 tmp_hits[1].momentum = tmp_hits[0].momentum * (1. - t) + tmp_hits[1].momentum * t;
1020 r1 = rout;
1021 }
1022
1023
1024 const auto dr_new = r1 - r0;
1025 tmp_hits[0].weight *= dr_new / dr_old;
1026 tmp_hits[1].weight *= dr_new / dr_old;
1027 }
1028
1029
1030 hits.push_back(std::move(tmp_hits[0]));
1031 hits.push_back(std::move(tmp_hits[1]));
1032 }
1033
1034
1035 cluster.truth_x = average<&interpolation_data_t::x>(hits);
1036 cluster.truth_y = average<&interpolation_data_t::y>(hits);
1037 cluster.truth_z = average<&interpolation_data_t::z>(hits);
1038
1039
1040 cluster.truth_px = average<&interpolation_data_t::px>(hits);
1041 cluster.truth_py = average<&interpolation_data_t::py>(hits);
1042 cluster.truth_pz = average<&interpolation_data_t::pz>(hits);
1043
1044 cluster.truth_r = get_r(cluster.truth_x, cluster.truth_y);
1045 cluster.truth_phi = std::atan2(cluster.truth_y, cluster.truth_x);
1046
1047
1048
1049
1050
1051 const auto cosphi(std::cos(cluster.truth_phi));
1052 const auto sinphi(std::sin(cluster.truth_phi));
1053 const auto truth_pphi = -cluster.truth_px * sinphi + cluster.truth_py * cosphi;
1054 const auto truth_pr = cluster.truth_px * cosphi + cluster.truth_py * sinphi;
1055
1056 cluster.truth_alpha = std::atan2(truth_pphi, truth_pr);
1057 cluster.truth_beta = std::atan2(cluster.truth_pz, truth_pr);
1058 }
1059
1060
1061 void TrackEvaluation::add_truth_information_micromegas(TrackEvaluationContainerv1::ClusterStruct& cluster, int tileid, const std::set<PHG4Hit*>& g4hits) const
1062 {
1063
1064 cluster.truth_size = g4hits.size();
1065
1066 const auto layer = cluster.layer;
1067 auto *const layergeom = dynamic_cast<CylinderGeomMicromegas*>(m_micromegas_geom_container->GetLayerGeom(layer));
1068 assert(layergeom);
1069
1070
1071 const TVector3 cluster_world(cluster.x, cluster.y, cluster.z);
1072 const auto cluster_local = layergeom->get_local_from_world_coords(tileid, m_tGeometry, cluster_world);
1073
1074
1075 interpolation_data_t::list hits;
1076 for (const auto& g4hit : g4hits)
1077 {
1078 const auto weight = g4hit->get_edep();
1079 for (int i = 0; i < 2; ++i)
1080 {
1081
1082 TVector3 g4hit_world(g4hit->get_x(i), g4hit->get_y(i), g4hit->get_z(i));
1083 auto g4hit_local = layergeom->get_local_from_world_coords(tileid, m_tGeometry, g4hit_world);
1084
1085
1086 TVector3 momentum_world(g4hit->get_px(i), g4hit->get_py(i), g4hit->get_pz(i));
1087 TVector3 momentum_local = layergeom->get_local_from_world_vect(tileid, m_tGeometry, momentum_world);
1088
1089 hits.push_back({.position = g4hit_local, .momentum = momentum_local, .weight = weight});
1090 }
1091 }
1092
1093
1094 const TVector3 interpolation_local(
1095 average<&interpolation_data_t::x>(hits),
1096 average<&interpolation_data_t::y>(hits),
1097 average<&interpolation_data_t::z>(hits));
1098
1099
1100 const TVector3 interpolation_world = layergeom->get_world_from_local_coords(tileid, m_tGeometry, interpolation_local);
1101
1102
1103 const TVector3 momentum_local(
1104 average<&interpolation_data_t::px>(hits),
1105 average<&interpolation_data_t::py>(hits),
1106 average<&interpolation_data_t::pz>(hits));
1107
1108
1109 const TVector3 momentum_world = layergeom->get_world_from_local_vect(tileid, m_tGeometry, momentum_local);
1110
1111 cluster.truth_x = interpolation_world.x();
1112 cluster.truth_y = interpolation_world.y();
1113 cluster.truth_z = interpolation_world.z();
1114 cluster.truth_r = get_r(cluster.truth_x, cluster.truth_y);
1115 cluster.truth_phi = std::atan2(cluster.truth_y, cluster.truth_x);
1116
1117
1118 cluster.truth_px = momentum_world.x();
1119 cluster.truth_py = momentum_world.y();
1120 cluster.truth_pz = momentum_world.z();
1121
1122
1123
1124
1125
1126 const auto cosphi(std::cos(cluster.truth_phi));
1127 const auto sinphi(std::sin(cluster.truth_phi));
1128 const auto truth_pphi = -cluster.truth_px * sinphi + cluster.truth_py * cosphi;
1129 const auto truth_pr = cluster.truth_px * cosphi + cluster.truth_py * sinphi;
1130
1131 cluster.truth_alpha = std::atan2(truth_pphi, truth_pr);
1132 cluster.truth_beta = std::atan2(cluster.truth_pz, truth_pr);
1133 }