Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:11:22

0001 #ifndef FASTSIMULATION_TRACKING_SEEDFINDER_H
0002 #define FASTSIMULATION_TRACKING_SEEDFINDER_H
0003 
0004 // system
0005 #include <vector>
0006 #include <functional>
0007 #include <array>
0008 
0009 // fastsim tracking
0010 #include "FastSimulation/Tracking/interface/SeedingTree.h"
0011 #include "FastSimulation/Tracking/interface/TrajectorySeedHitCandidate.h"
0012 #include "FastSimulation/Tracking/interface/SeedFinderSelector.h"
0013 
0014 class TrackerTopology;
0015 class FastTrackerRecHit;
0016 
0017 class SeedFinder {
0018 public:
0019 private:
0020   std::vector<std::vector<SeedFinderSelector*> > _selectorFunctionsByHits;
0021   const SeedingTree<TrackingLayer>& _seedingTree;
0022   const TrackerTopology* _trackerTopology;
0023 
0024 public:
0025   SeedFinder(const SeedingTree<TrackingLayer>& seedingTree, const TrackerTopology& trackerTopology)
0026       : _seedingTree(seedingTree), _trackerTopology(&trackerTopology) {}
0027 
0028   void addHitSelector(SeedFinderSelector* seedFinderSelector, unsigned int nHits) {
0029     if (_selectorFunctionsByHits.size() < nHits) {
0030       _selectorFunctionsByHits.resize(nHits);
0031       _selectorFunctionsByHits.reserve(nHits);
0032     }
0033     //shift indices by -1 so that _selectorFunctionsByHits[0] tests 1 hit
0034     _selectorFunctionsByHits[nHits - 1].push_back(seedFinderSelector);
0035   }
0036 
0037   std::vector<unsigned int> getSeed(const std::vector<const FastTrackerRecHit*>& trackerRecHits) const {
0038     std::vector<int> hitIndicesInTree(_seedingTree.numberOfNodes(), -1);
0039     //A SeedingNode is associated by its index to this list. The list stores the indices of the hits in 'trackerRecHits'
0040     /* example
0041            SeedingNode                     | hit index                 | hit
0042            -------------------------------------------------------------------------------
0043            index=  0:  [BPix1]             | hitIndicesInTree[0] (=1)  | trackerRecHits[1]
0044            index=  1:   -- [BPix2]         | hitIndicesInTree[1] (=3)  | trackerRecHits[3]
0045            index=  2:   --  -- [BPix3]     | hitIndicesInTree[2] (=4)  | trackerRecHits[4]
0046            index=  3:   --  -- [FPix1_pos] | hitIndicesInTree[3] (=6)  | trackerRecHits[6]
0047            index=  4:   --  -- [FPix1_neg] | hitIndicesInTree[4] (=7)  | trackerRecHits[7]
0048 
0049            The implementation has been chosen such that the tree only needs to be build once upon construction.
0050         */
0051     std::vector<TrajectorySeedHitCandidate> seedHitCandidates;
0052     for (const FastTrackerRecHit* trackerRecHit : trackerRecHits) {
0053       TrajectorySeedHitCandidate seedHitCandidate(trackerRecHit, _trackerTopology);
0054       seedHitCandidates.push_back(seedHitCandidate);
0055     }
0056     return iterateHits(0, seedHitCandidates, hitIndicesInTree, true);
0057 
0058     //TODO: create pairs of TrackingLayer -> remove TrajectorySeedHitCandidate class
0059   }
0060 
0061   //this method attempts to insert the hit at position 'trackerRecHit' in 'trackerRecHits'
0062   //into the seeding tree (stored as 'hitIndicesInTree')
0063   const SeedingNode<TrackingLayer>* insertHit(const std::vector<TrajectorySeedHitCandidate>& trackerRecHits,
0064                                               std::vector<int>& hitIndicesInTree,
0065                                               const SeedingNode<TrackingLayer>* node,
0066                                               unsigned int trackerHit) const {
0067     if (!node->getParent() || hitIndicesInTree[node->getParent()->getIndex()] >= 0) {
0068       if (hitIndicesInTree[node->getIndex()] < 0) {
0069         const TrajectorySeedHitCandidate& currentTrackerHit = trackerRecHits[trackerHit];
0070         if (currentTrackerHit.getTrackingLayer() != node->getData()) {
0071           return nullptr;
0072         }
0073 
0074         const unsigned int NHits = node->getDepth() + 1;
0075         if (_selectorFunctionsByHits.size() >= NHits) {
0076           //are there any selector functions stored for NHits?
0077           if (!_selectorFunctionsByHits[NHits - 1].empty()) {
0078             //fill vector of Hits from node to root to be passed to the selector function
0079             std::vector<const FastTrackerRecHit*> seedCandidateHitList(node->getDepth() + 1);
0080             seedCandidateHitList[node->getDepth()] = currentTrackerHit.hit();
0081             const SeedingNode<TrackingLayer>* parentNode = node->getParent();
0082             while (parentNode != nullptr) {
0083               seedCandidateHitList[parentNode->getDepth()] =
0084                   trackerRecHits[hitIndicesInTree[parentNode->getIndex()]].hit();
0085               parentNode = parentNode->getParent();
0086             }
0087 
0088             //loop over selector functions
0089             for (SeedFinderSelector* selectorFunction : _selectorFunctionsByHits[NHits - 1]) {
0090               if (!selectorFunction->pass(seedCandidateHitList)) {
0091                 return nullptr;
0092               }
0093             }
0094           }
0095         }
0096 
0097         //the hit was not rejected by all selector functions -> insert it into the tree
0098         hitIndicesInTree[node->getIndex()] = trackerHit;
0099         if (node->getChildrenSize() == 0) {
0100           return node;
0101         }
0102 
0103         return nullptr;
0104       } else {
0105         for (unsigned int ichild = 0; ichild < node->getChildrenSize(); ++ichild) {
0106           const SeedingNode<TrackingLayer>* seed =
0107               insertHit(trackerRecHits, hitIndicesInTree, node->getChild(ichild), trackerHit);
0108           if (seed) {
0109             return seed;
0110           }
0111         }
0112       }
0113     }
0114     return nullptr;
0115   }
0116 
0117   std::vector<unsigned int> iterateHits(unsigned int start,
0118                                         const std::vector<TrajectorySeedHitCandidate>& trackerRecHits,
0119                                         std::vector<int> hitIndicesInTree,
0120                                         bool processSkippedHits) const {
0121     for (unsigned int irecHit = start; irecHit < trackerRecHits.size(); ++irecHit) {
0122       // only accept hits that are on one of the requested layers
0123       if (_seedingTree.getSingleSet().find(trackerRecHits[irecHit].getTrackingLayer()) ==
0124           _seedingTree.getSingleSet().end()) {
0125         continue;
0126       }
0127 
0128       unsigned int currentHitIndex = irecHit;
0129 
0130       for (unsigned int inext = currentHitIndex + 1; inext < trackerRecHits.size(); ++inext) {
0131         //if multiple hits are on the same layer -> follow all possibilities by recusion
0132         if (trackerRecHits[currentHitIndex].getTrackingLayer() == trackerRecHits[inext].getTrackingLayer()) {
0133           if (processSkippedHits) {
0134             //recusively call the method again with hit 'inext' but skip all following on the same layer though 'processSkippedHits=false'
0135             std::vector<unsigned int> seedHits = iterateHits(inext, trackerRecHits, hitIndicesInTree, false);
0136             if (!seedHits.empty()) {
0137               return seedHits;
0138             }
0139           }
0140           irecHit += 1;
0141         } else {
0142           break;
0143         }
0144       }
0145 
0146       //processSkippedHits=true
0147 
0148       const SeedingNode<TrackingLayer>* seedNode = nullptr;
0149       for (unsigned int iroot = 0; seedNode == nullptr && iroot < _seedingTree.numberOfRoots(); ++iroot) {
0150         seedNode = insertHit(trackerRecHits, hitIndicesInTree, _seedingTree.getRoot(iroot), currentHitIndex);
0151       }
0152       if (seedNode) {
0153         std::vector<unsigned int> seedIndices(seedNode->getDepth() + 1);
0154         while (seedNode) {
0155           seedIndices[seedNode->getDepth()] = hitIndicesInTree[seedNode->getIndex()];
0156           seedNode = seedNode->getParent();
0157         }
0158         return seedIndices;
0159       }
0160     }
0161 
0162     return std::vector<unsigned int>();
0163   }
0164 };
0165 
0166 #endif