CMS 3D CMS Logo

EgAmbiguityTools.cc
Go to the documentation of this file.
7 
8 #include <algorithm>
9 
10 using namespace edm ;
11 using namespace std ;
12 using namespace reco ;
13 
14 namespace EgAmbiguityTools
15 {
16 
17 bool isBetter( reco::GsfElectron const& e1, reco::GsfElectron const& e2 )
18  { return (std::abs(e1.eSuperClusterOverP()-1) < std::abs(e2.eSuperClusterOverP() - 1)) ; }
19 
20 bool isInnerMost( reco::GsfElectron const& e1, reco::GsfElectron const& e2)
21 {
22  // retreive first valid hit
23  int gsfHitCounter1 = 0 ;
24  for(auto const& elHit : e1.gsfTrack()->recHits())
25  {
26  if (elHit->isValid()) break;
27  gsfHitCounter1++;
28  }
29 
30  int gsfHitCounter2 = 0 ;
31  for(auto const& elHit : e2.gsfTrack()->recHits())
32  {
33  if (elHit->isValid()) break;
34  gsfHitCounter2++;
35  }
36 
37  uint32_t gsfHit1 = e1.gsfTrack()->hitPattern().getHitPattern(HitPattern::TRACK_HITS, gsfHitCounter1);
38  uint32_t gsfHit2 = e2.gsfTrack()->hitPattern().getHitPattern(HitPattern::TRACK_HITS, gsfHitCounter2);
39 
40  if (HitPattern::getSubStructure(gsfHit1) != HitPattern::getSubStructure(gsfHit2)){
41  return (HitPattern::getSubStructure(gsfHit1) < HitPattern::getSubStructure(gsfHit2));
42  }else if (HitPattern::getLayer(gsfHit1) != HitPattern::getLayer(gsfHit2)){
43  return (HitPattern::getLayer(gsfHit1) < HitPattern::getLayer(gsfHit2));
44  }else{
45  return isBetter(e1, e2);
46  }
47 }
48 
49 int sharedHits(const GsfTrackRef & gsfTrackRef1, const GsfTrackRef & gsfTrackRef2 )
50 {
51  //get the Hit Pattern for the gsfTracks
52  HitPattern const& gsfHitPattern1 = gsfTrackRef1->hitPattern();
53  HitPattern const& gsfHitPattern2 = gsfTrackRef2->hitPattern();
54 
55  unsigned int shared = 0;
56 
57  int gsfHitCounter1 = 0;
58  for(trackingRecHit_iterator elHitsIt1 = gsfTrackRef1->recHitsBegin();
59  elHitsIt1 != gsfTrackRef1->recHitsEnd(); elHitsIt1++, gsfHitCounter1++) {
60  if(!(*elHitsIt1)->isValid()){
61  //count only valid Hits
62  continue;
63  }
64  //if (gsfHitCounter1>1) continue; // test only the first hit of the track 1
65  uint32_t gsfHit = gsfHitPattern1.getHitPattern(HitPattern::TRACK_HITS, gsfHitCounter1);
66  if(!(HitPattern::pixelHitFilter(gsfHit)
67  || HitPattern::stripTIBHitFilter(gsfHit)
68  || HitPattern::stripTOBHitFilter(gsfHit)
69  || HitPattern::stripTECHitFilter(gsfHit)
70  || HitPattern::stripTIDHitFilter(gsfHit))){
71  continue;
72  }
73 
74  int gsfHitsCounter2 = 0;
75  for(trackingRecHit_iterator gsfHitsIt2 = gsfTrackRef2->recHitsBegin();
76  gsfHitsIt2 != gsfTrackRef2->recHitsEnd(); gsfHitsIt2++, gsfHitsCounter2++) {
77  if(!(**gsfHitsIt2).isValid()){
78  //count only valid Hits
79  continue;
80  }
81  uint32_t gsfHit2 = gsfHitPattern2.getHitPattern(HitPattern::TRACK_HITS, gsfHitsCounter2);
82  if(!(HitPattern::pixelHitFilter(gsfHit2)
83  || HitPattern::stripTIBHitFilter(gsfHit2)
84  || HitPattern::stripTOBHitFilter(gsfHit2)
85  || HitPattern::stripTECHitFilter(gsfHit2)
86  || HitPattern::stripTIDHitFilter(gsfHit2))){
87  continue;
88  }
89  if((*elHitsIt1)->sharesInput(&(**gsfHitsIt2), TrackingRecHit::some)) {
90  //if (comp.equals(&(**elHitsIt1),&(**gsfHitsIt2))) {
92  shared++;
93  }
94  }//gsfHits2 iterator
95  }//gsfHits1 iterator
96 
97  //std::cout << "[sharedHits] number of shared hits " << shared << std::endl;
98  return shared;
99 }
100 
101 int sharedDets(const GsfTrackRef& gsfTrackRef1, const GsfTrackRef& gsfTrackRef2 )
102 {
103  //get the Hit Pattern for the gsfTracks
104  const HitPattern &gsfHitPattern1 = gsfTrackRef1->hitPattern();
105  const HitPattern &gsfHitPattern2 = gsfTrackRef2->hitPattern();
106 
107  unsigned int shared = 0;
108 
109  int gsfHitCounter1 = 0;
110  for(trackingRecHit_iterator elHitsIt1 = gsfTrackRef1->recHitsBegin();
111  elHitsIt1 != gsfTrackRef1->recHitsEnd(); elHitsIt1++, gsfHitCounter1++) {
112  if(!((**elHitsIt1).isValid())){
113  //count only valid Hits
114  continue;
115  }
116  //if (gsfHitCounter1>1) continue; // to test only the first hit of the track 1
117  uint32_t gsfHit = gsfHitPattern1.getHitPattern(HitPattern::TRACK_HITS, gsfHitCounter1);
118  if(!(HitPattern::pixelHitFilter(gsfHit)
119  || HitPattern::stripTIBHitFilter(gsfHit)
120  || HitPattern::stripTOBHitFilter(gsfHit)
121  || HitPattern::stripTECHitFilter(gsfHit)
122  || HitPattern::stripTIDHitFilter(gsfHit)))
123  {
124  continue;
125  }
126 
127  int gsfHitsCounter2 = 0;
128  for(trackingRecHit_iterator gsfHitsIt2 = gsfTrackRef2->recHitsBegin();
129  gsfHitsIt2 != gsfTrackRef2->recHitsEnd(); gsfHitsIt2++, gsfHitsCounter2++) {
130  if(!((**gsfHitsIt2).isValid())){
131  //count only valid Hits!
132  continue;
133  }
134 
135  uint32_t gsfHit2 = gsfHitPattern2.getHitPattern(HitPattern::TRACK_HITS, gsfHitsCounter2);
136  if(!(HitPattern::pixelHitFilter(gsfHit2)
137  || HitPattern::stripTIBHitFilter(gsfHit2)
138  || HitPattern::stripTOBHitFilter(gsfHit2)
139  || HitPattern::stripTECHitFilter(gsfHit2)
140  || HitPattern::stripTIDHitFilter(gsfHit2) )
141  ) continue;
142  if ((**elHitsIt1).geographicalId() == (**gsfHitsIt2).geographicalId()) shared++;
143  }//gsfHits2 iterator
144  }//gsfHits1 iterator
145 
146  //std::cout << "[sharedHits] number of shared dets " << shared << std::endl;
147  //return shared/min(gsfTrackRef1->numberOfValidHits(),gsfTrackRef2->numberOfValidHits());
148  return shared;
149 
150 }
151 
152 float sharedEnergy( CaloCluster const& clu1, CaloCluster const& clu2,
155 {
156 
157  double fractionShared = 0;
158 
159  for(auto const& h1 : clu1.hitsAndFractions()) {
160 
161  for(auto const& h2 : clu2.hitsAndFractions()) {
162 
163  if ( h1.first != h2.first ) continue;
164 
165  // here we have common Xtal id
167  if (h1.first.subdetId() == EcalBarrel) {
168  if ((itt=barrelRecHits.find(h1.first))!=barrelRecHits.end())
169  fractionShared += itt->energy();
170  } else if (h1.first.subdetId() == EcalEndcap) {
171  if ((itt=endcapRecHits.find(h1.first))!=endcapRecHits.end())
172  fractionShared += itt->energy();
173  }
174 
175  }
176  }
177 
178  //std::cout << "[sharedEnergy] shared energy /min(energy1,energy2) " << fractionShared << std::endl;
179  return fractionShared;
180 
181 }
182 
183 float sharedEnergy( SuperClusterRef const& sc1, SuperClusterRef const& sc2,
186 {
187  double energyShared = 0;
188  for(CaloCluster_iterator icl1=sc1->clustersBegin();icl1!=sc1->clustersEnd(); icl1++) {
189  for(CaloCluster_iterator icl2=sc2->clustersBegin();icl2!=sc2->clustersEnd(); icl2++) {
190  energyShared += sharedEnergy(**icl1,**icl2,barrelRecHits,endcapRecHits );
191  }
192  }
193  return energyShared;
194 
195 }
196 
197 }
GsfTrackRef gsfTrack() const override
reference to a GsfTrack
Definition: GsfElectron.h:186
int sharedDets(const GsfTrackRef &gsfTrackRef1, const GsfTrackRef &gsfTrackRef2)
int sharedHits(const GsfTrackRef &gsfTrackRef1, const GsfTrackRef &gsfTrackRef2)
float eSuperClusterOverP() const
Definition: GsfElectron.h:249
std::vector< EcalRecHit >::const_iterator const_iterator
bool isBetter(reco::GsfElectron const &, reco::GsfElectron const &)
const std::vector< std::pair< DetId, float > > & hitsAndFractions() const
Definition: CaloCluster.h:197
float sharedEnergy(reco::CaloCluster const &clu1, reco::CaloCluster const &clu2, EcalRecHitCollection const &barrelRecHits, EcalRecHitCollection const &endcapRecHits)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const_iterator end() const
iterator find(key_type k)
fixed size matrix
HLT enums.
bool isInnerMost(reco::GsfElectron const &, reco::GsfElectron const &)
uint16_t getHitPattern(HitCategory category, int position) const
Definition: HitPattern.h:553