CMS 3D CMS Logo

List of all members | Public Member Functions | Public Attributes
GsfElectronAlgo::ElectronData Struct Reference

Public Member Functions

void calculateMode (const MultiTrajectoryStateMode *mtsMode)
 
Candidate::LorentzVector calculateMomentum ()
 
bool calculateTSOS (const MultiTrajectoryStateTransform *, GsfConstraintAtVertex *)
 
void checkCtfTrack (edm::Handle< reco::TrackCollection > currentCtfTracks)
 
void computeCharge (int &charge, reco::GsfElectron::ChargeInfo &info)
 
 ElectronData (const reco::GsfElectronCoreRef &core, const reco::BeamSpot &bs)
 
CaloClusterPtr getEleBasicCluster (const MultiTrajectoryStateTransform *)
 
 ~ElectronData ()
 

Public Attributes

const reco::BeamSpot beamSpot
 
TrajectoryStateOnSurface constrainedVtxTSOS
 
const reco::GsfElectronCoreRef coreRef
 
reco::TrackRef ctfTrackRef
 
GlobalVector eleMom
 
GlobalPoint elePos
 
TrajectoryStateOnSurface eleTSOS
 
const reco::GsfTrackRef gsfTrackRef
 
GlobalVector innMom
 
GlobalPoint innPos
 
TrajectoryStateOnSurface innTSOS
 
GlobalVector outMom
 
GlobalPoint outPos
 
TrajectoryStateOnSurface outTSOS
 
GlobalVector sclMom
 
GlobalPoint sclPos
 
TrajectoryStateOnSurface sclTSOS
 
GlobalVector seedMom
 
GlobalPoint seedPos
 
TrajectoryStateOnSurface seedTSOS
 
float shFracInnerHits
 
const reco::SuperClusterRef superClusterRef
 
GlobalVector vtxMom
 
GlobalVector vtxMomWithConstraint
 
GlobalPoint vtxPos
 
TrajectoryStateOnSurface vtxTSOS
 

Detailed Description

Definition at line 285 of file GsfElectronAlgo.cc.

Constructor & Destructor Documentation

GsfElectronAlgo::ElectronData::ElectronData ( const reco::GsfElectronCoreRef core,
const reco::BeamSpot bs 
)

Definition at line 325 of file GsfElectronAlgo.cc.

327  : coreRef(core),
328  gsfTrackRef(coreRef->gsfTrack()),
329  superClusterRef(coreRef->superCluster()),
330  ctfTrackRef(coreRef->ctfTrack()), shFracInnerHits(coreRef->ctfGsfOverlap()),
331  beamSpot(bs)
332  {}
const reco::GsfTrackRef gsfTrackRef
const reco::GsfElectronCoreRef coreRef
const reco::SuperClusterRef superClusterRef
const reco::BeamSpot beamSpot
GsfElectronAlgo::ElectronData::~ElectronData ( )

Definition at line 334 of file GsfElectronAlgo.cc.

335  {}

Member Function Documentation

void GsfElectronAlgo::ElectronData::calculateMode ( const MultiTrajectoryStateMode mtsMode)

Definition at line 515 of file GsfElectronAlgo.cc.

References MultiTrajectoryStateMode::momentumFromModeCartesian(), and MultiTrajectoryStateMode::positionFromModeCartesian().

Referenced by GsfElectronAlgo::createElectron().

516  {
530  }
bool positionFromModeCartesian(const TrajectoryStateOnSurface tsos, GlobalPoint &position) const
TrajectoryStateOnSurface sclTSOS
bool momentumFromModeCartesian(const TrajectoryStateOnSurface tsos, GlobalVector &momentum) const
TrajectoryStateOnSurface innTSOS
TrajectoryStateOnSurface constrainedVtxTSOS
TrajectoryStateOnSurface seedTSOS
TrajectoryStateOnSurface eleTSOS
TrajectoryStateOnSurface outTSOS
TrajectoryStateOnSurface vtxTSOS
Candidate::LorentzVector GsfElectronAlgo::ElectronData::calculateMomentum ( )

Definition at line 532 of file GsfElectronAlgo.cc.

References Scenarios_cff::scale.

Referenced by GsfElectronAlgo::createElectron().

533  {
534  double scale = superClusterRef->energy()/vtxMom.mag() ;
536  ( vtxMom.x()*scale,vtxMom.y()*scale,vtxMom.z()*scale,
537  superClusterRef->energy() ) ;
538  }
T y() const
Definition: PV3DBase.h:63
T mag() const
Definition: PV3DBase.h:67
T z() const
Definition: PV3DBase.h:64
const reco::SuperClusterRef superClusterRef
T x() const
Definition: PV3DBase.h:62
math::PtEtaPhiELorentzVectorF LorentzVector
bool GsfElectronAlgo::ElectronData::calculateTSOS ( const MultiTrajectoryStateTransform mtsTransform,
GsfConstraintAtVertex constraintAtVtx 
)

Definition at line 481 of file GsfElectronAlgo.cc.

References ecalDrivenElectronSeedsParameters_cff::beamSpot, GsfConstraintAtVertex::constrainAtBeamSpot(), ele_convert(), MultiTrajectoryStateTransform::extrapolatedState(), MultiTrajectoryStateTransform::innerStateOnSurface(), MultiTrajectoryStateTransform::outerStateOnSurface(), and funct::true.

Referenced by GsfElectronAlgo::completeElectrons(), and getEleBasicCluster().

482  {
483  //at innermost point
484  innTSOS = mtsTransform->innerStateOnSurface(*gsfTrackRef);
485  if (!innTSOS.isValid()) return false;
486 
487  //at vertex
488  // innermost state propagation to the beam spot position
489  GlobalPoint bsPos ;
490  ele_convert(beamSpot.position(),bsPos) ;
491  vtxTSOS = mtsTransform->extrapolatedState(innTSOS,bsPos) ;
492  if (!vtxTSOS.isValid()) vtxTSOS=innTSOS;
493 
494  //at seed
495  outTSOS = mtsTransform->outerStateOnSurface(*gsfTrackRef);
496  if (!outTSOS.isValid()) return false;
497 
498  // TrajectoryStateOnSurface seedTSOS
499  seedTSOS = mtsTransform->extrapolatedState(outTSOS,
500  GlobalPoint(superClusterRef->seed()->position().x(),
501  superClusterRef->seed()->position().y(),
502  superClusterRef->seed()->position().z()));
504 
505  // at scl
507  if (!sclTSOS.isValid()) sclTSOS=outTSOS;
508 
509  // constrained momentum
511 
512  return true ;
513  }
TrajectoryStateOnSurface constrainAtBeamSpot(const reco::GsfTrack &, const reco::BeamSpot &) const
(multi)TSOS after including the beamspot
TrajectoryStateOnSurface sclTSOS
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
TrajectoryStateOnSurface innTSOS
TrajectoryStateOnSurface outerStateOnSurface(const reco::GsfTrack &tk) const
TrajectoryStateOnSurface constrainedVtxTSOS
TrajectoryStateOnSurface extrapolatedState(const TrajectoryStateOnSurface tsos, const GlobalPoint &point) const
const reco::GsfTrackRef gsfTrackRef
void ele_convert(const Type1 &obj1, Type2 &obj2)
TrajectoryStateOnSurface innerStateOnSurface(const reco::GsfTrack &tk) const
TrajectoryStateOnSurface seedTSOS
const reco::SuperClusterRef superClusterRef
const Point & position() const
position
Definition: BeamSpot.h:62
TrajectoryStateOnSurface outTSOS
TrajectoryStateOnSurface vtxTSOS
const reco::BeamSpot beamSpot
void GsfElectronAlgo::ElectronData::checkCtfTrack ( edm::Handle< reco::TrackCollection currentCtfTracks)

Definition at line 337 of file GsfElectronAlgo.cc.

References funct::abs(), TrackingRecHit::all, computeCharge(), counter, particleFlow_cfi::dEta, particleFlow_cfi::dPhi, reco::HitPattern::getHitPattern(), reco::HitPattern::hitPattern, min(), pi, edm::Handle< T >::product(), and mathSSE::sqrt().

Referenced by GsfElectronAlgo::createElectron().

338 {
339  if (!ctfTrackRef.isNull()) return ;
340 
341  // Code below from Puneeth Kalavase
342 
343  shFracInnerHits = 0 ;
344  const TrackCollection * ctfTrackCollection = currentCtfTracks.product() ;
345 
346  // get the Hit Pattern for the gsfTrack
347  const HitPattern &gsfHitPattern = gsfTrackRef->hitPattern() ;
348 
349  unsigned int counter ;
350  TrackCollection::const_iterator ctfTkIter ;
351  for (ctfTkIter = ctfTrackCollection->begin(), counter = 0;
352  ctfTkIter != ctfTrackCollection->end(); ctfTkIter++, counter++)
353  {
354  double dEta = gsfTrackRef->eta() - ctfTkIter->eta() ;
355  double dPhi = gsfTrackRef->phi() - ctfTkIter->phi() ;
356  double pi = acos(-1.);
357  if(std::abs(dPhi) > pi) dPhi = 2*pi - std::abs(dPhi) ;
358 
359  // dont want to look at every single track in the event!
360  if (sqrt(dEta*dEta + dPhi*dPhi) > 0.3) continue ;
361 
362  unsigned int shared = 0 ;
363  int gsfHitCounter = 0 ;
364  int numGsfInnerHits = 0 ;
365  int numCtfInnerHits = 0 ;
366 
367  // get the CTF Track Hit Pattern
368  const HitPattern &ctfHitPattern = ctfTkIter->hitPattern() ;
369 
370  trackingRecHit_iterator elHitsIt;
371  for (elHitsIt = gsfTrackRef->recHitsBegin();
372  elHitsIt != gsfTrackRef->recHitsEnd();
373  elHitsIt++, gsfHitCounter++)
374  {
375  if (!((**elHitsIt).isValid())) //count only valid Hits
376  { continue ; }
377 
378  // look only in the pixels/TIB/TID
379  uint32_t gsfHit = gsfHitPattern.getHitPattern(HitPattern::TRACK_HITS, gsfHitCounter) ;
380  if (!(HitPattern::pixelHitFilter(gsfHit)
381  || HitPattern::stripTIBHitFilter(gsfHit)
382  || HitPattern::stripTIDHitFilter(gsfHit))){
383  continue;
384  }
385 
386  numGsfInnerHits++ ;
387 
388  int ctfHitsCounter = 0 ;
389  numCtfInnerHits = 0 ;
390  trackingRecHit_iterator ctfHitsIt ;
391  for (ctfHitsIt = ctfTkIter->recHitsBegin();
392  ctfHitsIt != ctfTkIter->recHitsEnd();
393  ctfHitsIt++, ctfHitsCounter++ )
394  {
395  if(!((**ctfHitsIt).isValid())) //count only valid Hits!
396  { continue; }
397 
398  uint32_t ctfHit = ctfHitPattern.getHitPattern(HitPattern::TRACK_HITS, ctfHitsCounter);
399  if(!(HitPattern::pixelHitFilter(ctfHit)
400  || HitPattern::stripTIBHitFilter(ctfHit)
401  || HitPattern::stripTIDHitFilter(ctfHit)))
402  {
403  continue;
404  }
405 
406  numCtfInnerHits++ ;
407 
408  if((**elHitsIt).sharesInput(&(**ctfHitsIt), TrackingRecHit::all))
409  {
410  shared++ ;
411  break ;
412  }
413 
414  } //ctfHits iterator
415 
416  } //gsfHits iterator
417 
418  if ((numGsfInnerHits==0)||(numCtfInnerHits==0))
419  { continue ; }
420 
421  if ( static_cast<float>(shared)/min(numGsfInnerHits,numCtfInnerHits) > shFracInnerHits )
422  {
423  shFracInnerHits = static_cast<float>(shared)/min(numGsfInnerHits, numCtfInnerHits);
424  ctfTrackRef = TrackRef(currentCtfTracks,counter);
425  }
426  } //ctfTrack iterator
427 }
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
const Double_t pi
uint16_t hitPattern[ARRAY_LENGTH]
Definition: HitPattern.h:502
T sqrt(T t)
Definition: SSEVec.h:18
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
T min(T a, T b)
Definition: MathUtil.h:58
bool isNull() const
Checks for null.
Definition: Ref.h:248
T const * product() const
Definition: Handle.h:81
edm::Ref< TrackCollection > TrackRef
persistent reference to a Track
Definition: TrackFwd.h:20
const reco::GsfTrackRef gsfTrackRef
static std::atomic< unsigned int > counter
uint16_t getHitPattern(HitCategory category, int position) const
Definition: HitPattern.h:545
void GsfElectronAlgo::ElectronData::computeCharge ( int &  charge,
reco::GsfElectron::ChargeInfo info 
)

Definition at line 430 of file GsfElectronAlgo.cc.

References PV3DBase< T, PVType, FrameType >::barePhi(), ecalDrivenElectronSeedsParameters_cff::beamSpot, ele_convert(), getEleBasicCluster(), reco::GsfElectron::ChargeInfo::isGsfCtfConsistent, reco::GsfElectron::ChargeInfo::isGsfCtfScPixConsistent, reco::GsfElectron::ChargeInfo::isGsfScPixConsistent, normalized_phi(), and reco::GsfElectron::ChargeInfo::scPixCharge.

Referenced by checkCtfTrack(), and GsfElectronAlgo::createElectron().

431  {
432  // determine charge from SC
433  GlobalPoint orig, scpos ;
434  ele_convert(beamSpot.position(),orig) ;
435  ele_convert(superClusterRef->position(),scpos) ;
436  GlobalVector scvect(scpos-orig) ;
437  GlobalPoint inntkpos = innTSOS.globalPosition() ;
438  GlobalVector inntkvect = GlobalVector(inntkpos-orig) ;
439  float dPhiInnEle=normalized_phi(scvect.barePhi()-inntkvect.barePhi()) ;
440  if(dPhiInnEle>0) info.scPixCharge = -1 ;
441  else info.scPixCharge = 1 ;
442 
443  // flags
444  int chargeGsf = gsfTrackRef->charge() ;
445  info.isGsfScPixConsistent = ((chargeGsf*info.scPixCharge)>0) ;
446  info.isGsfCtfConsistent = (ctfTrackRef.isNonnull()&&((chargeGsf*ctfTrackRef->charge())>0)) ;
448 
449  // default charge
451  { charge = info.scPixCharge ; }
452  else
453  { charge = ctfTrackRef->charge() ; }
454  }
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:251
GlobalPoint globalPosition() const
TrajectoryStateOnSurface innTSOS
T barePhi() const
Definition: PV3DBase.h:68
RealType normalized_phi(RealType phi)
bool isNull() const
Checks for null.
Definition: Ref.h:248
const reco::GsfTrackRef gsfTrackRef
void ele_convert(const Type1 &obj1, Type2 &obj2)
const reco::SuperClusterRef superClusterRef
const Point & position() const
position
Definition: BeamSpot.h:62
Global3DVector GlobalVector
Definition: GlobalVector.h:10
const reco::BeamSpot beamSpot
CaloClusterPtr GsfElectronAlgo::ElectronData::getEleBasicCluster ( const MultiTrajectoryStateTransform mtsTransform)

Definition at line 457 of file GsfElectronAlgo.cc.

References funct::abs(), ecalDrivenElectronSeedsParameters_cff::beamSpot, calculateTSOS(), particleFlow_cfi::dPhi, MultiTrajectoryStateTransform::extrapolatedState(), and MultiTrajectoryStateTransform::outerStateOnSurface().

Referenced by computeCharge(), and GsfElectronAlgo::createElectron().

458  {
459  CaloClusterPtr eleRef ;
460  TrajectoryStateOnSurface tempTSOS ;
462  float dphimin = 1.e30 ;
463  for (CaloCluster_iterator bc=superClusterRef->clustersBegin(); bc!=superClusterRef->clustersEnd(); bc++)
464  {
465  GlobalPoint posclu((*bc)->position().x(),(*bc)->position().y(),(*bc)->position().z()) ;
466  tempTSOS = mtsTransform->extrapolatedState(outTSOS,posclu) ;
467  if (!tempTSOS.isValid()) tempTSOS=outTSOS ;
468  GlobalPoint extrap = tempTSOS.globalPosition() ;
469  float dphi = EleRelPointPair(posclu,extrap,beamSpot.position()).dPhi() ;
470  if (std::abs(dphi)<dphimin)
471  {
472  dphimin = std::abs(dphi) ;
473  eleRef = (*bc);
474  eleTSOS = tempTSOS ;
475  }
476  }
477  return eleRef ;
478  }
TrajectoryStateOnSurface outerStateOnSurface(const reco::GsfTrack &tk) const
TrajectoryStateOnSurface extrapolatedState(const TrajectoryStateOnSurface tsos, const GlobalPoint &point) const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const reco::GsfTrackRef gsfTrackRef
const reco::SuperClusterRef superClusterRef
TrajectoryStateOnSurface eleTSOS
const Point & position() const
position
Definition: BeamSpot.h:62
TrajectoryStateOnSurface outTSOS
const reco::BeamSpot beamSpot

Member Data Documentation

const reco::BeamSpot GsfElectronAlgo::ElectronData::beamSpot

Definition at line 293 of file GsfElectronAlgo.cc.

TrajectoryStateOnSurface GsfElectronAlgo::ElectronData::constrainedVtxTSOS

Definition at line 316 of file GsfElectronAlgo.cc.

const reco::GsfElectronCoreRef GsfElectronAlgo::ElectronData::coreRef

Definition at line 288 of file GsfElectronAlgo.cc.

Referenced by GsfElectronAlgo::createElectron().

reco::TrackRef GsfElectronAlgo::ElectronData::ctfTrackRef

Definition at line 291 of file GsfElectronAlgo.cc.

Referenced by GsfElectronAlgo::createElectron().

GlobalVector GsfElectronAlgo::ElectronData::eleMom

Definition at line 319 of file GsfElectronAlgo.cc.

Referenced by GsfElectronAlgo::createElectron().

GlobalPoint GsfElectronAlgo::ElectronData::elePos

Definition at line 320 of file GsfElectronAlgo.cc.

Referenced by GsfElectronAlgo::createElectron().

TrajectoryStateOnSurface GsfElectronAlgo::ElectronData::eleTSOS

Definition at line 315 of file GsfElectronAlgo.cc.

const reco::GsfTrackRef GsfElectronAlgo::ElectronData::gsfTrackRef

Definition at line 289 of file GsfElectronAlgo.cc.

GlobalVector GsfElectronAlgo::ElectronData::innMom

Definition at line 319 of file GsfElectronAlgo.cc.

Referenced by GsfElectronAlgo::createElectron().

GlobalPoint GsfElectronAlgo::ElectronData::innPos

Definition at line 320 of file GsfElectronAlgo.cc.

TrajectoryStateOnSurface GsfElectronAlgo::ElectronData::innTSOS

Definition at line 310 of file GsfElectronAlgo.cc.

GlobalVector GsfElectronAlgo::ElectronData::outMom

Definition at line 319 of file GsfElectronAlgo.cc.

Referenced by GsfElectronAlgo::createElectron().

GlobalPoint GsfElectronAlgo::ElectronData::outPos

Definition at line 320 of file GsfElectronAlgo.cc.

TrajectoryStateOnSurface GsfElectronAlgo::ElectronData::outTSOS

Definition at line 311 of file GsfElectronAlgo.cc.

GlobalVector GsfElectronAlgo::ElectronData::sclMom

Definition at line 319 of file GsfElectronAlgo.cc.

Referenced by GsfElectronAlgo::createElectron().

GlobalPoint GsfElectronAlgo::ElectronData::sclPos

Definition at line 320 of file GsfElectronAlgo.cc.

Referenced by GsfElectronAlgo::createElectron().

TrajectoryStateOnSurface GsfElectronAlgo::ElectronData::sclTSOS

Definition at line 313 of file GsfElectronAlgo.cc.

GlobalVector GsfElectronAlgo::ElectronData::seedMom

Definition at line 319 of file GsfElectronAlgo.cc.

Referenced by GsfElectronAlgo::createElectron().

GlobalPoint GsfElectronAlgo::ElectronData::seedPos

Definition at line 320 of file GsfElectronAlgo.cc.

Referenced by GsfElectronAlgo::createElectron().

TrajectoryStateOnSurface GsfElectronAlgo::ElectronData::seedTSOS

Definition at line 314 of file GsfElectronAlgo.cc.

float GsfElectronAlgo::ElectronData::shFracInnerHits

Definition at line 292 of file GsfElectronAlgo.cc.

Referenced by GsfElectronAlgo::createElectron().

const reco::SuperClusterRef GsfElectronAlgo::ElectronData::superClusterRef

Definition at line 290 of file GsfElectronAlgo.cc.

Referenced by GsfElectronAlgo::createElectron().

GlobalVector GsfElectronAlgo::ElectronData::vtxMom

Definition at line 319 of file GsfElectronAlgo.cc.

Referenced by GsfElectronAlgo::createElectron().

GlobalVector GsfElectronAlgo::ElectronData::vtxMomWithConstraint

Definition at line 321 of file GsfElectronAlgo.cc.

Referenced by GsfElectronAlgo::createElectron().

GlobalPoint GsfElectronAlgo::ElectronData::vtxPos

Definition at line 320 of file GsfElectronAlgo.cc.

Referenced by GsfElectronAlgo::createElectron().

TrajectoryStateOnSurface GsfElectronAlgo::ElectronData::vtxTSOS

Definition at line 312 of file GsfElectronAlgo.cc.

Referenced by GsfElectronAlgo::createElectron().