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