Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:26:08

0001 #ifndef DTSegment_DTCombinatorialExtendedPatternReco_h
0002 #define DTSegment_DTCombinatorialExtendedPatternReco_h
0003 
0004 /** \class DTCombinatorialExtendedPatternReco
0005  *
0006  * Algo for reconstructing 2d segment in DT using a combinatorial approach
0007  *  
0008  * \author Stefano Lacaprara - INFN Legnaro <stefano.lacaprara@pd.infn.it>
0009  * \author Riccardo Bellan - INFN TO <riccardo.bellan@cern.ch>
0010  *
0011  */
0012 
0013 /* Base Class Headers */
0014 #include "RecoLocalMuon/DTSegment/src/DTRecSegment2DBaseAlgo.h"
0015 
0016 /* Collaborating Class Declarations */
0017 namespace edm {
0018   class ParameterSet;
0019   class EventSetup;
0020   //  class ESHandle;
0021 }  // namespace edm
0022 class DTSegmentUpdator;
0023 class DTSegmentCleaner;
0024 class DTHitPairForFit;
0025 class DTSegmentCand;
0026 class DTSegmentExtendedCand;
0027 class MuonGeometryRecord;
0028 
0029 #include "DataFormats/DTRecHit/interface/DTSLRecCluster.h"
0030 #include "FWCore/Framework/interface/FrameworkfwdMostUsed.h"
0031 
0032 /* C++ Headers */
0033 #include <vector>
0034 #include <deque>
0035 #include <utility>
0036 
0037 #include "Geometry/DTGeometry/interface/DTGeometry.h"
0038 #include "FWCore/Framework/interface/ESHandle.h"
0039 #include "FWCore/Utilities/interface/ESGetToken.h"
0040 #include "RecoLocalMuon/DTSegment/src/DTSegmentCand.h"
0041 
0042 /* ====================================================================== */
0043 
0044 /* Class DTCombinatorialExtendedPatternReco Interface */
0045 
0046 class DTCombinatorialExtendedPatternReco : private DTRecSegment2DBaseAlgo {
0047 public:
0048   /// Constructor
0049   DTCombinatorialExtendedPatternReco(const edm::ParameterSet& pset, edm::ConsumesCollector);
0050 
0051   /// Destructor
0052   ~DTCombinatorialExtendedPatternReco() override;
0053 
0054   /* Operations */
0055 
0056   /// this function is called in the producer
0057   edm::OwnVector<DTSLRecSegment2D> reconstruct(const DTSuperLayer* sl,
0058                                                const std::vector<DTRecHit1DPair>& hits) override;
0059 
0060   /// return the algo name
0061   std::string algoName() const override { return theAlgoName; }
0062 
0063   /// Through this function the EventSetup is percolated to the
0064   /// objs which request it
0065   void setES(const edm::EventSetup& setup) override;
0066 
0067   // pass clusters to algo
0068   void setClusters(const std::vector<DTSLRecCluster>& clusters);
0069 
0070 protected:
0071 private:
0072   // create the DTHitPairForFit from the pairs for easy use
0073   std::vector<std::shared_ptr<DTHitPairForFit>> initHits(const DTSuperLayer* sl,
0074                                                          const std::vector<DTRecHit1DPair>& hits);
0075 
0076   // search for candidate, starting from pairs of hits in different layers
0077   std::vector<DTSegmentCand*> buildSegments(const DTSuperLayer* sl,
0078                                             const std::vector<std::shared_ptr<DTHitPairForFit>>& hits);
0079 
0080   // find all the hits compatible with the candidate
0081   std::vector<DTSegmentCand::AssPoint> findCompatibleHits(const LocalPoint& pos,
0082                                                           const LocalVector& dir,
0083                                                           const std::vector<std::shared_ptr<DTHitPairForFit>>& hits);
0084 
0085   // build segments from hits collection
0086   DTSegmentExtendedCand* buildBestSegment(std::vector<DTSegmentCand::AssPoint>& assHits, const DTSuperLayer* sl);
0087 
0088   bool checkDoubleCandidates(std::vector<DTSegmentCand*>& segs, DTSegmentCand* seg);
0089 
0090   /** build collection of compatible hits for L/R hits: the candidates is
0091      * updated with the segment candidates found */
0092   void buildPointsCollection(std::vector<DTSegmentCand::AssPoint>& points,
0093                              std::deque<std::shared_ptr<DTHitPairForFit>>& pointsNoLR,
0094                              std::vector<DTSegmentCand*>& candidates,
0095                              const DTSuperLayer* sl);
0096 
0097   /** extend the candidates with clusters from external SL */
0098   std::vector<DTSegmentExtendedCand*> extendCandidates(std::vector<DTSegmentCand*>& candidates, const DTSuperLayer* sl);
0099 
0100   bool closeSL(const DTSuperLayerId& id1, const DTSuperLayerId& id2);
0101 
0102 private:
0103   std::string theAlgoName;
0104   unsigned int theMaxAllowedHits;
0105   double theAlphaMaxTheta;
0106   double theAlphaMaxPhi;
0107   bool debug;
0108   bool usePairs;
0109   DTSegmentUpdator* theUpdator;  // the updator and fitter
0110   DTSegmentCleaner* theCleaner;  // the cleaner
0111 
0112   edm::ESHandle<DTGeometry> theDTGeometry;  // the DT geometry
0113   edm::ESGetToken<DTGeometry, MuonGeometryRecord> theDTGeometryToken;
0114 
0115 private:
0116   std::vector<std::vector<int>> theTriedPattern;
0117   std::vector<DTSLRecCluster> theClusters;
0118 };
0119 #endif  // DTSegment_DTCombinatorialExtendedPatternReco_h