CMS 3D CMS Logo

TTStubAlgorithm_official.h
Go to the documentation of this file.
1 
16 #ifndef L1_TRACK_TRIGGER_STUB_ALGO_official_H
17 #define L1_TRACK_TRIGGER_STUB_ALGO_official_H
18 
23 
26 
31 
32 #include <memory>
33 #include <string>
34 #include <map>
35 #include <typeinfo>
36 
37 template< typename T >
39 {
40  private:
44  bool m_tilted;
46 
47  std::vector< double > barrelCut;
48  std::vector< std::vector< double > > ringCut;
49  std::vector< std::vector< double > > tiltedCut;
50  std::vector< double > barrelNTilt;
51 
52  public:
54  TTStubAlgorithm_official( const TrackerGeometry* const theTrackerGeom, const TrackerTopology* const theTrackerTopo,
55  std::vector< double > setBarrelCut,
56  std::vector< std::vector< double > > setRingCut,
57  std::vector< std::vector< double > > setTiltedCut,
58  std::vector< double > setBarrelNTilt,
59  bool aPerformZMatchingPS, bool aPerformZMatching2S )
60  : TTStubAlgorithm< T >( theTrackerGeom, theTrackerTopo, __func__ )
61  {
62  barrelCut = setBarrelCut;
63  ringCut = setRingCut;
64  tiltedCut = setTiltedCut;
65  barrelNTilt = setBarrelNTilt;
66  mPerformZMatchingPS = aPerformZMatchingPS;
67  mPerformZMatching2S = aPerformZMatching2S;
68  }
69 
72 
74  void PatternHitCorrelation( bool &aConfirmation,
75  int &aDisplacement,
76  int &anOffset,
77  float &anROffset,
78  float &anHardBend,
79  const TTStub< T > &aTTStub ) const override;
80 
81  float degradeBend(bool psModule, int window, int bend) const;
82 
83 };
84 
92 template< >
95  int &aDisplacement,
96  int &anOffset,
97  float &anROffset,
98  float &anHardBend,
99  const TTStub< Ref_Phase2TrackerDigi_ > &aTTStub ) const;
100 
101 template< >
102 float TTStubAlgorithm_official< Ref_Phase2TrackerDigi_ >::degradeBend(bool psModule, int window, int bend) const;
103 
104 
113 template< typename T >
115 {
116  private:
118 
120  std::vector< double > setBarrelCut;
121  std::vector< std::vector< double > > setRingCut;
122  std::vector< std::vector< double > > setTiltedCut;
123 
124  std::vector< double > setBarrelNTilt;
125 
129 
130  public:
133  {
134  mPerformZMatchingPS = p.getParameter< bool >("zMatchingPS");
135  mPerformZMatching2S = p.getParameter< bool >("zMatching2S");
136  setBarrelCut = p.getParameter< std::vector< double > >("BarrelCut");
137  setBarrelNTilt = p.getParameter< std::vector< double > >("NTiltedRings");
138 
139  std::vector< edm::ParameterSet > vPSet = p.getParameter< std::vector< edm::ParameterSet > >("EndcapCutSet");
140  std::vector< edm::ParameterSet > vPSet2 = p.getParameter< std::vector< edm::ParameterSet > >("TiltedBarrelCutSet");
141 
142  std::vector< edm::ParameterSet >::const_iterator iPSet;
143  for ( iPSet = vPSet.begin(); iPSet != vPSet.end(); iPSet++ )
144  {
145  setRingCut.push_back( iPSet->getParameter< std::vector< double > >("EndcapCut") );
146  }
147 
148  for ( iPSet = vPSet2.begin(); iPSet != vPSet2.end(); iPSet++ )
149  {
150  setTiltedCut.push_back( iPSet->getParameter< std::vector< double > >("TiltedCut") );
151  }
152 
153 
154  setWhatProduced( this );
155  }
156 
159 
161  std::unique_ptr< TTStubAlgorithm< T > > produce( const TTStubAlgorithmRecord & record )
162  {
164  record.getRecord< TrackerDigiGeometryRecord >().get( tGeomHandle );
165  const TrackerGeometry* const theTrackerGeom = tGeomHandle.product();
166  edm::ESHandle<TrackerTopology> tTopoHandle;
167  record.getRecord< TrackerTopologyRcd >().get(tTopoHandle);
168  const TrackerTopology* const theTrackerTopo = tTopoHandle.product();
169 
170  TTStubAlgorithm< T >* TTStubAlgo = new TTStubAlgorithm_official< T >( theTrackerGeom, theTrackerTopo,
171  setBarrelCut, setRingCut, setTiltedCut, setBarrelNTilt,
173 
174  return std::unique_ptr< TTStubAlgorithm< T > >( TTStubAlgo );
175  }
176 
177 };
178 
179 #endif
180 
181 
Class for "official" algorithm to be used in TTStubBuilder.
void PatternHitCorrelation(bool &aConfirmation, int &aDisplacement, int &anOffset, float &anROffset, float &anHardBend, const TTStub< T > &aTTStub) const override
Matching operations.
T getParameter(std::string const &) const
Base class for any algorithm to be used in TTStubBuilder.
JetCorrectorParameters::Record record
Definition: classes.h:7
Class to declare the algorithm to the framework.
float degradeBend(bool psModule, int window, int bend) const
~TTStubAlgorithm_official() override
Destructor.
~ES_TTStubAlgorithm_official() override
Destructor.
ES_TTStubAlgorithm_official(const edm::ParameterSet &p)
Constructor.
std::vector< std::vector< double > > setRingCut
def window(xmin, xmax, ymin, ymax, x=0, y=0, width=100, height=100, xlogbase=None, ylogbase=None, minusInfinity=-1000, flipx=False, flipy=True)
Definition: svgfig.py:642
TTStubAlgorithm_official(const TrackerGeometry *const theTrackerGeom, const TrackerTopology *const theTrackerTopo, std::vector< double > setBarrelCut, std::vector< std::vector< double > > setRingCut, std::vector< std::vector< double > > setTiltedCut, std::vector< double > setBarrelNTilt, bool aPerformZMatchingPS, bool aPerformZMatching2S)
Constructor.
std::vector< double > barrelNTilt
Class to store the TTStubAlgorithm used in TTStubBuilder.
Class to store the L1 Track Trigger stubs.
Definition: TTStub.h:22
bool mPerformZMatchingPS
Data members.
std::vector< std::vector< double > > ringCut
std::vector< std::vector< double > > setTiltedCut
std::vector< std::vector< double > > tiltedCut
std::unique_ptr< TTStubAlgorithm< T > > produce(const TTStubAlgorithmRecord &record)
Implement the producer.
std::vector< double > barrelCut
long double T
T const * product() const
Definition: ESHandle.h:86
std::vector< double > setBarrelCut
Data members.