File indexing completed on 2025-08-06 08:18:32
0001
0002
0003
0004
0005
0006
0007 #include "PHSimpleKFProp.h"
0008 #include "ALICEKF.h"
0009 #include "GPUTPCTrackLinearisation.h"
0010 #include "GPUTPCTrackParam.h"
0011 #include "PHGhostRejection.h"
0012 #include "nanoflann.hpp"
0013
0014 #include <ffamodules/CDBInterface.h>
0015
0016 #include <g4detectors/PHG4TpcCylinderGeom.h>
0017 #include <g4detectors/PHG4TpcCylinderGeomContainer.h>
0018
0019 #include <phfield/PHField.h>
0020 #include <phfield/PHFieldConfig.h>
0021 #include <phfield/PHFieldConfigv1.h>
0022 #include <phfield/PHFieldUtility.h>
0023
0024
0025 #include <tpc/TpcDistortionCorrectionContainer.h>
0026
0027 #include <trackbase/ActsGeometry.h>
0028 #include <trackbase/TrackFitUtils.h>
0029 #include <trackbase/TrkrCluster.h>
0030 #include <trackbase/TrkrClusterContainer.h>
0031 #include <trackbase/TrkrClusterIterationMapv1.h>
0032 #include <trackbase/TrkrDefs.h>
0033
0034 #include <trackbase_historic/ActsTransformations.h>
0035 #include <trackbase_historic/TrackSeedContainer.h>
0036 #include <trackbase_historic/TrackSeed_v2.h>
0037 #include <trackbase_historic/TrackSeedHelper.h>
0038
0039 #include <fun4all/Fun4AllReturnCodes.h>
0040
0041 #include <phool/PHTimer.h>
0042 #include <phool/getClass.h>
0043 #include <phool/phool.h> // for PHWHERE
0044
0045 #include <Acts/MagneticField/InterpolatedBFieldMap.hpp>
0046 #include <Acts/MagneticField/MagneticFieldProvider.hpp>
0047
0048 #include <TSystem.h>
0049
0050 #include <Geant4/G4SystemOfUnits.hh>
0051
0052 #include <Eigen/Core>
0053 #include <Eigen/Dense>
0054
0055 #include <omp.h>
0056
0057 #include <filesystem>
0058 #include <iostream>
0059 #include <syncstream>
0060 #include <vector>
0061
0062
0063 namespace
0064 {
0065
0066 template <class T>
0067 inline constexpr T square(const T& x)
0068 {
0069 return x * x;
0070 }
0071 }
0072
0073 using keylist = std::vector<TrkrDefs::cluskey>;
0074
0075 PHSimpleKFProp::PHSimpleKFProp(const std::string& name)
0076 : SubsysReco(name)
0077 {}
0078
0079
0080 PHSimpleKFProp::~PHSimpleKFProp()
0081 {
0082 if( m_own_fieldmap )
0083 { delete _field_map; }
0084 }
0085
0086
0087 int PHSimpleKFProp::End(PHCompositeNode* )
0088 {
0089 return Fun4AllReturnCodes::EVENT_OK;
0090 }
0091
0092 int PHSimpleKFProp::InitRun(PHCompositeNode* topNode)
0093 {
0094 int ret = get_nodes(topNode);
0095 if (ret != Fun4AllReturnCodes::EVENT_OK)
0096 {
0097 return ret;
0098 }
0099
0100 PHFieldConfigv1 fcfg;
0101 fcfg.set_field_config(PHFieldConfig::FieldConfigTypes::Field3DCartesian);
0102
0103 if (std::filesystem::path(m_magField).extension() != ".root")
0104 { m_magField = CDBInterface::instance()->getUrl(m_magField); }
0105
0106 if (!_use_const_field)
0107 {
0108 if (!std::filesystem::exists(m_magField))
0109 {
0110 if (m_magField.empty())
0111 { m_magField = "empty string"; }
0112 std::cout << PHWHERE << "Fieldmap " << m_magField << " does not exist" << std::endl;
0113 gSystem->Exit(1);
0114 }
0115
0116 fcfg.set_filename(m_magField);
0117 }
0118 else
0119 {
0120 fcfg.set_field_config(PHFieldConfig::FieldConfigTypes::kFieldUniform);
0121 fcfg.set_magfield_rescale(_const_field);
0122 }
0123
0124
0125
0126
0127
0128
0129
0130 const auto node_fcfg = PHFieldUtility::GetFieldConfigNode(&fcfg, topNode);
0131 if( fcfg == *node_fcfg )
0132 {
0133
0134 std::cout << "PHSimpleKFProp::InitRun - using node tree field map" << std::endl;
0135 _field_map = PHFieldUtility::GetFieldMapNode(&fcfg, topNode);
0136 m_own_fieldmap = false;
0137 } else {
0138
0139 std::cout << "PHSimpleKFProp::InitRun - using own field map" << std::endl;
0140 _field_map = PHFieldUtility::BuildFieldMap(&fcfg);
0141 m_own_fieldmap = true;
0142 }
0143
0144
0145 fitter = std::make_unique<ALICEKF>(topNode, _cluster_map, _field_map, _fieldDir, _min_clusters_per_track, _max_sin_phi, Verbosity());
0146 fitter->setNeonFraction(Ne_frac);
0147 fitter->setArgonFraction(Ar_frac);
0148 fitter->setCF4Fraction(CF4_frac);
0149 fitter->setNitrogenFraction(N2_frac);
0150 fitter->setIsobutaneFraction(isobutane_frac);
0151 fitter->useConstBField(_use_const_field);
0152 fitter->setConstBField(_const_field);
0153 fitter->useFixedClusterError(_use_fixed_clus_err);
0154 fitter->setFixedClusterError(0, _fixed_clus_err.at(0));
0155 fitter->setFixedClusterError(1, _fixed_clus_err.at(1));
0156 fitter->setFixedClusterError(2, _fixed_clus_err.at(2));
0157
0158
0159
0160
0161 std::cout << "PHSimpleKFProp::InitRun - m_num_threads: " << m_num_threads << std::endl;
0162 if( m_num_threads >= 1 ) { omp_set_num_threads( m_num_threads ); }
0163
0164 return Fun4AllReturnCodes::EVENT_OK;
0165 }
0166
0167 double PHSimpleKFProp::get_Bz(double x, double y, double z) const
0168 {
0169 if (_use_const_field || std::abs(z) > 105.5)
0170 {
0171 return _const_field;
0172 }
0173 double p[4] = {x * cm, y * cm, z * cm, 0. * cm};
0174 double bfield[3];
0175
0176
0177 if( omp_get_thread_num() == 0 )
0178 {
0179 _field_map->GetFieldValue(p, bfield);
0180 } else {
0181 _field_map->GetFieldValue_nocache(p, bfield);
0182 }
0183
0184
0185
0186
0187
0188
0189
0190
0191
0192
0193
0194
0195
0196 return bfield[2] / tesla;
0197 }
0198
0199 int PHSimpleKFProp::get_nodes(PHCompositeNode* topNode)
0200 {
0201
0202
0203 m_tgeometry = findNode::getClass<ActsGeometry>(topNode, "ActsGeometry");
0204 if (!m_tgeometry)
0205 {
0206 std::cout << "No Acts tracking geometry, exiting." << std::endl;
0207 return Fun4AllReturnCodes::ABORTEVENT;
0208 }
0209
0210
0211 m_globalPositionWrapper.loadNodes(topNode);
0212
0213
0214 if (_use_truth_clusters)
0215 {
0216 _cluster_map = findNode::getClass<TrkrClusterContainer>(topNode, "TRKR_CLUSTER_TRUTH");
0217 }
0218 else
0219 {
0220 _cluster_map = findNode::getClass<TrkrClusterContainer>(topNode, "TRKR_CLUSTER");
0221 }
0222 if (!_cluster_map)
0223 {
0224 std::cerr << PHWHERE << " ERROR: Can't find node TRKR_CLUSTER" << std::endl;
0225 return Fun4AllReturnCodes::ABORTEVENT;
0226 }
0227
0228
0229 _track_map = findNode::getClass<TrackSeedContainer>(topNode, "TpcTrackSeedContainer");
0230 if (!_track_map)
0231 {
0232 std::cerr << PHWHERE << " ERROR: Can't find TrackSeedContainer " << std::endl;
0233 return Fun4AllReturnCodes::ABORTEVENT;
0234 }
0235
0236
0237 auto geom_container = findNode::getClass<PHG4TpcCylinderGeomContainer>(topNode, "CYLINDERCELLGEOM_SVTX");
0238 if (!geom_container)
0239 {
0240 std::cerr << PHWHERE << "ERROR: Can't find node CYLINDERCELLGEOM_SVTX" << std::endl;
0241 return Fun4AllReturnCodes::ABORTRUN;
0242 }
0243
0244 for (int i = 7; i <= 54; i++)
0245 {
0246 radii.push_back(geom_container->GetLayerCellGeom(i)->get_radius());
0247 }
0248
0249 return Fun4AllReturnCodes::EVENT_OK;
0250 }
0251
0252 int PHSimpleKFProp::process_event(PHCompositeNode* topNode)
0253 {
0254 if (_n_iteration != 0)
0255 {
0256 _iteration_map = findNode::getClass<TrkrClusterIterationMapv1>(topNode, "CLUSTER_ITERATION_MAP");
0257 if (!_iteration_map)
0258 {
0259 std::cerr << PHWHERE << "Cluster Iteration Map missing, aborting." << std::endl;
0260 return Fun4AllReturnCodes::ABORTEVENT;
0261 }
0262 }
0263
0264 PHTimer timer("KFPropTimer");
0265 timer.restart();
0266
0267
0268 if(_max_seeds > 0 && _track_map->size() > _max_seeds)
0269 {
0270 std::cout << "number of TPC seeds: " << _track_map->size() << std::endl;
0271 std::cout << PHWHERE << "number of TPC seeds > " << _max_seeds << " aborting event." << std::endl;
0272 return Fun4AllReturnCodes::ABORTEVENT;
0273 }
0274
0275 const auto globalPositions = PrepareKDTrees();
0276 if (Verbosity())
0277 { std::cout << "PHSimpleKFProp::process_event - PrepareKDTrees time: " << timer.elapsed() << " ms" << std::endl; }
0278
0279
0280 std::vector<std::vector<TrkrDefs::cluskey>> new_chains;
0281 std::vector<TrackSeed_v2> unused_tracks;
0282
0283 timer.restart();
0284 #pragma omp parallel
0285 {
0286 if (Verbosity())
0287 {
0288 std::osyncstream(std::cout)
0289 << "PHSimpleKFProp -"
0290 << " num_threads: " << omp_get_num_threads()
0291 << " this thread: " << omp_get_thread_num()
0292 << std::endl;
0293 }
0294
0295 PHTimer timer_mp("KFPropTimer_parallel");
0296
0297 std::vector<std::vector<TrkrDefs::cluskey>> local_chains;
0298 std::vector<TrackSeed_v2> local_unused;
0299
0300 #pragma omp for schedule(static)
0301 for (size_t track_it = 0; track_it != _track_map->size(); ++track_it)
0302 {
0303 if (Verbosity())
0304 {
0305 std::osyncstream(std::cout)
0306 << "PHSimpleKFProp -"
0307 << " num_threads: " << omp_get_num_threads()
0308 << " this thread: " << omp_get_thread_num()
0309 << " processing seed " << track_it << std::endl;
0310 }
0311
0312
0313 auto track = _track_map->get(track_it);
0314 const bool is_tpc = std::any_of(
0315 track->begin_cluster_keys(),
0316 track->end_cluster_keys(),
0317 [](const TrkrDefs::cluskey& key)
0318 { return TrkrDefs::getTrkrId(key) == TrkrDefs::tpcId; });
0319
0320 if (is_tpc)
0321 {
0322
0323
0324 std::vector<std::vector<TrkrDefs::cluskey>> keylist_A(1);
0325 std::copy(track->begin_cluster_keys(), track->end_cluster_keys(), std::back_inserter(keylist_A[0]));
0326
0327
0328 std::map<TrkrDefs::cluskey, Acts::Vector3> trackClusPositions;
0329 std::transform(track->begin_cluster_keys(), track->end_cluster_keys(), std::inserter(trackClusPositions, trackClusPositions.end()),
0330 [globalPositions](const auto& key)
0331 { return std::make_pair(key, globalPositions.at(key)); });
0332
0333
0334 if (keylist_A[0].size() < 3)
0335 {
0336 continue;
0337 }
0338
0339
0340
0341 std::vector<float> trackChi2;
0342
0343 timer_mp.restart();
0344 auto seedpair = fitter->ALICEKalmanFilter(keylist_A, false, trackClusPositions, trackChi2);
0345
0346 if (Verbosity() > 3)
0347 {
0348 std::cout << "PHSimpleKFProp::process_event - single track ALICEKF time " << timer_mp.elapsed() << " ms" << std::endl;
0349 }
0350
0351 timer_mp.restart();
0352
0353
0354 TrackSeedHelper::circleFitByTaubin(track, trackClusPositions, 7, 55);
0355 TrackSeedHelper::lineFit(track, trackClusPositions, 7, 55);
0356 track->set_phi(TrackSeedHelper::get_phi(track, trackClusPositions));
0357 if (Verbosity() > 3)
0358 {
0359 std::cout << "PHSimpleKFProp::process_event - single track circle fit time " << timer_mp.elapsed() << " ms" << std::endl;
0360 }
0361
0362 if (seedpair.first.empty()|| seedpair.second.empty())
0363 {
0364 continue;
0365 }
0366
0367 if (Verbosity())
0368 {
0369 std::cout << "is tpc track" << std::endl;
0370 }
0371
0372 timer_mp.restart();
0373
0374 if (Verbosity())
0375 {
0376 std::cout << "propagate first round" << std::endl;
0377 }
0378
0379 auto preseed = PropagateTrack(track, PropagationDirection::Inward, seedpair.second.at(0), globalPositions);
0380 if (Verbosity())
0381 {
0382 std::cout << "preseed size " << preseed.size() << std::endl;
0383 }
0384
0385 std::vector<std::vector<TrkrDefs::cluskey>> kl = {preseed};
0386 if (Verbosity())
0387 {
0388 std::cout << "kl size " << kl.size() << std::endl;
0389 }
0390 std::vector<float> pretrackChi2;
0391
0392 auto prepair = fitter->ALICEKalmanFilter(kl, false, globalPositions, pretrackChi2);
0393 if (prepair.first.empty() || prepair.second.empty())
0394 {
0395 continue;
0396 }
0397
0398 std::reverse(kl.at(0).begin(), kl.at(0).end());
0399
0400 auto pretrack = prepair.first.at(0);
0401
0402
0403 std::map<TrkrDefs::cluskey, Acts::Vector3> pretrackClusPositions;
0404 std::transform(pretrack.begin_cluster_keys(), pretrack.end_cluster_keys(), std::inserter(pretrackClusPositions, pretrackClusPositions.end()),
0405 [globalPositions](const auto& key)
0406 { return std::make_pair(key, globalPositions.at(key)); });
0407
0408
0409 TrackSeedHelper::circleFitByTaubin(&pretrack,pretrackClusPositions, 7, 55);
0410 TrackSeedHelper::lineFit(&pretrack, pretrackClusPositions, 7, 55);
0411 pretrack.set_phi(TrackSeedHelper::get_phi(&pretrack, pretrackClusPositions));
0412
0413 prepair.second.at(0).SetDzDs(-prepair.second.at(0).GetDzDs());
0414 const auto finalchain = PropagateTrack(&pretrack, kl.at(0), PropagationDirection::Outward, prepair.second.at(0), globalPositions);
0415
0416 if (finalchain.size() > kl.at(0).size())
0417 {
0418 local_chains.push_back(std::move(finalchain));
0419 }
0420 else
0421 {
0422 local_chains.push_back(std::move(kl.at(0)));
0423 }
0424
0425 if (Verbosity() > 3)
0426 {
0427 std::cout << "PHSimpleKFProp::process_event - propagate track time " << timer_mp.elapsed() << " ms" << std::endl;
0428 }
0429 }
0430 else
0431 {
0432 if (Verbosity())
0433 {
0434 std::cout << "is NOT tpc track" << std::endl;
0435 }
0436 local_unused.emplace_back(*track);
0437 }
0438 }
0439
0440
0441 std::sort(local_chains.begin(),local_chains.end());
0442 local_chains.erase(std::unique(local_chains.begin(),local_chains.end()),local_chains.end());
0443
0444
0445 #pragma omp critical
0446 {
0447 new_chains.reserve(new_chains.size()+local_chains.size());
0448 new_chains.insert(new_chains.end(), std::make_move_iterator(local_chains.begin()), std::make_move_iterator(local_chains.end()));
0449
0450 unused_tracks.reserve(unused_tracks.size()+local_unused.size());
0451 unused_tracks.insert(unused_tracks.end(), std::make_move_iterator(local_unused.begin()), std::make_move_iterator(local_unused.end()));
0452 }
0453 }
0454 if (Verbosity())
0455 { std::cout << "PHSimpleKFProp::process_event - first seed loop time: " << timer.elapsed() << " ms" << std::endl; }
0456
0457
0458 timer.restart();
0459 std::sort(new_chains.begin(),new_chains.end());
0460 new_chains.erase(std::unique(new_chains.begin(),new_chains.end()),new_chains.end());
0461
0462 if (Verbosity())
0463 { std::cout << "PHSimpleKFProp::process_event - first cleanup time: " << timer.elapsed() << " ms" << std::endl; }
0464
0465 if( Verbosity() )
0466 {
0467
0468 std::cout << "PHSimpleKFProp::process_event - new_chains size: " << new_chains.size() << std::endl;
0469 for( const auto& chain:new_chains )
0470 {
0471 std::cout << "PHSimpleKFProp::process_event - { ";
0472 for( const auto& key:chain )
0473 {
0474 std::cout << key << " ";
0475 }
0476 std::cout << "}" << std::endl << std::endl;
0477 }
0478 }
0479
0480
0481 new_chains.erase( std::remove_if( new_chains.begin(), new_chains.end(),
0482 [](const auto& chain) { return chain.size()<3; } ),
0483 new_chains.end() );
0484
0485
0486 timer.restart();
0487 std::vector<float> trackChi2;
0488 auto seeds = fitter->ALICEKalmanFilter(new_chains, true, globalPositions, trackChi2);
0489 if (Verbosity())
0490 { std::cout << "PHSimpleKFProp::process_event - ALICEKalmanFilter time: " << timer.elapsed() << " ms" << std::endl; }
0491
0492
0493 _track_map->Reset();
0494
0495
0496 timer.restart();
0497 if (m_ghostrejection)
0498 {
0499 rejectAndPublishSeeds(seeds.first, globalPositions, trackChi2);
0500 }
0501 else
0502 {
0503 publishSeeds(seeds.first);
0504 }
0505
0506
0507 publishSeeds(unused_tracks);
0508
0509 if (Verbosity())
0510 { std::cout << "PHSimpleKFProp::process_event - publishSeeds time: " << timer.elapsed() << " ms" << std::endl; }
0511
0512 return Fun4AllReturnCodes::EVENT_OK;
0513 }
0514
0515 Acts::Vector3 PHSimpleKFProp::getGlobalPosition(TrkrDefs::cluskey key, TrkrCluster* cluster) const
0516 {
0517
0518 return _pp_mode ?
0519 m_tgeometry->getGlobalPosition(key, cluster):
0520 m_globalPositionWrapper.getGlobalPositionDistortionCorrected( key, cluster, 0 );
0521 }
0522
0523 PositionMap PHSimpleKFProp::PrepareKDTrees()
0524 {
0525 PositionMap globalPositions;
0526
0527 std::vector<std::vector<std::vector<double>>> kdhits;
0528 kdhits.resize(58);
0529 if (!_cluster_map)
0530 {
0531 std::cout << "WARNING: (tracking.PHTpcTrackerUtil.convert_clusters_to_hits) cluster map is not provided" << std::endl;
0532 return globalPositions;
0533 }
0534
0535 for (const auto& hitsetkey : _cluster_map->getHitSetKeys(TrkrDefs::TrkrId::tpcId))
0536 {
0537 auto range = _cluster_map->getClusters(hitsetkey);
0538 for (TrkrClusterContainer::ConstIterator it = range.first; it != range.second; ++it)
0539 {
0540 const auto& [cluskey,cluster] = *it;
0541 if (!cluster)
0542 {
0543 continue;
0544 }
0545
0546
0547 if (_n_iteration && _iteration_map && _iteration_map->getIteration(cluskey) > 0)
0548 {
0549 continue;
0550 }
0551
0552 const auto globalpos = getGlobalPosition(cluskey, cluster);
0553 globalPositions.emplace(cluskey, globalpos);
0554
0555 const int layer = TrkrDefs::getLayer(cluskey);
0556 std::vector<double> kdhit{ globalpos.x(), globalpos.y(), globalpos.z(), 0 };
0557 const uint64_t key = cluskey;
0558 std::memcpy(&kdhit[3], &key, sizeof(key));
0559
0560
0561
0562
0563
0564
0565 kdhits[layer].push_back(std::move(kdhit));
0566 }
0567 }
0568 _ptclouds.resize(kdhits.size());
0569 _kdtrees.resize(kdhits.size());
0570 for (size_t l = 0; l < kdhits.size(); ++l)
0571 {
0572 if (Verbosity() > 1)
0573 {
0574 std::cout << "l: " << l << std::endl;
0575 }
0576 _ptclouds[l] = std::make_shared<KDPointCloud<double>>();
0577 _ptclouds[l]->pts = std::move(kdhits[l]);
0578
0579 _kdtrees[l] = std::make_shared<nanoflann::KDTreeSingleIndexAdaptor<nanoflann::L2_Simple_Adaptor<double, KDPointCloud<double>>, KDPointCloud<double>, 3>>(3, *(_ptclouds[l]), nanoflann::KDTreeSingleIndexAdaptorParams(10));
0580 _kdtrees[l]->buildIndex();
0581 }
0582
0583 return globalPositions;
0584 }
0585
0586 bool PHSimpleKFProp::TransportAndRotate(
0587 double old_radius,
0588 double new_radius,
0589 double& phi,
0590 GPUTPCTrackParam& kftrack,
0591 GPUTPCTrackParam::GPUTPCTrackFitParam& fp) const
0592 {
0593 if (Verbosity() > 2)
0594 {
0595 std::cout << "old_radius " << old_radius << ", new_radius " << new_radius << std::endl;
0596 }
0597 const float transport_spacing = .05;
0598 const int Ndivisions = floor(std::abs(new_radius - old_radius) / transport_spacing);
0599 if (Verbosity() > 2)
0600 {
0601 std::cout << "Ndivisions: " << Ndivisions << std::endl;
0602 }
0603 for (int i = 1; i <= Ndivisions + 1; i++)
0604 {
0605 if (std::isnan(kftrack.GetX()) ||
0606 std::isnan(kftrack.GetY()) ||
0607 std::isnan(kftrack.GetZ()))
0608 {
0609 if (Verbosity() > 1)
0610 {
0611 std::cout << "position is NaN, exiting" << std::endl;
0612 }
0613 return false;
0614 }
0615
0616 if (new_radius > 78.)
0617 {
0618 if (Verbosity() > 1)
0619 {
0620 std::cout << "outside TPC, exiting" << std::endl;
0621 }
0622 return false;
0623 }
0624
0625 double r_div;
0626 if (i == Ndivisions + 1)
0627 {
0628 r_div = new_radius;
0629 }
0630 else if (old_radius < new_radius)
0631 {
0632 r_div = old_radius + transport_spacing * i;
0633 }
0634 else
0635 {
0636 r_div = old_radius - transport_spacing * i;
0637 }
0638
0639 if (Verbosity() > 2)
0640 {
0641 std::cout << "transporting to intermediate radius " << r_div << std::endl;
0642 }
0643
0644
0645 const double tX = kftrack.GetX();
0646 const double tY = kftrack.GetY();
0647 const double tz = kftrack.GetZ();
0648
0649
0650 const double tx = tX * cos(phi) - tY * sin(phi);
0651 const double ty = tX * sin(phi) + tY * cos(phi);
0652
0653 const double Bz = _Bzconst * get_Bz(tx, ty, tz);
0654
0655 kftrack.CalculateFitParameters(fp);
0656
0657
0658 if (!kftrack.TransportToXWithMaterial(r_div, fp, Bz, 1.))
0659 {
0660 if (Verbosity() > 1)
0661 {
0662 std::cout << "transport failed" << std::endl;
0663 }
0664 return false;
0665 }
0666
0667 const double new_tX = kftrack.GetX();
0668 const double new_tY = kftrack.GetY();
0669
0670 const double new_tx = new_tX * cos(phi) - new_tY * sin(phi);
0671 const double new_ty = new_tX * sin(phi) + new_tY * cos(phi);
0672 const double new_phi = atan2(new_ty, new_tx);
0673 const double alpha = new_phi - phi;
0674
0675 if (!kftrack.Rotate(alpha, 1.))
0676 {
0677 if (Verbosity() > 1)
0678 {
0679 std::cout << "rotate failed" << std::endl;
0680 }
0681 return false;
0682 }
0683 phi = new_phi;
0684 }
0685
0686 if (Verbosity() > 1)
0687 {
0688 const double final_tX = kftrack.GetX();
0689 const double final_tY = kftrack.GetY();
0690 const double final_tx = final_tX * cos(phi) - final_tY * sin(phi);
0691 const double final_ty = final_tX * sin(phi) + final_tY * cos(phi);
0692 const double final_tz = kftrack.GetZ();
0693 std::cout << "track position after transport: (" << final_tx << ", " << final_ty << ", " << final_tz << ")" << std::endl;
0694 }
0695 return true;
0696 }
0697
0698
0699 bool PHSimpleKFProp::PropagateStep(
0700 unsigned int& current_layer,
0701 double& current_phi,
0702 PropagationDirection& direction,
0703 std::vector<TrkrDefs::cluskey>& propagated_track,
0704 const std::vector<TrkrDefs::cluskey>& ckeys,
0705 GPUTPCTrackParam& kftrack,
0706 GPUTPCTrackParam::GPUTPCTrackFitParam& fp,
0707 const PositionMap& globalPositions) const
0708 {
0709
0710 if (std::isnan(kftrack.GetX()) ||
0711 std::isnan(kftrack.GetY()) ||
0712 std::isnan(kftrack.GetZ()))
0713 {
0714 if (Verbosity() > 1)
0715 {
0716 std::cout << "position is NaN, exiting loop" << std::endl;
0717 }
0718 return false;
0719 }
0720
0721 if (Verbosity() > 1)
0722 {
0723 std::cout << "current layer: " << current_layer << std::endl;
0724 }
0725 if (Verbosity() > 5)
0726 {
0727 std::cout << "original seed size: " << ckeys.size() << std::endl;
0728 }
0729
0730
0731
0732 int next_layer;
0733 if (direction == PropagationDirection::Outward)
0734 {
0735 next_layer = current_layer + 1;
0736 }
0737 else
0738 {
0739 next_layer = current_layer - 1;
0740 }
0741
0742 if (Verbosity() > 1)
0743 {
0744 std::cout << "next layer: " << next_layer << std::endl;
0745 }
0746
0747
0748 if (next_layer < 7 || next_layer > 54)
0749 {
0750 if (Verbosity() > 1)
0751 {
0752 std::cout << "reached end of TPC, exiting loop" << std::endl;
0753 }
0754 return false;
0755 }
0756
0757
0758 const double tX = kftrack.GetX();
0759 const double tY = kftrack.GetY();
0760 const double tz = kftrack.GetZ();
0761
0762
0763 const double tx = tX * cos(current_phi) - tY * sin(current_phi);
0764 const double ty = tX * sin(current_phi) + tY * cos(current_phi);
0765
0766 if (Verbosity() > 1)
0767 {
0768 std::cout << std::endl;
0769 std::cout << "track position: (" << tx << ", " << ty << ", " << tz << ")" << std::endl;
0770 }
0771
0772 if (Verbosity() > 1)
0773 {
0774 std::cout << "transporting to radius: " << radii[next_layer - 7] << std::endl;
0775 }
0776
0777 if (!TransportAndRotate(kftrack.GetX(), radii[next_layer - 7], current_phi, kftrack, fp))
0778 {
0779 if (std::isnan(kftrack.GetX()) ||
0780 std::isnan(kftrack.GetY()) ||
0781 std::isnan(kftrack.GetZ()))
0782 {
0783 return false;
0784 }
0785
0786 if (Verbosity() > 1)
0787 {
0788 std::cout << "track turned around near X=" << kftrack.GetX() << std::endl;
0789 }
0790
0791
0792
0793
0794
0795
0796 double pt = 1. / std::abs(kftrack.GetQPt());
0797 double end_tx = kftrack.GetX() * cos(current_phi) - kftrack.GetY() * sin(current_phi);
0798 double end_ty = kftrack.GetX() * sin(current_phi) + kftrack.GetY() * cos(current_phi);
0799 double end_tz = kftrack.GetZ();
0800 if (Verbosity() > 1)
0801 {
0802 std::cout << "current parameters: pt=" << pt << ", (x, y, z) = (" << end_tx << ", " << end_ty << ", " << end_tz << ")" << std::endl;
0803 }
0804
0805 double R = 100. * pt / (0.3 * get_Bz(end_tx, end_ty, end_tz));
0806 if (Verbosity() > 1)
0807 {
0808 std::cout << "R=" << R << std::endl;
0809 }
0810 double pX = pt * kftrack.GetCosPhi();
0811 double pY = pt * kftrack.GetSinPhi();
0812 if (Verbosity() > 1)
0813 {
0814 std::cout << "(pX, pY) = (" << pX << ", " << pY << ")" << std::endl;
0815 }
0816 double px = pX * cos(current_phi) - pY * sin(current_phi);
0817 double py = pX * sin(current_phi) + pY * cos(current_phi);
0818 double tangent_phi = std::atan2(py, px);
0819 double center_phi = 0;
0820 if (kftrack.GetQPt() > 0)
0821 {
0822 center_phi = tangent_phi + M_PI / 2.;
0823 }
0824 else
0825 {
0826 center_phi = tangent_phi - M_PI / 2.;
0827 }
0828 if (center_phi > M_PI)
0829 {
0830 center_phi -= 2. * M_PI;
0831 }
0832 if (center_phi < -M_PI)
0833 {
0834 center_phi += 2. * M_PI;
0835 }
0836 double xc = end_tx - R * cos(center_phi);
0837 double yc = end_ty - R * sin(center_phi);
0838 if (Verbosity() > 1)
0839 {
0840 std::cout << "(xc, yc) = (" << xc << ", " << yc << ")" << std::endl;
0841 }
0842 auto circle_output = TrackFitUtils::circle_circle_intersection(
0843 std::sqrt(square(kftrack.GetX()) + square(kftrack.GetY())), R, xc, yc);
0844
0845 double new_tx;
0846 double new_ty;
0847 double xplus = std::get<0>(circle_output);
0848 double yplus = std::get<1>(circle_output);
0849 double xminus = std::get<2>(circle_output);
0850 double yminus = std::get<3>(circle_output);
0851 if (Verbosity() > 1)
0852 {
0853 std::cout << "circle-circle intersection: (" << xplus << ", " << yplus << "), (" << xminus << ", " << yminus << ")" << std::endl;
0854 }
0855
0856 if (std::sqrt(square(end_tx - xplus) + square(end_ty - yplus)) > std::sqrt(square(end_tx - xminus) + square(end_ty - yminus)))
0857 {
0858 new_tx = xplus;
0859 new_ty = yplus;
0860 }
0861 else
0862 {
0863 new_tx = xminus;
0864 new_ty = yminus;
0865 }
0866
0867 if (Verbosity() > 1)
0868 {
0869 std::cout << "track now at (" << new_tx << ", " << new_ty << ")" << std::endl;
0870 }
0871 const double rot_phi = atan2(new_ty, new_tx);
0872
0873
0874
0875 const double new_tX = new_tx * cos(rot_phi) + new_ty * sin(rot_phi);
0876 const double new_tY = -new_tx * sin(rot_phi) + new_ty * cos(rot_phi);
0877 const double new_centerphi = atan2(new_ty - yc, new_tx - xc);
0878 double dcenterphi = new_centerphi - center_phi;
0879 if (dcenterphi > M_PI)
0880 {
0881 dcenterphi = 2. * M_PI - dcenterphi;
0882 }
0883 if (dcenterphi < -M_PI)
0884 {
0885 dcenterphi = 2. * M_PI + dcenterphi;
0886 }
0887 const double ds = R * std::abs(dcenterphi);
0888 const double dz = kftrack.GetDzDs() * ds;
0889
0890 current_phi = rot_phi;
0891 kftrack.SetX(new_tX);
0892 kftrack.SetY(new_tY);
0893 kftrack.SetZ(end_tz + dz);
0894
0895 kftrack.SetSignCosPhi(-kftrack.GetSignCosPhi());
0896
0897
0898 if (!TransportAndRotate(kftrack.GetX(), radii[current_layer - 7], current_phi, kftrack, fp))
0899 {
0900 return false;
0901 }
0902
0903 if (direction == PropagationDirection::Outward)
0904 {
0905 direction = PropagationDirection::Inward;
0906 }
0907 else
0908 {
0909 direction = PropagationDirection::Outward;
0910 }
0911
0912
0913 next_layer = current_layer;
0914 }
0915
0916
0917 if(!_pp_mode)
0918 {
0919 if (Verbosity() > 1)
0920 {
0921 std::cout << "doing distortion correction" << std::endl;
0922 }
0923
0924
0925
0926
0927
0928 const double uncorr_tX = kftrack.GetX();
0929 const double uncorr_tY = kftrack.GetY();
0930 const double uncorr_tx = uncorr_tX * cos(current_phi) - uncorr_tY * sin(current_phi);
0931 const double uncorr_ty = uncorr_tX * sin(current_phi) + uncorr_tY * cos(current_phi);
0932 const double uncorr_tz = kftrack.GetZ();
0933
0934 Acts::Vector3 proj_pt(uncorr_tx, uncorr_ty, uncorr_tz);
0935
0936 const double proj_radius = sqrt(square(uncorr_tx) + square(uncorr_ty));
0937 if (proj_radius > 78.0)
0938 {
0939
0940 return false;
0941 }
0942
0943 if (Verbosity() > 2)
0944 {
0945 std::cout << "distortion correction for (" << tx << ", " << ty << ", " << tz << "), layer " << next_layer << ", radius " << proj_radius << std::endl;
0946 }
0947
0948
0949 proj_pt = m_globalPositionWrapper.applyDistortionCorrections(proj_pt);
0950
0951
0952 const double dist_radius = sqrt(square(proj_pt[0]) + square(proj_pt[1]));
0953
0954 if (Verbosity() > 2)
0955 {
0956 std::cout << "after correction: (" << proj_pt[0] << ", " << proj_pt[1] << ", " << proj_pt[2] << "), radius " << dist_radius << std::endl;
0957 }
0958
0959 if (!TransportAndRotate(kftrack.GetX(), dist_radius, current_phi, kftrack, fp))
0960 {
0961 return false;
0962 }
0963
0964 if (std::isnan(kftrack.GetX()) ||
0965 std::isnan(kftrack.GetY()) ||
0966 std::isnan(kftrack.GetZ()))
0967 {
0968 return false;
0969 }
0970 }
0971
0972 const double new_tX = kftrack.GetX();
0973 const double new_tY = kftrack.GetY();
0974 const double new_tx = new_tX * cos(current_phi) - new_tY * sin(current_phi);
0975 const double new_ty = new_tX * sin(current_phi) + new_tY * cos(current_phi);
0976 const double new_tz = kftrack.GetZ();
0977
0978
0979 double query_pt[3] = {new_tx, new_ty, new_tz};
0980 std::vector<long unsigned int> index_out(1);
0981 std::vector<double> distance_out(1);
0982 int n_results = _kdtrees[next_layer]->knnSearch(&query_pt[0], 1, &index_out[0], &distance_out[0]);
0983
0984
0985 if (!n_results)
0986 {
0987 if (Verbosity() > 1)
0988 {
0989 std::cout << "no clusters found in search window, moving on" << std::endl;
0990 }
0991 current_layer = next_layer;
0992 return true;
0993 }
0994 const std::vector<double>& point = _ptclouds[next_layer]->pts[index_out[0]];
0995 TrkrDefs::cluskey closest_ckey = (*((int64_t*) &point[3]));
0996 TrkrCluster* clusterCandidate = _cluster_map->findCluster(closest_ckey);
0997 const auto candidate_globalpos = globalPositions.at(closest_ckey);
0998 const double cand_x = candidate_globalpos(0);
0999 const double cand_y = candidate_globalpos(1);
1000 const double cand_z = candidate_globalpos(2);
1001
1002 if (Verbosity() > 1)
1003 {
1004 std::cout << "found closest cluster candidate at (" << cand_x << ", " << cand_y << ", " << cand_z << ")" << std::endl;
1005 }
1006
1007
1008 const double tYerr = std::sqrt(kftrack.GetCov(0));
1009 const double txerr = std::abs(tYerr * sin(current_phi));
1010 const double tyerr = std::abs(tYerr * cos(current_phi));
1011 const double tzerr = sqrt(kftrack.GetCov(5));
1012
1013 if (Verbosity() > 1)
1014 {
1015 std::cout << "track state position errors: (" << txerr << ", " << tyerr << ", " << tzerr << ")" << std::endl;
1016 }
1017
1018 const double cand_xerr = sqrt(fitter->getClusterError(clusterCandidate, closest_ckey, candidate_globalpos, 0, 0));
1019 const double cand_yerr = sqrt(fitter->getClusterError(clusterCandidate, closest_ckey, candidate_globalpos, 1, 1));
1020 const double cand_zerr = sqrt(fitter->getClusterError(clusterCandidate, closest_ckey, candidate_globalpos, 2, 2));
1021
1022 if (Verbosity() > 1)
1023 {
1024 std::cout << "cluster position errors: (" << cand_xerr << ", " << cand_yerr << ", " << cand_zerr << ")" << std::endl;
1025 }
1026
1027
1028 if (std::abs(new_tx - cand_x) < _max_dist * std::sqrt(square(txerr) + square(cand_xerr)) &&
1029 std::abs(new_ty - cand_y) < _max_dist * std::sqrt(square(tyerr) + square(cand_yerr)) &&
1030 std::abs(new_tz - cand_z) < _max_dist * std::sqrt(square(tzerr) + square(cand_zerr)))
1031 {
1032 if (Verbosity() > 1)
1033 {
1034 std::cout << "added cluster" << std::endl;
1035 }
1036 propagated_track.push_back(closest_ckey);
1037
1038
1039 if (std::find(ckeys.begin(), ckeys.end(), closest_ckey) == ckeys.end())
1040 {
1041 const double cand_Y = -cand_x * std::sin(current_phi) + cand_y * std::cos(current_phi);
1042 const double cand_xycov2 = fitter->getClusterError(clusterCandidate, closest_ckey, candidate_globalpos, 0, 1);
1043 const double cand_Yerr2 = cand_xerr * cand_xerr * sin(current_phi) * sin(current_phi) + cand_xycov2 * sin(current_phi) * cos(current_phi) + cand_yerr * cand_yerr * cos(current_phi) * cos(current_phi);
1044 const double cand_zerr2 = cand_zerr * cand_zerr;
1045
1046 if (Verbosity() > 1)
1047 {
1048 std::cout << "Filtering cluster with Y=" << cand_Y << ", z=" << cand_z << ", Yerr2=" << cand_Yerr2 << ", zerr2=" << cand_zerr2 << std::endl;
1049 }
1050
1051 kftrack.Filter(cand_Y, cand_z, cand_Yerr2, cand_zerr2, _max_sin_phi);
1052 }
1053 }
1054
1055
1056 current_layer = next_layer;
1057
1058 if (Verbosity() > 1)
1059 {
1060 std::cout << "track current parameters:" << std::endl;
1061 std::cout << "(X, Y, Z) = (" << kftrack.GetX() << ", " << kftrack.GetY() << ", " << kftrack.GetZ() << ")" << std::endl;
1062 std::cout << "pt = " << 1. / fabs(kftrack.GetQPt()) << std::endl;
1063 std::cout << "sinPhi = " << kftrack.GetSinPhi() << std::endl;
1064 std::cout << "cosPhi = " << kftrack.GetCosPhi() << std::endl;
1065 std::cout << "dzds = " << kftrack.GetDzDs() << std::endl;
1066 }
1067
1068
1069 return true;
1070 }
1071
1072
1073 std::vector<TrkrDefs::cluskey> PHSimpleKFProp::PropagateTrack(TrackSeed* track, PropagationDirection direction, GPUTPCTrackParam& aliceSeed, const PositionMap& globalPositions) const
1074 {
1075
1076 std::vector<TrkrDefs::cluskey> ckeys;
1077 if (direction == PropagationDirection::Inward)
1078 {
1079 std::reverse_copy(track->begin_cluster_keys(), track->end_cluster_keys(), std::back_inserter(ckeys));
1080 } else {
1081 std::copy(track->begin_cluster_keys(), track->end_cluster_keys(), std::back_inserter(ckeys));
1082 }
1083 return PropagateTrack(track, ckeys, direction, aliceSeed, globalPositions);
1084 }
1085
1086
1087 std::vector<TrkrDefs::cluskey> PHSimpleKFProp::PropagateTrack(TrackSeed* track, std::vector<TrkrDefs::cluskey>& ckeys, PropagationDirection direction, GPUTPCTrackParam& aliceSeed, const PositionMap& globalPositions) const
1088 {
1089 if (direction == PropagationDirection::Inward)
1090 {
1091 aliceSeed.SetDzDs(-aliceSeed.GetDzDs());
1092 }
1093
1094 const auto track_pos = TrackSeedHelper::get_xyz(track);
1095 double track_x = track_pos.x();
1096 double track_y = track_pos.y();
1097 double track_z = track_pos.z();
1098
1099 if (Verbosity() > 1)
1100 {
1101 std::cout << "layers:" << std::endl;
1102 for (auto& c : ckeys)
1103 {
1104 std::cout << (int) TrkrDefs::getLayer(c) << ", ";
1105 }
1106 std::cout << std::endl;
1107 }
1108
1109 if (Verbosity() > 1)
1110 {
1111 std::cout << "track (x,y,z) = (" << track_x << ", " << track_y << ", " << track_z << ")" << std::endl;
1112 }
1113
1114 double track_px = NAN;
1115 double track_py = NAN;
1116 double track_pz = NAN;
1117 if (_use_const_field)
1118 {
1119 float pt = fabs(1. / track->get_qOverR()) * (0.3 / 100) * _const_field;
1120 float phi = track->get_phi();
1121 track_px = pt * std::cos(phi);
1122 track_py = pt * std::sin(phi);
1123 track_pz = pt * std::cosh(track->get_eta()) * std::cos(track->get_theta());
1124 }
1125 else
1126 {
1127 track_px = track->get_px();
1128 track_py = track->get_py();
1129 track_pz = track->get_pz();
1130 }
1131
1132 if (Verbosity() > 1)
1133 {
1134 std::cout << "track (px,py,pz) = (" << track_px << ", " << track_py << ", " << track_pz << ")" << std::endl;
1135 }
1136
1137 std::vector<Acts::Vector3> trkGlobPos;
1138 for (const auto& ckey : ckeys)
1139 {
1140 if (TrkrDefs::getTrkrId(ckey) == TrkrDefs::tpcId)
1141 {
1142 trkGlobPos.push_back(globalPositions.at(ckey));
1143 }
1144 }
1145
1146
1147 size_t inner_index;
1148 size_t second_index;
1149 if (TrkrDefs::getLayer(ckeys[0]) > TrkrDefs::getLayer(ckeys.back()))
1150 {
1151 inner_index = ckeys.size() - 1;
1152 second_index = ckeys.size() - 2;
1153 }
1154 else
1155 {
1156 inner_index = 0;
1157 second_index = 1;
1158 }
1159
1160 double xc = track->get_X0();
1161 double yc = track->get_Y0();
1162 double cluster_x = trkGlobPos.at(inner_index)(0);
1163 double cluster_y = trkGlobPos.at(inner_index)(1);
1164 double dy = cluster_y - yc;
1165 double dx = cluster_x - xc;
1166 double phi = atan2(dy, dx);
1167 double second_dx = trkGlobPos.at(second_index)(0) - xc;
1168 double second_dy = trkGlobPos.at(second_index)(1) - yc;
1169 double second_phi = atan2(second_dy, second_dx);
1170 double dphi = second_phi - phi;
1171
1172 if (dphi > 0)
1173 {
1174 phi += M_PI / 2.0;
1175 }
1176 else
1177 {
1178 phi -= M_PI / 2.0;
1179 }
1180
1181 double pt = sqrt(track_px * track_px + track_py * track_py);
1182
1183 track_px = pt * cos(phi);
1184 track_py = pt * sin(phi);
1185 track_x = trkGlobPos.at(0)(0);
1186 track_y = trkGlobPos.at(0)(1);
1187 track_z = trkGlobPos.at(0)(2);
1188
1189 if (Verbosity() > 1)
1190 {
1191 std::cout << "new track (x,y,z) = " << track_x << ", " << track_y << ", " << track_z << ")" << std::endl;
1192 }
1193 if (Verbosity() > 1)
1194 {
1195 std::cout << "new track (px,py,pz) = " << track_px << ", " << track_py << ", " << track_pz << ")" << std::endl;
1196 }
1197
1198
1199
1200 double track_pt = sqrt(track_px * track_px + track_py * track_py);
1201 if (Verbosity() > 1)
1202 {
1203 std::cout << "track pt: " << track_pt << std::endl;
1204 }
1205
1206
1207 GPUTPCTrackParam kftrack{};
1208 kftrack.InitParam();
1209 kftrack.setNeonFraction(Ne_frac);
1210 kftrack.setArgonFraction(Ar_frac);
1211 kftrack.setCF4Fraction(CF4_frac);
1212 kftrack.setNitrogenFraction(N2_frac);
1213 kftrack.setIsobutaneFraction(isobutane_frac);
1214
1215 float track_phi = atan2(track_y, track_x);
1216 if (Verbosity() > 1)
1217 {
1218 std::cout << "track phi: " << track_phi << std::endl;
1219 std::cout << "track charge: " << track->get_charge() << std::endl;
1220 }
1221 kftrack.SetQPt(track->get_charge() / track_pt);
1222
1223 float track_pX = track_px * cos(track_phi) + track_py * sin(track_phi);
1224 float track_pY = -track_px * sin(track_phi) + track_py * cos(track_phi);
1225
1226 kftrack.SetSignCosPhi(track_pX / track_pt);
1227 kftrack.SetSinPhi(track_pY / track_pt);
1228 kftrack.SetDzDs(track_pz / sqrt(track_pt * track_pt + track_pz * track_pz));
1229
1230 if (Verbosity() > 1)
1231 {
1232 std::cout << "track pX = " << track_pX << ", pY = " << track_pY << ", CosPhi = " << track_pX / track_pt << ", signCosPhi = " << kftrack.GetSignCosPhi() << std::endl;
1233 }
1234
1235
1236
1237
1238
1239
1240
1241
1242
1243
1244
1245
1246
1247
1248
1249
1250
1251
1252
1253
1254
1255
1256
1257
1258
1259
1260
1261
1262
1263
1264
1265
1266
1267
1268
1269
1270
1271
1272
1273
1274
1275
1276
1277
1278
1279
1280
1281
1282
1283
1284
1285
1286
1287
1288
1289
1290
1291
1292
1293
1294 std::vector<TrkrDefs::cluskey> propagated_track;
1295
1296 kftrack.SetX(track_x * cos(track_phi) + track_y * sin(track_phi));
1297 kftrack.SetY(-track_x * sin(track_phi) + track_y * cos(track_phi));
1298 kftrack.SetZ(track_z);
1299
1300 if (kftrack.GetSignCosPhi() < 0)
1301 {
1302 kftrack.SetSinPhi(-kftrack.GetSinPhi());
1303 kftrack.SetDzDs(-kftrack.GetDzDs());
1304
1305 kftrack.SetCov(3, -kftrack.GetCov(3));
1306 kftrack.SetCov(4, -kftrack.GetCov(4));
1307 kftrack.SetCov(6, -kftrack.GetCov(6));
1308 kftrack.SetCov(7, -kftrack.GetCov(7));
1309 kftrack.SetCov(10, -kftrack.GetCov(10));
1310 kftrack.SetCov(11, -kftrack.GetCov(11));
1311 }
1312
1313 if (Verbosity() > 1)
1314 {
1315 std::cout << "initial track params:" << std::endl;
1316 std::cout << "X: " << kftrack.GetX() << std::endl;
1317 std::cout << "Y: " << kftrack.GetY() << std::endl;
1318 std::cout << "Z: " << kftrack.GetZ() << std::endl;
1319 std::cout << "SinPhi: " << kftrack.GetSinPhi() << std::endl;
1320 std::cout << "DzDs: " << kftrack.GetDzDs() << std::endl;
1321 std::cout << "QPt: " << kftrack.GetQPt() << std::endl;
1322 std::cout << "cov: " << std::endl;
1323 for (int i = 0; i < 15; i++)
1324 {
1325 std::cout << kftrack.GetCov(i) << ", ";
1326 }
1327 std::cout << std::endl;
1328 }
1329
1330
1331 std::vector<unsigned int> layers;
1332 std::transform(ckeys.begin(), ckeys.end(), std::back_inserter(layers), [](const TrkrDefs::cluskey& key)
1333 { return TrkrDefs::getLayer(key); });
1334
1335 double old_phi = track_phi;
1336 unsigned int old_layer = TrkrDefs::getLayer(ckeys[0]);
1337 if (Verbosity() > 1)
1338 {
1339 std::cout << "first layer: " << old_layer << std::endl;
1340 std::cout << "cluster (x,y,z) = (" << globalPositions.at(ckeys[0])(0) << ", " << globalPositions.at(ckeys[0])(1) << ", " << globalPositions.at(ckeys[0])(2) << ")" << std::endl;
1341 }
1342
1343 propagated_track.push_back(ckeys[0]);
1344
1345 GPUTPCTrackLinearisation kfline(aliceSeed);
1346 GPUTPCTrackParam::GPUTPCTrackFitParam fp{};
1347 aliceSeed.CalculateFitParameters(fp);
1348
1349 for (int step = 0; step <= _max_propagation_steps; ++step)
1350 {
1351 if (Verbosity() > 1)
1352 {
1353 std::cout << std::endl
1354 << "------------------------" << std::endl
1355 << "step " << step << std::endl;
1356 }
1357 if (!PropagateStep(old_layer, old_phi, direction, propagated_track, ckeys, aliceSeed, fp, globalPositions))
1358 {
1359 break;
1360 }
1361 }
1362
1363 if (Verbosity() > 1)
1364 {
1365 std::cout << "propagated track has " << propagated_track.size() << " clusters, at layers: ";
1366 for (auto c : propagated_track)
1367 {
1368 std::cout << (int) TrkrDefs::getLayer(c) << ", ";
1369 }
1370 std::cout << std::endl;
1371 }
1372 return propagated_track;
1373 }
1374
1375 std::vector<keylist> PHSimpleKFProp::RemoveBadClusters(const std::vector<keylist>& chains, const PositionMap& ) const
1376 {
1377 if (Verbosity())
1378 {
1379 std::cout << "removing bad clusters" << std::endl;
1380 }
1381 std::vector<keylist> clean_chains;
1382
1383
1384
1385
1386
1387
1388
1389
1390
1391
1392
1393
1394
1395
1396
1397
1398
1399
1400
1401
1402
1403
1404
1405
1406
1407
1408
1409
1410
1411
1412
1413
1414
1415
1416
1417
1418
1419
1420
1421
1422
1423
1424 std::copy_if(chains.begin(), chains.end(), std::back_inserter(clean_chains),
1425 [](const keylist& chain)
1426 { return chain.size() >= 3; });
1427
1428 return clean_chains;
1429 }
1430
1431 void PHSimpleKFProp::rejectAndPublishSeeds(std::vector<TrackSeed_v2>& seeds, const PositionMap& positions, std::vector<float>& trackChi2)
1432 {
1433
1434 PHTimer timer("KFPropTimer");
1435
1436
1437 timer.restart();
1438
1439
1440 PHGhostRejection rejector(Verbosity(), seeds);
1441 rejector.set_phi_cut(_ghost_phi_cut);
1442 rejector.set_eta_cut(_ghost_eta_cut);
1443 rejector.set_x_cut(_ghost_x_cut);
1444 rejector.set_y_cut(_ghost_y_cut);
1445 rejector.set_z_cut(_ghost_z_cut);
1446
1447
1448
1449
1450
1451
1452 for (unsigned int itrack = 0; itrack < seeds.size(); ++itrack)
1453 { rejector.cut_from_clusters(itrack); }
1454
1455 #pragma omp parallel
1456 {
1457
1458 #pragma omp for schedule(static)
1459 for (unsigned int itrack = 0; itrack < seeds.size(); ++itrack)
1460 {
1461
1462 if (rejector.is_rejected(itrack))
1463 { continue; }
1464
1465 auto& seed = seeds[itrack];
1466
1467 const int q = seed.get_charge();
1468
1469 PositionMap local;
1470 std::transform(seed.begin_cluster_keys(), seed.end_cluster_keys(), std::inserter(local, local.end()),
1471 [positions](const auto& key)
1472 { return std::make_pair(key, positions.at(key)); });
1473 TrackSeedHelper::circleFitByTaubin(&seed,local, 7, 55);
1474 TrackSeedHelper::lineFit(&seed,local, 7, 55);
1475 seed.set_phi(TrackSeedHelper::get_phi(&seed,local));
1476 seed.set_qOverR(fabs(seed.get_qOverR()) * q);
1477 }
1478
1479 }
1480
1481 if (Verbosity())
1482 { std::cout << "PHSimpleKFProp::rejectAndPublishSeeds - circle fit: " << timer.elapsed() << " ms" << std::endl; }
1483
1484
1485 timer.restart();
1486 rejector.find_ghosts(trackChi2);
1487 if (Verbosity())
1488 { std::cout << "PHSimpleKFProp::rejectAndPublishSeeds - ghost rejection: " << timer.elapsed() << " ms" << std::endl; }
1489
1490 timer.restart();
1491 for (unsigned int itrack = 0; itrack < seeds.size(); ++itrack)
1492 {
1493 if (rejector.is_rejected(itrack))
1494 {
1495 if (Verbosity() > 0)
1496 {
1497 std::cout << " Seed " << ((int) itrack) << " rejected. Not getting published." << std::endl;
1498 }
1499 continue;
1500 }
1501 auto& seed = seeds[itrack];
1502 _track_map->insert(&seed);
1503
1504 int q = seed.get_charge();
1505 if (Verbosity() > 0)
1506 {
1507 std::cout << "Publishing seed " << ((int) itrack)
1508 << " q " << q
1509 << " qOverR " << fabs(seed.get_qOverR()) * q
1510 << " x " << TrackSeedHelper::get_x(&seed)
1511 << " y " << TrackSeedHelper::get_y(&seed)
1512 << " z " << TrackSeedHelper::get_z(&seed)
1513 << " pT " << seed.get_pt()
1514 << " eta " << seed.get_eta()
1515 << " phi " << seed.get_phi()
1516 << std::endl;
1517 }
1518 }
1519 if (Verbosity())
1520 { std::cout << "PHSimpleKFProp::rejectAndPublishSeeds - publication: " << timer.elapsed() << " ms" << std::endl; }
1521
1522 }
1523
1524 void PHSimpleKFProp::publishSeeds(const std::vector<TrackSeed_v2>& seeds)
1525 {
1526 for (const auto& seed : seeds)
1527 {
1528 _track_map->insert(&seed);
1529 }
1530 }