CMS 3D CMS Logo

SiStripElectronSeedGenerator.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: EgammaElectronAlgos
4 // Class: SiStripElectronSeedGenerator.
5 //
11 //
12 //
13 
14 #include <vector>
15 #include <utility>
16 
19 
21 
22 
23 
25 
28 
33 
36 
38 
39 // files for retrieving hits using measurement tracker
40 
46 
47 
49 
51  : beamSpotTag_(tokens.token_bs),
52  theUpdator(nullptr),thePropagator(nullptr),theMeasurementTracker(nullptr),
53  theMeasurementTrackerEventTag(tokens.token_mte),
54  theSetup(nullptr), theMatcher_(nullptr),
55  cacheIDMagField_(0),cacheIDCkfComp_(0),cacheIDTrkGeom_(0),
56  tibOriginZCut_(pset.getParameter<double>("tibOriginZCut")),
57  tidOriginZCut_(pset.getParameter<double>("tidOriginZCut")),
58  tecOriginZCut_(pset.getParameter<double>("tecOriginZCut")),
59  monoOriginZCut_(pset.getParameter<double>("monoOriginZCut")),
60  tibDeltaPsiCut_(pset.getParameter<double>("tibDeltaPsiCut")),
61  tidDeltaPsiCut_(pset.getParameter<double>("tidDeltaPsiCut")),
62  tecDeltaPsiCut_(pset.getParameter<double>("tecDeltaPsiCut")),
63  monoDeltaPsiCut_(pset.getParameter<double>("monoDeltaPsiCut")),
64  tibPhiMissHit2Cut_(pset.getParameter<double>("tibPhiMissHit2Cut")),
65  tidPhiMissHit2Cut_(pset.getParameter<double>("tidPhiMissHit2Cut")),
66  tecPhiMissHit2Cut_(pset.getParameter<double>("tecPhiMissHit2Cut")),
67  monoPhiMissHit2Cut_(pset.getParameter<double>("monoPhiMissHit2Cut")),
68  tibZMissHit2Cut_(pset.getParameter<double>("tibZMissHit2Cut")),
69  tidRMissHit2Cut_(pset.getParameter<double>("tidRMissHit2Cut")),
70  tecRMissHit2Cut_(pset.getParameter<double>("tecRMissHit2Cut")),
71  tidEtaUsage_(pset.getParameter<double>("tidEtaUsage")),
72  tidMaxHits_(pset.getParameter<int>("tidMaxHits")),
73  tecMaxHits_(pset.getParameter<int>("tecMaxHits")),
74  monoMaxHits_(pset.getParameter<int>("monoMaxHits")),
75  maxSeeds_(pset.getParameter<int>("maxSeeds"))
76 {
77  // use of a theMeasurementTrackerName
78  if (pset.exists("measurementTrackerName"))
79  { theMeasurementTrackerName = pset.getParameter<std::string>("measurementTrackerName") ; }
80 
81 
82  // new beamSpot tag
83  /*
84  if (pset.exists("beamSpot"))
85  { beamSpotTag_ = pset.getParameter<edm::InputTag>("beamSpot") ; }
86  */
87 
88  theUpdator = new KFUpdator();
90 }
91 
92 
94  delete thePropagator;
95  delete theUpdator;
96 }
97 
98 
100 
103  cacheIDMagField_=setup.get<IdealMagneticFieldRecord>().cacheIdentifier();
104  if (thePropagator) delete thePropagator;
106  }
107 
110  cacheIDCkfComp_=setup.get<CkfComponentsRecord>().cacheIdentifier();
112  }
113 
116  cacheIDTrkGeom_=setup.get<TrackerDigiGeometryRecord>().cacheIdentifier();
117  }
118 
119 }
120 
124  theSetup= &setup;
125 
129 
130 
131  for (unsigned int i=0;i<clusters->size();++i) {
132  edm::Ref<reco::SuperClusterCollection> theClusB(clusters,i);
133  // Find the seeds
134  LogDebug ("run") << "new cluster, calling findSeedsFromCluster";
135  findSeedsFromCluster(theClusB,theBeamSpot,*data,out);
136  }
137 
138  LogDebug ("run") << ": For event "<<e.id();
139  LogDebug ("run") <<"Nr of superclusters: "<<clusters->size()
140  <<", no. of ElectronSeeds found = " << out.size();
141 }
142 
143 
144 // Find seeds using a supercluster
148  const MeasurementTrackerEvent & trackerData,
150  {
151  // clear the member vectors of good hits
152  layer1Hits_.clear() ;
153  layer2Hits_.clear() ;
154  backupLayer2Hits_.clear() ;
155 
156  using namespace std;
157 
158  double sCenergy = seedCluster->energy();
159  math::XYZPoint sCposition = seedCluster->position();
160  double scEta = seedCluster->eta();
161 
162  double scz = sCposition.z();
163  double scr = sqrt(pow(sCposition.x(),2)+pow(sCposition.y(),2));
164 
165  double pT = sCenergy * seedCluster->position().rho()/sqrt(seedCluster->x()*seedCluster->x()+seedCluster->y()*seedCluster->y()+seedCluster->z()*seedCluster->z());
166 
167  // FIXME use nominal field (see below)
168  double magneticField = 3.8;
169 
170  // cf Jackson p. 581-2, a little geometry
171  double phiVsRSlope = -3.00e-3 * magneticField / pT / 2.;
172 
173 
174  //Need to create TSOS to feed MeasurementTracker
175  GlobalPoint beamSpot(bs->x0(),bs->y0(),bs->z0());
176  GlobalPoint superCluster(sCposition.x(),sCposition.y(),sCposition.z());
177  double r0 = beamSpot.perp();
178  double z0 = beamSpot.z();
179 
180  //We need to pick a charge for the particle we want to reconstruct before hits can be retrieved
181  //Choosing both charges improves seeding efficiency by less than 0.5% for signal events
182  //If we pick a single charge, this reduces fake rate and CPU time
183  //So we pick a charge that is equally likely to be positive or negative
184 
185  int chargeHypothesis;
186  double chargeSelector = sCenergy - (int)sCenergy;
187  if(chargeSelector >= 0.5) chargeHypothesis = -1;
188  else {chargeHypothesis = 1;}
189 
190  //Use BeamSpot and SC position to estimate 3rd point
191  double rFake = 25.;
192  double phiFake = phiDiff(superCluster.phi(),chargeHypothesis * phiVsRSlope * (scr - rFake));
193  double zFake = (rFake*(scz-z0)-r0*scz+scr*z0)/(scr-r0);
194  double xFake = rFake * cos(phiFake);
195  double yFake = rFake * sin(phiFake);
196  GlobalPoint fakePoint(xFake,yFake,zFake);
197 
198  // FIXME optmize, move outside loop
200  theSetup->get<IdealMagneticFieldRecord>().get(bfield);
201  float nomField = bfield->nominalValue();
202 
203  //Use 3 points to make helix
204  FastHelix initialHelix(superCluster,fakePoint,beamSpot,nomField,&*bfield);
205 
206  //Use helix to get FTS
207  FreeTrajectoryState initialFTS(initialHelix.stateAtVertex());
208 
209  //Use FTS and BeamSpot to create TSOS
211  TrajectoryStateOnSurface initialTSOS = tipe->extrapolate(initialFTS,beamSpot);
212 
213  //Use GST to retrieve hits from various DetLayers using layerMeasurements class
214  const GeometricSearchTracker* gst = theMeasurementTracker->geometricSearchTracker();
215 
216  std::vector<const BarrelDetLayer*> tibLayers = gst->tibLayers();
217  const DetLayer* tib1 = tibLayers.at(0);
218  const DetLayer* tib2 = tibLayers.at(1);
219 
220  std::vector<const ForwardDetLayer*> tecLayers;
221  std::vector<const ForwardDetLayer*> tidLayers;
222  if(scEta < 0){
223  tecLayers = gst->negTecLayers();
224  tidLayers = gst->negTidLayers();
225  }
226  if(scEta > 0){
227  tecLayers = gst->posTecLayers();
228  tidLayers = gst->posTidLayers();
229  }
230 
231  const DetLayer* tid1 = tidLayers.at(0);
232  const DetLayer* tid2 = tidLayers.at(1);
233  const DetLayer* tid3 = tidLayers.at(2);
234  const DetLayer* tec1 = tecLayers.at(0);
235  const DetLayer* tec2 = tecLayers.at(1);
236  const DetLayer* tec3 = tecLayers.at(2);
237 
238  //Figure out which DetLayers to use based on SC Eta
239  std::vector<bool> useDL = useDetLayer(scEta);
240  bool useTID = false;
241 
242  //Use counters to restrict the number of hits in TID and TEC layers
243  //This reduces seed multiplicity
244  int tid1MHC = 0;
245  int tid2MHC = 0;
246  int tid3MHC = 0;
247  int tid1BHC = 0;
248  int tid2BHC = 0;
249  int tid3BHC = 0;
250  int tec1MHC = 0;
251  int tec2MHC = 0;
252  int tec3MHC = 0;
253 
254  //Use counter to limit the allowed number of seeds
255  int seedCounter = 0;
256 
257  bool hasLay1Hit = false;
258  bool hasLay2Hit = false;
259  bool hasBackupHit = false;
260 
261  LayerMeasurements layerMeasurements(*theMeasurementTracker, trackerData);
262 
263  std::vector<TrajectoryMeasurement> tib1measurements;
264  if(useDL.at(0)) tib1measurements = layerMeasurements.measurements(*tib1,initialTSOS,*thePropagator,*theEstimator);
265  std::vector<TrajectoryMeasurement> tib2measurements;
266  if(useDL.at(1)) tib2measurements = layerMeasurements.measurements(*tib2,initialTSOS,*thePropagator,*theEstimator);
267 
268  //Basic idea: Retrieve hits from a given DetLayer
269  //Check if it is a Matched Hit and satisfies some cuts
270  //If yes, accept hit for seed making
271 
272  for(std::vector<TrajectoryMeasurement>::const_iterator tmIter = tib1measurements.begin(); tmIter != tib1measurements.end(); ++ tmIter){
273  ConstRecHitPointer hit = tmIter->recHit();
274  const SiStripMatchedRecHit2D* matchedHit = matchedHitConverter(hit);
275  if(matchedHit){
276  GlobalPoint position = trackerGeometryHandle->idToDet(matchedHit->geographicalId())->surface().toGlobal(matchedHit->localPosition());
277  if(preselection(position, superCluster, phiVsRSlope, 1)){
278  hasLay1Hit = true;
279  layer1Hits_.push_back(matchedHit);
280  }
281  }
282  }
283 
284  for(std::vector<TrajectoryMeasurement>::const_iterator tmIter = tib2measurements.begin(); tmIter != tib2measurements.end(); ++ tmIter){
285  ConstRecHitPointer hit = tmIter->recHit();
286  const SiStripMatchedRecHit2D* matchedHit = matchedHitConverter(hit);
287  if(matchedHit){
288  GlobalPoint position = trackerGeometryHandle->idToDet(matchedHit->geographicalId())->surface().toGlobal(matchedHit->localPosition());
289  if(preselection(position, superCluster, phiVsRSlope, 1)){
290  hasLay2Hit = true;
291  layer2Hits_.push_back(matchedHit);
292  }
293  }
294  }
295 
296  if(!(hasLay1Hit && hasLay2Hit)) useTID = true;
297  if(std::abs(scEta) > tidEtaUsage_) useTID = true;
298  std::vector<TrajectoryMeasurement> tid1measurements;
299  if(useDL.at(2) && useTID) tid1measurements = layerMeasurements.measurements(*tid1,initialTSOS,*thePropagator,*theEstimator);
300  std::vector<TrajectoryMeasurement> tid2measurements;
301  if(useDL.at(3) && useTID) tid2measurements = layerMeasurements.measurements(*tid2,initialTSOS,*thePropagator,*theEstimator);
302  std::vector<TrajectoryMeasurement> tid3measurements;
303  if(useDL.at(4) && useTID) tid3measurements = layerMeasurements.measurements(*tid3,initialTSOS,*thePropagator,*theEstimator);
304 
305  for(std::vector<TrajectoryMeasurement>::const_iterator tmIter = tid1measurements.begin(); tmIter != tid1measurements.end(); ++ tmIter){
306  if(tid1MHC < tidMaxHits_){
307  ConstRecHitPointer hit = tmIter->recHit();
308  const SiStripMatchedRecHit2D* matchedHit = matchedHitConverter(hit);
309  if(matchedHit){
310  GlobalPoint position = trackerGeometryHandle->idToDet(matchedHit->geographicalId())->surface().toGlobal(matchedHit->localPosition());
311  if(preselection(position, superCluster, phiVsRSlope, 2)){
312  tid1MHC++;
313  hasLay1Hit = true;
314  layer1Hits_.push_back(matchedHit);
315  hasLay2Hit = true;
316  layer2Hits_.push_back(matchedHit);
317  }
318  }else if(useDL.at(8) && tid1BHC < monoMaxHits_){
319  const SiStripRecHit2D* backupHit = backupHitConverter(hit);
320  if(backupHit){
321  GlobalPoint position = trackerGeometryHandle->idToDet(backupHit->geographicalId())->surface().toGlobal(backupHit->localPosition());
322  if(preselection(position, superCluster, phiVsRSlope, 4) && position.perp() > 37.){
323  tid1BHC++;
324  hasBackupHit = true;
325  backupLayer2Hits_.push_back(backupHit);
326  }
327  }
328  }
329  }
330  }
331 
332  for(std::vector<TrajectoryMeasurement>::const_iterator tmIter = tid2measurements.begin(); tmIter != tid2measurements.end(); ++ tmIter){
333  if(tid2MHC < tidMaxHits_){
334  ConstRecHitPointer hit = tmIter->recHit();
335  const SiStripMatchedRecHit2D* matchedHit = matchedHitConverter(hit);
336  if(matchedHit){
337  GlobalPoint position = trackerGeometryHandle->idToDet(matchedHit->geographicalId())->surface().toGlobal(matchedHit->localPosition());
338  if(preselection(position, superCluster, phiVsRSlope, 2)){
339  tid2MHC++;
340  hasLay1Hit = true;
341  layer1Hits_.push_back(matchedHit);
342  hasLay2Hit = true;
343  layer2Hits_.push_back(matchedHit);
344  }
345  }else if(useDL.at(8) && tid2BHC < monoMaxHits_){
346  const SiStripRecHit2D* backupHit = backupHitConverter(hit);
347  if(backupHit){
348  GlobalPoint position = trackerGeometryHandle->idToDet(backupHit->geographicalId())->surface().toGlobal(backupHit->localPosition());
349  if(preselection(position, superCluster, phiVsRSlope, 4) && position.perp() > 37.){
350  tid2BHC++;
351  hasBackupHit = true;
352  backupLayer2Hits_.push_back(backupHit);
353  }
354  }
355  }
356  }
357  }
358 
359  for(std::vector<TrajectoryMeasurement>::const_iterator tmIter = tid3measurements.begin(); tmIter != tid3measurements.end(); ++ tmIter){
360  if(tid3MHC < tidMaxHits_){
361  ConstRecHitPointer hit = tmIter->recHit();
362  const SiStripMatchedRecHit2D* matchedHit = matchedHitConverter(hit);
363  if(matchedHit){
364  GlobalPoint position = trackerGeometryHandle->idToDet(matchedHit->geographicalId())->surface().toGlobal(matchedHit->localPosition());
365  if(preselection(position, superCluster, phiVsRSlope, 2)){
366  tid3MHC++;
367  hasLay1Hit = true;
368  layer1Hits_.push_back(matchedHit);
369  hasLay2Hit = true;
370  layer2Hits_.push_back(matchedHit);
371  }
372  }else if(useDL.at(8) && tid3BHC < monoMaxHits_){
373  const SiStripRecHit2D* backupHit = backupHitConverter(hit);
374  if(backupHit){
375  GlobalPoint position = trackerGeometryHandle->idToDet(backupHit->geographicalId())->surface().toGlobal(backupHit->localPosition());
376  if(preselection(position, superCluster, phiVsRSlope, 4) && position.perp() > 37.){
377  tid3BHC++;
378  hasBackupHit = true;
379  backupLayer2Hits_.push_back(backupHit);
380  }
381  }
382  }
383  }
384  }
385 
386  std::vector<TrajectoryMeasurement> tec1measurements;
387  if(useDL.at(5)) tec1measurements = layerMeasurements.measurements(*tec1,initialTSOS,*thePropagator,*theEstimator);
388  std::vector<TrajectoryMeasurement> tec2measurements;
389  if(useDL.at(6)) tec2measurements = layerMeasurements.measurements(*tec2,initialTSOS,*thePropagator,*theEstimator);
390  std::vector<TrajectoryMeasurement> tec3measurements;
391  if(useDL.at(7)) tec3measurements = layerMeasurements.measurements(*tec3,initialTSOS,*thePropagator,*theEstimator);
392 
393  for(std::vector<TrajectoryMeasurement>::const_iterator tmIter = tec1measurements.begin(); tmIter != tec1measurements.end(); ++ tmIter){
394  if(tec1MHC < tecMaxHits_){
395  ConstRecHitPointer hit = tmIter->recHit();
396  const SiStripMatchedRecHit2D* matchedHit = matchedHitConverter(hit);
397  if(matchedHit){
398  GlobalPoint position = trackerGeometryHandle->idToDet(matchedHit->geographicalId())->surface().toGlobal(matchedHit->localPosition());
399  if(preselection(position, superCluster, phiVsRSlope, 3)){
400  tec1MHC++;
401  hasLay1Hit = true;
402  layer1Hits_.push_back(matchedHit);
403  hasLay2Hit = true;
404  layer2Hits_.push_back(matchedHit);
405  }
406  }
407  }
408  }
409 
410  for(std::vector<TrajectoryMeasurement>::const_iterator tmIter = tec2measurements.begin(); tmIter != tec2measurements.end(); ++ tmIter){
411  if(tec2MHC < tecMaxHits_){
412  ConstRecHitPointer hit = tmIter->recHit();
413  const SiStripMatchedRecHit2D* matchedHit = matchedHitConverter(hit);
414  if(matchedHit){
415  GlobalPoint position = trackerGeometryHandle->idToDet(matchedHit->geographicalId())->surface().toGlobal(matchedHit->localPosition());
416  if(preselection(position, superCluster, phiVsRSlope, 3)){
417  tec2MHC++;
418  hasLay1Hit = true;
419  layer1Hits_.push_back(matchedHit);
420  hasLay2Hit = true;
421  layer2Hits_.push_back(matchedHit);
422  }
423  }
424  }
425  }
426 
427  for(std::vector<TrajectoryMeasurement>::const_iterator tmIter = tec3measurements.begin(); tmIter != tec3measurements.end(); ++ tmIter){
428  if(tec3MHC < tecMaxHits_){
429  ConstRecHitPointer hit = tmIter->recHit();
430  const SiStripMatchedRecHit2D* matchedHit = matchedHitConverter(hit);
431  if(matchedHit){
432  GlobalPoint position = trackerGeometryHandle->idToDet(matchedHit->geographicalId())->surface().toGlobal(matchedHit->localPosition());
433  if(preselection(position, superCluster, phiVsRSlope, 3)){
434  tec3MHC++;
435  hasLay2Hit = true;
436  layer2Hits_.push_back(matchedHit);
437  }
438  }
439  }
440  }
441 
442  // We have 2 arrays of hits, combine them to form seeds
443  if( hasLay1Hit && hasLay2Hit ){
444 
445  for (std::vector<const SiStripMatchedRecHit2D*>::const_iterator hit1 = layer1Hits_.begin() ; hit1!= layer1Hits_.end(); ++hit1) {
446  for (std::vector<const SiStripMatchedRecHit2D*>::const_iterator hit2 = layer2Hits_.begin() ; hit2!= layer2Hits_.end(); ++hit2) {
447 
448  if(seedCounter < maxSeeds_){
449 
450  if(checkHitsAndTSOS(hit1,hit2,scr,scz,pT,scEta)) {
451 
452  seedCounter++;
453 
454  recHits_.clear();
455 
457  hit=new SiStripMatchedRecHit2D(*(dynamic_cast <const SiStripMatchedRecHit2D *> (*hit1) ) );
458  recHits_.push_back(hit);
459  hit=new SiStripMatchedRecHit2D(*(dynamic_cast <const SiStripMatchedRecHit2D *> (*hit2) ) );
460  recHits_.push_back(hit);
461 
464  reco::ElectronSeed::CaloClusterRef caloCluster(seedCluster) ;
465  seed.setCaloCluster(caloCluster) ;
466  result.push_back(seed);
467 
468  }
469 
470  }
471 
472  }// end of hit 2 loop
473 
474  }// end of hit 1 loop
475 
476  }//end of seed making
477 
478  //Make seeds using TID Ring 3 if necessary
479 
480  if(hasLay1Hit && hasBackupHit && seedCounter == 0){
481 
482  for (std::vector<const SiStripMatchedRecHit2D*>::const_iterator hit1 = layer1Hits_.begin() ; hit1!= layer1Hits_.end(); ++hit1) {
483  for (std::vector<const SiStripRecHit2D*>::const_iterator hit2 = backupLayer2Hits_.begin() ; hit2!= backupLayer2Hits_.end(); ++hit2) {
484 
485  if(seedCounter < maxSeeds_){
486 
487  if(altCheckHitsAndTSOS(hit1,hit2,scr,scz,pT,scEta)) {
488 
489  seedCounter++;
490 
491  recHits_.clear();
492 
493  SiStripMatchedRecHit2D *innerHit;
494  innerHit=new SiStripMatchedRecHit2D(*(dynamic_cast <const SiStripMatchedRecHit2D *> (*hit1) ) );
495  recHits_.push_back(innerHit);
496  SiStripRecHit2D *outerHit;
497  outerHit=new SiStripRecHit2D(*(dynamic_cast <const SiStripRecHit2D *> (*hit2) ) );
498  recHits_.push_back(outerHit);
499 
502  reco::ElectronSeed::CaloClusterRef caloCluster(seedCluster) ;
503  seed.setCaloCluster(caloCluster) ;
504  result.push_back(seed);
505 
506  }
507 
508  }
509 
510  }// end of hit 2 loop
511 
512  }// end of hit 1 loop
513 
514  }// end of backup seed making
515 
516 } // end of findSeedsFromCluster
517 
518 
519 
520 bool SiStripElectronSeedGenerator::checkHitsAndTSOS(std::vector<const SiStripMatchedRecHit2D*>::const_iterator hit1,
521  std::vector<const SiStripMatchedRecHit2D*>::const_iterator hit2,
522  double rc,double zc,double pT,double scEta) {
523 
524  bool seedCutsSatisfied = false;
525 
526  using namespace std;
527 
528  GlobalPoint hit1Pos = trackerGeometryHandle->idToDet((*hit1)->geographicalId())->surface().toGlobal((*hit1)->localPosition());
529  double r1 = sqrt(hit1Pos.x()*hit1Pos.x() + hit1Pos.y()*hit1Pos.y());
530  double phi1 = hit1Pos.phi();
531  double z1=hit1Pos.z();
532 
533  GlobalPoint hit2Pos = trackerGeometryHandle->idToDet((*hit2)->geographicalId())->surface().toGlobal((*hit2)->localPosition());
534  double r2 = sqrt(hit2Pos.x()*hit2Pos.x() + hit2Pos.y()*hit2Pos.y());
535  double phi2 = hit2Pos.phi();
536  double z2 = hit2Pos.z();
537 
538  if(r2 > r1 && (std::abs(z2) > std::abs(z1) || std::abs(scEta) < 0.25)) {
539 
540  //Consider the circle made of IP and Hit 1; Calculate it's radius using pT
541 
542  double curv = pT*100*.877;
543 
544  //Predict phi of hit 2
545  double a = (r2-r1)/(2*curv);
546  double b = phiDiff(phi2,phi1);
547  //UB added '=0' to avoid compiler warning
548  double phiMissHit2=0;
549  if(std::abs(b - a)<std::abs(b + a)) phiMissHit2 = b - a;
550  if(std::abs(b - a)>std::abs(b + a)) phiMissHit2 = b + a;
551 
552  double zMissHit2 = z2 - (r2*(zc-z1)-r1*zc+rc*z1)/(rc-r1);
553 
554  double rPredHit2 = r1 + (rc-r1)/(zc-z1)*(z2-z1);
555  double rMissHit2 = r2 - rPredHit2;
556 
557  int subdetector = whichSubdetector(hit2);
558 
559  bool zDiff = true;
560  double zVar1 = std::abs(z1);
561  double zVar2 = std::abs(z2 - z1);
562  if(zVar1 > 75 && zVar1 < 95 && (zVar2 > 18 || zVar2 < 5)) zDiff = false;
563  if(zVar1 > 100 && zVar1 < 110 && (zVar2 > 35 || zVar2 < 5)) zDiff = false;
564  if(zVar1 > 125 && zVar1 < 150 && (zVar2 > 18 || zVar2 < 5)) zDiff = false;
565 
566  if(subdetector == 1){
567  int tibExtraCut = 0;
568  if(r1 > 23 && r1 < 28 && r2 > 31 && r2 < 37) tibExtraCut = 1;
569  if(std::abs(phiMissHit2) < tibPhiMissHit2Cut_ && std::abs(zMissHit2) < tibZMissHit2Cut_ && tibExtraCut == 1) seedCutsSatisfied = true;
570  }else if(subdetector == 2){
571  int tidExtraCut = 0;
572  if(r1 > 23 && r1 < 34 && r2 > 26 && r2 < 42) tidExtraCut = 1;
573  if(std::abs(phiMissHit2) < tidPhiMissHit2Cut_ && std::abs(rMissHit2) < tidRMissHit2Cut_ && tidExtraCut == 1 && zDiff) seedCutsSatisfied = true;
574  }else if(subdetector == 3){
575  int tecExtraCut = 0;
576  if(r1 > 23 && r1 < 32 && r2 > 26 && r2 < 42) tecExtraCut = 1;
577  if(std::abs(phiMissHit2) < tecPhiMissHit2Cut_ && std::abs(rMissHit2) < tecRMissHit2Cut_ && tecExtraCut == 1 && zDiff) seedCutsSatisfied = true;
578  }
579 
580  }
581 
582  if(!seedCutsSatisfied) return false;
583 
584  // seed checks borrowed from pixel-based algoritm
585 
586 
587 
589 
590  double vertexZ = z1 - (r1 * (zc - z1) ) / (rc - r1);
591  GlobalPoint eleVertex(0.,0.,vertexZ);
592 
593  // FIXME optimize: move outside loop
595  theSetup->get<IdealMagneticFieldRecord>().get(bfield);
596  float nomField = bfield->nominalValue();
597 
598  // make a spiral
599  FastHelix helix(hit2Pos,hit1Pos,eleVertex,nomField,&*bfield);
600  if (!helix.isValid()) return false;
601 
602  FreeTrajectoryState fts(helix.stateAtVertex());
603  TSOS propagatedState = thePropagator->propagate(fts,(*hit1)->det()->surface());
604 
605  if (!propagatedState.isValid()) return false;
606 
607  TSOS updatedState = theUpdator->update(propagatedState, **hit1);
608  TSOS propagatedState_out = thePropagator->propagate(fts,(*hit2)->det()->surface()) ;
609 
610  if (!propagatedState_out.isValid()) return false;
611 
612  // the seed has now passed all the cuts
613 
614  TSOS updatedState_out = theUpdator->update(propagatedState_out, **hit2);
615 
616  pts_ = trajectoryStateTransform::persistentState(updatedState_out, (*hit2)->geographicalId().rawId());
617 
618  return true;
619 }
620 
621 bool SiStripElectronSeedGenerator::altCheckHitsAndTSOS(std::vector<const SiStripMatchedRecHit2D*>::const_iterator hit1,
622  std::vector<const SiStripRecHit2D*>::const_iterator hit2,
623  double rc,double zc,double pT,double scEta) {
624 
625  bool seedCutSatisfied = false;
626 
627  using namespace std;
628 
629  GlobalPoint hit1Pos = trackerGeometryHandle->idToDet((*hit1)->geographicalId())->surface().toGlobal((*hit1)->localPosition());
630  double r1 = sqrt(hit1Pos.x()*hit1Pos.x() + hit1Pos.y()*hit1Pos.y());
631  double phi1 = hit1Pos.phi();
632  double z1=hit1Pos.z();
633 
634  GlobalPoint hit2Pos = trackerGeometryHandle->idToDet((*hit2)->geographicalId())->surface().toGlobal((*hit2)->localPosition());
635  double r2 = sqrt(hit2Pos.x()*hit2Pos.x() + hit2Pos.y()*hit2Pos.y());
636  double phi2 = hit2Pos.phi();
637  double z2 = hit2Pos.z();
638 
639  if(r2 > r1 && std::abs(z2) > std::abs(z1)) {
640 
641  //Consider the circle made of IP and Hit 1; Calculate it's radius using pT
642 
643  double curv = pT*100*.877;
644 
645  //Predict phi of hit 2
646  double a = (r2-r1)/(2*curv);
647  double b = phiDiff(phi2,phi1);
648  double phiMissHit2 = 0;
649  if(std::abs(b - a)<std::abs(b + a)) phiMissHit2 = b - a;
650  if(std::abs(b - a)>std::abs(b + a)) phiMissHit2 = b + a;
651 
652  if(std::abs(phiMissHit2) < monoPhiMissHit2Cut_) seedCutSatisfied = true;
653 
654  }
655 
656  if(!seedCutSatisfied) return false;
657 
658  // seed checks borrowed from pixel-based algoritm
659 
660 
661 
662 
664 
665  double vertexZ = z1 - (r1 * (zc - z1) ) / (rc - r1);
666  GlobalPoint eleVertex(0.,0.,vertexZ);
667 
668  // FIXME optimize: move outside loop
670  theSetup->get<IdealMagneticFieldRecord>().get(bfield);
671  float nomField = bfield->nominalValue();
672 
673  // make a spiral
674  FastHelix helix(hit2Pos,hit1Pos,eleVertex,nomField,&*bfield);
675  if (!helix.isValid()) return false;
676 
677  FreeTrajectoryState fts(helix.stateAtVertex());
678  TSOS propagatedState = thePropagator->propagate(fts,(*hit1)->det()->surface());
679 
680  if (!propagatedState.isValid()) return false;
681 
682  TSOS updatedState = theUpdator->update(propagatedState, **hit1);
683  TSOS propagatedState_out = thePropagator->propagate(fts,(*hit2)->det()->surface()) ;
684 
685  if (!propagatedState_out.isValid()) return false;
686 
687  // the seed has now passed all the cuts
688 
689  TSOS updatedState_out = theUpdator->update(propagatedState_out, **hit2);
690 
691  pts_ = trajectoryStateTransform::persistentState(updatedState_out, (*hit2)->geographicalId().rawId());
692 
693  return true;
694 }
695 
696 
697 bool SiStripElectronSeedGenerator::preselection(GlobalPoint position,GlobalPoint superCluster,double phiVsRSlope,int hitLayer){
698  double r = position.perp();
699  double phi = position.phi();
700  double z = position.z();
701  double scr = superCluster.perp();
702  double scphi = superCluster.phi();
703  double scz = superCluster.z();
704  double psi = phiDiff(phi,scphi);
705  double deltaPsi = psi - (scr-r)*phiVsRSlope;
706  double antiDeltaPsi = psi - (r-scr)*phiVsRSlope;
707  double dP;
708  if (std::abs(deltaPsi)<std::abs(antiDeltaPsi)){
709  dP = deltaPsi;
710  }else{
711  dP = antiDeltaPsi;
712  }
713  double originZ = (scr*z - r*scz)/(scr-r);
714 
715  bool result = false;
716 
717  if(hitLayer == 1){
718  if(std::abs(originZ) < tibOriginZCut_ && std::abs(dP) < tibDeltaPsiCut_) result = true;
719  }else if(hitLayer == 2){
720  if(std::abs(originZ) < tidOriginZCut_ && std::abs(dP) < tidDeltaPsiCut_) result = true;
721  }else if(hitLayer == 3){
722  if(std::abs(originZ) < tecOriginZCut_ && std::abs(dP) < tecDeltaPsiCut_) result = true;
723  }else if(hitLayer == 4){
724  if(std::abs(originZ) < monoOriginZCut_ && std::abs(dP) < monoDeltaPsiCut_) result = true;
725  }
726 
727  return result;
728 }
729 
730 // Helper algorithms
731 
732 int SiStripElectronSeedGenerator::whichSubdetector(std::vector<const SiStripMatchedRecHit2D*>::const_iterator hit){
733  int result = 0;
734  if(((*hit)->geographicalId()).subdetId() == StripSubdetector::TIB){
735  result = 1;
736  }else if(((*hit)->geographicalId()).subdetId() == StripSubdetector::TID){
737  result = 2;
738  }else if(((*hit)->geographicalId()).subdetId() == StripSubdetector::TEC){
739  result = 3;
740  }
741  return result;
742 }
743 
745  const TrackingRecHit* trh = crhp->hit();
746  const SiStripMatchedRecHit2D* matchedHit = dynamic_cast<const SiStripMatchedRecHit2D*>(trh);
747  return matchedHit;
748 }
749 
751  const TrackingRecHit* trh = crhp->hit();
752  const SiStripRecHit2D* backupHit = dynamic_cast<const SiStripRecHit2D*>(trh);
753  return backupHit;
754 }
755 
756 std::vector<bool> SiStripElectronSeedGenerator::useDetLayer(double scEta){
757  std::vector<bool> useDetLayer;
758  double variable = std::abs(scEta);
759  if(variable > 0 && variable < 1.8){
760  useDetLayer.push_back(true);
761  }else{
762  useDetLayer.push_back(false);
763  }
764  if(variable > 0 && variable < 1.5){
765  useDetLayer.push_back(true);
766  }else{
767  useDetLayer.push_back(false);
768  }
769  if(variable > 1 && variable < 2.1){
770  useDetLayer.push_back(true);
771  }else{
772  useDetLayer.push_back(false);
773  }
774  if(variable > 1 && variable < 2.2){
775  useDetLayer.push_back(true);
776  }else{
777  useDetLayer.push_back(false);
778  }
779  if(variable > 1 && variable < 2.3){
780  useDetLayer.push_back(true);
781  }else{
782  useDetLayer.push_back(false);
783  }
784  if(variable > 1.8 && variable < 2.5){
785  useDetLayer.push_back(true);
786  }else{
787  useDetLayer.push_back(false);
788  }
789  if(variable > 1.8 && variable < 2.5){
790  useDetLayer.push_back(true);
791  }else{
792  useDetLayer.push_back(false);
793  }
794  if(variable > 1.8 && variable < 2.5){
795  useDetLayer.push_back(true);
796  }else{
797  useDetLayer.push_back(false);
798  }
799  if(variable > 1.2 && variable < 1.6){
800  useDetLayer.push_back(true);
801  }else{
802  useDetLayer.push_back(false);
803  }
804  return useDetLayer;
805 }
806 
807 
808 
#define LogDebug(id)
T getParameter(std::string const &) const
unsigned long long cacheIdentifier() const
double z0() const
z coordinate
Definition: BeamSpot.h:68
std::vector< TrajectoryMeasurement > measurements(const DetLayer &layer, const TrajectoryStateOnSurface &startingState, const Propagator &prop, const MeasurementEstimator &est) const
bool isValid() const
Definition: FastHelix.h:63
bool checkHitsAndTSOS(std::vector< const SiStripMatchedRecHit2D * >::const_iterator hit1, std::vector< const SiStripMatchedRecHit2D * >::const_iterator hit2, double scr, double scz, double pT, double scEta)
T perp() const
Definition: PV3DBase.h:72
std::vector< bool > useDetLayer(double scEta)
int nominalValue() const
The nominal field value for this map in kGauss.
Definition: MagneticField.h:58
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
TString subdetector
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:54
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
const MeasurementTracker * theMeasurementTracker
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
PTrajectoryStateOnDet persistentState(const TrajectoryStateOnSurface &ts, unsigned int detid)
T y() const
Definition: PV3DBase.h:63
edm::ESHandle< MagneticField > theMagField
bool exists(std::string const &parameterName) const
checks if a parameter exists
TransientTrackingRecHit::ConstRecHitPointer ConstRecHitPointer
Chi2MeasurementEstimator * theEstimator
def setup(process, global_tag, zero_tesla=False)
Definition: GeneralSetup.py:1
const SiStripMatchedRecHit2D * matchedHitConverter(ConstRecHitPointer crhp)
PropagationDirection
#define nullptr
edm::EDGetTokenT< reco::BeamSpot > beamSpotTag_
edm::EDGetTokenT< MeasurementTrackerEvent > theMeasurementTrackerEventTag
std::map< std::string, int, std::less< std::string > > psi
void push_back(D *&d)
Definition: OwnVector.h:290
TrajectoryStateOnSurface update(const TrajectoryStateOnSurface &, const TrackingRecHit &) const override
Definition: KFUpdator.cc:169
T sqrt(T t)
Definition: SSEVec.h:18
T z() const
Definition: PV3DBase.h:64
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
void clear()
Definition: OwnVector.h:445
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
SiStripElectronSeedGenerator(const edm::ParameterSet &, const Tokens &)
std::vector< ElectronSeed > ElectronSeedCollection
collection of ElectronSeed objects
edm::ESHandle< TrackerGeometry > trackerGeometryHandle
edm::Handle< reco::BeamSpot > theBeamSpot
int whichSubdetector(std::vector< const SiStripMatchedRecHit2D * >::const_iterator hit)
void findSeedsFromCluster(edm::Ref< reco::SuperClusterCollection >, edm::Handle< reco::BeamSpot >, const MeasurementTrackerEvent &trackerData, reco::ElectronSeedCollection &)
bool preselection(GlobalPoint position, GlobalPoint superCluster, double phiVsRSlope, int hitLayer)
TrajectoryStateOnSurface TSOS
Definition: TestHits.cc:19
virtual TrackingRecHit const * hit() const
std::vector< BarrelDetLayer const * > const & tibLayers() const
const SiStripRecHit2D * backupHitConverter(ConstRecHitPointer crhp)
bool altCheckHitsAndTSOS(std::vector< const SiStripMatchedRecHit2D * >::const_iterator hit1, std::vector< const SiStripRecHit2D * >::const_iterator hit2, double scr, double scz, double pT, double scEta)
void setupES(const edm::EventSetup &setup)
std::vector< ForwardDetLayer const * > const & posTecLayers() const
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
double b
Definition: hdecay.h:120
std::vector< const SiStripMatchedRecHit2D * > layer2Hits_
TrajectoryStateOnSurface propagate(STA const &state, SUR const &surface) const
Definition: Propagator.h:53
std::vector< ForwardDetLayer const * > const & negTidLayers() const
double phiDiff(double phi1, double phi2)
void setCaloCluster(const CaloClusterRef &clus)
Definition: ElectronSeed.h:89
std::vector< ForwardDetLayer const * > const & posTidLayers() const
edm::EventID id() const
Definition: EventBase.h:60
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
double a
Definition: hdecay.h:121
static int position[264][3]
Definition: ReadPGInfo.cc:509
T get() const
Definition: EventSetup.h:63
const TrackerGeomDet * idToDet(DetId) const override
double y0() const
y coordinate
Definition: BeamSpot.h:66
LocalPoint localPosition() const final
edm::ESHandle< MeasurementTracker > measurementTrackerHandle
std::vector< const SiStripRecHit2D * > backupLayer2Hits_
void run(edm::Event &, const edm::EventSetup &setup, const edm::Handle< reco::SuperClusterCollection > &, reco::ElectronSeedCollection &)
GlobalTrajectoryParameters stateAtVertex() const
Definition: FastHelix.h:65
DetId geographicalId() const
dbl *** dir
Definition: mlp_gen.cc:35
std::vector< ForwardDetLayer const * > const & negTecLayers() const
std::vector< const SiStripMatchedRecHit2D * > layer1Hits_
T x() const
Definition: PV3DBase.h:62
T const * product() const
Definition: ESHandle.h:86
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
double x0() const
x coordinate
Definition: BeamSpot.h:64