Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 /*! \class   TTStubAlgorithm_official
0002  *  \brief   Class for "official" algorithm to be used
0003  *           in TTStubBuilder
0004  *  \details HW-friendly algorithm: layer-wise LUT.
0005  *           After moving from SimDataFormats to DataFormats,
0006  *           the template structure of the class was maintained
0007  *           in order to accomodate any types other than PixelDigis
0008  *           in case there is such a need in the future.
0009  *
0010  *  \author Nicola Pozzobon
0011  *  \author Sebastien Viret
0012  *  \date   2013, Jul 18
0013  *
0014  */
0015 
0016 #ifndef L1_TRACK_TRIGGER_STUB_ALGO_official_H
0017 #define L1_TRACK_TRIGGER_STUB_ALGO_official_H
0018 
0019 #include "FWCore/Framework/interface/EventSetup.h"
0020 #include "FWCore/Framework/interface/ModuleFactory.h"
0021 #include "FWCore/Framework/interface/ESProducer.h"
0022 
0023 #include "L1Trigger/TrackTrigger/interface/TTStubAlgorithm.h"
0024 #include "L1Trigger/TrackTrigger/interface/TTStubAlgorithmRecord.h"
0025 
0026 #include "DataFormats/GeometryCommonDetAlgo/interface/MeasurementPoint.h"
0027 #include "Geometry/CommonTopologies/interface/Topology.h"
0028 #include "Geometry/CommonDetUnit/interface/PixelGeomDetUnit.h"
0029 #include "Geometry/CommonTopologies/interface/PixelTopology.h"
0030 
0031 #include <memory>
0032 #include <string>
0033 #include <map>
0034 #include <typeinfo>
0035 
0036 template <typename T>
0037 class TTStubAlgorithm_official : public TTStubAlgorithm<T> {
0038 private:
0039   /// Data members
0040   bool mPerformZMatchingPS;
0041   bool mPerformZMatching2S;
0042   bool m_tilted;
0043   std::string className_;
0044 
0045   std::vector<double> barrelCut;
0046   std::vector<std::vector<double>> ringCut;
0047   std::vector<std::vector<double>> tiltedCut;
0048   std::vector<double> barrelNTilt;
0049 
0050 public:
0051   /// Constructor
0052   TTStubAlgorithm_official(const TrackerGeometry *const theTrackerGeom,
0053                            const TrackerTopology *const theTrackerTopo,
0054                            std::vector<double> setBarrelCut,
0055                            std::vector<std::vector<double>> setRingCut,
0056                            std::vector<std::vector<double>> setTiltedCut,
0057                            std::vector<double> setBarrelNTilt,
0058                            bool aPerformZMatchingPS,
0059                            bool aPerformZMatching2S)
0060       : TTStubAlgorithm<T>(theTrackerGeom, theTrackerTopo, __func__) {
0061     barrelCut = setBarrelCut;
0062     ringCut = setRingCut;
0063     tiltedCut = setTiltedCut;
0064     barrelNTilt = setBarrelNTilt;
0065     mPerformZMatchingPS = aPerformZMatchingPS;
0066     mPerformZMatching2S = aPerformZMatching2S;
0067   }
0068 
0069   /// Destructor
0070   ~TTStubAlgorithm_official() override {}
0071 
0072   /// Matching operations
0073   void PatternHitCorrelation(bool &aConfirmation,
0074                              int &aDisplacement,
0075                              int &anOffset,
0076                              float &anHardBend,
0077                              const TTStub<T> &aTTStub) const override;
0078 
0079   float degradeBend(bool psModule, int window, int bend) const;
0080 
0081 };  /// Close class
0082 
0083 /*! \brief   Implementation of methods
0084  *  \details Here, in the header file, the methods which do not depend
0085  *           on the specific type <T> that can fit the template.
0086  *           Other methods, with type-specific features, are implemented
0087  *           in the source file.
0088  */
0089 
0090 /// Matching operations
0091 template <>
0092 void TTStubAlgorithm_official<Ref_Phase2TrackerDigi_>::PatternHitCorrelation(
0093     bool &aConfirmation,
0094     int &aDisplacement,
0095     int &anOffset,
0096     float &anHardBend,
0097     const TTStub<Ref_Phase2TrackerDigi_> &aTTStub) const;
0098 
0099 template <>
0100 float TTStubAlgorithm_official<Ref_Phase2TrackerDigi_>::degradeBend(bool psModule, int window, int bend) const;
0101 
0102 /*! \class   ES_TTStubAlgorithm_official
0103  *  \brief   Class to declare the algorithm to the framework
0104  *
0105  *  \author Nicola Pozzobon
0106  *  \date   2013, Jul 18
0107  *
0108  */
0109 
0110 template <typename T>
0111 class ES_TTStubAlgorithm_official : public edm::ESProducer {
0112 private:
0113   /// Data members
0114   edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> mGeomToken;
0115   edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> mTopoToken;
0116 
0117   /// Windows
0118   std::vector<double> setBarrelCut;
0119   std::vector<std::vector<double>> setRingCut;
0120   std::vector<std::vector<double>> setTiltedCut;
0121 
0122   std::vector<double> setBarrelNTilt;
0123 
0124   /// Z-matching
0125   bool mPerformZMatchingPS;
0126   bool mPerformZMatching2S;
0127 
0128 public:
0129   /// Constructor
0130   ES_TTStubAlgorithm_official(const edm::ParameterSet &p) {
0131     mPerformZMatchingPS = p.getParameter<bool>("zMatchingPS");
0132     mPerformZMatching2S = p.getParameter<bool>("zMatching2S");
0133     setBarrelCut = p.getParameter<std::vector<double>>("BarrelCut");
0134     setBarrelNTilt = p.getParameter<std::vector<double>>("NTiltedRings");
0135 
0136     std::vector<edm::ParameterSet> vPSet = p.getParameter<std::vector<edm::ParameterSet>>("EndcapCutSet");
0137     std::vector<edm::ParameterSet> vPSet2 = p.getParameter<std::vector<edm::ParameterSet>>("TiltedBarrelCutSet");
0138 
0139     std::vector<edm::ParameterSet>::const_iterator iPSet;
0140     for (iPSet = vPSet.begin(); iPSet != vPSet.end(); iPSet++) {
0141       setRingCut.push_back(iPSet->getParameter<std::vector<double>>("EndcapCut"));
0142     }
0143 
0144     for (iPSet = vPSet2.begin(); iPSet != vPSet2.end(); iPSet++) {
0145       setTiltedCut.push_back(iPSet->getParameter<std::vector<double>>("TiltedCut"));
0146     }
0147 
0148     auto cc = setWhatProduced(this);
0149     mGeomToken = cc.consumes();
0150     mTopoToken = cc.consumes();
0151   }
0152 
0153   /// Destructor
0154   ~ES_TTStubAlgorithm_official() override {}
0155 
0156   /// Implement the producer
0157   std::unique_ptr<TTStubAlgorithm<T>> produce(const TTStubAlgorithmRecord &record) {
0158     return std::make_unique<TTStubAlgorithm_official<T>>(&record.get(mGeomToken),
0159                                                          &record.get(mTopoToken),
0160                                                          setBarrelCut,
0161                                                          setRingCut,
0162                                                          setTiltedCut,
0163                                                          setBarrelNTilt,
0164                                                          mPerformZMatchingPS,
0165                                                          mPerformZMatching2S);
0166   }
0167 };
0168 
0169 #endif