CMS 3D CMS Logo

PixelTripletLowPtGenerator.cc
Go to the documentation of this file.
5 
8 
10 
13 
19 
20 #undef Debug
21 
22 using namespace std;
23 
24 /*****************************************************************************/
26  : HitTripletGeneratorFromPairAndLayers(), // no theMaxElement used in this class
28  m_topoToken(iC.esConsumes<TrackerTopology, TrackerTopologyRcd>()),
29  m_magfieldToken(iC.esConsumes()),
30  m_ttrhBuilderToken(iC.esConsumes(edm::ESInputTag("", cfg.getParameter<string>("TTRHBuilder")))),
31  m_msmakerToken(iC.esConsumes()),
32  m_clusterFilterToken(iC.esConsumes(edm::ESInputTag("", "ClusterShapeHitFilter"))),
33  theTracker(nullptr),
34  theClusterShapeCacheToken(
35  iC.consumes<SiPixelClusterShapeCache>(cfg.getParameter<edm::InputTag>("clusterShapeCacheSrc"))) {
36  checkMultipleScattering = cfg.getParameter<bool>("checkMultipleScattering");
37  nSigMultipleScattering = cfg.getParameter<double>("nSigMultipleScattering");
38  checkClusterShape = cfg.getParameter<bool>("checkClusterShape");
39  rzTolerance = cfg.getParameter<double>("rzTolerance");
40  maxAngleRatio = cfg.getParameter<double>("maxAngleRatio");
41 }
42 
43 /*****************************************************************************/
45 
46 /*****************************************************************************/
47 
48 /*****************************************************************************/
50  if (theTracker == nullptr) {
51  // Get tracker geometry
53  }
54 
55  if (!theFilter) {
56  theFilter = std::make_unique<TripletFilter>(&es.getData(m_clusterFilterToken));
57  }
58 }
59 
60 /*****************************************************************************/
62  DetId detId = recHit->geographicalId();
63 
64  return theTracker->idToDet(detId)->toGlobal(recHit->localPosition());
65 }
66 
67 /*****************************************************************************/
70  const edm::Event& ev,
71  const edm::EventSetup& es,
72  const SeedingLayerSetsHits::SeedingLayerSet& pairLayers,
73  const std::vector<SeedingLayerSetsHits::SeedingLayer>& thirdLayers) {
74  //Retrieve tracker topology from geometry
75  const TrackerTopology* tTopo = &es.getData(m_topoToken);
76  const auto& magfield = es.getData(m_magfieldToken);
77  const auto& ttrhBuilder = es.getData(m_ttrhBuilderToken);
78  const auto& msmaker = es.getData(m_msmakerToken);
79 
82 
83  // Generate pairs
84  OrderedHitPairs pairs;
85  pairs.reserve(30000);
86  thePairGenerator->hitPairs(region, pairs, ev, es, pairLayers);
87 
88  if (pairs.empty())
89  return;
90 
91  int size = thirdLayers.size();
92 
93  // Set aliases
94  const RecHitsSortedInPhi** thirdHitMap = new const RecHitsSortedInPhi*[size];
95  for (int il = 0; il < size; il++)
96  thirdHitMap[il] = &(*theLayerCache)(thirdLayers[il], region);
97 
98  // Get tracker
99  getTracker(es);
100 
101  // Look at all generated pairs
102  for (OrderedHitPairs::const_iterator ip = pairs.begin(); ip != pairs.end(); ip++) {
103  // Fill rechits and points
104  vector<const TrackingRecHit*> recHits(3);
105  vector<GlobalPoint> points(3);
106 
107  recHits[0] = (*ip).inner()->hit();
108  recHits[1] = (*ip).outer()->hit();
109 
110 #ifdef Debug
111  cerr << " RecHits " + HitInfo::getInfo(*recHits[0]) + HitInfo::getInfo(*recHits[1]) << endl;
112 #endif
113 
114  for (int i = 0; i < 2; i++)
115  points[i] = getGlobalPosition(recHits[i]);
116 
117  // Initialize helix prediction
118  ThirdHitPrediction thePrediction(
119  region, points[0], points[1], magfield, ttrhBuilder, nSigMultipleScattering, maxAngleRatio);
120 
121  // Look at all layers
122  for (int il = 0; il < size; il++) {
123  const DetLayer* layer = thirdLayers[il].detLayer();
124 
125 #ifdef Debug
126  cerr << " check layer " << layer->subDetector() << " " << layer->location() << endl;
127 #endif
128 
129  // Get ranges for the third hit
130  float phi[2], rz[2];
131  thePrediction.getRanges(layer, phi, rz);
132 
133  PixelRecoRange<float> phiRange(phi[0], phi[1]);
134  PixelRecoRange<float> rzRange(rz[0] - rzTolerance, rz[1] + rzTolerance);
135 
136  // Get third hit candidates from cache
138  vector<Hit> thirdHits = thirdHitMap[il]->hits(phiRange.min(), phiRange.max());
139  typedef vector<Hit>::const_iterator IH;
140 
141  for (IH th = thirdHits.begin(), eh = thirdHits.end(); th < eh; ++th) {
142  // Fill rechit and point
143  recHits[2] = (*th)->hit();
144  points[2] = getGlobalPosition(recHits[2]);
145 
146 #ifdef Debug
147  cerr << " third hit " + HitInfo::getInfo(*recHits[2]) << endl;
148 #endif
149 
150  // Check if third hit is compatible with multiple scattering
151  vector<GlobalVector> globalDirs;
152  if (thePrediction.isCompatibleWithMultipleScattering(points[2], recHits, globalDirs, msmaker) == false) {
153 #ifdef Debug
154  cerr << " not compatible: multiple scattering" << endl;
155 #endif
157  continue;
158  }
159 
160  // Convert to localDirs
161  /*
162  vector<LocalVector> localDirs;
163  vector<GlobalVector>::const_iterator globalDir = globalDirs.begin();
164  for(vector<const TrackingRecHit *>::const_iterator
165  recHit = recHits.begin();
166  recHit != recHits.end(); recHit++)
167  {
168  localDirs.push_back(theTracker->idToDet(
169  (*recHit)->geographicalId())->toLocal(*globalDir));
170  globalDir++;
171  }
172 */
173 
174  // Check if the cluster shapes are compatible with thrusts
175  if (checkClusterShape) {
176  if (!theFilter->checkTrack(recHits, globalDirs, tTopo, *clusterShapeCache)) {
177 #ifdef Debug
178  cerr << " not compatible: cluster shape" << endl;
179 #endif
180  continue;
181  }
182  }
183 
184  // All checks passed, put triplet back
185  result.push_back(OrderedHitTriplet((*ip).inner(), (*ip).outer(), *th));
186  }
187  }
188  }
189  delete[] thirdHitMap;
190 
191  return;
192 }
195  const edm::EventSetup& es,
196  const HitDoublets& doublets,
197  const RecHitsSortedInPhi** thirdHitMap,
198  const std::vector<const DetLayer*>& thirdLayerDetLayer,
199  const int nThirdLayers) {
200  throw cms::Exception("Error") << "PixelTripletLowPtGenerator::hitTriplets is not implemented \n";
201 }
size
Write out results.
PixelTripletLowPtGenerator(const edm::ParameterSet &cfg, edm::ConsumesCollector &iC)
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
void getTracker(const edm::EventSetup &es)
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > m_geomToken
void hitTriplets(const TrackingRegion &region, OrderedHitTriplets &trs, const edm::Event &ev, const edm::EventSetup &es, const SeedingLayerSetsHits::SeedingLayerSet &pairLayers, const std::vector< SeedingLayerSetsHits::SeedingLayer > &thirdLayers) override
edm::ESGetToken< ClusterShapeHitFilter, CkfComponentsRecord > m_clusterFilterToken
std::unique_ptr< TripletFilter > theFilter
static std::string getInfo(const DetId &id, const TrackerTopology *tTopo)
Definition: HitInfo.cc:19
edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > m_magfieldToken
const TrackerGeomDet * idToDet(DetId) const override
std::vector< Hit > hits(float phiMin, float phiMax) const
edm::EDGetTokenT< SiPixelClusterShapeCache > theClusterShapeCacheToken
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:49
SeedingHitSet::ConstRecHitPointer Hit
BaseTrackerRecHit const * Hit
Definition: DetId.h:17
std::unique_ptr< HitPairGeneratorFromLayerPair > thePairGenerator
const TrackerGeometry * theTracker
HLT enums.
GlobalPoint getGlobalPosition(const TrackingRecHit *recHit)
void getRanges(const DetLayer *layer, float phi[], float rz[])
edm::ESGetToken< MultipleScatteringParametrisationMaker, TrackerMultipleScatteringRecord > m_msmakerToken
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > m_topoToken
edm::ESGetToken< TransientTrackingRecHitBuilder, TransientRecHitRecord > m_ttrhBuilderToken
bool isCompatibleWithMultipleScattering(GlobalPoint g3, const std::vector< const TrackingRecHit *> &h, std::vector< GlobalVector > &localDirs, const MultipleScatteringParametrisationMaker &msmaker)