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  const TTStub< T > &aTTStub ) const;
78 
79 };
80 
88 template< >
91  int &aDisplacement,
92  int &anOffset,
93  const TTStub< Ref_Phase2TrackerDigi_ > &aTTStub ) const;
94 
103 template< typename T >
105 {
106  private:
108  std::shared_ptr< TTStubAlgorithm< T > > _theAlgo;
109 
111  std::vector< double > setBarrelCut;
112  std::vector< std::vector< double > > setRingCut;
113  std::vector< std::vector< double > > setTiltedCut;
114 
115  std::vector< double > setBarrelNTilt;
116 
120 
121  public:
124  {
125  mPerformZMatchingPS = p.getParameter< bool >("zMatchingPS");
126  mPerformZMatching2S = p.getParameter< bool >("zMatching2S");
127  setBarrelCut = p.getParameter< std::vector< double > >("BarrelCut");
128  setBarrelNTilt = p.getParameter< std::vector< double > >("NTiltedRings");
129 
130  std::vector< edm::ParameterSet > vPSet = p.getParameter< std::vector< edm::ParameterSet > >("EndcapCutSet");
131  std::vector< edm::ParameterSet > vPSet2 = p.getParameter< std::vector< edm::ParameterSet > >("TiltedBarrelCutSet");
132 
133  std::vector< edm::ParameterSet >::const_iterator iPSet;
134  for ( iPSet = vPSet.begin(); iPSet != vPSet.end(); iPSet++ )
135  {
136  setRingCut.push_back( iPSet->getParameter< std::vector< double > >("EndcapCut") );
137  }
138 
139  for ( iPSet = vPSet2.begin(); iPSet != vPSet2.end(); iPSet++ )
140  {
141  setTiltedCut.push_back( iPSet->getParameter< std::vector< double > >("TiltedCut") );
142  }
143 
144 
145  setWhatProduced( this );
146  }
147 
150 
152  std::shared_ptr< TTStubAlgorithm< T > > produce( const TTStubAlgorithmRecord & record )
153  {
155  record.getRecord< TrackerDigiGeometryRecord >().get( tGeomHandle );
156  const TrackerGeometry* const theTrackerGeom = tGeomHandle.product();
157  edm::ESHandle<TrackerTopology> tTopoHandle;
158  record.getRecord< TrackerTopologyRcd >().get(tTopoHandle);
159  const TrackerTopology* const theTrackerTopo = tTopoHandle.product();
160 
161  TTStubAlgorithm< T >* TTStubAlgo = new TTStubAlgorithm_official< T >( theTrackerGeom, theTrackerTopo,
162  setBarrelCut, setRingCut, setTiltedCut, setBarrelNTilt,
164 
165  _theAlgo = std::shared_ptr< TTStubAlgorithm< T > >( TTStubAlgo );
166  return _theAlgo;
167  }
168 
169 };
170 
171 #endif
172 
173 
Class for "official" algorithm to be used in TTStubBuilder.
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.
std::shared_ptr< TTStubAlgorithm< T > > _theAlgo
Data members.
ES_TTStubAlgorithm_official(const edm::ParameterSet &p)
Constructor.
std::vector< std::vector< double > > setRingCut
void PatternHitCorrelation(bool &aConfirmation, int &aDisplacement, int &anOffset, const TTStub< T > &aTTStub) const
Matching operations.
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
std::shared_ptr< TTStubAlgorithm< T > > produce(const TTStubAlgorithmRecord &record)
Implement the producer.
bool mPerformZMatchingPS
Data members.
std::vector< std::vector< double > > ringCut
std::vector< std::vector< double > > setTiltedCut
std::vector< std::vector< double > > tiltedCut
virtual ~ES_TTStubAlgorithm_official()
Destructor.
std::vector< double > barrelCut
long double T
T const * product() const
Definition: ESHandle.h:86
std::vector< double > setBarrelCut
Windows.