CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
GsfElectronCoreBaseProducer.cc
Go to the documentation of this file.
1 
3 
7 
14 //#include "DataFormats/Common/interface/ValueMap.h"
15 
16 //#include <map>
17 
18 using namespace reco ;
19 
21  {
22  desc.add<edm::InputTag>("gsfPfRecTracks",edm::InputTag("pfTrackElec")) ;
23  desc.add<edm::InputTag>("gsfTracks",edm::InputTag("electronGsfTracks")) ;
24  desc.add<edm::InputTag>("ctfTracks",edm::InputTag("generalTracks")) ;
25  desc.add<bool>("useGsfPfRecTracks",true) ;
26  }
27 
29  {
30  produces<GsfElectronCoreCollection>() ;
31  gsfPfRecTracksTag_ = mayConsume<reco::GsfPFRecTrackCollection>(config.getParameter<edm::InputTag>("gsfPfRecTracks")) ;
32  gsfTracksTag_ = consumes<reco::GsfTrackCollection>(config.getParameter<edm::InputTag>("gsfTracks"));
33  ctfTracksTag_ = consumes<reco::TrackCollection>(config.getParameter<edm::InputTag>("ctfTracks"));
34  useGsfPfRecTracks_ = config.getParameter<bool>("useGsfPfRecTracks") ;
35  }
36 
38  {}
39 
40 
41 //=======================================================================================
42 // For derived producers
43 //=======================================================================================
44 
45 // to be called at the beginning of each new event
47  {
48  if (useGsfPfRecTracks_)
49  { event.getByToken(gsfPfRecTracksTag_,gsfPfRecTracksH_) ; }
50  event.getByToken(gsfTracksTag_,gsfTracksH_) ;
51  event.getByToken(ctfTracksTag_,ctfTracksH_) ;
52  }
53 
55  {
56  const GsfTrackRef & gsfTrackRef = eleCore->gsfTrack() ;
57 
58  std::pair<TrackRef,float> ctfpair = getCtfTrackRef(gsfTrackRef) ;
59  eleCore->setCtfTrack(ctfpair.first,ctfpair.second) ;
60  }
61 
62 
63 //=======================================================================================
64 // Code from Puneeth Kalavase
65 //=======================================================================================
66 
67 std::pair<TrackRef,float> GsfElectronCoreBaseProducer::getCtfTrackRef
68  ( const GsfTrackRef & gsfTrackRef )
69  {
70  float maxFracShared = 0;
71  TrackRef ctfTrackRef = TrackRef() ;
72  const TrackCollection * ctfTrackCollection = ctfTracksH_.product() ;
73 
74  // get the Hit Pattern for the gsfTrack
75  const HitPattern& gsfHitPattern = gsfTrackRef->hitPattern();
76 
77  unsigned int counter ;
78  TrackCollection::const_iterator ctfTkIter ;
79  for ( ctfTkIter = ctfTrackCollection->begin() , counter = 0 ;
80  ctfTkIter != ctfTrackCollection->end() ; ctfTkIter++, counter++ )
81  {
82 
83  double dEta = gsfTrackRef->eta() - ctfTkIter->eta();
84  double dPhi = gsfTrackRef->phi() - ctfTkIter->phi();
85  double pi = acos(-1.);
86  if(std::abs(dPhi) > pi) dPhi = 2*pi - std::abs(dPhi);
87 
88  // dont want to look at every single track in the event!
89  if(sqrt(dEta*dEta + dPhi*dPhi) > 0.3) continue;
90 
91  unsigned int shared = 0 ;
92  int gsfHitCounter = 0 ;
93  int numGsfInnerHits = 0 ;
94  int numCtfInnerHits = 0 ;
95  // get the CTF Track Hit Pattern
96  const HitPattern& ctfHitPattern = ctfTkIter->hitPattern();
97 
98  trackingRecHit_iterator elHitsIt ;
99  for ( elHitsIt = gsfTrackRef->recHitsBegin() ;
100  elHitsIt != gsfTrackRef->recHitsEnd() ;
101  elHitsIt++, gsfHitCounter++ )
102  {
103  if(!((**elHitsIt).isValid())) //count only valid Hits
104  { continue ; }
105 
106  // look only in the pixels/TIB/TID
107  uint32_t gsfHit = gsfHitPattern.getHitPattern(HitPattern::TRACK_HITS, gsfHitCounter);
108  if (!(HitPattern::pixelHitFilter(gsfHit) ||
111  { continue ; }
112 
113  numGsfInnerHits++ ;
114 
115  int ctfHitsCounter = 0 ;
116  numCtfInnerHits = 0 ;
117  trackingRecHit_iterator ctfHitsIt ;
118  for ( ctfHitsIt = ctfTkIter->recHitsBegin() ;
119  ctfHitsIt != ctfTkIter->recHitsEnd() ;
120  ctfHitsIt++, ctfHitsCounter++ )
121  {
122  if(!((**ctfHitsIt).isValid())) //count only valid Hits!
123  { continue ; }
124 
125  uint32_t ctfHit = ctfHitPattern.getHitPattern(HitPattern::TRACK_HITS, ctfHitsCounter);
126  if(!(HitPattern::pixelHitFilter(ctfHit) ||
129  { continue ; }
130 
131  numCtfInnerHits++ ;
132 
133  if( (**elHitsIt).sharesInput(&(**ctfHitsIt),TrackingRecHit::all) )
134  {
135  shared++ ;
136  break ;
137  }
138 
139  } //ctfHits iterator
140 
141  } //gsfHits iterator
142 
143  if ((numGsfInnerHits==0)||(numCtfInnerHits==0))
144  { continue ; }
145 
146  if ( static_cast<float>(shared)/std::min(numGsfInnerHits,numCtfInnerHits) > maxFracShared )
147  {
148  maxFracShared = static_cast<float>(shared)/std::min(numGsfInnerHits, numCtfInnerHits);
149  ctfTrackRef = TrackRef(ctfTracksH_,counter);
150  }
151 
152  } //ctfTrack iterator
153 
154  return make_pair(ctfTrackRef,maxFracShared) ;
155  }
156 
157 
158 
void initEvent(edm::Event &event, const edm::EventSetup &setup)
T getParameter(std::string const &) const
void setCtfTrack(const TrackRef &closestCtfTrack, float ctfGsfOverlap)
static bool pixelHitFilter(uint16_t pattern)
Definition: HitPattern.h:551
const GsfTrackRef & gsfTrack() const
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
GsfElectronCoreBaseProducer(const edm::ParameterSet &conf)
const Double_t pi
uint16_t hitPattern[ARRAY_LENGTH]
Definition: HitPattern.h:459
double dPhi(double phi1, double phi2)
Definition: JetUtil.h:30
T sqrt(T t)
Definition: SSEVec.h:48
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
T min(T a, T b)
Definition: MathUtil.h:58
ParameterDescriptionBase * add(U const &iLabel, T const &value)
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
std::pair< reco::TrackRef, float > getCtfTrackRef(const reco::GsfTrackRef &)
static void fillDescription(edm::ParameterSetDescription &)
edm::Ref< TrackCollection > TrackRef
persistent reference to a Track
Definition: TrackFwd.h:20
static bool stripTIBHitFilter(uint16_t pattern)
Definition: HitPattern.h:597
static std::atomic< unsigned int > counter
uint16_t getHitPattern(HitCategory category, int position) const
Definition: HitPattern.h:502
static bool stripTIDHitFilter(uint16_t pattern)
Definition: HitPattern.h:602
void fillElectronCore(reco::GsfElectronCore *)
void setup(std::vector< TH2F > &depth, std::string name, std::string units="")
TrackingRecHitCollection::base::const_iterator trackingRecHit_iterator
iterator over a vector of reference to TrackingRecHit in the same collection