File indexing completed on 2025-08-06 08:18:27
0001
0002
0003
0004
0005
0006
0007
0008
0009 #include "PHActsTrkFitter.h"
0010 #include "MakeSourceLinks.h"
0011
0012 #include <tpc/TpcDistortionCorrectionContainer.h>
0013
0014
0015 #include <trackbase/Calibrator.h>
0016 #include <trackbase/ClusterErrorPara.h>
0017 #include <trackbase/InttDefs.h>
0018 #include <trackbase/MvtxDefs.h>
0019 #include <trackbase/TpcDefs.h>
0020 #include <trackbase/TrkrCluster.h>
0021 #include <trackbase/TrkrClusterContainer.h>
0022
0023 #include <trackbase_historic/ActsTransformations.h>
0024 #include <trackbase_historic/SvtxAlignmentStateMap_v1.h>
0025 #include <trackbase_historic/SvtxTrackMap_v2.h>
0026
0027 #include <trackbase_historic/SvtxTrackState_v3.h>
0028 #include <trackbase_historic/SvtxTrack_v4.h>
0029 #include <trackbase_historic/TrackSeed.h>
0030 #include <trackbase_historic/TrackSeedContainer.h>
0031 #include <trackbase_historic/TrackSeedHelper.h>
0032
0033 #include <g4detectors/PHG4TpcCylinderGeomContainer.h>
0034
0035 #include <micromegas/MicromegasDefs.h>
0036
0037 #include <ffamodules/CDBInterface.h>
0038
0039 #include <fun4all/Fun4AllReturnCodes.h>
0040
0041 #include <phool/PHCompositeNode.h>
0042 #include <phool/PHDataNode.h>
0043 #include <phool/PHNode.h>
0044 #include <phool/PHNodeIterator.h>
0045 #include <phool/PHObject.h>
0046 #include <phool/PHTimer.h>
0047 #include <phool/getClass.h>
0048 #include <phool/phool.h>
0049
0050 #include <Acts/EventData/MultiTrajectory.hpp>
0051 #include <Acts/EventData/MultiTrajectoryHelpers.hpp>
0052 #include <Acts/EventData/SourceLink.hpp>
0053 #include <Acts/EventData/TrackParameters.hpp>
0054 #include <Acts/Surfaces/PerigeeSurface.hpp>
0055 #include <Acts/Surfaces/Surface.hpp>
0056 #include <Acts/TrackFitting/GainMatrixSmoother.hpp>
0057 #include <Acts/TrackFitting/GainMatrixUpdater.hpp>
0058
0059 #include <cmath>
0060 #include <filesystem>
0061 #include <iostream>
0062 #include <vector>
0063
0064 namespace
0065 {
0066
0067 inline bool is_valid(const Acts::Vector3& vec)
0068 {
0069 return !(std::isnan(vec.x()) || std::isnan(vec.y()) || std::isnan(vec.z()));
0070 }
0071 template <class T>
0072 inline T square(const T& x)
0073 {
0074 return x * x;
0075 }
0076 }
0077
0078 #include <trackbase/alignmentTransformationContainer.h>
0079 #include <Eigen/Dense>
0080 #include <Eigen/Geometry>
0081
0082 PHActsTrkFitter::PHActsTrkFitter(const std::string& name)
0083 : SubsysReco(name)
0084 , m_trajectories(nullptr)
0085 {
0086 }
0087
0088 int PHActsTrkFitter::InitRun(PHCompositeNode* topNode)
0089 {
0090 if (Verbosity() > 1)
0091 {
0092 std::cout << "Setup PHActsTrkFitter" << std::endl;
0093 }
0094
0095 if (createNodes(topNode) != Fun4AllReturnCodes::EVENT_OK)
0096 {
0097 return Fun4AllReturnCodes::ABORTEVENT;
0098 }
0099
0100 if (getNodes(topNode) != Fun4AllReturnCodes::EVENT_OK)
0101 {
0102 return Fun4AllReturnCodes::ABORTEVENT;
0103 }
0104
0105
0106 m_alignStates.loadNodes(topNode);
0107 m_alignStates.verbosity(Verbosity());
0108 m_alignStates.fieldMap(m_fieldMap);
0109
0110
0111 std::istringstream stringline(m_fieldMap);
0112 stringline >> fieldstrength;
0113 if (!stringline.fail())
0114 {
0115 m_ConstField = true;
0116 }
0117
0118 auto level = Acts::Logging::FATAL;
0119 if (Verbosity() > 5)
0120 {
0121 level = Acts::Logging::VERBOSE;
0122 }
0123
0124 m_fitCfg.fit = ActsTrackFittingAlgorithm::makeKalmanFitterFunction(
0125 m_tGeometry->geometry().tGeometry,
0126 m_tGeometry->geometry().magField,
0127 true, true, 0.0, Acts::FreeToBoundCorrection(), *Acts::getDefaultLogger("Kalman", level));
0128
0129 m_fitCfg.dFit = ActsTrackFittingAlgorithm::makeDirectedKalmanFitterFunction(
0130 m_tGeometry->geometry().tGeometry,
0131 m_tGeometry->geometry().magField);
0132
0133 MaterialSurfaceSelector selector;
0134 if (m_fitSiliconMMs || m_directNavigation)
0135 {
0136 m_tGeometry->geometry().tGeometry->visitSurfaces(selector,false);
0137
0138 m_materialSurfaces = selector.surfaces;
0139 }
0140
0141 m_outlierFinder.verbosity = Verbosity();
0142 std::map<long unsigned int, float> chi2Cuts;
0143 chi2Cuts.insert(std::make_pair(10, 4));
0144 chi2Cuts.insert(std::make_pair(12, 4));
0145 chi2Cuts.insert(std::make_pair(14, 9));
0146 chi2Cuts.insert(std::make_pair(16, 4));
0147 m_outlierFinder.chi2Cuts = chi2Cuts;
0148 if (m_useOutlierFinder)
0149 {
0150 m_outlierFinder.m_tGeometry = m_tGeometry;
0151 m_fitCfg.fit->outlierFinder(m_outlierFinder);
0152 }
0153
0154 if (m_timeAnalysis)
0155 {
0156 m_timeFile = new TFile(std::string(Name() + ".root").c_str(),
0157 "RECREATE");
0158 h_eventTime = new TH1F("h_eventTime", ";time [ms]",
0159 100000, 0, 10000);
0160 h_fitTime = new TH2F("h_fitTime", ";p_{T} [GeV];time [ms]",
0161 80, 0, 40, 100000, 0, 1000);
0162 h_updateTime = new TH1F("h_updateTime", ";time [ms]",
0163 100000, 0, 1000);
0164
0165 h_rotTime = new TH1F("h_rotTime", ";time [ms]",
0166 100000, 0, 1000);
0167 h_stateTime = new TH1F("h_stateTime", ";time [ms]",
0168 100000, 0, 1000);
0169 }
0170
0171 if (m_actsEvaluator)
0172 {
0173 m_evaluator = std::make_unique<ActsEvaluator>(m_evalname);
0174 m_evaluator->Init(topNode);
0175 if(m_actsEvaluator && !m_simActsEvaluator)
0176 {
0177 m_evaluator->isData();
0178 }
0179 m_evaluator->verbosity(Verbosity());
0180 }
0181
0182 _tpccellgeo = findNode::getClass<PHG4TpcCylinderGeomContainer>(topNode, "CYLINDERCELLGEOM_SVTX");
0183
0184 if (Verbosity() > 1)
0185 {
0186 std::cout << "Finish PHActsTrkFitter Setup" << std::endl;
0187 }
0188
0189 return Fun4AllReturnCodes::EVENT_OK;
0190 }
0191
0192 int PHActsTrkFitter::process_event(PHCompositeNode* topNode)
0193 {
0194 PHTimer eventTimer("eventTimer");
0195 eventTimer.stop();
0196 eventTimer.restart();
0197 m_nBadFits = 0;
0198 m_event++;
0199
0200 auto logLevel = Acts::Logging::FATAL;
0201
0202 if (m_actsEvaluator)
0203 {
0204 m_evaluator->next_event(topNode);
0205 }
0206
0207 if (Verbosity() > 1)
0208 {
0209 std::cout << PHWHERE << "Events processed: " << m_event << std::endl;
0210 std::cout << "Start PHActsTrkFitter::process_event" << std::endl;
0211 if (Verbosity() > 4)
0212 {
0213 logLevel = Acts::Logging::VERBOSE;
0214 }
0215 }
0216
0217
0218 m_trackMap->Reset();
0219
0220 loopTracks(logLevel);
0221
0222 eventTimer.stop();
0223 auto eventTime = eventTimer.get_accumulated_time();
0224
0225 if (Verbosity() > 1)
0226 {
0227 std::cout << "PHActsTrkFitter total event time "
0228 << eventTime << std::endl;
0229 }
0230
0231 if (m_timeAnalysis)
0232 {
0233 h_eventTime->Fill(eventTime);
0234 }
0235
0236 if (Verbosity() > 1)
0237 {
0238 std::cout << "PHActsTrkFitter::process_event finished"
0239 << std::endl;
0240 }
0241
0242
0243 if (Verbosity() > 0)
0244 {
0245 std::cout << "The Acts track fitter had " << m_nBadFits
0246 << " fits return an error" << std::endl;
0247 std::cout << " seed map size " << m_seedMap->size() << std::endl;
0248
0249 std::cout << " SvtxTrackMap size is now " << m_trackMap->size()
0250 << std::endl;
0251 }
0252
0253 return Fun4AllReturnCodes::EVENT_OK;
0254 }
0255
0256 int PHActsTrkFitter::ResetEvent(PHCompositeNode* )
0257 {
0258 if (Verbosity() > 1)
0259 {
0260 std::cout << "Reset PHActsTrkFitter" << std::endl;
0261 }
0262
0263 m_trajectories->clear();
0264
0265 return Fun4AllReturnCodes::EVENT_OK;
0266 }
0267
0268 int PHActsTrkFitter::End(PHCompositeNode* )
0269 {
0270 if (m_timeAnalysis)
0271 {
0272 m_timeFile->cd();
0273 h_fitTime->Write();
0274 h_eventTime->Write();
0275 h_rotTime->Write();
0276 h_stateTime->Write();
0277 h_updateTime->Write();
0278 m_timeFile->Write();
0279 m_timeFile->Close();
0280 }
0281
0282 if (m_actsEvaluator)
0283 {
0284 m_evaluator->End();
0285 }
0286 if(m_useOutlierFinder)
0287 {
0288 m_outlierFinder.Write();
0289 }
0290 if (Verbosity() > 0)
0291 {
0292 std::cout << "Finished PHActsTrkFitter" << std::endl;
0293 }
0294 return Fun4AllReturnCodes::EVENT_OK;
0295 }
0296
0297 void PHActsTrkFitter::loopTracks(Acts::Logging::Level logLevel)
0298 {
0299 auto logger = Acts::getDefaultLogger("PHActsTrkFitter", logLevel);
0300
0301 for (auto track : *m_seedMap)
0302 {
0303 if (!track)
0304 {
0305 continue;
0306 }
0307
0308 unsigned int tpcid = track->get_tpc_seed_index();
0309 unsigned int siid = track->get_silicon_seed_index();
0310
0311
0312
0313 short silicon_crossing = SHRT_MAX;
0314 auto siseed = m_siliconSeeds->get(siid);
0315 if(siseed)
0316 {
0317 silicon_crossing = siseed->get_crossing();
0318 }
0319 short crossing = silicon_crossing;
0320 short int crossing_estimate = crossing;
0321
0322 if(m_enable_crossing_estimate)
0323 {
0324 crossing_estimate = track->get_crossing_estimate();
0325 }
0326
0327
0328
0329
0330 if (m_fitSiliconMMs)
0331 {
0332 if( (siid == std::numeric_limits<unsigned int>::max()) || (silicon_crossing == SHRT_MAX))
0333 {
0334 continue;
0335 }
0336 }
0337
0338
0339 if(!siseed)
0340 {
0341 crossing = 0;
0342 }
0343
0344 if (Verbosity() > 1)
0345 {
0346 if(siseed)
0347 {
0348 std::cout << "tpc and si id " << tpcid << ", " << siid << " silicon_crossing " << silicon_crossing
0349 << " crossing " << crossing << " crossing estimate " << crossing_estimate << std::endl;
0350 }
0351 }
0352
0353 auto tpcseed = m_tpcSeeds->get(tpcid);
0354
0355
0356 if (!tpcseed)
0357 {
0358 std::cout << "no tpc seed" << std::endl;
0359 continue;
0360 }
0361
0362 if (Verbosity() > 0)
0363 {
0364 if (siseed)
0365 {
0366 const auto si_position = TrackSeedHelper::get_xyz(siseed);
0367 const auto tpc_position = TrackSeedHelper::get_xyz(tpcseed);
0368 std::cout << " silicon seed position is (x,y,z) = " << si_position.x() << " " << si_position.y() << " " << si_position.z() << std::endl;
0369 std::cout << " tpc seed position is (x,y,z) = " << tpc_position.x() << " " << tpc_position.y() << " " << tpc_position.z() << std::endl;
0370 }
0371 }
0372
0373 PHTimer trackTimer("TrackTimer");
0374 trackTimer.stop();
0375 trackTimer.restart();
0376
0377 if (Verbosity() > 1 && siseed)
0378 {
0379 std::cout << " m_pp_mode " << m_pp_mode << " m_enable_crossing_estimate " << m_enable_crossing_estimate
0380 << " INTT crossing " << crossing << " crossing_estimate " << crossing_estimate << std::endl;
0381 }
0382
0383 short int this_crossing = crossing;
0384 bool use_estimate = false;
0385 short int nvary = 0;
0386 std::vector<float> chisq_ndf;
0387 std::vector<SvtxTrack_v4> svtx_vec;
0388
0389 if(m_pp_mode)
0390 {
0391 if (m_enable_crossing_estimate && crossing == SHRT_MAX)
0392 {
0393
0394
0395 use_estimate = true;
0396 nvary = max_bunch_search;
0397 if (Verbosity() > 1)
0398 {
0399 std::cout << " No INTT crossing: use crossing_estimate " << crossing_estimate << " with nvary " << nvary << std::endl;
0400 }
0401 }
0402 else
0403 {
0404
0405 crossing_estimate = crossing;
0406 }
0407 }
0408 else
0409 {
0410
0411 if(siseed && silicon_crossing != 0)
0412 {
0413 crossing = 0;
0414
0415 }
0416 crossing_estimate = crossing;
0417 }
0418
0419
0420
0421
0422
0423 for (short int ivary = -nvary; ivary <= nvary; ++ivary)
0424 {
0425 this_crossing = crossing_estimate + ivary;
0426
0427 if (Verbosity() > 1)
0428 {
0429 std::cout << " nvary " << nvary << " trial fit with ivary " << ivary << " this_crossing = " << this_crossing << std::endl;
0430 }
0431
0432 ActsTrackFittingAlgorithm::MeasurementContainer measurements;
0433
0434 SourceLinkVec sourceLinks;
0435
0436 MakeSourceLinks makeSourceLinks;
0437 makeSourceLinks.initialize(_tpccellgeo);
0438 makeSourceLinks.setVerbosity(Verbosity());
0439 makeSourceLinks.set_pp_mode(m_pp_mode);
0440 for(const auto& layer : m_ignoreLayer)
0441 {
0442 makeSourceLinks.ignoreLayer(layer);
0443 }
0444
0445
0446 makeSourceLinks.resetTransientTransformMap(
0447 m_alignmentTransformationMapTransient,
0448 m_transient_id_set,
0449 m_tGeometry);
0450
0451 if (m_use_clustermover)
0452 {
0453
0454
0455 if (siseed && !m_ignoreSilicon)
0456 {
0457
0458 sourceLinks = makeSourceLinks.getSourceLinksClusterMover(
0459 siseed,
0460 measurements,
0461 m_clusterContainer,
0462 m_tGeometry,
0463 m_globalPositionWrapper,
0464 this_crossing);
0465 }
0466
0467
0468 const auto tpcSourceLinks = makeSourceLinks.getSourceLinksClusterMover(
0469 tpcseed,
0470 measurements,
0471 m_clusterContainer,
0472 m_tGeometry,
0473 m_globalPositionWrapper,
0474 this_crossing);
0475
0476
0477 sourceLinks.insert(sourceLinks.end(), tpcSourceLinks.begin(), tpcSourceLinks.end());
0478 }
0479 else
0480 {
0481
0482
0483 if(Verbosity() > 1)
0484 { std::cout << "Calling getSourceLinks for si seed, siid " << siid << " and tpcid " << tpcid << std::endl; }
0485
0486 if (siseed && !m_ignoreSilicon)
0487 {
0488
0489 sourceLinks = makeSourceLinks.getSourceLinks(
0490 siseed,
0491 measurements,
0492 m_clusterContainer,
0493 m_tGeometry,
0494 m_globalPositionWrapper,
0495 m_alignmentTransformationMapTransient,
0496 m_transient_id_set,
0497 this_crossing);
0498 }
0499
0500 if(Verbosity() > 1)
0501 { std::cout << "Calling getSourceLinks for tpc seed, siid " << siid << " and tpcid " << tpcid << std::endl; }
0502
0503
0504 const auto tpcSourceLinks = makeSourceLinks.getSourceLinks(
0505 tpcseed,
0506 measurements,
0507 m_clusterContainer,
0508 m_tGeometry,
0509 m_globalPositionWrapper,
0510 m_alignmentTransformationMapTransient,
0511 m_transient_id_set,
0512 this_crossing);
0513
0514
0515 sourceLinks.insert(sourceLinks.end(), tpcSourceLinks.begin(), tpcSourceLinks.end());
0516 }
0517
0518
0519 m_transient_geocontext = m_alignmentTransformationMapTransient;
0520
0521
0522 Acts::Vector3 position(0, 0, 0);
0523 if (siseed)
0524 {
0525 position = TrackSeedHelper::get_xyz(siseed)*Acts::UnitConstants::cm;
0526 }
0527 if(!siseed || !is_valid(position) || m_ignoreSilicon)
0528 {
0529 position = TrackSeedHelper::get_xyz(tpcseed)*Acts::UnitConstants::cm;
0530 }
0531 if (!is_valid(position))
0532 {
0533 if(Verbosity() > 4)
0534 {
0535 std::cout << "Invalid position of " << position.transpose() << std::endl;
0536 }
0537 continue;
0538 }
0539
0540 if (sourceLinks.empty())
0541 {
0542 continue;
0543 }
0544
0545
0546 SurfacePtrVec surfaces_tmp;
0547 SurfacePtrVec surfaces;
0548 if (m_fitSiliconMMs || m_directNavigation)
0549 {
0550 sourceLinks = getSurfaceVector(sourceLinks, surfaces_tmp);
0551
0552
0553 if (surfaces_tmp.empty())
0554 {
0555 continue;
0556 }
0557 for (const auto& surface_apr : m_materialSurfaces)
0558 {
0559 if(m_forceSiOnlyFit)
0560 {
0561 if(surface_apr->geometryId().volume() >12)
0562 {
0563 continue;
0564 }
0565 }
0566 bool pop_flag = false;
0567 if(surface_apr->geometryId().approach() == 1)
0568 {
0569 surfaces.push_back(surface_apr);
0570 }
0571 else
0572 {
0573 pop_flag = true;
0574 for (const auto& surface_sns: surfaces_tmp)
0575 {
0576 if (surface_apr->geometryId().volume() == surface_sns->geometryId().volume())
0577 {
0578 if ( surface_apr->geometryId().layer()==surface_sns->geometryId().layer())
0579 {
0580 pop_flag = false;
0581 surfaces.push_back(surface_sns);
0582 }
0583 }
0584 }
0585 if (!pop_flag)
0586 {
0587 surfaces.push_back(surface_apr);
0588 }
0589 else
0590 {
0591 surfaces.pop_back();
0592 pop_flag = false;
0593 }
0594 if (surface_apr->geometryId().volume() == 12&& surface_apr->geometryId().layer()==8)
0595 {
0596 for (const auto& surface_sns: surfaces_tmp)
0597 {
0598 if (14 == surface_sns->geometryId().volume())
0599 {
0600 surfaces.push_back(surface_sns);
0601 }
0602 }
0603 }
0604 }
0605 }
0606 checkSurfaceVec(surfaces);
0607 if (Verbosity() > 1)
0608 {
0609 for (const auto& surf : surfaces)
0610 {
0611 std::cout << "Surface vector : " << surf->geometryId() << std::endl;
0612 }
0613 }
0614 if (m_fitSiliconMMs)
0615 {
0616
0617 if (m_useMicromegas &&
0618 std::none_of(surfaces.begin(), surfaces.end(), [this](const auto& surface)
0619 { return m_tGeometry->maps().isMicromegasSurface(surface); }))
0620 {
0621 continue;
0622 }
0623 }
0624 }
0625
0626 float px = std::numeric_limits<float>::quiet_NaN();
0627 float py = std::numeric_limits<float>::quiet_NaN();
0628 float pz = std::numeric_limits<float>::quiet_NaN();
0629
0630
0631 float seedphi = 0;
0632 float seedtheta = 0;
0633 float seedeta = 0;
0634 if(siseed)
0635 {
0636 seedphi = siseed->get_phi();
0637 seedtheta = siseed->get_theta();
0638 seedeta = siseed->get_eta();
0639 }
0640 else
0641 {
0642 seedphi = tpcseed->get_phi();
0643 seedtheta = tpcseed->get_theta();
0644 seedeta = tpcseed->get_eta();
0645 }
0646
0647 float seedpt = tpcseed->get_pt();
0648
0649 if (m_ConstField)
0650 {
0651 float pt = fabs(1. / tpcseed->get_qOverR()) * (0.3 / 100) * fieldstrength;
0652 float phi = seedphi;
0653 float eta = seedeta;
0654 float theta = seedtheta;
0655 px = pt * std::cos(phi);
0656 py = pt * std::sin(phi);
0657 pz = pt * std::cosh(eta) * std::cos(theta);
0658 }
0659 else
0660 {
0661 px = seedpt * cos(seedphi);
0662 py = seedpt * sin(seedphi);
0663 pz = seedpt * std::cosh(seedeta) * std::cos(seedtheta);
0664 }
0665
0666 Acts::Vector3 momentum(px, py, pz);
0667 if (!is_valid(momentum))
0668 {
0669 if(Verbosity() > 4)
0670 {
0671 std::cout << "Invalid momentum of " << momentum.transpose() << std::endl;
0672 }
0673 continue;
0674 }
0675
0676 auto pSurface = Acts::Surface::makeShared<Acts::PerigeeSurface>(
0677 position);
0678
0679 auto actsFourPos = Acts::Vector4(position(0), position(1),
0680 position(2),
0681 10 * Acts::UnitConstants::ns);
0682 Acts::BoundSquareMatrix cov = setDefaultCovariance();
0683
0684 int charge = tpcseed->get_charge();
0685
0686
0687 auto seed = ActsTrackFittingAlgorithm::TrackParameters::create(
0688 pSurface,
0689 m_transient_geocontext,
0690 actsFourPos,
0691 momentum,
0692 charge / momentum.norm(),
0693 cov,
0694 Acts::ParticleHypothesis::pion())
0695 .value();
0696
0697 if (Verbosity() > 2)
0698 {
0699 printTrackSeed(seed);
0700 }
0701
0702
0703 Acts::PropagatorPlainOptions ppPlainOptions;
0704
0705 auto calibptr = std::make_unique<Calibrator>();
0706 CalibratorAdapter calibrator{*calibptr, measurements};
0707
0708 auto magcontext = m_tGeometry->geometry().magFieldContext;
0709 auto calibcontext = m_tGeometry->geometry().calibContext;
0710
0711 ActsTrackFittingAlgorithm::GeneralFitterOptions
0712 kfOptions{
0713 m_transient_geocontext,
0714 magcontext,
0715 calibcontext,
0716 pSurface.get(),
0717 ppPlainOptions};
0718
0719 PHTimer fitTimer("FitTimer");
0720 fitTimer.stop();
0721 fitTimer.restart();
0722
0723 auto trackContainer =
0724 std::make_shared<Acts::VectorTrackContainer>();
0725 auto trackStateContainer =
0726 std::make_shared<Acts::VectorMultiTrajectory>();
0727 ActsTrackFittingAlgorithm::TrackContainer
0728 tracks(trackContainer, trackStateContainer);
0729
0730 if(Verbosity() > 1)
0731 { std::cout << "Calling fitTrack for track with siid " << siid << " tpcid " << tpcid << " crossing " << crossing << std::endl; }
0732
0733 auto result = fitTrack(sourceLinks, seed, kfOptions,
0734 surfaces, calibrator, tracks);
0735 fitTimer.stop();
0736 auto fitTime = fitTimer.get_accumulated_time();
0737
0738 if (Verbosity() > 1)
0739 {
0740 std::cout << "PHActsTrkFitter Acts fit time " << fitTime << std::endl;
0741 }
0742
0743
0744 if (result.ok())
0745 {
0746 if (use_estimate)
0747 {
0748
0749
0750
0751 SvtxTrack_v4 newTrack;
0752 newTrack.set_tpc_seed(tpcseed);
0753 newTrack.set_crossing(this_crossing);
0754 newTrack.set_silicon_seed(siseed);
0755
0756 if (getTrackFitResult(result, track, &newTrack, tracks, measurements))
0757 {
0758 float chi2ndf = newTrack.get_quality();
0759 chisq_ndf.push_back(chi2ndf);
0760 svtx_vec.push_back(newTrack);
0761 if (Verbosity() > 1)
0762 {
0763 std::cout << " tpcid " << tpcid << " siid " << siid << " ivary " << ivary << " this_crossing " << this_crossing << " chi2ndf " << chi2ndf << std::endl;
0764 }
0765 }
0766
0767 if (ivary != nvary)
0768 {
0769 if(Verbosity() > 3)
0770 {
0771 std::cout << "Skipping track fit for trial variation" << std::endl;
0772 }
0773 continue;
0774 }
0775
0776
0777 if (Verbosity() > 1)
0778 {
0779 std::cout << "Finished with trial fits, chisq_ndf size is " << chisq_ndf.size() << " chisq_ndf values are:" << std::endl;
0780 }
0781 float best_chisq = 1000.0;
0782 short int best_ivary = 0;
0783 for (unsigned int i = 0; i < chisq_ndf.size(); ++i)
0784 {
0785 if (chisq_ndf[i] < best_chisq)
0786 {
0787 best_chisq = chisq_ndf[i];
0788 best_ivary = i;
0789 }
0790 if (Verbosity() > 1)
0791 {
0792 std::cout << " trial " << i << " chisq_ndf " << chisq_ndf[i] << " best_chisq " << best_chisq << " best_ivary " << best_ivary << std::endl;
0793 }
0794 }
0795 unsigned int trid = m_trackMap->size();
0796 svtx_vec[best_ivary].set_id(trid);
0797
0798 m_trackMap->insertWithKey(&svtx_vec[best_ivary], trid);
0799 }
0800 else
0801 {
0802 SvtxTrack_v4 newTrack;
0803 newTrack.set_tpc_seed(tpcseed);
0804 newTrack.set_crossing(this_crossing);
0805 newTrack.set_silicon_seed(siseed);
0806
0807 if (m_fitSiliconMMs)
0808 {
0809 unsigned int trid = m_directedTrackMap->size();
0810 newTrack.set_id(trid);
0811
0812 if (getTrackFitResult(result, track, &newTrack, tracks, measurements))
0813 {
0814 m_directedTrackMap->insertWithKey(&newTrack, trid);
0815 }
0816 }
0817 else
0818 {
0819 unsigned int trid = m_trackMap->size();
0820 newTrack.set_id(trid);
0821
0822 if (getTrackFitResult(result, track, &newTrack, tracks, measurements))
0823 {
0824 m_trackMap->insertWithKey(&newTrack, trid);
0825 }
0826 }
0827 }
0828 }
0829 else if (!m_fitSiliconMMs)
0830 {
0831
0832 m_nBadFits++;
0833 if (Verbosity() > 1)
0834 {
0835 std::cout << "Track fit failed for track " << m_seedMap->find(track)
0836 << " with Acts error message "
0837 << result.error() << ", " << result.error().message()
0838 << std::endl;
0839 }
0840 }
0841 }
0842
0843 trackTimer.stop();
0844 auto trackTime = trackTimer.get_accumulated_time();
0845
0846 if (Verbosity() > 1)
0847 {
0848 std::cout << "PHActsTrkFitter total single track time " << trackTime << std::endl;
0849 }
0850 }
0851
0852 return;
0853 }
0854
0855 bool PHActsTrkFitter::getTrackFitResult(FitResult& fitOutput,
0856 TrackSeed* seed, SvtxTrack* track,
0857 ActsTrackFittingAlgorithm::TrackContainer& tracks,
0858 const ActsTrackFittingAlgorithm::MeasurementContainer& measurements)
0859 {
0860
0861
0862 std::vector<Acts::MultiTrajectoryTraits::IndexType> trackTips;
0863 trackTips.reserve(1);
0864 auto& outtrack = fitOutput.value();
0865 if (outtrack.hasReferenceSurface())
0866 {
0867 trackTips.emplace_back(outtrack.tipIndex());
0868 Trajectory::IndexedParameters indexedParams;
0869 indexedParams.emplace(std::pair{outtrack.tipIndex(),
0870 ActsExamples::TrackParameters{outtrack.referenceSurface().getSharedPtr(),
0871 outtrack.parameters(), outtrack.covariance(), outtrack.particleHypothesis()}});
0872
0873 if (Verbosity() > 2)
0874 {
0875 std::cout << "Fitted parameters for track" << std::endl;
0876 std::cout << " position : " << outtrack.referenceSurface().localToGlobal(m_transient_geocontext, Acts::Vector2(outtrack.loc0(), outtrack.loc1()), Acts::Vector3(1, 1, 1)).transpose()
0877
0878 << std::endl;
0879 int otcharge = outtrack.qOverP() > 0 ? 1 : -1;
0880 std::cout << "charge: " << otcharge << std::endl;
0881 std::cout << " momentum : " << outtrack.momentum().transpose()
0882 << std::endl;
0883 std::cout << "For trackTip == " << outtrack.tipIndex() << std::endl;
0884 }
0885
0886
0887
0888 PHTimer updateTrackTimer("UpdateTrackTimer");
0889 updateTrackTimer.stop();
0890 updateTrackTimer.restart();
0891 updateSvtxTrack(trackTips, indexedParams, tracks, track);
0892
0893 if (m_commissioning)
0894 {
0895 if (track->get_silicon_seed() && track->get_tpc_seed())
0896 {
0897 m_alignStates.fillAlignmentStateMap(tracks, trackTips,
0898 track, measurements);
0899 }
0900 }
0901
0902 updateTrackTimer.stop();
0903 auto updateTime = updateTrackTimer.get_accumulated_time();
0904
0905 if (Verbosity() > 1)
0906 {
0907 std::cout << "PHActsTrkFitter update SvtxTrack time "
0908 << updateTime << std::endl;
0909 }
0910
0911 if (m_timeAnalysis)
0912 {
0913 h_updateTime->Fill(updateTime);
0914 }
0915
0916 Trajectory trajectory(tracks.trackStateContainer(),
0917 trackTips, indexedParams);
0918
0919 m_trajectories->insert(std::make_pair(track->get_id(), trajectory));
0920
0921 if (m_actsEvaluator)
0922 {
0923 m_evaluator->evaluateTrackFit(tracks, trackTips, indexedParams, track,
0924 seed, measurements);
0925 }
0926
0927 return true;
0928 }
0929
0930 return false;
0931 }
0932
0933 ActsTrackFittingAlgorithm::TrackFitterResult PHActsTrkFitter::fitTrack(
0934 const std::vector<Acts::SourceLink>& sourceLinks,
0935 const ActsTrackFittingAlgorithm::TrackParameters& seed,
0936 const ActsTrackFittingAlgorithm::GeneralFitterOptions& kfOptions,
0937 const SurfacePtrVec& surfSequence,
0938 const CalibratorAdapter& calibrator,
0939 ActsTrackFittingAlgorithm::TrackContainer& tracks)
0940 {
0941 if (m_fitSiliconMMs)
0942 {
0943 return (*m_fitCfg.dFit)(sourceLinks, seed, kfOptions,
0944 surfSequence, calibrator, tracks);
0945 }
0946 else
0947 {
0948 if(m_directNavigation)
0949 {
0950 return (*m_fitCfg.dFit)(sourceLinks, seed, kfOptions,
0951 surfSequence, calibrator, tracks);
0952 }
0953 else
0954 {
0955 return (*m_fitCfg.fit)(sourceLinks, seed, kfOptions,
0956 calibrator, tracks);
0957 }
0958 }
0959 }
0960
0961 SourceLinkVec PHActsTrkFitter::getSurfaceVector(const SourceLinkVec& sourceLinks,
0962 SurfacePtrVec& surfaces) const
0963 {
0964 SourceLinkVec siliconMMSls;
0965
0966
0967
0968
0969 for (const auto& sl : sourceLinks)
0970 {
0971 const ActsSourceLink asl = sl.get<ActsSourceLink>();
0972 if (Verbosity() > 1)
0973 {
0974 std::cout << "SL available on : " << asl.geometryId() << std::endl;
0975 }
0976
0977 const auto surf = m_tGeometry->geometry().tGeometry->findSurface(asl.geometryId());
0978 if (m_fitSiliconMMs)
0979 {
0980
0981 if (m_tGeometry->maps().isTpcSurface(surf))
0982 {
0983 continue;
0984 }
0985
0986
0987 if (m_tGeometry->maps().isMicromegasSurface(surf) && !m_useMicromegas)
0988 {
0989 continue;
0990 }
0991 }
0992
0993 if(m_forceSiOnlyFit)
0994 {
0995 if(m_tGeometry->maps().isMicromegasSurface(surf)||m_tGeometry->maps().isTpcSurface(surf))
0996 {
0997 continue;
0998 }
0999 }
1000
1001 siliconMMSls.push_back(sl);
1002 surfaces.push_back(surf);
1003 }
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013 if (Verbosity() > 10)
1014 {
1015 for (const auto& surf : surfaces)
1016 {
1017 std::cout << "Surface vector : " << surf->geometryId() << std::endl;
1018 }
1019 }
1020
1021 return siliconMMSls;
1022 }
1023
1024 void PHActsTrkFitter::checkSurfaceVec(SurfacePtrVec& surfaces) const
1025 {
1026 for (unsigned int i = 0; i < surfaces.size() - 1; i++)
1027 {
1028 const auto& surface = surfaces.at(i);
1029 const auto thisVolume = surface->geometryId().volume();
1030 const auto thisLayer = surface->geometryId().layer();
1031
1032 const auto nextSurface = surfaces.at(i + 1);
1033 const auto nextVolume = nextSurface->geometryId().volume();
1034 const auto nextLayer = nextSurface->geometryId().layer();
1035
1036
1037 if (nextVolume == thisVolume)
1038 {
1039 if (nextLayer < thisLayer)
1040 {
1041 std::cout
1042 << "PHActsTrkFitter::checkSurfaceVec - "
1043 << "Surface not in order... removing surface"
1044 << surface->geometryId() << std::endl;
1045
1046 surfaces.erase(surfaces.begin() + i);
1047
1048
1049 --i;
1050 continue;
1051 }
1052 }
1053 else
1054 {
1055 if (nextVolume < thisVolume)
1056 {
1057 std::cout
1058 << "PHActsTrkFitter::checkSurfaceVec - "
1059 << "Volume not in order... removing surface"
1060 << surface->geometryId() << std::endl;
1061
1062 surfaces.erase(surfaces.begin() + i);
1063
1064
1065 --i;
1066 continue;
1067 }
1068 }
1069 }
1070 }
1071
1072 void PHActsTrkFitter::updateSvtxTrack(std::vector<Acts::MultiTrajectoryTraits::IndexType>& tips,
1073 Trajectory::IndexedParameters& paramsMap,
1074 ActsTrackFittingAlgorithm::TrackContainer& tracks,
1075 SvtxTrack* track)
1076 {
1077 const auto& mj = tracks.trackStateContainer();
1078
1079
1080 auto& trackTip = tips.front();
1081
1082 if (Verbosity() > 2)
1083 {
1084 std::cout << "Identify (proto) track before updating with acts results " << std::endl;
1085 track->identify();
1086 }
1087
1088 if (!m_fitSiliconMMs && !m_forceSiOnlyFit)
1089 {
1090 track->clear_states();
1091 }
1092
1093
1094
1095 float pathlength = 0.0;
1096
1097 SvtxTrackState_v3 out(pathlength);
1098 out.set_localX(0.0);
1099 out.set_localY(0.0);
1100 out.set_x(0.0);
1101 out.set_y(0.0);
1102 out.set_z(0.0);
1103 track->insert_state(&out);
1104
1105 auto trajState =
1106 Acts::MultiTrajectoryHelpers::trajectoryState(mj, trackTip);
1107
1108 const auto& params = paramsMap.find(trackTip)->second;
1109
1110
1111 track->set_x(params.position(m_transient_geocontext)(0) / Acts::UnitConstants::cm);
1112 track->set_y(params.position(m_transient_geocontext)(1) / Acts::UnitConstants::cm);
1113 track->set_z(params.position(m_transient_geocontext)(2) / Acts::UnitConstants::cm);
1114
1115 track->set_px(params.momentum()(0));
1116 track->set_py(params.momentum()(1));
1117 track->set_pz(params.momentum()(2));
1118
1119 track->set_charge(params.charge());
1120 track->set_chisq(trajState.chi2Sum);
1121 track->set_ndf(trajState.NDF);
1122
1123 ActsTransformations rotater;
1124 rotater.setVerbosity(Verbosity());
1125
1126 if (params.covariance())
1127 {
1128 auto rotatedCov = rotater.rotateActsCovToSvtxTrack(params);
1129
1130 for (int i = 0; i < 6; i++)
1131 {
1132 for (int j = 0; j < 6; j++)
1133 {
1134 track->set_error(i, j, rotatedCov(i, j));
1135 }
1136 }
1137 }
1138
1139
1140
1141 PHTimer trackStateTimer("TrackStateTimer");
1142 trackStateTimer.stop();
1143 trackStateTimer.restart();
1144
1145 if (m_fillSvtxTrackStates)
1146 {
1147 rotater.fillSvtxTrackStates(mj, trackTip, track,
1148 m_transient_geocontext);
1149 }
1150
1151 trackStateTimer.stop();
1152 auto stateTime = trackStateTimer.get_accumulated_time();
1153
1154 if (Verbosity() > 1)
1155 {
1156 std::cout << "PHActsTrkFitter update SvtxTrackStates time "
1157 << stateTime << std::endl;
1158 }
1159
1160 if (m_timeAnalysis)
1161 {
1162 h_stateTime->Fill(stateTime);
1163 }
1164
1165 if (Verbosity() > 2)
1166 {
1167 std::cout << " Identify fitted track after updating track states:"
1168 << std::endl;
1169 track->identify();
1170 }
1171
1172 return;
1173 }
1174
1175 Acts::BoundSquareMatrix PHActsTrkFitter::setDefaultCovariance() const
1176 {
1177 Acts::BoundSquareMatrix cov = Acts::BoundSquareMatrix::Zero();
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189 if (m_fitSiliconMMs)
1190 {
1191 cov << 1000 * Acts::UnitConstants::um, 0., 0., 0., 0., 0.,
1192 0., 1000 * Acts::UnitConstants::um, 0., 0., 0., 0.,
1193 0., 0., 0.1, 0., 0., 0.,
1194 0., 0., 0., 0.1, 0., 0.,
1195 0., 0., 0., 0., 0.005, 0.,
1196 0., 0., 0., 0., 0., 1.;
1197 }
1198 else
1199 {
1200
1201 double sigmaD0 = 50 * Acts::UnitConstants::um;
1202 double sigmaZ0 = 50 * Acts::UnitConstants::um;
1203
1204 double sigmaPhi = 1. * Acts::UnitConstants::degree;
1205 double sigmaTheta = 1. * Acts::UnitConstants::degree;
1206 double sigmaT = 1. * Acts::UnitConstants::ns;
1207
1208 cov(Acts::eBoundLoc0, Acts::eBoundLoc0) = sigmaD0 * sigmaD0;
1209 cov(Acts::eBoundLoc1, Acts::eBoundLoc1) = sigmaZ0 * sigmaZ0;
1210 cov(Acts::eBoundTime, Acts::eBoundTime) = sigmaT * sigmaT;
1211 cov(Acts::eBoundPhi, Acts::eBoundPhi) = sigmaPhi * sigmaPhi;
1212 cov(Acts::eBoundTheta, Acts::eBoundTheta) = sigmaTheta * sigmaTheta;
1213
1214
1215 cov(Acts::eBoundQOverP, Acts::eBoundQOverP) = 0.025;
1216 }
1217
1218 return cov;
1219 }
1220
1221 void PHActsTrkFitter::printTrackSeed(const ActsTrackFittingAlgorithm::TrackParameters& seed) const
1222 {
1223 std::cout
1224 << PHWHERE
1225 << " Processing proto track:"
1226 << std::endl;
1227
1228 std::cout
1229 << "position: " << seed.position(m_transient_geocontext).transpose()
1230 << std::endl
1231 << "momentum: " << seed.momentum().transpose()
1232 << std::endl;
1233
1234 std::cout << "charge : " << seed.charge() << std::endl;
1235 std::cout << "absolutemom : " << seed.absoluteMomentum() << std::endl;
1236 }
1237
1238 int PHActsTrkFitter::createNodes(PHCompositeNode* topNode)
1239 {
1240 PHNodeIterator iter(topNode);
1241
1242 PHCompositeNode* dstNode = dynamic_cast<PHCompositeNode*>(iter.findFirst("PHCompositeNode", "DST"));
1243
1244 if (!dstNode)
1245 {
1246 std::cerr << "DST node is missing, quitting" << std::endl;
1247 throw std::runtime_error("Failed to find DST node in PHActsTrkFitter::createNodes");
1248 }
1249
1250 PHNodeIterator dstIter(topNode);
1251 PHCompositeNode* svtxNode = dynamic_cast<PHCompositeNode*>(dstIter.findFirst("PHCompositeNode", "SVTX"));
1252
1253 if (!svtxNode)
1254 {
1255 svtxNode = new PHCompositeNode("SVTX");
1256 dstNode->addNode(svtxNode);
1257 }
1258
1259 if (m_fitSiliconMMs)
1260 {
1261 m_directedTrackMap = findNode::getClass<SvtxTrackMap>(topNode,
1262 "SvtxSiliconMMTrackMap");
1263 if (!m_directedTrackMap)
1264 {
1265
1266 m_directedTrackMap = new SvtxTrackMap_v2;
1267
1268 PHIODataNode<PHObject>* trackNode =
1269 new PHIODataNode<PHObject>(m_directedTrackMap, "SvtxSiliconMMTrackMap", "PHObject");
1270 svtxNode->addNode(trackNode);
1271 }
1272 }
1273
1274 m_trajectories = findNode::getClass<std::map<const unsigned int, Trajectory>>(topNode, m_trajectories_name);
1275 if (!m_trajectories)
1276 {
1277 m_trajectories = new std::map<const unsigned int, Trajectory>;
1278 auto node =
1279 new PHDataNode<std::map<const unsigned int, Trajectory>>(m_trajectories, m_trajectories_name);
1280 svtxNode->addNode(node);
1281 }
1282
1283 m_trackMap = findNode::getClass<SvtxTrackMap>(topNode, _track_map_name);
1284
1285 if (!m_trackMap)
1286 {
1287 m_trackMap = new SvtxTrackMap_v2;
1288 PHIODataNode<PHObject>* node = new PHIODataNode<PHObject>(m_trackMap, _track_map_name, "PHObject");
1289 svtxNode->addNode(node);
1290 }
1291
1292 m_alignmentStateMap = findNode::getClass<SvtxAlignmentStateMap>(topNode, _svtx_alignment_state_map_name);
1293 if (!m_alignmentStateMap)
1294 {
1295 m_alignmentStateMap = new SvtxAlignmentStateMap_v1;
1296 auto node = new PHDataNode<SvtxAlignmentStateMap>(m_alignmentStateMap, _svtx_alignment_state_map_name, "PHObject");
1297 svtxNode->addNode(node);
1298 }
1299
1300 return Fun4AllReturnCodes::EVENT_OK;
1301 }
1302
1303 int PHActsTrkFitter::getNodes(PHCompositeNode* topNode)
1304 {
1305 m_alignmentTransformationMap = findNode::getClass<alignmentTransformationContainer>(topNode, "alignmentTransformationContainer");
1306 if (!m_alignmentTransformationMap)
1307 {
1308 std::cout << PHWHERE << "alignmentTransformationContainer not on node tree. Bailing"
1309 << std::endl;
1310 return Fun4AllReturnCodes::ABORTEVENT;
1311 }
1312
1313 m_alignmentTransformationMapTransient = findNode::getClass<alignmentTransformationContainer>(topNode, "alignmentTransformationContainerTransient");
1314 if (!m_alignmentTransformationMapTransient)
1315 {
1316 std::cout << PHWHERE << "alignmentTransformationContainerTransient not on node tree. Bailing"
1317 << std::endl;
1318 return Fun4AllReturnCodes::ABORTEVENT;
1319 }
1320
1321
1322 m_tpcSeeds = findNode::getClass<TrackSeedContainer>(topNode, "TpcTrackSeedContainer");
1323 if (!m_tpcSeeds)
1324 {
1325 std::cout << PHWHERE << "TpcTrackSeedContainer not on node tree. Bailing"
1326 << std::endl;
1327 return Fun4AllReturnCodes::ABORTEVENT;
1328 }
1329
1330
1331 m_siliconSeeds = findNode::getClass<TrackSeedContainer>(topNode, "SiliconTrackSeedContainer");
1332 if (!m_siliconSeeds)
1333 {
1334 std::cout << PHWHERE << "SiliconTrackSeedContainer not on node tree. Bailing"
1335 << std::endl;
1336 return Fun4AllReturnCodes::ABORTEVENT;
1337 }
1338
1339
1340 m_clusterContainer = findNode::getClass<TrkrClusterContainer>(topNode, m_clusterContainerName);
1341 if (!m_clusterContainer)
1342 {
1343 std::cout << PHWHERE
1344 << "No trkr cluster container, exiting." << std::endl;
1345 return Fun4AllReturnCodes::ABORTEVENT;
1346 }
1347
1348
1349 m_tGeometry = findNode::getClass<ActsGeometry>(topNode, "ActsGeometry");
1350 if (!m_tGeometry)
1351 {
1352 std::cout << "ActsGeometry not on node tree. Exiting."
1353 << std::endl;
1354
1355 return Fun4AllReturnCodes::ABORTEVENT;
1356 }
1357
1358
1359 m_seedMap = findNode::getClass<TrackSeedContainer>(topNode, _svtx_seed_map_name);
1360 if (!m_seedMap)
1361 {
1362 std::cout << "No Svtx seed map on node tree. Exiting."
1363 << std::endl;
1364 return Fun4AllReturnCodes::ABORTEVENT;
1365 }
1366
1367
1368 m_globalPositionWrapper.loadNodes(topNode);
1369 m_globalPositionWrapper.set_suppressCrossing(true);
1370
1371 return Fun4AllReturnCodes::EVENT_OK;
1372 }