Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #ifndef DTSegment_DTMeantimerPatternReco4D_h
0002 #define DTSegment_DTMeantimerPatternReco4D_h
0003 
0004 /** \class DTMeantimerPatternReco4D
0005  *
0006  * Algo for reconstructing 4d segment in DT using a Meantimer 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/DTRecSegment4DBaseAlgo.h"
0015 
0016 class DTRecSegment2DBaseAlgo;
0017 
0018 // Collaborating Class Declarations
0019 namespace edm {
0020   class ParameterSet;
0021   class EventSetup;
0022 }  // namespace edm
0023 class DTSegmentUpdator;
0024 class MuonGeometryRecord;
0025 //class DTSegmentCleaner;
0026 
0027 // C++ Headers
0028 #include <vector>
0029 //#include <utility>
0030 
0031 #include "Geometry/DTGeometry/interface/DTGeometry.h"
0032 #include "FWCore/Framework/interface/ESHandle.h"
0033 #include "FWCore/Utilities/interface/ESGetToken.h"
0034 #include "FWCore/Framework/interface/FrameworkfwdMostUsed.h"
0035 
0036 // ======================================================================
0037 class DTSegmentCand;
0038 class DTMeantimerPatternReco;
0039 class DTHitPairForFit;
0040 
0041 // Class DTMeantimerPatternReco4D Interface
0042 
0043 class DTMeantimerPatternReco4D : public DTRecSegment4DBaseAlgo {
0044 public:
0045   /// Constructor
0046   DTMeantimerPatternReco4D(const edm::ParameterSet& pset, edm::ConsumesCollector);
0047 
0048   /// Destructor
0049   ~DTMeantimerPatternReco4D() override;
0050 
0051   /// Operations
0052   edm::OwnVector<DTRecSegment4D> reconstruct() override;
0053 
0054   std::string algoName() const override { return theAlgoName; }
0055 
0056   void setES(const edm::EventSetup& setup) override;
0057   void setDTRecHit1DContainer(edm::Handle<DTRecHitCollection> all1DHits) override;
0058   void setDTRecSegment2DContainer(edm::Handle<DTRecSegment2DCollection> all2DSegments) override;
0059   void setChamber(const DTChamberId& chId) override;
0060   bool wants2DSegments() override { return !allDTRecHits; }
0061 
0062 protected:
0063 private:
0064   std::vector<DTSegmentCand*> buildPhiSuperSegmentsCandidates(
0065       std::vector<std::shared_ptr<DTHitPairForFit>>& pairPhiOwned);
0066   DTRecSegment4D* segmentSpecialZed(DTRecSegment4D* seg);
0067 
0068   std::string theAlgoName;
0069 
0070   bool debug;
0071   // DTSegmentUpdator* theUpdator; // the updator and fitter
0072   // DTSegmentCleaner* theCleaner; // the cleaner
0073 
0074   edm::ESHandle<DTGeometry> theDTGeometry;  // the DT geometry
0075   const edm::ESGetToken<DTGeometry, MuonGeometryRecord> theDTGeometryToken;
0076 
0077   // The reconstruction 2D algorithm
0078   // For the 2D reco I use thei reconstructor!
0079   DTMeantimerPatternReco* the2DAlgo;
0080 
0081   // the updator
0082   DTSegmentUpdator* theUpdator;
0083 
0084   const DTChamber* theChamber;
0085 
0086   //the input type
0087   bool allDTRecHits;
0088   bool applyT0corr;
0089   bool computeT0corr;
0090 
0091   //  std::vector<DTRecHit1D> the1DPhiHits;
0092   std::vector<DTSLRecSegment2D> theSegments2DTheta;
0093   std::vector<DTRecHit1DPair> theHitsFromPhi1;
0094   std::vector<DTRecHit1DPair> theHitsFromTheta;
0095   std::vector<DTRecHit1DPair> theHitsFromPhi2;
0096 };
0097 #endif