CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SiStripElectronAlgo.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: EgammaElectronAlgos
4 // Class : SiStripElectronAlgo
5 //
6 // Implementation:
7 // <Notes on implementation>
8 //
9 // Original Author: Jim Pivarski
10 // Created: Fri May 26 16:12:04 EDT 2006
11 // $Id: SiStripElectronAlgo.cc,v 1.41 2013/01/02 18:59:12 dlange Exp $
12 //
13 
14 // system include files
15 
16 // user include files
20 
21 
24 
28 
29 
30 //
31 // constants, enums and typedefs
32 //
33 
34 //
35 // static data member definitions
36 //
37 namespace {
38  struct CompareHits {
39  CompareHits( const TrackingRecHitLessFromGlobalPosition& iLess ) :
40  less_(iLess) {}
41  bool operator()(const TrackingRecHit* iLHS, const TrackingRecHit* iRHS) {
42  return less_(*iLHS,*iRHS);
43  }
44 
46  };
47 }
48 
49 
50 //
51 // constructors and destructor
52 //
53 SiStripElectronAlgo::SiStripElectronAlgo(unsigned int maxHitsOnDetId,
54  double originUncertainty,
55  double phiBandWidth,
56  double maxNormResid,
57  unsigned int minHits,
58  double maxReducedChi2)
59  : maxHitsOnDetId_(maxHitsOnDetId)
60  , originUncertainty_(originUncertainty)
61  , phiBandWidth_(phiBandWidth)
62  , maxNormResid_(maxNormResid)
63  , minHits_(minHits)
64  , maxReducedChi2_(maxReducedChi2)
65 {
66 }
67 
68 // SiStripElectronAlgo::SiStripElectronAlgo(const SiStripElectronAlgo& rhs)
69 // {
70 // // do actual copying here;
71 // }
72 
74 {
75 }
76 
77 //
78 // assignment operators
79 //
80 // const SiStripElectronAlgo& SiStripElectronAlgo::operator=(const SiStripElectronAlgo& rhs)
81 // {
82 // //An exception safe implementation is
83 // SiStripElectronAlgo temp(rhs);
84 // swap(rhs);
85 //
86 // return *this;
87 // }
88 
89 //
90 // member functions
91 //
92 
97  const edm::ESHandle<MagneticField>& magneticField)
98 {
99  LogDebug("") << " In prepareEvent " ;
100 
101  tracker_p_ = tracker.product();
102  rphiHits_p_ = rphiHits.product();
103  stereoHits_p_ = stereoHits.product();
104  matchedHits_p_ = matchedHits.product();
105  magneticField_p_ = magneticField.product();
106 
107  rphiHits_hp_ = &rphiHits;
108  stereoHits_hp_ = &stereoHits;
109 
110  // Keep a table that relates hit pointers to their index (key) in the collections
111  rphiKey_.clear();
112  stereoKey_.clear();
113  // Keep track of which hits have been used already (so a hit is assigned to only one electron)
114  hitUsed_.clear();
115  matchedHitUsed_.clear();
116 
117  unsigned int counter = 0;
119  rphiKey_[&(*it)] = counter;
120  hitUsed_[&(*it)] = false;
121  counter++;
122  }
123 
124  counter = 0;
126  stereoKey_[&(*it)] = counter;
127  hitUsed_[&(*it)] = false;
128  counter++;
129  }
130 
131  counter = 0;
133  matchedKey_[&(*it)] = counter;
134  matchedHitUsed_[&(*it)] = false;
135  counter++;
136  }
137 
138  LogDebug("") << " Leaving prepareEvent " ;
139 }
140 
141 
142 
143 // returns true iff an electron was found
144 // inserts electrons and trackcandiates into electronOut and trackCandidateOut
146  TrackCandidateCollection& trackCandidateOut,
147  const reco::SuperClusterRef& superclusterIn,
148  const TrackerTopology *tTopo)
149 {
150  // Try each of the two charge hypotheses, but only take one
151  bool electronSuccess = projectPhiBand(-1., superclusterIn, tTopo);
152  bool positronSuccess = projectPhiBand( 1., superclusterIn, tTopo);
153 
154  // electron hypothesis did better than electron
155  if ((electronSuccess && !positronSuccess) ||
156  (electronSuccess && positronSuccess && redchi2_neg_ <= redchi2_pos_)) {
157 
158  // Initial uncertainty for tracking
159  AlgebraicSymMatrix55 errors; // makes 5x5 matrix, indexed from (0,0) to (44)
160  errors(0,0) = 3.; // uncertainty**2 in 1/momentum
161  errors(1,1) = 0.01; // uncertainty**2 in lambda (lambda == pi/2 - polar angle theta)
162  errors(2,2) = 0.0001; // uncertainty**2 in phi
163  errors(3,3) = 0.01; // uncertainty**2 in x_transverse (where x is in cm)
164  errors(4,4) = 0.01; // uncertainty**2 in y_transverse (where y is in cm)
165 
166 #ifdef EDM_ML_DEBUG
167  // JED Debugging possible double hit sort problem
168  std::ostringstream debugstr6;
169  debugstr6 << " HERE BEFORE SORT electron case " << " \n" ;
170  for (std::vector<const TrackingRecHit*>::iterator itHit=outputHits_neg_.begin();
171  itHit != outputHits_neg_.end(); ++itHit) {
172  debugstr6 <<" HIT "<<((*itHit)->geographicalId()).rawId()<<" \n"
173  <<" Local Position: "<<(*itHit)->localPosition()<<" \n"
174  <<" Global Rho: "
175  <<(((TrackingGeometry*)(tracker_p_))->idToDet((*itHit)->geographicalId())->surface().toGlobal((*itHit)->localPosition()).perp())
176  <<" Phi "
177  <<(((TrackingGeometry*)(tracker_p_))->idToDet((*itHit)->geographicalId())->surface().toGlobal((*itHit)->localPosition()).phi())
178  << " Z "
179  <<(((TrackingGeometry*)(tracker_p_))->idToDet((*itHit)->geographicalId())->surface().toGlobal((*itHit)->localPosition()).z())<< " \n";
180  }
181  // end of dump
182 
183 
184 
185  // JED call dump
186  debugstr6 << " Calling sort alongMomentum " << " \n";
187 #endif
188 
191 
192 #ifdef EDM_ML_DEBUG
193  debugstr6 << " Done with sort " << " \n";
194 
195  debugstr6 << " HERE AFTER SORT electron case " << " \n";
196  for (std::vector<const TrackingRecHit*>::iterator itHit=outputHits_neg_.begin();
197  itHit != outputHits_neg_.end(); ++itHit) {
198  debugstr6 <<" HIT "<<((*itHit)->geographicalId()).rawId()<<" \n"
199  <<" Local Position: "<<(*itHit)->localPosition()<<" \n"
200  <<" Global Rho: "
201  <<(((TrackingGeometry*)(tracker_p_))->idToDet((*itHit)->geographicalId())->surface().toGlobal((*itHit)->localPosition()).perp())
202  <<" Phi "
203  <<(((TrackingGeometry*)(tracker_p_))->idToDet((*itHit)->geographicalId())->surface().toGlobal((*itHit)->localPosition()).phi())
204  << " Z "
205  <<(((TrackingGeometry*)(tracker_p_))->idToDet((*itHit)->geographicalId())->surface().toGlobal((*itHit)->localPosition()).z())<< " \n";;
206  }
207  // end of dump
208 
209  LogDebug("")<< debugstr6.str();
210 #endif
211 
212  //create an OwnVector needed by the classes which will be stored to the Event
214  for(std::vector<const TrackingRecHit*>::iterator itHit =outputHits_neg_.begin();
215  itHit != outputHits_neg_.end();
216  ++itHit) {
217  hits.push_back( (*itHit)->clone());
218  if( !(hitUsed_.find(*itHit) != hitUsed_.end()) ) {
219  LogDebug("") << " Assert failure " ;
220  assert(hitUsed_.find(*itHit) != hitUsed_.end());
221  }
222  hitUsed_[*itHit] = true;
223  }
224 
228  tracker_p_->idToDet(innerhit_neg_->geographicalId())->surface());
229 
230 
231  PTrajectoryStateOnDet const & PTraj = trajectoryStateTransform::persistentState(state, innerhit_neg_->geographicalId().rawId());
232  TrajectorySeed trajectorySeed(PTraj, hits, alongMomentum);
233  trackCandidateOut.push_back(TrackCandidate(hits, trajectorySeed, PTraj));
234 
235  electronOut.push_back(reco::SiStripElectron(superclusterIn,
236  -1,
240  slope_neg_,
241  intercept_neg_ + superclusterIn->position().phi(),
242  chi2_neg_,
243  ndof_neg_,
245  pZ_neg_,
250 
251  return true;
252  }
253 
254  // positron hypothesis did better than electron
255  if ((!electronSuccess && positronSuccess) ||
256  (electronSuccess && positronSuccess && redchi2_neg_ > redchi2_pos_)) {
257 
258  // Initial uncertainty for tracking
259  AlgebraicSymMatrix55 errors; // same as abovr
260  errors(0,0) = 3.; // uncertainty**2 in 1/momentum
261  errors(1,1) = 0.01; // uncertainty**2 in lambda (lambda == pi/2 - polar angle theta)
262  errors(2,2) = 0.0001; // uncertainty**2 in phi
263  errors(3,3) = 0.01; // uncertainty**2 in x_transverse (where x is in cm)
264  errors(4,4) = 0.01; // uncertainty**2 in y_transverse (where y is in cm)
265 
266 #ifdef EDM_ML_DEBUG
267  // JED Debugging possible double hit sort problem
268  std::ostringstream debugstr7;
269  debugstr7 << " HERE BEFORE SORT Positron case " << " \n";
270  for (std::vector<const TrackingRecHit*>::iterator itHit = outputHits_pos_.begin();
271  itHit != outputHits_pos_.end(); ++itHit) {
272  debugstr7 <<" HIT "<<((*itHit)->geographicalId()).rawId()<<" \n"
273  <<" Local Position: "<<(*itHit)->localPosition()<<" \n"
274  <<" Global Rho: "
275  <<(((TrackingGeometry*)(tracker_p_))->idToDet((*itHit)->geographicalId())->surface().toGlobal((*itHit)->localPosition()).perp())
276  <<" Phi "
277  <<(((TrackingGeometry*)(tracker_p_))->idToDet((*itHit)->geographicalId())->surface().toGlobal((*itHit)->localPosition()).phi())
278  << " Z "
279  <<(((TrackingGeometry*)(tracker_p_))->idToDet((*itHit)->geographicalId())->surface().toGlobal((*itHit)->localPosition()).z())<< " \n";
280  }
281  // end of dump
282 
283  debugstr7 << " Calling sort alongMomentum " << " \n";
284 #endif
285 
288 
289 #ifdef EDM_ML_DEBUG
290  debugstr7 << " Done with sort " << " \n";
291 
292  // JED Debugging possible double hit sort problem
293  debugstr7 << " HERE AFTER SORT Positron case " << " \n";
294  for (std::vector<const TrackingRecHit*>::iterator itHit = outputHits_pos_.begin();
295  itHit != outputHits_pos_.end(); ++itHit) {
296  debugstr7 <<" HIT "<<((*itHit)->geographicalId()).rawId()<<" \n"
297  <<" Local Position: "<<(*itHit)->localPosition()<<" \n"
298  <<" Global Rho: "
299  <<(((TrackingGeometry*)(tracker_p_))->idToDet((*itHit)->geographicalId())->surface().toGlobal((*itHit)->localPosition()).perp())
300  <<" Phi "
301  <<(((TrackingGeometry*)(tracker_p_))->idToDet((*itHit)->geographicalId())->surface().toGlobal((*itHit)->localPosition()).phi())
302  << " Z "
303  <<(((TrackingGeometry*)(tracker_p_))->idToDet((*itHit)->geographicalId())->surface().toGlobal((*itHit)->localPosition()).z())<< " \n";
304  }
305  // end of dump
306  LogDebug("") << debugstr7.str();
307 #endif
308 
309  //create an OwnVector needed by the classes which will be stored to the Event
311  for(std::vector<const TrackingRecHit*>::iterator itHit =outputHits_pos_.begin();
312  itHit != outputHits_pos_.end();
313  ++itHit) {
314  hits.push_back( (*itHit)->clone());
315  assert(hitUsed_.find(*itHit) != hitUsed_.end());
316  hitUsed_[*itHit] = true;
317  }
318 
322  tracker_p_->idToDet(innerhit_pos_->geographicalId())->surface());
323 
324 
325  PTrajectoryStateOnDet const & PTraj = trajectoryStateTransform::persistentState(state, innerhit_pos_->geographicalId().rawId());
326  TrajectorySeed trajectorySeed(PTraj, hits, alongMomentum);
327  trackCandidateOut.push_back(TrackCandidate(hits, trajectorySeed, PTraj));
328 
329  electronOut.push_back(reco::SiStripElectron(superclusterIn,
330  1,
334  slope_pos_,
335  intercept_pos_ + superclusterIn->position().phi(),
336  chi2_pos_,
337  ndof_pos_,
339  pZ_pos_,
344 
345  return true;
346  }
347 
348  return false;
349 }
350 
351 // inserts pointers to good hits into hitPointersOut
352 // selects hits on DetIds that have no more than maxHitsOnDetId_
353 // selects from stereo if stereo == true, rphi otherwise
354 // selects from TID or TEC if endcap == true, TIB or TOB otherwise
355 void SiStripElectronAlgo::coarseHitSelection(std::vector<const SiStripRecHit2D*>& hitPointersOut,
356  const TrackerTopology *tTopo,
357  bool stereo, bool endcap)
358 {
359  // This function is not time-efficient. If you want to improve the
360  // speed of the algorithm, you'll probably want to change this
361  // function. There may be a more efficienct way to extract hits,
362  // and it would definitely help to put a geographical cut on the
363  // DetIds. (How does one determine the global position of a given
364  // DetId? Is tracker_p_->idToDet(id)->surface().toGlobal(LocalPosition(0,0,0))
365  // expensive?)
366 
367  // Loop over the detector ids
370  for (; itdet != eddet; ++itdet) {
371  // Get the hits on this detector id
372  SiStripRecHit2DCollection::DetSet hits = *itdet;
373  DetId id(hits.detId());
374 
375  // Count the number of hits on this detector id
376  unsigned int numberOfHits = hits.size();
377 
378  // Only take the hits if there aren't too many
379  // (Would it be better to loop only once, fill a temporary list,
380  // and copy that if numberOfHits <= maxHitsOnDetId_?)
381  if (numberOfHits <= maxHitsOnDetId_) {
383  hit != hits.end(); ++hit) {
384  // check that hit is valid first !
385  if(!(*hit).isValid()) {
386  LogDebug("") << " InValid hit skipped in coarseHitSelection " << std::endl ;
387  continue ;
388  }
389  std::string theDet = "null";
390  int theLayer = -999;
391  bool isStereoDet = false ;
392  if(tracker_p_->idToDetUnit(hit->geographicalId())->type().subDetector() == GeomDetEnumerators::TIB) {
393  theDet = "TIB" ;
394  theLayer = tTopo->tibLayer(id);
395  if(tTopo->tibStereo(id)==1) { isStereoDet = true ; }
396  } else if
397  (tracker_p_->idToDetUnit(hit->geographicalId())->type().subDetector() == GeomDetEnumerators::TOB) {
398  theDet = "TOB" ;
399  theLayer = tTopo->tobLayer(id);
400  if(tTopo->tobStereo(id)==1) { isStereoDet = true ; }
401  }else if
402  (tracker_p_->idToDetUnit(hit->geographicalId())->type().subDetector() == GeomDetEnumerators::TID) {
403  theDet = "TID" ;
404  theLayer = tTopo->tidWheel(id); // or ring ?
405  if(tTopo->tidStereo(id)==1) { isStereoDet = true ; }
406  }else if
407  (tracker_p_->idToDetUnit(hit->geographicalId())->type().subDetector() == GeomDetEnumerators::TEC) {
408  theDet = "TEC" ;
409  theLayer = tTopo->tecWheel(id); // or ring or petal ?
410  if(tTopo->tecStereo(id)==1) { isStereoDet = true ; }
411  } else {
412  LogDebug("") << " UHOH BIG PROBLEM - Unrecognized SI Layer" ;
413  LogDebug("") << " Det "<< theDet << " Lay " << theLayer ;
414  assert(1!=1) ;
415  }
416 
417  if ((endcap && stereo && (theDet=="TID" || theDet== "TEC") && isStereoDet ) ||
418  (endcap && !stereo && (theDet=="TID" || theDet== "TEC") && !isStereoDet ) ||
419  (!endcap && stereo && (theDet=="TIB" || theDet=="TOB") && isStereoDet ) ||
420  (!endcap && !stereo && (theDet=="TIB" || theDet=="TOB" )&& !isStereoDet )
421  ) {
422 
423 
424  hitPointersOut.push_back(&(*hit));
425 
426  } // end if this is the right subdetector
427  } // end loop over hits
428  } // end if this detector id doesn't have too many hits on it
429  } // end loop over detector ids
430 }
431 
432 // select all matched hits for now
433 
434 void SiStripElectronAlgo::coarseMatchedHitSelection(std::vector<const SiStripMatchedRecHit2D*>& coarseMatchedHitPointersOut)
435 {
436 
437  // Loop over the detector ids
439  for (; itdet != eddet; ++itdet) {
440 
441  // Get the hits on this detector id
443 
444  // Count the number of hits on this detector id
445  unsigned int numberOfHits = 0;
447  if ( !((hit->geographicalId()).subdetId() == StripSubdetector::TIB) &&
448  !( (hit->geographicalId()).subdetId() == StripSubdetector::TOB )) { break;}
449  numberOfHits++;
450  if (numberOfHits > maxHitsOnDetId_) { break; }
451  }
452 
453  // Only take the hits if there aren't too many
454  if (numberOfHits <= maxHitsOnDetId_) {
456  if(!(*hit).isValid()) {
457  LogDebug("") << " InValid hit skipped in coarseMatchedHitSelection " << std::endl ;
458  continue ;
459  }
460  if ( !((hit->geographicalId()).subdetId() == StripSubdetector::TIB) &&
461  !( (hit->geographicalId()).subdetId() == StripSubdetector::TOB )) { break;}
462 
463  coarseMatchedHitPointersOut.push_back(&(*hit));
464  } // end loop over hits
465 
466  } // end if this detector id doesn't have too many hits on it
467  } // end loop over detector ids
468 
469 
470 }// end of matchedHitSelection
471 
472 
473 
474 
475 // projects a phi band of width phiBandWidth_ from supercluster into tracker (given a chargeHypothesis)
476 // fills *_pos_ or *_neg_ member data with the results
477 // returns true iff the electron/positron passes cuts
478 bool SiStripElectronAlgo::projectPhiBand(float chargeHypothesis, const reco::SuperClusterRef& superclusterIn,
479  const TrackerTopology *tTopo)
480 {
481  // This algorithm projects a phi band into the tracker three times:
482  // (a) for all stereo hits, (b) for barrel rphi hits, and (c) for
483  // endcap zphi hits. While accumulating stereo hits in step (a),
484  // we fit r vs z to a line. This resolves the ambiguity in z for
485  // rphi hits and the ambiguity in r for zphi hits. We can then cut
486  // on the z of rphi hits (a little wider than one strip length),
487  // and we can convert the z of zphi hits into r to apply the phi
488  // band cut. (We don't take advantage of the endcap strips'
489  // segmentation in r.)
490  //
491  // As we project a phi band into the tracker, we count hits within
492  // that band and performs a linear fit for phi vs r. The number of
493  // hits and reduced chi^2 from the fit are used to select a good
494  // candidate.
495 
496  // set limits on Si hits - smallest error that makes sense = 1 micron ?
497  static const double rphiHitSmallestError = 0.0001 ;
498  static const double stereoHitSmallestError = 0.0001 ;
499  //not used since endcap is not implemented
500  // static const double zphiHitSmallestError = 0.0001 ;
501 
502 
503  static const double stereoPitchAngle = 0.100 ; // stereo angle in rad
504  static const double cos2SiPitchAngle = cos(stereoPitchAngle)*cos(stereoPitchAngle) ;
505  static const double sin2SiPitchAngle = sin(stereoPitchAngle)*sin(stereoPitchAngle) ;
506  // overall misalignment fudge to be added in quad to position errors.
507  // this is a rough approx to values reported in tracking meet 5/16/2007
508  static const double rphiErrFudge = 0.0200 ;
509  static const double stereoErrFudge = 0.0200 ;
510 
511  // max chi2 of a hit on an SiDet relative to the prediction
512  static const double chi2HitMax = 25.0 ;
513 
514  // Minimum number of hits to consider on a candidate
515  static const unsigned int nHitsLeftMinimum = 3 ;
516 
517  // Create and fill vectors of pointers to hits
518  std::vector<const SiStripRecHit2D*> stereoHits;
519  std::vector<const SiStripRecHit2D*> rphiBarrelHits;
520  std::vector<const SiStripRecHit2D*> zphiEndcapHits;
521 
522  // stereo? endcap?
523  coarseHitSelection(stereoHits, tTopo, true, false);
524 
525  // skip endcap stereo for now
526  // LogDebug("") << " Getting endcap stereo hits " ;
527  // coarseHitSelection(stereoHits, true, true);
528 
529  std::ostringstream debugstr1;
530  debugstr1 << " Getting barrel rphi hits " << " \n" ;
531 
532  coarseHitSelection(rphiBarrelHits, tTopo, false, false);
533 
534  // LogDebug("") << " Getting endcap zphi hits " ;
535  // coarseHitSelection(zphiEndcapHits, false, true);
536 
537  debugstr1 << " Getting matched hits " << " \n" ;
538  std::vector<const SiStripMatchedRecHit2D*> matchedHits;
539  coarseMatchedHitSelection(matchedHits);
540 
541 
542  // Determine how to project from the supercluster into the tracker
543  double energy = superclusterIn->energy();
544  double pT = energy * superclusterIn->position().rho()/sqrt(superclusterIn->x()*superclusterIn->x() +
545  superclusterIn->y()*superclusterIn->y() +
546  superclusterIn->z()*superclusterIn->z());
547  // cf Jackson p. 581-2, a little geometry
548  double phiVsRSlope = -3.00e-3 * chargeHypothesis * magneticField_p_->inTesla(GlobalPoint(superclusterIn->x(), superclusterIn->y(), 0.)).z() / pT / 2.;
549 
550  // Shorthand for supercluster radius, z
551  const double scr = superclusterIn->position().rho();
552  const double scz = superclusterIn->position().z();
553 
554  // These are used to fit all hits to a line in phi(r)
555  std::vector<bool> uselist;
556  std::vector<double> rlist, philist, w2list;
557  std::vector<int> typelist; // stereo = 0, rphi barrel = 1, and zphi disk = 2 (only used in this function)
558  std::vector<const SiStripRecHit2D*> hitlist;
559 
560  std::vector<bool> matcheduselist;
561  std::vector<const SiStripMatchedRecHit2D*> matchedhitlist;
562 
563  // These are used to fit the stereo hits to a line in z(r), constrained to pass through supercluster
564  double zSlopeFitNumer = 0.;
565  double zSlopeFitDenom = 0.;
566 
567 
568  debugstr1 << " There are a total of " << stereoHits.size() << " stereoHits in this event " << " \n"
569  << " There are a total of " << rphiBarrelHits.size() << " rphiBarrelHits in this event " << " \n"
570  << " There are a total of " << zphiEndcapHits.size() << " zphiEndcapHits in this event " << " \n\n";
571 
572 
573  LogDebug("") << debugstr1.str() ;
574 
576 
577 
578  // Loop over all matched hits
579  // make a list of good matched rechits.
580  // in the stereo and rphi loops check to see if the hit is associated with a matchedhit
581  LogDebug("") << " Loop over matched hits " << " \n";
582 
583  for (std::vector<const SiStripMatchedRecHit2D*>::const_iterator hit = matchedHits.begin() ;
584  hit != matchedHits.end() ; ++ hit) {
585 
586  assert(matchedHitUsed_.find(*hit) != matchedHitUsed_.end());
587 
588  if (!matchedHitUsed_[*hit]) {
589 
590  // Calculate the 3-D position of this hit
591  GlobalPoint position = tracker_p_->idToDet((*hit)->geographicalId())->surface().toGlobal((*hit)->localPosition());
592  double r = sqrt(position.x()*position.x() + position.y()*position.y());
593  double phi = unwrapPhi(position.phi() - superclusterIn->position().phi()); // phi is relative to supercluster
594  double z = position.z();
595 
596  // Cut a triangle in the z-r plane using the supercluster's eta/dip angle information
597  // and the fact that the electron originated *near* the origin
598  if ((-originUncertainty_ + (scz + originUncertainty_)*(r/scr)) < z && z < (originUncertainty_ + (scz - originUncertainty_)*(r/scr))) {
599 
600  // Cut a narrow band around the supercluster's projection in phi
601  if (unwrapPhi((r-scr)*phiVsRSlope - phiBandWidth_) < phi && phi < unwrapPhi((r-scr)*phiVsRSlope + phiBandWidth_)) {
602 
603 
604  matcheduselist.push_back(true);
605  matchedhitlist.push_back(*hit);
606 
607 
608 
609  // Use this hit to fit z(r)
610  zSlopeFitNumer += -(scr - r) * (z - scz);
611  zSlopeFitDenom += (scr - r) * (scr - r);
612 
613  } // end cut on phi band
614  } // end cut on electron originating *near* the origin
615  } // end assign disjoint sets of hits to electrons
616  } // end loop over matched hits
617 
618  // Calculate the linear fit for z(r)
619  double zVsRSlope;
620  if (zSlopeFitDenom > 0.) {
621  zVsRSlope = zSlopeFitNumer / zSlopeFitDenom;
622  }
623  else {
624  // zVsRSlope assumes electron is from origin if there were no stereo hits
625  zVsRSlope = scz/scr;
626  }
627 
628  // // Loop over all stereo hits
629  LogDebug("") << " Loop over stereo hits" << " \n";
630 
631  // check if the stereo hit is matched to one of the matched hit
632  unsigned int numberOfStereoHits = 0;
633  for (std::vector<const SiStripRecHit2D*>::const_iterator hit = stereoHits.begin(); hit != stereoHits.end(); ++hit) {
634  assert(hitUsed_.find(*hit) != hitUsed_.end());
635 
636  // Calculate the 3-D position of this hit
637  GlobalPoint position = tracker_p_->idToDet((*hit)->geographicalId())->surface().toGlobal((*hit)->localPosition());
638  double r_stereo = sqrt(position.x()*position.x() + position.y()*position.y());
639  double phi_stereo = unwrapPhi(position.phi() - superclusterIn->position().phi()); // phi is relative to supercluster
640  // double z_stereo = position.z();
641 
642  // stereo is has a pitch of 100 mrad - so consider both components
643  double r_stereo_err = sqrt((*hit)->localPositionError().xx()*cos2SiPitchAngle +
644  (*hit)->localPositionError().yy()*sin2SiPitchAngle ) ;
645 
646 
647 
648  // make sure that the error for this hit is sensible, ie > 1 micron
649  // otherwise skip this hit
650  if(r_stereo_err > stereoHitSmallestError ) {
651  r_stereo_err = sqrt(r_stereo_err*r_stereo_err+stereoErrFudge*stereoErrFudge);
652 
653  OmniClusterRef const & stereocluster=(*hit)->omniClusterRef();
654 
655  bool thisHitIsMatched = false ;
656 
657  if (!hitUsed_[*hit]) {
658 
659  unsigned int matcheduselist_size = matcheduselist.size();
660  for (unsigned int i = 0; i < matcheduselist_size; i++) {
661  if (matcheduselist[i]) {
662  OmniClusterRef const & mystereocluster = matchedhitlist[i]->stereoClusterRef();
663  if( stereocluster == mystereocluster ) {
664  thisHitIsMatched = true ;
665  // LogDebug("")<< " This hit is matched " << tracker_p_->idToDet(matchedhitlist[i]->stereoHit()->geographicalId())->surface().toGlobal(matchedhitlist[i]->stereoHit()->localPosition()) << std::endl;
666  // break ;
667  }
668  } // check if matcheduselist okay
669  }// loop over matched hits
670 
671  if(thisHitIsMatched) {
672  // Use this hit to fit phi(r)
673  uselist.push_back(true);
674  rlist.push_back(r_stereo);
675  philist.push_back(phi_stereo);
676  w2list.push_back(1./(r_stereo_err/r_stereo)/(r_stereo_err/r_stereo)); // weight**2 == 1./uncertainty**2
677  typelist.push_back(0);
678  hitlist.push_back(*hit);
679  } // thisHitIsMatched
680  } // if(!hitUsed)
681 
682  } // end of check on hit position error size
683 
684  } // end loop over stereo hits
685 
686  LogDebug("") << " There are " << uselist.size() << " good hits after stereo loop " ;
687 
688 
689  // Loop over barrel rphi hits
690  LogDebug("") << " Looping over barrel rphi hits " ;
691  unsigned int rphiMatchedCounter = 0 ;
692  unsigned int rphiUnMatchedCounter = 0 ;
693  unsigned int numberOfBarrelRphiHits = 0;
694  for (std::vector<const SiStripRecHit2D*>::const_iterator hit = rphiBarrelHits.begin(); hit != rphiBarrelHits.end(); ++hit) {
695  assert(hitUsed_.find(*hit) != hitUsed_.end());
696  // Calculate the 2.5-D position of this hit
697  GlobalPoint position = tracker_p_->idToDet((*hit)->geographicalId())->surface().toGlobal((*hit)->localPosition());
698  double r = sqrt(position.x()*position.x() + position.y()*position.y());
699  double phi = unwrapPhi(position.phi() - superclusterIn->position().phi()); // phi is relative to supercluster
700  double z = position.z();
701  double r_mono_err = sqrt((*hit)->localPositionError().xx()) ;
702 
703  // only consider hits with errors that make sense
704  if( r_mono_err > rphiHitSmallestError) {
705  // inflate the reported error
706  r_mono_err=sqrt(r_mono_err*r_mono_err+rphiErrFudge*rphiErrFudge);
707 
708  OmniClusterRef const & monocluster=(*hit)->omniClusterRef();
709 
710 
711  if (!hitUsed_[*hit]) {
712  if( (tracker_p_->idToDetUnit((*hit)->geographicalId())->type().subDetector() == GeomDetEnumerators::TIB &&
713  (tTopo->tibLayer((*hit)->geographicalId())==1 || tTopo->tibLayer((*hit)->geographicalId())==2)) ||
714  (tracker_p_->idToDetUnit((*hit)->geographicalId())->type().subDetector() == GeomDetEnumerators::TOB
715  && (tTopo->tobLayer((*hit)->geographicalId())==1 || tTopo->tobLayer((*hit)->geographicalId())==2)) ) {
716  bool thisHitIsMatched = false ;
717  unsigned int matcheduselist_size = matcheduselist.size();
718  for (unsigned int i = 0; i < matcheduselist_size; i++) {
719  if (matcheduselist[i]) {
720  OmniClusterRef const & mymonocluster = matchedhitlist[i]->monoClusterRef();
721  if( monocluster == mymonocluster ) {
722  thisHitIsMatched = true ;
723  }
724  } // check if matcheduselist okay
725  }// loop over matched hits
726 
727 
728  if( thisHitIsMatched ) {
729  // Use this hit to fit phi(r)
730  uselist.push_back(true);
731  rlist.push_back(r);
732  philist.push_back(phi);
733  w2list.push_back(1./(r_mono_err/r)/(r_mono_err/r)); // weight**2 == 1./uncertainty**2
734  typelist.push_back(1);
735  hitlist.push_back(*hit);
736  rphiMatchedCounter++;
737  } // end of matched hit check
738 
739  } else {
740 
741 
742  // The expected z position of this hit, according to the z(r) fit
743  double zFit = zVsRSlope * (r - scr) + scz;
744 
745  // Cut on the Z of the strip
746  // TIB strips are 11 cm long, TOB strips are 19 cm long (can I get these from a function?)
747  if ((tracker_p_->idToDetUnit((*hit)->geographicalId())->type().subDetector() == GeomDetEnumerators::TIB &&
748  std::abs(z - zFit) < 12.) ||
749  (tracker_p_->idToDetUnit((*hit)->geographicalId())->type().subDetector() == GeomDetEnumerators::TOB &&
750  std::abs(z - zFit) < 20.) ) {
751 
752  // Cut a narrow band around the supercluster's projection in phi
753  if (unwrapPhi((r-scr)*phiVsRSlope - phiBandWidth_) < phi && phi < unwrapPhi((r-scr)*phiVsRSlope + phiBandWidth_)) {
754 
755  // Use this hit to fit phi(r)
756  uselist.push_back(true);
757  rlist.push_back(r);
758  philist.push_back(phi);
759  w2list.push_back(1./(r_mono_err/r)/(r_mono_err/r)); // weight**2 == 1./uncertainty**2
760  typelist.push_back(1);
761  hitlist.push_back(*hit);
762  rphiUnMatchedCounter++;
763 
764  } // end cut on phi band
765  } // end cut on strip z
766  } // loop over TIB/TOB layer 1,2
767  } // end assign disjoint sets of hits to electrons
768  } // end of check on rphi hit position error size
769  } // end loop over barrel rphi hits
770 
771  LogDebug("") << " There are " << rphiMatchedCounter <<" matched rphi hits";
772  LogDebug("") << " There are " << rphiUnMatchedCounter <<" unmatched rphi hits";
773  LogDebug("") << " There are " << uselist.size() << " good stereo+rphi hits " ;
774 
775 
776 
777 
778 
779 
781 
782  // Loop over endcap zphi hits
783  LogDebug("") << " Looping over barrel zphi hits " ;
784 
785 
786  unsigned int numberOfEndcapZphiHits = 0;
787  for (std::vector<const SiStripRecHit2D*>::const_iterator hit = zphiEndcapHits.begin();
788  hit != zphiEndcapHits.end(); ++hit) {
789  assert(hitUsed_.find(*hit) != hitUsed_.end());
790  if (!hitUsed_[*hit]) {
791  // Calculate the 2.5-D position of this hit
792  GlobalPoint position = tracker_p_->idToDet((*hit)->geographicalId())->surface().toGlobal((*hit)->localPosition());
793  double z = position.z();
794  double phi = unwrapPhi(position.phi() - superclusterIn->position().phi()); // phi is relative to supercluster
795  // double r=sqrt(position.x()*position.x()+position.y()*position.y()) ;
796 
797  // The expected r position of this hit, according to the z(r) fit
798  double rFit = (z - scz)/zVsRSlope + scr;
799 
800  // I don't know the r widths of the endcap strips, otherwise I
801  // would apply a cut on r similar to the rphi z cut
802 
803  // Cut a narrow band around the supercluster's projection in phi,
804  // and be sure the hit isn't in a reflected band (extrapolation of large z differences into negative r)
805  if (rFit > 0. &&
806  unwrapPhi((rFit-scr)*phiVsRSlope - phiBandWidth_) < phi && phi < unwrapPhi((rFit-scr)*phiVsRSlope + phiBandWidth_)) {
807 
808  // Use this hit to fit phi(r)
809  uselist.push_back(true);
810  rlist.push_back(rFit);
811  philist.push_back(phi);
812  w2list.push_back(1./(0.05/rFit)/(0.05/rFit)); // weight**2 == 1./uncertainty**2
813  typelist.push_back(2);
814  hitlist.push_back(*hit);
815 
816  } // end cut on phi band
817  } // end assign disjoint sets of hits to electrons
818  } // end loop over endcap zphi hits
819 
820  LogDebug("") << " There are " << uselist.size() << " good stereo+rphi+zphi hits " ;
822 
823 #ifdef EDM_ML_DEBUG
824  std::ostringstream debugstr5;
825  debugstr5 << " List of hits before uniqification " << " \n" ;
826  for (unsigned int i = 0; i < uselist.size(); i++) {
827  if ( uselist[i] ) {
828  const SiStripRecHit2D* hit = hitlist[i];
829  debugstr5 << " DetID " << ((hit)->geographicalId()).rawId()
830  << " R " << rlist[i]
831  << " Phi " << philist[i]
832  << " Weight " << w2list[i]
833  << " PhiPred " << (rlist[i]-scr)*phiVsRSlope
834  << " Chi2 " << (philist[i]-(rlist[i]-scr)*phiVsRSlope)*(philist[i]-(rlist[i]-scr)*phiVsRSlope)*w2list[i]
835  << " \n" ;
836  }
837  }
838  debugstr5 << " \n\n\n" ;
839 
840  debugstr5 << " Counting number of unique detectors " << " \n" ;
841 
842  debugstr5 << " These are all the detectors with hits " << " \n";
843 #endif
844  // Loop over hits, and find the best hit for each SiDetUnit - keep only those
845  // get a list of DetIds in uselist
846  std::vector<uint32_t> detIdList ;
847 
848  for (unsigned int i = 0; i < uselist.size(); i++) {
849  if (uselist[i]) {
850  const SiStripRecHit2D* hit = hitlist[i];
851  uint32_t detIDUsed = ((hit)->geographicalId()).rawId() ;
852 #ifdef EDM_ML_DEBUG
853  debugstr5<< " DetId " << detIDUsed << "\n";
854 #endif
855  detIdList.push_back(detIDUsed) ;
856  }
857  }
858 
859 #ifdef EDM_ML_DEBUG
860  debugstr5 << " There are " << detIdList.size() << " hits on detectors \n";
861 #endif
862  // now sort and then uniq this list of detId
863  std::sort(detIdList.begin(), detIdList.end());
864  detIdList.erase(
865  std::unique(detIdList.begin(), detIdList.end()), detIdList.end());
866 #ifdef EDM_ML_DEBUG
867  debugstr5 << " There are " << detIdList.size() << " unique detectors \n";
868 #endif
869  //now we have a list of unique detectors used
870 
871 
872 #ifdef EDM_ML_DEBUG
873  debugstr5 << " Looping over detectors to choose best hit " << " \n";
874 #endif
875  // Loop over detectors ID and hits to create list of hits on same DetId
876  for (unsigned int idet = 0 ; idet < detIdList.size() ; idet++ ) {
877  for (unsigned int i = 0; i+1 < uselist.size(); i++) {
878  if (uselist[i]) {
879  // Get Chi2 of this hit relative to predicted hit
880  const SiStripRecHit2D* hit1 = hitlist[i];
881  double r_hit1 = rlist[i];
882  double phi_hit1 = philist[i];
883  double w2_hit1 = w2list[i] ;
884  double phi_pred1 = (r_hit1-scr)*phiVsRSlope ;
885  double chi1 = (phi_hit1-phi_pred1)*(phi_hit1-phi_pred1)*w2_hit1;
886  if(detIdList[idet]== ((hit1)->geographicalId()).rawId()) {
887  for (unsigned int j = i+1; j < uselist.size(); j++) {
888  if (uselist[j]) {
889  const SiStripRecHit2D* hit2 = hitlist[j];
890  if(detIdList[idet]== ((hit2)->geographicalId()).rawId()) {
891 #ifdef EDM_ML_DEBUG
892  debugstr5 << " Found 2 hits on same Si Detector "
893  << ((hit2)->geographicalId()).rawId() << "\n";
894 #endif
895  // Get Chi2 of this hit relative to predicted hit
896  double r_hit2 = rlist[j];
897  double phi_hit2 = philist[j];
898  double w2_hit2 = w2list[j] ;
899  double phi_pred2 = (r_hit2-scr)*phiVsRSlope ;
900  double chi2 = (phi_hit2-phi_pred2)*(phi_hit2-phi_pred2)*w2_hit2;
901 #ifdef EDM_ML_DEBUG
902  debugstr5 << " Chi1 " << chi1 << " Chi2 " << chi2 <<"\n";
903 #endif
904  if( chi1< chi2 ){
905  uselist[j] = 0;
906  }else{
907  uselist[i] = 0;
908  }
909 
910  } // end of Det check
911  } // end of j- usehit check
912  } // end of j hit list loop
913 
914  } // end of detector check
915  } // end of uselist check
916  } // end of i-hit list loop
917 
918 
919  } // end of DetId loop
920 
921 
922 
923  // now let's through out hits with a predicted chi > chi2HitMax
924  for ( unsigned int i = 0; i < uselist.size(); i++ ) {
925  if ( uselist[i] ) {
926  double localchi2 = (philist[i]-(rlist[i]-scr)*phiVsRSlope)*(philist[i]-(rlist[i]-scr)*phiVsRSlope)*w2list[i] ;
927  if(localchi2 > chi2HitMax ) {
928 #ifdef EDM_ML_DEBUG
929  const SiStripRecHit2D* hit = hitlist[i];
930  debugstr5 << " Throwing out "
931  <<" DetID " << ((hit)->geographicalId()).rawId()
932  << " R " << rlist[i]
933  << " Phi " << philist[i]
934  << " Weight " << w2list[i]
935  << " PhiPred " << (rlist[i]-scr)*phiVsRSlope
936  << " Chi2 " << (philist[i]-(rlist[i]-scr)*phiVsRSlope)*(philist[i]-(rlist[i]-scr)*phiVsRSlope)*w2list[i]
937  << " \n" ;
938 #endif
939  uselist[i]=0 ;
940  }
941  }
942  }
943 
944 #ifdef EDM_ML_DEBUG
945  debugstr5 << " List of hits after uniqification " << " \n" ;
946  for (unsigned int i = 0; i < uselist.size(); i++) {
947  if ( uselist[i] ) {
948  const SiStripRecHit2D* hit = hitlist[i];
949  debugstr5 << " DetID " << ((hit)->geographicalId()).rawId()
950  << " R " << rlist[i]
951  << " Phi " << philist[i]
952  << " Weight " << w2list[i]
953  << " PhiPred " << (rlist[i]-scr)*phiVsRSlope
954  << " Chi2 " << (philist[i]-(rlist[i]-scr)*phiVsRSlope)*(philist[i]-(rlist[i]-scr)*phiVsRSlope)*w2list[i]
955  << " \n" ;
956  }
957  }
958  debugstr5 << " \n\n\n" ;
959 #endif
960 
961  // need to check that there are more than 2 hits left here!
962  unsigned int nHitsLeft =0;
963  for (unsigned int i = 0; i < uselist.size(); i++) {
964  if ( uselist[i] ) {
965  nHitsLeft++;
966  }
967  }
968  if(nHitsLeft < nHitsLeftMinimum ) {
969 #ifdef EDM_ML_DEBUG
970  debugstr5 << " Too few hits "<< nHitsLeft
971  << " left - return false " << " \n";
972 #endif
973  return false ;
974  }
975 #ifdef EDM_ML_DEBUG
976  LogDebug("") << debugstr5.str();
977 #endif
978 
980  // Calculate a linear phi(r) fit and drop hits until the biggest contributor to chi^2 is less than maxNormResid_
981  bool done = false;
982  double intercept = 0 ;
983  double slope = 0 ;
984  double chi2 = 0;
985 
986  std::ostringstream debugstr4;
987  debugstr4 <<" Calc of phi(r) "<<" \n";
988 
989  while (!done) {
990  // Use an iterative update of <psi>, <r> etc instead in an
991  // attempt to minize divisions of large numbers
992  // The linear fit psi = slope*r + intercept
993  // slope is (<psi*r>-<psi>*<r>) / (<r^2>-<r>^2>)
994  // intercept is <psi> - slope*<r>
995 
996  double phiBar = 0.;
997  double phiBarOld = 0.;
998  double rBar = 0.;
999  double rBarOld = 0.;
1000  double r2Bar = 0.;
1001  double r2BarOld = 0.;
1002  double phirBar = 0.;
1003  double phirBarOld = 0.;
1004  double totalWeight = 0.;
1005  double totalWeightOld = 0.;
1006  unsigned int uselist_size = uselist.size();
1007  for (unsigned int i = 0; i < uselist_size; i++) {
1008  if (uselist[i]) {
1009  double r = rlist[i];
1010  double phi = philist[i];
1011  double weight2 = w2list[i];
1012 
1013  totalWeight = totalWeightOld + weight2 ;
1014 
1015  //weight2 is 1/sigma^2. Typically sigma is 100micron pitch
1016  // over root(12) = 30 microns so weight2 is about 10^5-10^6
1017 
1018  double totalWeightRatio = totalWeightOld/totalWeight ;
1019  double localWeightRatio = weight2/totalWeight ;
1020 
1021  phiBar = phiBarOld*totalWeightRatio + phi*localWeightRatio ;
1022  rBar = rBarOld*totalWeightRatio + r*localWeightRatio ;
1023  r2Bar= r2BarOld*totalWeightRatio + r*r*localWeightRatio ;
1024  phirBar = phirBarOld*totalWeightRatio + phi*r*localWeightRatio ;
1025 
1026  totalWeightOld = totalWeight ;
1027  phiBarOld = phiBar ;
1028  rBarOld = rBar ;
1029  r2BarOld = r2Bar ;
1030  phirBarOld = phirBar ;
1031 #ifdef EDM_ML_DEBUG
1032  debugstr4 << " totalWeight " << totalWeight
1033  << " totalWeightRatio " << totalWeightRatio
1034  << " localWeightRatio "<< localWeightRatio
1035  << " phiBar " << phiBar
1036  << " rBar " << rBar
1037  << " r2Bar " << r2Bar
1038  << " phirBar " << phirBar
1039  << " \n ";
1040 #endif
1041 
1042  } // end of use hit loop
1043  } // end of hit loop to calculate a linear fit
1044  slope = (phirBar-phiBar*rBar)/(r2Bar-rBar*rBar);
1045  intercept = phiBar-slope*rBar ;
1046 
1047  debugstr4 << " end of loop slope " << slope
1048  << " intercept " << intercept << " \n";
1049 
1050  // Calculate chi^2
1051  chi2 = 0.;
1052  double biggest_normresid = -1.;
1053  unsigned int biggest_index = 0;
1054  for (unsigned int i = 0; i < uselist_size; i++) {
1055  if (uselist[i]) {
1056  double r = rlist[i];
1057  double phi = philist[i];
1058  double weight2 = w2list[i];
1059 
1060  double normresid = weight2 * (phi - intercept - slope*r)*(phi - intercept - slope*r);
1061  chi2 += normresid;
1062 
1063  if (normresid > biggest_normresid) {
1064  biggest_normresid = normresid;
1065  biggest_index = i;
1066  }
1067  }
1068  } // end loop over hits to calculate chi^2 and find its biggest contributer
1069 
1070  if (biggest_normresid > maxNormResid_) {
1071 #ifdef EDM_ML_DEBUG
1072  debugstr4 << "Dropping hit from fit due to Chi2 " << " \n" ;
1073  const SiStripRecHit2D* hit = hitlist[biggest_index];
1074  debugstr4 << " DetID " << ((hit)->geographicalId()).rawId()
1075  << " R " << rlist[biggest_index]
1076  << " Phi " << philist[biggest_index]
1077  << " Weight " << w2list[biggest_index]
1078  << " PhiPred " << (rlist[biggest_index]-scr)*phiVsRSlope
1079  << " Chi2 " << (philist[biggest_index]-(rlist[biggest_index]-scr)*phiVsRSlope)*(philist[biggest_index]-(rlist[biggest_index]-scr)*phiVsRSlope)*w2list[biggest_index]
1080  << " normresid " << biggest_normresid
1081  << " \n\n";
1082 #endif
1083  uselist[biggest_index] = false;
1084  }
1085  else {
1086  done = true;
1087  }
1088  } // end loop over trial fits
1089 #ifdef EDM_ML_DEBUG
1090  debugstr4 <<" Fit "
1091  << " Intercept " << intercept
1092  << " slope " << slope
1093  << " chi2 " << chi2
1094  << " \n" ;
1095 
1096  LogDebug("") << debugstr4.str() ;
1097 #endif
1098 
1099  // Now we have intercept, slope, and chi2; uselist to tell us which hits are used, and hitlist for the hits
1100 
1101  // Identify the innermost hit
1102  const SiStripRecHit2D* innerhit = (SiStripRecHit2D*)(0);
1103  double innerhitRadius = -1.; // meaningless until innerhit is defined
1104 
1105  // Copy hits into an OwnVector, which we put in the TrackCandidate
1106  std::vector<const TrackingRecHit*> outputHits;
1107  // Reference rphi and stereo hits into RefVectors, which we put in the SiStripElectron
1108  std::vector<SiStripRecHit2D> outputRphiHits;
1109  std::vector<SiStripRecHit2D> outputStereoHits;
1110 
1111  typedef edm::Ref<SiStripRecHit2DCollection,SiStripRecHit2D> SiStripRecHit2DRef;
1112 
1113 
1114  for (unsigned int i = 0; i < uselist.size(); i++) {
1115  if (uselist[i]) {
1116  double r = rlist[i];
1117  const SiStripRecHit2D* hit = hitlist[i];
1118 
1119  // Keep track of the innermost hit
1120  if (innerhit == (SiStripRecHit2D*)(0) || r < innerhitRadius) {
1121  innerhit = hit;
1122  innerhitRadius = r;
1123  }
1124 
1125  if (typelist[i] == 0) {
1126  numberOfStereoHits++;
1127 
1128  // Copy this hit for the TrajectorySeed
1129  outputHits.push_back(hit);
1130  outputStereoHits.push_back(*hit);
1131  }
1132  else if (typelist[i] == 1) {
1133  numberOfBarrelRphiHits++;
1134 
1135  // Copy this hit for the TrajectorySeed
1136  outputHits.push_back(hit);
1137  outputRphiHits.push_back(*hit);
1138  }
1139  else if (typelist[i] == 2) {
1140  numberOfEndcapZphiHits++;
1141 
1142  // Copy this hit for the TrajectorySeed
1143  outputHits.push_back(hit);
1144  outputRphiHits.push_back(*hit);
1145  }
1146  }
1147  } // end loop over all hits, after having culled the ones with big residuals
1148 
1149  unsigned int totalNumberOfHits = numberOfStereoHits + numberOfBarrelRphiHits + numberOfEndcapZphiHits;
1150  double reducedChi2 = (totalNumberOfHits > 2 ? chi2 / (totalNumberOfHits - 2) : 1e10);
1151 
1152  // Select this candidate if it passes minHits_ and maxReducedChi2_ cuts
1153  if (totalNumberOfHits >= minHits_ && reducedChi2 <= maxReducedChi2_) {
1154  // GlobalTrajectoryParameters evaluated at the position of the innerhit
1155  GlobalPoint position = tracker_p_->idToDet(innerhit->geographicalId())->surface().toGlobal(innerhit->localPosition());
1156 
1157  // Use our phi(r) linear fit to correct pT (pT is inversely proportional to slope)
1158  // (By applying a correction instead of going back to first
1159  // principles, any correction to the phiVsRSlope definition
1160  // above will be automatically propagated here.)
1161  double correct_pT = pT * phiVsRSlope / slope;
1162 
1163  // Our phi(r) linear fit returns phi relative to the supercluster phi for a given radius
1164  double phi = intercept + slope*sqrt(position.x()*position.x() + position.y()*position.y()) + superclusterIn->position().phi();
1165 
1166  // Our z(r) linear fit gives us a better measurement of eta/dip angle
1167  double pZ = correct_pT * zVsRSlope;
1168 
1169  // Now we can construct a trajectory momentum out of linear fits to hits
1170  GlobalVector momentum = GlobalVector(correct_pT*cos(phi), correct_pT*sin(phi), pZ);
1171 
1172  if (chargeHypothesis > 0.) {
1173  redchi2_pos_ = chi2 / double(totalNumberOfHits - 2);
1175  momentum_pos_ = momentum;
1176  innerhit_pos_ = innerhit;
1177  outputHits_pos_ = outputHits;
1178  outputRphiHits_pos_ = outputRphiHits;
1179  outputStereoHits_pos_ = outputStereoHits;
1180  phiVsRSlope_pos_ = phiVsRSlope;
1181  slope_pos_ = slope;
1182  intercept_pos_ = intercept;
1183  chi2_pos_ = chi2;
1184  ndof_pos_ = totalNumberOfHits - 2;
1185  correct_pT_pos_ = correct_pT;
1186  pZ_pos_ = pZ;
1187  zVsRSlope_pos_ = zVsRSlope;
1188  numberOfStereoHits_pos_ = numberOfStereoHits;
1189  numberOfBarrelRphiHits_pos_ = numberOfBarrelRphiHits;
1190  numberOfEndcapZphiHits_pos_ = numberOfEndcapZphiHits;
1191  }
1192  else {
1193  redchi2_neg_ = chi2 / double(totalNumberOfHits - 2);
1195  momentum_neg_ = momentum;
1196  innerhit_neg_ = innerhit;
1197  outputHits_neg_ = outputHits;
1198  outputRphiHits_neg_ = outputRphiHits;
1199  outputStereoHits_neg_ = outputStereoHits;
1200  phiVsRSlope_neg_ = phiVsRSlope;
1201  slope_neg_ = slope;
1202  intercept_neg_ = intercept;
1203  chi2_neg_ = chi2;
1204  ndof_neg_ = totalNumberOfHits - 2;
1205  correct_pT_neg_ = correct_pT;
1206  pZ_neg_ = pZ;
1207  zVsRSlope_neg_ = zVsRSlope;
1208  numberOfStereoHits_neg_ = numberOfStereoHits;
1209  numberOfBarrelRphiHits_neg_ = numberOfBarrelRphiHits;
1210  numberOfEndcapZphiHits_neg_ = numberOfEndcapZphiHits;
1211  }
1212 
1213  LogDebug("") << " return from projectFindBand with True \n" << std::endl ;
1214  return true;
1215  } // end if this is a good electron candidate
1216 
1217  // Signal for a failed electron candidate
1218  LogDebug("") << " return from projectFindBand with False \n" << std::endl ;
1219  return false;
1220 }
1221 
1222 //
1223 // const member functions
1224 //
1225 
1226 //
1227 // static member functions
1228 //
1229 
#define LogDebug(id)
type
Definition: HCALResponse.h:21
int i
Definition: DBlmapReader.cc:9
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
const_iterator begin() const
std::vector< const TrackingRecHit * > outputHits_pos_
unsigned int tibLayer(const DetId &id) const
std::vector< SiStripRecHit2D > outputRphiHits_pos_
const edm::Handle< SiStripRecHit2DCollection > * stereoHits_hp_
const SiStripRecHit2D * innerhit_neg_
virtual GlobalVector inTesla(const GlobalPoint &gp) const =0
Field value ad specified global point, in Tesla.
static const double slope[3]
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:47
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
std::vector< TrackCandidate > TrackCandidateCollection
uint32_t tobStereo(const DetId &id) const
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
PTrajectoryStateOnDet persistentState(const TrajectoryStateOnSurface &ts, unsigned int detid)
T y() const
Definition: PV3DBase.h:63
#define abs(x)
Definition: mlp_lapack.h:159
unsigned int tidWheel(const DetId &id) const
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55
const SiStripRecHit2D * innerhit_pos_
float float float z
static int position[TOTALCHAMBERS][3]
Definition: ReadPGInfo.cc:509
void push_back(D *&d)
Definition: OwnVector.h:273
std::vector< SiStripRecHit2D > outputStereoHits_pos_
uint32_t tidStereo(const DetId &id) const
std::map< const TrackingRecHit *, bool > matchedHitUsed_
std::vector< SiStripElectron > SiStripElectronCollection
collectin of SiStripElectron objects
T sqrt(T t)
Definition: SSEVec.h:48
T z() const
Definition: PV3DBase.h:64
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
data_type const * data(size_t cell) const
const SiStripRecHit2DCollection * rphiHits_p_
int j
Definition: DBlmapReader.cc:9
double unwrapPhi(double phi) const
std::map< const SiStripRecHit2D *, unsigned int > stereoKey_
const_iterator end() const
virtual const GeomDet * idToDet(DetId) const
std::vector< SiStripRecHit2D > outputRphiHits_neg_
const SiStripRecHit2DCollection * stereoHits_p_
bool findElectron(reco::SiStripElectronCollection &electronOut, TrackCandidateCollection &trackCandidateOut, const reco::SuperClusterRef &superclusterIn, const TrackerTopology *tTopo)
unsigned int numberOfStereoHits_pos_
Definition: DetId.h:20
bool projectPhiBand(float chargeHypothesis, const reco::SuperClusterRef &superclusterIn, const TrackerTopology *tTopo)
unsigned int numberOfBarrelRphiHits_pos_
std::vector< SiStripRecHit2D > outputStereoHits_neg_
void coarseHitSelection(std::vector< const SiStripRecHit2D * > &hitPointersOut, const TrackerTopology *tTopo, bool stereo, bool endcap)
const TrackerGeometry * tracker_p_
T const * product() const
Definition: ESHandle.h:62
const edm::Handle< SiStripRecHit2DCollection > * rphiHits_hp_
id_type detId() const
Definition: DetSetNew.h:72
const SiStripMatchedRecHit2DCollection * matchedHits_p_
T const * product() const
Definition: Handle.h:74
char state
Definition: procUtils.cc:75
virtual const GeomDetUnit * idToDetUnit(DetId) const
Return the pointer to the GeomDetUnit corresponding to a given DetId.
std::map< const SiStripMatchedRecHit2D *, unsigned int > matchedKey_
T perp() const
Magnitude of transverse component.
iterator end()
Definition: DetSetNew.h:59
uint32_t tecStereo(const DetId &id) const
std::map< const TrackingRecHit *, bool > hitUsed_
void coarseMatchedHitSelection(std::vector< const SiStripMatchedRecHit2D * > &coarseMatchedHitPointersOut)
const MagneticField * magneticField_p_
std::map< const SiStripRecHit2D *, unsigned int > rphiKey_
uint32_t tibStereo(const DetId &id) const
size_type size() const
Definition: DetSetNew.h:75
SiStripElectronAlgo(unsigned int maxHitsOnDetId, double originUncertainty, double phiBandWidth, double maxNormResid, unsigned int minHits, double maxReducedChi2)
T x() const
Definition: PV3DBase.h:62
unsigned int tecWheel(const DetId &id) const
unsigned int numberOfEndcapZphiHits_pos_
unsigned int tobLayer(const DetId &id) const
void prepareEvent(const edm::ESHandle< TrackerGeometry > &tracker, const edm::Handle< SiStripRecHit2DCollection > &rphiHits, const edm::Handle< SiStripRecHit2DCollection > &stereoHits, const edm::Handle< SiStripMatchedRecHit2DCollection > &matchedHits, const edm::ESHandle< MagneticField > &magneticField)
Global3DVector GlobalVector
Definition: GlobalVector.h:10
std::vector< const TrackingRecHit * > outputHits_neg_
Definition: DDAxes.h:10
iterator begin()
Definition: DetSetNew.h:56