CMS 3D CMS Logo

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