Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-08-06 08:18:33

0001 #include "PHSimpleVertexFinder.h"
0002 
0003 /// Tracking includes
0004 
0005 #include <trackbase/TrackVertexCrossingAssoc_v1.h>
0006 #include <trackbase/TrkrCluster.h>  // for TrkrCluster
0007 #include <trackbase/TrkrClusterContainer.h>
0008 #include <trackbase/TrkrDefs.h>            // for cluskey, getLayer, TrkrId
0009 #include <trackbase/TrackFitUtils.h>
0010 #include <trackbase_historic/SvtxTrack.h>  // for SvtxTrack, SvtxTrack::C...
0011 #include <trackbase_historic/SvtxTrackMap_v2.h>
0012 
0013 #include <globalvertex/SvtxVertexMap_v1.h>
0014 #include <globalvertex/SvtxVertex_v2.h>
0015 
0016 #include <fun4all/Fun4AllReturnCodes.h>
0017 #include <phool/PHCompositeNode.h>
0018 #include <phool/PHNode.h>
0019 #include <phool/PHNodeIterator.h>
0020 
0021 #include <phool/getClass.h>
0022 #include <phool/phool.h>
0023 
0024 #include <Acts/Surfaces/PerigeeSurface.hpp>
0025 
0026 #include <cmath>  // for sqrt, fabs, atan2, cos
0027 #include <iomanip>
0028 #include <iostream>  // for operator<<, basic_ostream
0029 #include <map>       // for map
0030 #include <set>       // for _Rb_tree_const_iterator
0031 #include <utility>   // for pair, make_pair
0032 
0033 #include <algorithm>
0034 #include <cassert>
0035 #include <functional>
0036 #include <iostream>
0037 #include <numeric>
0038 #include <vector>
0039 
0040 #include <Eigen/Dense>
0041 
0042 //____________________________________________________________________________..
0043 PHSimpleVertexFinder::PHSimpleVertexFinder(const std::string &name)
0044   : SubsysReco(name)
0045 {
0046 }
0047 
0048 //____________________________________________________________________________..
0049 int PHSimpleVertexFinder::InitRun(PHCompositeNode *topNode)
0050 {
0051   if (Verbosity() > 0)
0052   {
0053     std::cout << __PRETTY_FUNCTION__ << " is pp mode? " << _pp_mode << std::endl;
0054   }
0055 
0056   int ret = GetNodes(topNode);
0057   if (ret != Fun4AllReturnCodes::EVENT_OK)
0058   {
0059     return ret;
0060   }
0061   ret = CreateNodes(topNode);
0062   if (ret != Fun4AllReturnCodes::EVENT_OK)
0063   {
0064     return ret;
0065   }
0066   return ret;
0067 }
0068 
0069 //____________________________________________________________________________..
0070 int PHSimpleVertexFinder::process_event(PHCompositeNode * /*topNode*/)
0071 {
0072   if (Verbosity() > 0)
0073   {
0074     std::cout << PHWHERE << " track map size " << _track_map->size() << std::endl;
0075   }
0076 
0077   _active_dcacut = _base_dcacut;
0078 
0079   // in case these objects are in the input file, we clear the nodes and replace them
0080     _svtx_vertex_map->Reset();
0081     _track_vertex_crossing_map->Reset();
0082 
0083   // Write to a new map on the node tree that contains (crossing, trackid) pairs for all tracks
0084   // Later, will add to it a map  containing (crossing, vertexid)
0085 
0086   std::set<short int> crossings;
0087   for (const auto &[trackkey, track] : *_track_map)
0088   {
0089     auto crossing = track->get_crossing();
0090 
0091     if (Verbosity() > 0)
0092       {
0093     std::cout << "track id " << trackkey << " crossing " << crossing 
0094           << " x y z " << track->get_x() << "  " << track->get_y() << "  " << track->get_z() << std::endl;
0095       }
0096 
0097     // crossing zero contains unmatched TPC tracks
0098     // Here we skip those crossing = zero tracks that do not have silicon seeds
0099     if (_pp_mode)
0100     {
0101       auto siseed = track->get_silicon_seed();
0102       if (crossing == 0 && !siseed)
0103       {
0104       continue;
0105       }
0106     }
0107     
0108     crossings.insert(crossing);
0109     _track_vertex_crossing_map->addTrackAssoc(crossing, trackkey);    
0110   }
0111   
0112   unsigned int vertex_id = 0;
0113 
0114   for (auto cross : crossings)
0115   {
0116     // reset maps for each crossing
0117     _vertex_track_map.clear();
0118     _track_pair_map.clear();
0119     _track_pair_pca_map.clear();
0120     _vertex_position_map.clear();
0121     _vertex_covariance_map.clear();
0122     _vertex_set.clear();
0123 
0124     if (Verbosity() > 0)
0125     {
0126       std::cout << "process tracks for beam crossing " << cross << std::endl;
0127     }
0128 
0129     // get the subset of tracks for this crossing
0130     auto crossing_track_index = _track_vertex_crossing_map->getTracks(cross);
0131     SvtxTrackMap *crossing_tracks = new SvtxTrackMap_v2;
0132     for (auto iter = crossing_track_index.first; iter != crossing_track_index.second; ++iter)
0133     {
0134       unsigned int trackkey = (*iter).second;
0135       SvtxTrack *track = _track_map->get(trackkey);
0136       if(!track)
0137       {
0138         continue;
0139       }
0140       crossing_tracks->insertWithKey(track, trackkey);
0141     }
0142 
0143     // Find all instances where two tracks have a dca of < _dcacut,  and capture the pair details
0144     // Fills _track_pair_map and _track_pair_pca_map
0145     if(_zero_field)
0146       {
0147     checkDCAsZF(crossing_tracks);
0148       }
0149     else
0150       {
0151     checkDCAs(crossing_tracks);
0152       }
0153 
0154     /// If we didn't find any matches, try again with a slightly larger DCA cut
0155     if (_track_pair_map.size() == 0)
0156     {
0157       _active_dcacut = 3.0 * _base_dcacut;
0158       if(_zero_field)
0159     {
0160       checkDCAsZF(crossing_tracks);
0161     }
0162       else
0163     {
0164       checkDCAs(crossing_tracks);
0165     }
0166     }
0167     
0168     if (Verbosity() > 0)
0169     {
0170       std::cout << "crossing " << cross << " track pair map size " << _track_pair_map.size() << std::endl;
0171     }
0172 
0173     // get all connected pairs of tracks by looping over the track_pair map
0174     std::vector<std::set<unsigned int>> connected_tracks = findConnectedTracks();
0175 
0176     // make vertices - each set of connected tracks is a vertex
0177     for (unsigned int ivtx = 0; ivtx < connected_tracks.size(); ++ivtx)
0178     {
0179       if (Verbosity() > 0)
0180       {
0181         std::cout << "process vertex " << ivtx + vertex_id << std::endl;
0182       }
0183 
0184       for (auto it : connected_tracks[ivtx])
0185       {
0186         unsigned int id = it;
0187         _vertex_track_map.insert(std::make_pair(ivtx, id));
0188         if (Verbosity() > 0)
0189         {
0190           std::cout << "  adding track " << id << " to vertex " << ivtx + vertex_id << std::endl;
0191         }
0192       }
0193     }
0194 
0195     // make a list of vertices
0196     for (auto it : _vertex_track_map)
0197     {
0198       if (Verbosity() > 1)
0199       {
0200         std::cout << " vertex " << it.first + vertex_id << " track " << it.second << std::endl;
0201       }
0202       _vertex_set.insert(it.first);
0203     }
0204 
0205     // this finds average vertex positions after removal of outlying track pairs
0206     removeOutlierTrackPairs();
0207 
0208     // average covariance for accepted tracks
0209     for (auto it : _vertex_set)
0210     {
0211       matrix_t avgCov = matrix_t::Zero();
0212       double cov_wt = 0.0;
0213 
0214       auto ret = _vertex_track_map.equal_range(it);
0215       for (auto cit = ret.first; cit != ret.second; ++cit)
0216       {
0217         unsigned int trid = cit->second;
0218         matrix_t cov;
0219         auto track = _track_map->get(trid);
0220         for (int i = 0; i < 3; ++i)
0221         {
0222           for (int j = 0; j < 3; ++j)
0223           {
0224             cov(i, j) = track->get_error(i, j);
0225           }
0226         }
0227 
0228         avgCov += cov;
0229         cov_wt++;
0230       }
0231 
0232       avgCov /= sqrt(cov_wt);
0233       if (Verbosity() > 2)
0234       {
0235         std::cout << "Average covariance for vertex " << it << " is:" << std::endl;
0236         std::cout << std::setprecision(8) << avgCov << std::endl;
0237       }
0238       _vertex_covariance_map.insert(std::make_pair(it, avgCov));
0239     }
0240 
0241     // Write the vertices to the vertex map on the node tree
0242     //==============================================
0243 
0244     for (auto it : _vertex_set)
0245     {
0246       unsigned int thisid = it + vertex_id;  // the address of the vertex in the event
0247 
0248       auto svtxVertex = std::make_unique<SvtxVertex_v2>();
0249 
0250       svtxVertex->set_chisq(0.0);
0251       svtxVertex->set_ndof(0);
0252       svtxVertex->set_t0(0);
0253       svtxVertex->set_id(thisid);
0254       svtxVertex->set_beam_crossing(cross);
0255 
0256       auto ret = _vertex_track_map.equal_range(it);
0257       for (auto cit = ret.first; cit != ret.second; ++cit)
0258       {
0259         unsigned int trid = cit->second;
0260 
0261         if (Verbosity() > 1)
0262         {
0263           std::cout << "   vertex " << thisid << " insert track " << trid << std::endl;
0264         }
0265         svtxVertex->insert_track(trid);
0266         _track_map->get(trid)->set_vertex_id(thisid);
0267       }
0268 
0269       Eigen::Vector3d pos = _vertex_position_map.find(it)->second;
0270       svtxVertex->set_x(pos.x());
0271       svtxVertex->set_y(pos.y());
0272       svtxVertex->set_z(pos.z());
0273       if (Verbosity() > 1)
0274       {
0275         std::cout << "   vertex " << thisid << " insert pos.x " << pos.x() << " pos.y " << pos.y() << " pos.z " << pos.z() << std::endl;
0276       }
0277 
0278       auto vtxCov = _vertex_covariance_map.find(it)->second;
0279       for (int i = 0; i < 3; ++i)
0280       {
0281         for (int j = 0; j < 3; ++j)
0282         {
0283           svtxVertex->set_error(i, j, vtxCov(i, j));
0284         }
0285       }
0286 
0287       _svtx_vertex_map->insert(svtxVertex.release());
0288     }
0289 
0290     vertex_id += _vertex_set.size();
0291 
0292     /// Iterate through the tracks and assign the closest vtx id to
0293     /// the track position for propagating back to the vtx. Catches any
0294     /// tracks that were missed or were not  compatible with any of the
0295     /// identified vertices
0296     //=================================================
0297     for (const auto &[trackkey, track] : *crossing_tracks)
0298     {
0299       auto thistrack = _track_map->get(trackkey);  // get the original, not the copy
0300       auto vtxid = thistrack->get_vertex_id();
0301       if (Verbosity() > 1)
0302       {
0303         std::cout << "        track " << trackkey << " track vtxid " << vtxid << std::endl;
0304       }
0305 
0306       /// If there is a vertex already assigned, keep going
0307       if (vtxid != std::numeric_limits<unsigned int>::max())
0308       {
0309         continue;
0310       }
0311 
0312       float maxdz = std::numeric_limits<float>::max();
0313       unsigned int newvtxid = std::numeric_limits<unsigned int>::max();
0314 
0315       for (auto it : _vertex_set)
0316       {
0317         unsigned int thisid = it + vertex_id - _vertex_set.size();
0318 
0319         if (Verbosity() > 1)
0320         {
0321           std::cout << "                test vertex " << thisid << std::endl;
0322         }
0323 
0324         auto thisvertex = _svtx_vertex_map->get(thisid);
0325         float dz = thistrack->get_z() - thisvertex->get_z();
0326         if (std::fabs(dz) < maxdz)
0327         {
0328           maxdz = dz;
0329           newvtxid = thisid;
0330         }
0331       }
0332 
0333       // this updates the track, but does not add it to the vertex primary track list
0334       if (newvtxid != std::numeric_limits<unsigned int>::max())
0335       {
0336         thistrack->set_vertex_id(newvtxid);
0337         if (Verbosity() > 1)
0338         {
0339           std::cout << "                assign vertex " << newvtxid << " to additional track " << trackkey << std::endl;
0340         }
0341       }
0342     }
0343 
0344     delete crossing_tracks;
0345 
0346   }  // end loop over crossings
0347 
0348   // update the crossing vertex map with the results
0349   for (const auto &[vtxkey, vertex] : *_svtx_vertex_map)
0350   {
0351     short int crossing = vertex->get_beam_crossing();
0352     _track_vertex_crossing_map->addVertexAssoc(crossing, vtxkey);
0353 
0354     if (Verbosity() > 1)
0355     {
0356       std::cout << "Vertex ID: " << vtxkey << " vertex crossing " << crossing << " list of tracks: " << std::endl;
0357       for (auto trackiter = vertex->begin_tracks(); trackiter != vertex->end_tracks(); ++trackiter)
0358       {
0359         SvtxTrack *track = _track_map->get(*trackiter);
0360         if (!track)
0361         {
0362           continue;
0363         }
0364 
0365         auto siseed = track->get_silicon_seed();
0366         short int intt_crossing = siseed->get_crossing();
0367 
0368         // the track crossing may be from the INTT clusters or from geometric matching if there are no INTT clusters
0369         short int track_crossing = track->get_crossing();
0370         std::cout << " vtxid " << vtxkey << " vertex crossing " << crossing
0371                   << " track crossing " << track_crossing
0372                   << " intt crossing " << intt_crossing
0373                   << " trackID " << *trackiter
0374                   << " track Z " << track->get_z()
0375                   << " X " << track->get_x()
0376                   << " Y " << track->get_y()
0377                   << " quality " << track->get_quality()
0378                   << " pt " << track->get_pt()
0379                   << std::endl;
0380         if (Verbosity() > 3)
0381         {
0382           siseed->identify();
0383         }
0384       }
0385     }
0386   }
0387 
0388   if (Verbosity() > 2)
0389   {
0390     _track_vertex_crossing_map->identify();
0391   }
0392 
0393   return Fun4AllReturnCodes::EVENT_OK;
0394 }
0395 
0396 int PHSimpleVertexFinder::End(PHCompositeNode * /*topNode*/)
0397 {
0398   return Fun4AllReturnCodes::EVENT_OK;
0399 }
0400 
0401 int PHSimpleVertexFinder::CreateNodes(PHCompositeNode *topNode)
0402 {
0403   PHNodeIterator iter(topNode);
0404 
0405   PHCompositeNode *dstNode = dynamic_cast<PHCompositeNode *>(iter.findFirst("PHCompositeNode", "DST"));
0406 
0407   /// Check that it is there
0408   if (!dstNode)
0409   {
0410     std::cerr << "DST Node missing, quitting" << std::endl;
0411     throw std::runtime_error("failed to find DST node in PHActsInitialVertexFinder::createNodes");
0412   }
0413 
0414   /// Get the tracking subnode
0415   PHCompositeNode *svtxNode = dynamic_cast<PHCompositeNode *>(iter.findFirst("PHCompositeNode", "SVTX"));
0416 
0417   if (!svtxNode)
0418   {
0419     svtxNode = new PHCompositeNode("SVTX");
0420     dstNode->addNode(svtxNode);
0421   }
0422 
0423   _svtx_vertex_map = findNode::getClass<SvtxVertexMap>(topNode, _vertex_map_name);
0424 
0425   if (!_svtx_vertex_map)
0426   {
0427     _svtx_vertex_map = new SvtxVertexMap_v1;
0428     PHIODataNode<PHObject> *vertexNode = new PHIODataNode<PHObject>(
0429         _svtx_vertex_map, _vertex_map_name, "PHObject");
0430 
0431     svtxNode->addNode(vertexNode);
0432   }
0433 
0434   _track_vertex_crossing_map = findNode::getClass<TrackVertexCrossingAssoc>(topNode, "TrackVertexCrossingAssocMap");
0435   if (!_track_vertex_crossing_map)
0436   {
0437     _track_vertex_crossing_map = new TrackVertexCrossingAssoc_v1;
0438     PHIODataNode<PHObject> *trackvertexNode = new PHIODataNode<PHObject>(
0439         _track_vertex_crossing_map, "TrackVertexCrossingAssocMap", "PHObject");
0440 
0441     svtxNode->addNode(trackvertexNode);
0442   }
0443 
0444   return Fun4AllReturnCodes::EVENT_OK;
0445 }
0446 int PHSimpleVertexFinder::GetNodes(PHCompositeNode *topNode)
0447 {
0448   _track_map = findNode::getClass<SvtxTrackMap>(topNode, _track_map_name);
0449   if (!_track_map)
0450   {
0451     std::cout << PHWHERE << " ERROR: Can't find SvtxTrackMap: " << _track_map_name << std::endl;
0452     return Fun4AllReturnCodes::ABORTEVENT;
0453   }
0454 
0455 
0456   _cluster_map = findNode::getClass<TrkrClusterContainer>(topNode, m_clusterContainerName);
0457   if (!_cluster_map)
0458   {
0459     std::cout << PHWHERE << " ERROR: Can't find node TRKR_CLUSTER container" << std::endl;
0460     return Fun4AllReturnCodes::ABORTEVENT;
0461   }
0462 
0463   _tGeometry = findNode::getClass<ActsGeometry>(topNode, "ActsGeometry");
0464   if (!_tGeometry)
0465   {
0466     std::cout << PHWHERE << "Error, can't find acts tracking geometry" << std::endl;
0467     return Fun4AllReturnCodes::ABORTEVENT;
0468   }
0469 
0470   return Fun4AllReturnCodes::EVENT_OK;
0471 }
0472 
0473 void PHSimpleVertexFinder::checkDCAs(SvtxTrackMap *track_map)
0474 {
0475   // Loop over tracks and check for close DCA match with all other tracks
0476   for (auto tr1_it = track_map->begin(); tr1_it != track_map->end(); ++tr1_it)
0477   {
0478     auto id1 = tr1_it->first;
0479     auto tr1 = tr1_it->second;
0480     if (tr1->get_quality() > _qual_cut)
0481     {
0482       continue;
0483     }
0484     if (_require_mvtx)
0485     {
0486       unsigned int nmvtx = 0;
0487       TrackSeed *siliconseed = tr1->get_silicon_seed();
0488       if (!siliconseed)
0489       {
0490         continue;
0491       }
0492 
0493       for (auto clusit = siliconseed->begin_cluster_keys(); clusit != siliconseed->end_cluster_keys(); ++clusit)
0494       {
0495         if (TrkrDefs::getTrkrId(*clusit) == TrkrDefs::mvtxId)
0496         {
0497           nmvtx++;
0498         }
0499         if (nmvtx >= _nmvtx_required)
0500         {
0501           break;
0502         }
0503       }
0504       if (nmvtx < _nmvtx_required)
0505       {
0506         continue;
0507       }
0508       if (Verbosity() > 3)
0509       {
0510         std::cout << " tr1 id " << id1 << " has nmvtx at least " << nmvtx << std::endl;
0511       }
0512     }
0513 
0514     // look for close DCA matches with all other such tracks
0515     for (auto tr2_it = std::next(tr1_it); tr2_it != track_map->end(); ++tr2_it)
0516     {
0517       auto id2 = tr2_it->first;
0518       auto tr2 = tr2_it->second;
0519       if (tr2->get_quality() > _qual_cut)
0520       {
0521         continue;
0522       }
0523       if (_require_mvtx)
0524       {
0525         unsigned int nmvtx = 0;
0526         TrackSeed *siliconseed = tr2->get_silicon_seed();
0527         if (!siliconseed)
0528         {
0529           continue;
0530         }
0531 
0532         for (auto clusit = siliconseed->begin_cluster_keys(); clusit != siliconseed->end_cluster_keys(); ++clusit)
0533         {
0534           if (TrkrDefs::getTrkrId(*clusit) == TrkrDefs::mvtxId)
0535           {
0536             nmvtx++;
0537           }
0538           if (nmvtx >= _nmvtx_required)
0539           {
0540             break;
0541           }
0542         }
0543         if (nmvtx < _nmvtx_required)
0544         {
0545           continue;
0546         }
0547         if (Verbosity() > 3)
0548         {
0549           std::cout << " tr2 id " << id2 << " has nmvtx at least " << nmvtx << std::endl;
0550         }
0551       }
0552 
0553       // find DCA of these two tracks
0554       if (Verbosity() > 3)
0555       {
0556         std::cout << "Check DCA for tracks " << id1 << " and  " << id2 << std::endl;
0557       }
0558 
0559       findDcaTwoTracks(tr1, tr2);
0560     }
0561   }
0562 }
0563 
0564 void PHSimpleVertexFinder::checkDCAsZF(SvtxTrackMap *track_map)
0565 {
0566   // ZF tracks do not have an Acts fit, and the seeding does not give
0567   // reliable track parameters - refit clusters with straight lines
0568   // No distortion corrections applied in TPC at present
0569 
0570   std::vector<unsigned int> cumulative_trackid_vec;
0571   std::vector<unsigned int> cumulative_nmvtx_vec;
0572   std::vector<unsigned int> cumulative_nintt_vec;
0573   std::vector<std::vector<Acts::Vector3>> cumulative_global_vec;
0574   std::vector<std::vector<TrkrDefs::cluskey>> cumulative_cluskey_vec;
0575   std::vector<std::vector<float>> cumulative_fitpars_vec;
0576 
0577   for (auto & tr1_it : *track_map)
0578   {
0579     auto id1 = tr1_it.first;
0580     auto tr1 = tr1_it.second;
0581 
0582     //    tr1->identify();
0583  
0584     TrackSeed *siliconseed = tr1->get_silicon_seed();
0585     if (_require_mvtx)
0586       {
0587     if (!siliconseed)
0588       {
0589         continue;
0590       }
0591       }
0592     TrackSeed *tpcseed = tr1->get_tpc_seed();
0593 
0594     std::vector<Acts::Vector3> global_vec;
0595     std::vector<TrkrDefs::cluskey> cluskey_vec;
0596     
0597     // Get a vector of cluster keys from the silicon seed and TPC seed
0598     if(siliconseed)
0599       {
0600     getTrackletClusterList(siliconseed, cluskey_vec);
0601     if(Verbosity() > 0) 
0602       { 
0603         std::cout << "  after silicon: silicon cluskey_vec size " << cluskey_vec.size() << std::endl; 
0604         for(unsigned long i : cluskey_vec)
0605           {  
0606         std::cout << i << std::endl;
0607           }
0608       }
0609       }
0610     if(tpcseed)
0611       {
0612 
0613     getTrackletClusterList(tpcseed, cluskey_vec);
0614     if(Verbosity() > 0) 
0615       { 
0616         std::cout << "  after tpc: cluskey_vec size " << cluskey_vec.size() << std::endl; 
0617         for(unsigned long i : cluskey_vec)
0618           {  
0619         std::cout << i << std::endl;
0620           }
0621       }
0622       }
0623 
0624     unsigned int nmvtx = 0;
0625     unsigned int nintt = 0;
0626     for (auto& key : cluskey_vec)
0627       {
0628         if (TrkrDefs::getTrkrId(key) == TrkrDefs::mvtxId)
0629       {
0630         nmvtx++;
0631       }
0632     if(TrkrDefs::getTrkrId(key) == TrkrDefs::inttId)
0633       {
0634         nintt++;
0635       }
0636       }
0637     
0638     // store cluster global positions in a vector
0639     TrackFitUtils::getTrackletClusters(_tGeometry, _cluster_map, global_vec, cluskey_vec);
0640     
0641     std::vector<float> fitpars = TrackFitUtils::fitClustersZeroField(global_vec, cluskey_vec, true);
0642     
0643     if (Verbosity() > 1)
0644       {
0645     if(fitpars.size() == 4)
0646       {
0647         std::cout << " Track " << id1 << " dy/dx " << fitpars[0] << " y intercept " << fitpars[1] 
0648               << " dx/dz " << fitpars[2] << " Z0 " << fitpars[3] << std::endl;
0649       }
0650     else
0651       {
0652         std::cout << " Track " << id1 << " ZF line fits failed, fitpars is empty" << std::endl;
0653       }
0654       }
0655 
0656     cumulative_trackid_vec.push_back(id1);   
0657     cumulative_nmvtx_vec.push_back(nmvtx);   
0658     cumulative_nintt_vec.push_back(nintt);   
0659     cumulative_global_vec.push_back(global_vec);
0660     cumulative_cluskey_vec.push_back(cluskey_vec);
0661     cumulative_fitpars_vec.push_back(fitpars);
0662   }
0663 
0664   for(unsigned int i1 = 0; i1 < cumulative_trackid_vec.size(); ++i1)
0665     {
0666       if(cumulative_fitpars_vec[i1].size() == 0) { continue; }
0667 
0668       for(unsigned int i2 = i1; i2 < cumulative_trackid_vec.size(); ++i2)
0669     {
0670       if(cumulative_fitpars_vec[i2].size() == 0) { continue; }
0671 
0672       //  For straight line: fitpars[4] = { xyslope, y0, xzslope, z0 }
0673       Eigen::Vector3d a1(0.0, cumulative_fitpars_vec[i1][1],cumulative_fitpars_vec[i1][3]);      // point on track 1 at x = 0
0674       Eigen::Vector3d a2(0.0, cumulative_fitpars_vec[i2][1],cumulative_fitpars_vec[i2][3]);      // point on track 2 at x = 0
0675       // direction vectors made from dy/dx = xyslope and dz/dx = xzslope
0676       Eigen::Vector3d b1(1.0, cumulative_fitpars_vec[i1][0],cumulative_fitpars_vec[i1][2]);      // direction vector of track 1
0677       Eigen::Vector3d b2(1.0, cumulative_fitpars_vec[i2][0],cumulative_fitpars_vec[i2][2]);      // direction vector of track 2
0678           
0679       Eigen::Vector3d PCA1(0, 0, 0);
0680       Eigen::Vector3d PCA2(0, 0, 0);
0681       double dca = dcaTwoLines(a1, b1, a2, b2, PCA1, PCA2);
0682 
0683 
0684       // check dca cut is satisfied, and that PCA is close to beam line
0685       if (fabs(dca) < _active_dcacut && (fabs(PCA1.x()) < _beamline_xy_cut && fabs(PCA1.y()) < _beamline_xy_cut))
0686         {
0687           int id1 = cumulative_trackid_vec[i1];
0688           int id2 = cumulative_trackid_vec[i2];
0689 
0690           if (Verbosity() > 3)
0691         {
0692           std::cout << " good match for tracks " << id1 << " and " << id2 << std::endl;
0693           std::cout << "    a1.x " << a1.x() << " a1.y " << a1.y() << " a1.z " << a1.z() << std::endl;
0694           std::cout << "    a2.x  " << a2.x() << " a2.y " << a2.y() << " a2.z " << a2.z() << std::endl;
0695           std::cout << "    PCA1.x() " << PCA1.x() << " PCA1.y " << PCA1.y() << " PCA1.z " << PCA1.z() << std::endl;
0696           std::cout << "    PCA2.x() " << PCA2.x() << " PCA2.y " << PCA2.y() << " PCA2.z " << PCA2.z() << std::endl;
0697           std::cout << "    dca " << dca << std::endl;
0698         }
0699           
0700           // capture the results for successful matches
0701           _track_pair_map.insert(std::make_pair(id1, std::make_pair(id2, dca)));
0702           _track_pair_pca_map.insert(std::make_pair(id1, std::make_pair(id2, std::make_pair(PCA1, PCA2))));
0703         }         
0704     }
0705     }
0706 
0707   return; 
0708 }
0709 
0710 void PHSimpleVertexFinder::getTrackletClusterList(TrackSeed* tracklet, std::vector<TrkrDefs::cluskey>& cluskey_vec)
0711 {
0712   for (auto clusIter = tracklet->begin_cluster_keys();
0713        clusIter != tracklet->end_cluster_keys();
0714        ++clusIter)
0715   {
0716     auto key = *clusIter;
0717     auto cluster = _cluster_map->findCluster(key);
0718     if (!cluster)
0719     {
0720       std::cout << PHWHERE << "Failed to get cluster with key " << key << std::endl;
0721       continue;
0722     }
0723 
0724     /// Make a safety check for clusters that couldn't be attached to a surface
0725     auto surf = _tGeometry->maps().getSurface(key, cluster);
0726     if (!surf)
0727     {
0728       continue;
0729     }
0730 
0731     // drop some bad layers in the TPC completely
0732     unsigned int layer = TrkrDefs::getLayer(key);
0733     if (layer == 7 || layer == 22 || layer == 23 || layer == 38 || layer == 39)
0734     {
0735       continue;
0736     }
0737 
0738     // drop INTT clusters for now  -- TEMPORARY!
0739     if (layer > 2 && layer < 7)
0740     {
0741       continue;
0742     }
0743 
0744 
0745     cluskey_vec.push_back(key);
0746 
0747   }  // end loop over clusters for this track
0748 }
0749 
0750 void PHSimpleVertexFinder::checkDCAs()
0751 {
0752   // Loop over tracks and check for close DCA match with all other tracks
0753   for (auto tr1_it = _track_map->begin(); tr1_it != _track_map->end(); ++tr1_it)
0754   {
0755     auto id1 = tr1_it->first;
0756     auto tr1 = tr1_it->second;
0757     if (tr1->get_quality() > _qual_cut)
0758     {
0759       continue;
0760     }
0761     if (_require_mvtx)
0762     {
0763       unsigned int nmvtx = 0;
0764       TrackSeed *siliconseed = tr1->get_silicon_seed();
0765       if (!siliconseed)
0766       {
0767         continue;
0768       }
0769 
0770       for (auto clusit = siliconseed->begin_cluster_keys(); clusit != siliconseed->end_cluster_keys(); ++clusit)
0771       {
0772         if (TrkrDefs::getTrkrId(*clusit) == TrkrDefs::mvtxId)
0773         {
0774           nmvtx++;
0775         }
0776         if (nmvtx >= _nmvtx_required)
0777         {
0778           break;
0779         }
0780       }
0781       if (nmvtx < _nmvtx_required)
0782       {
0783         continue;
0784       }
0785       if (Verbosity() > 3)
0786       {
0787         std::cout << " tr1 id " << id1 << " has nmvtx at least " << nmvtx << std::endl;
0788       }
0789     }
0790 
0791     // look for close DCA matches with all other such tracks
0792     for (auto tr2_it = std::next(tr1_it); tr2_it != _track_map->end(); ++tr2_it)
0793     {
0794       auto id2 = tr2_it->first;
0795       auto tr2 = tr2_it->second;
0796       if (tr2->get_quality() > _qual_cut)
0797       {
0798         continue;
0799       }
0800       if (_require_mvtx)
0801       {
0802         unsigned int nmvtx = 0;
0803         TrackSeed *siliconseed = tr2->get_silicon_seed();
0804         if (!siliconseed)
0805         {
0806           continue;
0807         }
0808 
0809         for (auto clusit = siliconseed->begin_cluster_keys(); clusit != siliconseed->end_cluster_keys(); ++clusit)
0810         {
0811           if (TrkrDefs::getTrkrId(*clusit) == TrkrDefs::mvtxId)
0812           {
0813             nmvtx++;
0814           }
0815           if (nmvtx >= _nmvtx_required)
0816           {
0817             break;
0818           }
0819         }
0820         if (nmvtx < _nmvtx_required)
0821         {
0822           continue;
0823         }
0824         if (Verbosity() > 3)
0825         {
0826           std::cout << " tr2 id " << id2 << " has nmvtx at least " << nmvtx << std::endl;
0827         }
0828       }
0829 
0830       // find DCA of these two tracks
0831       if (Verbosity() > 3)
0832       {
0833         std::cout << "Check DCA for tracks " << id1 << " and  " << id2 << std::endl;
0834       }
0835 
0836       findDcaTwoTracks(tr1, tr2);
0837     }
0838   }
0839 }
0840 
0841 void PHSimpleVertexFinder::findDcaTwoTracks(SvtxTrack *tr1, SvtxTrack *tr2)
0842 {
0843   if (tr1->get_pt() < _track_pt_cut)
0844   {
0845     return;
0846   }
0847   if (tr2->get_pt() < _track_pt_cut)
0848   {
0849     return;
0850   }
0851 
0852   unsigned int id1 = tr1->get_id();
0853   unsigned int id2 = tr2->get_id();
0854 
0855   // get the line equations for the tracks
0856 
0857   Eigen::Vector3d a1(tr1->get_x(), tr1->get_y(), tr1->get_z());
0858   Eigen::Vector3d b1(tr1->get_px() / tr1->get_p(), tr1->get_py() / tr1->get_p(), tr1->get_pz() / tr1->get_p());
0859   Eigen::Vector3d a2(tr2->get_x(), tr2->get_y(), tr2->get_z());
0860   Eigen::Vector3d b2(tr2->get_px() / tr2->get_p(), tr2->get_py() / tr2->get_p(), tr2->get_pz() / tr2->get_p());
0861 
0862   Eigen::Vector3d PCA1(0, 0, 0);
0863   Eigen::Vector3d PCA2(0, 0, 0);
0864   double dca = dcaTwoLines(a1, b1, a2, b2, PCA1, PCA2);
0865 
0866   if (Verbosity() > 3)
0867   {
0868     std::cout << " pair dca is " << dca << " _active_dcacut is " << _active_dcacut
0869               << " PCA1.x " << PCA1.x() << " PCA1.y " << PCA1.y()
0870               << " PCA2.x " << PCA2.x() << " PCA2.y " << PCA2.y() << std::endl;
0871   }
0872 
0873   // check dca cut is satisfied, and that PCA is close to beam line
0874   if (fabs(dca) < _active_dcacut && (fabs(PCA1.x()) < _beamline_xy_cut && fabs(PCA1.y()) < _beamline_xy_cut))
0875   {
0876     if (Verbosity() > 3)
0877     {
0878       std::cout << " good match for tracks " << tr1->get_id() << " and " << tr2->get_id() << " with pT " << tr1->get_pt() << " and " << tr2->get_pt() << std::endl;
0879       std::cout << "    a1.x " << a1.x() << " a1.y " << a1.y() << " a1.z " << a1.z() << std::endl;
0880       std::cout << "    a2.x  " << a2.x() << " a2.y " << a2.y() << " a2.z " << a2.z() << std::endl;
0881       std::cout << "    PCA1.x() " << PCA1.x() << " PCA1.y " << PCA1.y() << " PCA1.z " << PCA1.z() << std::endl;
0882       std::cout << "    PCA2.x() " << PCA2.x() << " PCA2.y " << PCA2.y() << " PCA2.z " << PCA2.z() << std::endl;
0883       std::cout << "    dca " << dca << std::endl;
0884     }
0885 
0886     // capture the results for successful matches
0887     _track_pair_map.insert(std::make_pair(id1, std::make_pair(id2, dca)));
0888     _track_pair_pca_map.insert(std::make_pair(id1, std::make_pair(id2, std::make_pair(PCA1, PCA2))));
0889   }
0890 
0891   return;
0892 }
0893 
0894 double PHSimpleVertexFinder::dcaTwoLines(const Eigen::Vector3d &a1, const Eigen::Vector3d &b1,
0895                                          const Eigen::Vector3d &a2, const Eigen::Vector3d &b2,
0896                                          Eigen::Vector3d &PCA1, Eigen::Vector3d &PCA2)
0897 {
0898   // The shortest distance between two skew lines described by
0899   //  a1 + c * b1
0900   //  a2 + d * b2
0901   // where a1, a2, are vectors representing points on the lines, b1, b2 are direction vectors, and c and d are scalars
0902   // is:
0903   // dca = (b1 x b2) .(a2-a1) / |b1 x b2|
0904 
0905   // bcrossb/mag_bcrossb is a unit vector perpendicular to both direction vectors b1 and b2
0906   auto bcrossb = b1.cross(b2);
0907   auto mag_bcrossb = bcrossb.norm();
0908   // a2-a1 is the vector joining any arbitrary points on the two lines
0909   auto aminusa = a2 - a1;
0910 
0911   // The DCA of these two lines is the projection of a2-a1 along the direction of the perpendicular to both
0912   // remember that a2-a1 is longer than (or equal to) the dca by definition
0913   double dca = 999;
0914   if (mag_bcrossb != 0)
0915   {
0916     dca = bcrossb.dot(aminusa) / mag_bcrossb;
0917   }
0918   else
0919   {
0920     return dca;  // same track, skip combination
0921   }
0922 
0923   // get the points at which the normal to the lines intersect the lines, where the lines are perpendicular
0924   // Assume the shortest distance occurs at points A on line 1, and B on line 2, call the line AB
0925   //    AB = a2+d*b2 - (a1+c*b1)
0926   // we need to find c and d where AB is perpendicular to both lines. so AB.b1 = 0 and AB.b2 = 0
0927   //    ( (a2 -a1) + d*b2 -c*b1 ).b1 = 0
0928   //    ( (a2 -a1) + d*b2 -c*b1 ).b2 = 0
0929   // so we have two simultaneous equations in 2 unknowns
0930   //    (a2-a1).b1 + d*b2.b1 -c*b1.b1 = 0 => d*b2.b1 = c*b1.b1 - (a2-a1).b1 => d = (1/b2.b1) * (c*b1.b1 - (a2-a1).b1)
0931   //    (a2-a1).b2 + d*b2.b2 - c*b1.b2 = 0 => c*b1.b2 =  (a2-a1).b2 + [(1/b2.b1) * (c*b1*b1 -(a2-a1).b1)}*b2.b2
0932   //    c*b1.b2 - (1/b2.b1) * c*b1.b1*b2.b2 = (a2-a1).b2 - (1/b2.b1)*(a2-a1).b1*b2.b2
0933   //    c *(b1.b2 - (1/b2.b1)*b1.b1*b2.b2)  = (a2-a1).b2 - (1/b2.b1)*(a2-a1).b1*b2.b2
0934   // call this: c*X = Y
0935   // plug back into the d equation
0936   //     d = c*b1.b1 / b2.b1 - (a2-a1).b1 / b2.b1
0937   // and call the d equation: d = c * F - G
0938 
0939   double X = b1.dot(b2) - b1.dot(b1) * b2.dot(b2) / b2.dot(b1);
0940   double Y = (a2.dot(b2) - a1.dot(b2)) - (a2.dot(b1) - a1.dot(b1)) * b2.dot(b2) / b2.dot(b1);
0941   double c = Y / X;
0942 
0943   double F = b1.dot(b1) / b2.dot(b1);
0944   double G = -(a2.dot(b1) - a1.dot(b1)) / b2.dot(b1);
0945   double d = c * F + G;
0946 
0947   // then the points of closest approach are:
0948   PCA1 = a1 + c * b1;
0949   PCA2 = a2 + d * b2;
0950 
0951   return dca;
0952 }
0953 
0954 std::vector<std::set<unsigned int>> PHSimpleVertexFinder::findConnectedTracks()
0955 {
0956   std::vector<std::set<unsigned int>> connected_tracks;
0957   std::set<unsigned int> connected;
0958   std::set<unsigned int> used;
0959   for (auto it : _track_pair_map)
0960   {
0961     unsigned int id1 = it.first;
0962     unsigned int id2 = it.second.first;
0963 
0964     if ((used.find(id1) != used.end()) && (used.find(id2) != used.end()))
0965     {
0966       if (Verbosity() > 3)
0967       {
0968         std::cout << " tracks " << id1 << " and " << id2 << " are both in used , skip them" << std::endl;
0969       }
0970       continue;
0971     }
0972     else if ((used.find(id1) == used.end()) && (used.find(id2) == used.end()))
0973     {
0974       if (Verbosity() > 3)
0975       {
0976         std::cout << " tracks " << id1 << " and " << id2 << " are both not in used , start a new connected set" << std::endl;
0977       }
0978       // close out and start a new connections set
0979       if (connected.size() > 0)
0980       {
0981         connected_tracks.push_back(connected);
0982         connected.clear();
0983         if (Verbosity() > 3)
0984         {
0985           std::cout << "           closing out set " << std::endl;
0986         }
0987       }
0988     }
0989 
0990     // get everything connected to id1 and id2
0991     connected.insert(id1);
0992     used.insert(id1);
0993     connected.insert(id2);
0994     used.insert(id2);
0995     for (auto cit : _track_pair_map)
0996     {
0997       unsigned int id3 = cit.first;
0998       unsigned int id4 = cit.second.first;
0999       if ((connected.find(id3) != connected.end()) || (connected.find(id4) != connected.end()))
1000       {
1001         if (Verbosity() > 3)
1002         {
1003           std::cout << " found connection to " << id3 << " and " << id4 << std::endl;
1004         }
1005         connected.insert(id3);
1006         used.insert(id3);
1007         connected.insert(id4);
1008         used.insert(id4);
1009       }
1010     }
1011   }
1012 
1013   // close out the last set
1014   if (connected.size() > 0)
1015   {
1016     connected_tracks.push_back(connected);
1017     connected.clear();
1018     if (Verbosity() > 3)
1019     {
1020       std::cout << "           closing out last set " << std::endl;
1021     }
1022   }
1023 
1024   if (Verbosity() > 3)
1025   {
1026     std::cout << "connected_tracks size " << connected_tracks.size() << std::endl;
1027   }
1028 
1029   return connected_tracks;
1030 }
1031 
1032 void PHSimpleVertexFinder::removeOutlierTrackPairs()
1033 {
1034   //  Note: std::multimap<unsigned int, std::pair<unsigned int, std::pair<Eigen::Vector3d,  Eigen::Vector3d>>>  _track_pair_pca_map
1035 
1036   for (auto it : _vertex_set)
1037   {
1038     unsigned int vtxid = it;
1039     if (Verbosity() > 1)
1040     {
1041       std::cout << "calculate average position for vertex " << vtxid << std::endl;
1042     }
1043 
1044     // we need the median values of the x and y positions
1045     std::vector<double> vx;
1046     std::vector<double> vy;
1047     std::vector<double> vz;
1048 
1049     double pca_median_x = 0.;
1050     double pca_median_y = 0.;
1051     double pca_median_z = 0.;
1052 
1053     Eigen::Vector3d new_pca_avge(0., 0., 0.);
1054     double new_wt = 0.0;
1055 
1056     auto ret = _vertex_track_map.equal_range(vtxid);
1057 
1058     // Start by getting the positions for this vertex into vectors for the median calculation
1059     for (auto cit = ret.first; cit != ret.second; ++cit)
1060     {
1061       unsigned int tr1id = cit->second;
1062       if (Verbosity() > 2)
1063       {
1064         std::cout << "   vectors: get entries for track " << tr1id << " for vertex " << vtxid << std::endl;
1065       }
1066 
1067       // find all pairs for this vertex with tr1id
1068       auto pca_range = _track_pair_pca_map.equal_range(tr1id);
1069       for (auto pit = pca_range.first; pit != pca_range.second; ++pit)
1070       {
1071         unsigned int tr2id = pit->second.first;
1072 
1073         Eigen::Vector3d PCA1 = pit->second.second.first;
1074         Eigen::Vector3d PCA2 = pit->second.second.second;
1075 
1076         if (Verbosity() > 2)
1077         {
1078           std::cout << " vectors: tr1id " << tr1id << " tr2id " << tr2id
1079                     << " PCA1 " << PCA1.x() << "  " << PCA1.y() << "  " << PCA1.z()
1080                     << " PCA2 " << PCA2.x() << "  " << PCA2.y() << "  " << PCA2.z()
1081                     << std::endl;
1082         }
1083 
1084         vx.push_back(PCA1.x());
1085         vx.push_back(PCA2.x());
1086         vy.push_back(PCA1.y());
1087         vy.push_back(PCA2.y());
1088         vz.push_back(PCA1.z());
1089         vz.push_back(PCA2.z());
1090       }
1091     }
1092 
1093     // Get the medians for this vertex
1094     // Using the median as a reference for rejecting outliers only makes sense for more than 2 tracks
1095     if (vx.size() < 3)
1096     {
1097       new_pca_avge.x() = getAverage(vx);
1098       new_pca_avge.y() = getAverage(vy);
1099       new_pca_avge.z() = getAverage(vz);
1100       _vertex_position_map.insert(std::make_pair(vtxid, new_pca_avge));
1101       if (Verbosity() > 1)
1102       {
1103         std::cout << " Vertex has only 2 tracks, use average for PCA: " << new_pca_avge.x() << "  " << new_pca_avge.y() << "  " << new_pca_avge.z() << std::endl;
1104       }
1105 
1106       // done with this vertex
1107       continue;
1108     }
1109 
1110     pca_median_x = getMedian(vx);
1111     pca_median_y = getMedian(vy);
1112     pca_median_z = getMedian(vz);
1113     if (Verbosity() > 1)
1114     {
1115       std::cout << "Median values: x " << pca_median_x << " y " << pca_median_y << " z : " << pca_median_z << std::endl;
1116     }
1117 
1118     // Make the average vertex position with outlier rejection wrt the median
1119     for (auto cit = ret.first; cit != ret.second; ++cit)
1120     {
1121       unsigned int tr1id = cit->second;
1122       if (Verbosity() > 2)
1123       {
1124         std::cout << "   average: get entries for track " << tr1id << " for vertex " << vtxid << std::endl;
1125       }
1126 
1127       // find all pairs for this vertex with tr1id
1128       auto pca_range = _track_pair_pca_map.equal_range(tr1id);
1129       for (auto pit = pca_range.first; pit != pca_range.second; ++pit)
1130       {
1131         unsigned int tr2id = pit->second.first;
1132 
1133         Eigen::Vector3d PCA1 = pit->second.second.first;
1134         Eigen::Vector3d PCA2 = pit->second.second.second;
1135 
1136         if (
1137             fabs(PCA1.x() - pca_median_x) < _outlier_cut &&
1138             fabs(PCA1.y() - pca_median_y) < _outlier_cut &&
1139             fabs(PCA2.x() - pca_median_x) < _outlier_cut &&
1140             fabs(PCA2.y() - pca_median_y) < _outlier_cut)
1141         {
1142           // good track pair, add to new average
1143 
1144           new_pca_avge += PCA1;
1145           new_wt++;
1146           new_pca_avge += PCA2;
1147           new_wt++;
1148         }
1149         else
1150         {
1151           if (Verbosity() > 1)
1152           {
1153             std::cout << "Reject pair with tr1id " << tr1id << " tr2id " << tr2id << std::endl;
1154           }
1155         }
1156       }
1157     }
1158     if (new_wt > 0.0)
1159     {
1160       new_pca_avge = new_pca_avge / new_wt;
1161     }
1162     else
1163     {
1164       // There were no pairs that survived the track cuts, use the median values
1165       new_pca_avge.x() = pca_median_x;
1166       new_pca_avge.y() = pca_median_y;
1167       new_pca_avge.z() = pca_median_z;
1168     }
1169 
1170     _vertex_position_map.insert(std::make_pair(vtxid, new_pca_avge));
1171   }
1172 
1173   return;
1174 }
1175 
1176 double PHSimpleVertexFinder::getMedian(std::vector<double> &v)
1177 {
1178   double median = 0.0;
1179 
1180   if ((v.size() % 2) == 0)
1181   {
1182     // even number of entries
1183     // we want the average of the middle two numbers, v.size()/2 and v.size()/2-1
1184     auto m1 = v.begin() + v.size() / 2;
1185     std::nth_element(v.begin(), m1, v.end());
1186     double median1 = v[v.size() / 2];
1187 
1188     auto m2 = v.begin() + v.size() / 2 - 1;
1189     std::nth_element(v.begin(), m2, v.end());
1190     double median2 = v[v.size() / 2 - 1];
1191 
1192     median = (median1 + median2) / 2.0;
1193     if (Verbosity() > 2)
1194     {
1195       std::cout << "The vector size is " << v.size()
1196                 << " element m is " << v.size() / 2 << " = " << v[v.size() / 2]
1197                 << " element m-1 is " << v.size() / 2 - 1 << " = " << v[v.size() / 2 - 1]
1198                 << std::endl;
1199     }
1200   }
1201   else
1202   {
1203     // odd number of entries
1204     auto m = v.begin() + v.size() / 2;
1205     std::nth_element(v.begin(), m, v.end());
1206     median = v[v.size() / 2];
1207     if (Verbosity() > 2)
1208     {
1209       std::cout << "The vector size is " << v.size() << " element m is " << v.size() / 2 << " = " << v[v.size() / 2] << std::endl;
1210     }
1211   }
1212 
1213   return median;
1214 }
1215 double PHSimpleVertexFinder::getAverage(std::vector<double> &v)
1216 {
1217   double avge = 0.0;
1218   double wt = 0.0;
1219   for (auto it : v)
1220   {
1221     avge += it;
1222     wt++;
1223   }
1224 
1225   avge /= wt;
1226   if (Verbosity() > 2)
1227   {
1228     std::cout << " average = " << avge << std::endl;
1229   }
1230 
1231   return avge;
1232 }