CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SiStripElectronAssociator.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: SiStripElectronAssociator
4 // Class: SiStripElectronAssociator
5 //
13 //
14 // Original Author: Jim Pivarski
15 // Created: Tue Aug 1 15:24:02 EDT 2006
16 // $Id: SiStripElectronAssociator.cc,v 1.5 2008/10/31 15:30:56 nancy Exp $
17 //
18 //
19 
20 
21 #include <map>
22 #include <sstream>
23 
25 
32 
40 
41 //
42 // constants, enums and typedefs
43 //
44 
45 //
46 // static data member definitions
47 //
48 
49 //
50 // constructors and destructor
51 //
53 {
54  //register your products
55  electronsLabel_ = iConfig.getParameter<edm::InputTag>("electronsLabel");
56  produces<reco::ElectronCollection>(electronsLabel_.label());
57 
58  //now do what ever other initialization is needed
59  //siStripElectronProducer_ = iConfig.getParameter<edm::InputTag>("siStripElectronProducer");
60  siStripElectronCollection_ = iConfig.getParameter<edm::InputTag>("siStripElectronCollection");
61  //trackProducer_ = iConfig.getParameter<edm::InputTag>("trackProducer");
62  trackCollection_ = iConfig.getParameter<edm::InputTag>("trackCollection");
63 }
64 
65 
67 {}
68 
69 
70 void
72 {
73 
74  // position tolerance for equality of 2 hits set to 10 microns
75  static const double positionTol = 1e-3 ;
76 
78  iEvent.getByLabel(siStripElectronCollection_, siStripElectrons);
79 
81  iEvent.getByLabel(trackCollection_, tracks);
82 
83  std::map<const reco::SiStripElectron*, bool> alreadySeen;
84  for (reco::SiStripElectronCollection::const_iterator strippyIter = siStripElectrons->begin(); strippyIter != siStripElectrons->end(); ++strippyIter) {
85  alreadySeen[&(*strippyIter)] = false;
86  }
87 
88  // Output the high-level Electrons
89  std::auto_ptr<reco::ElectronCollection> output(new reco::ElectronCollection);
90 
91  LogDebug("SiStripElectronAssociator") << " About to loop over tracks " << std::endl ;
92  LogDebug("SiStripElectronAssociator") << " Number of tracks in Associator " << tracks.product()->size() ;
93  LogDebug("SiStripElectronAssociator") << " Number of SiStripElectron Candidates " << siStripElectrons->size();
94 
95  // The reco::Track's hits are a (improper?) subset of the reco::SiStripElectron's
96  // countSiElFit counts electron candidates that return w/ a good track fit
97  int countSiElFit = 0 ;
98  for (unsigned int i = 0; i < tracks.product()->size(); i++) {
99  const reco::Track* trackPtr = &(*reco::TrackRef(tracks, i));
100 
101  // If the reco::Track and the reco::SiStripElectron share even
102  // one hit in common, they belong to each other. (Disjoint sets
103  // of hits are assigned to electrons.) So let's look at one hit.
104 
105  // But first, make sure the track's hit list is not empty.
106  if (trackPtr->recHitsBegin() == trackPtr->recHitsEnd()) { continue; }
107 
108  // Detector id is not enough to completely specify a hit
109  uint32_t id = (*trackPtr->recHitsBegin())->geographicalId().rawId();
110  LocalPoint pos = (*trackPtr->recHitsBegin())->localPosition();
111 
112  LogDebug("SiStripElectronAssociator") << " New Track Candidate " << i
113  << " DetId " << id
114  << " pos " << pos << "\n";
115 
116  // Find the electron with that hit!
117  bool foundElectron = false;
118  for (reco::SiStripElectronCollection::const_iterator strippyIter = siStripElectrons->begin(); strippyIter != siStripElectrons->end(); ++strippyIter) {
119  if (!alreadySeen[&(*strippyIter)]) {
120 
121  bool hitInCommon = false;
122  LogDebug("SiStripElectronAssociator") << " Looping over Mono hits " << "\n" ;
123 
124  for (std::vector<SiStripRecHit2D>::const_iterator hitIter = strippyIter->rphiRecHits().begin(); hitIter != strippyIter->rphiRecHits().end(); ++hitIter) {
125 
126  LogDebug("SiStripElectronAssociator") << " SiStripCand "
127  << " DetId " << hitIter->geographicalId().rawId()
128  << " localPos " << hitIter->localPosition()
129  << " deltasPos " << (hitIter->localPosition() - pos).mag() ;
130 
131  if (hitIter->geographicalId().rawId() == id &&
132  (hitIter->localPosition() - pos).mag() < positionTol ) {
133  hitInCommon = true;
134  LogDebug("SiStripElectronAssociator") << " hitInCommon True " << "\n" ;
135  break;
136  } else {
137  LogDebug("SiStripElectronAssociator") << " hitInCommon False " << "\n" ;
138  }
139  } // end loop over rphi hits
140 
141  if(!hitInCommon) {
142  LogDebug("SiStripElectronAssociator") << " Looping over Stereo hits " << "\n" ;
143 
144  for (std::vector<SiStripRecHit2D>::const_iterator hitIter = strippyIter->stereoRecHits().begin(); hitIter != strippyIter->stereoRecHits().end(); ++hitIter) {
145 
146  LogDebug("SiStripElectronAssociator") << " SiStripCand "
147  << " DetId " << hitIter->geographicalId().rawId()
148  << " localPos " << hitIter->localPosition()
149  << " deltasPos " << (hitIter->localPosition() - pos).mag() ;
150 
151  if (hitIter->geographicalId().rawId() == id &&
152  (hitIter->localPosition() - pos).mag() < positionTol) {
153  hitInCommon = true;
154  LogDebug("SiStripElectronAssociator") << " hitInCommon True " << "\n" ;
155  break;
156  } else {
157  LogDebug("SiStripElectronAssociator") << " hitInCommon False " << "\n" ;
158  }
159 
160  } // end loop over stereo hits
161  } // end of hitInCommon check for loop over stereo hits
162 
163  if (hitInCommon) {
164  LogDebug("SiStripElectronAssociator") << " Hit in Common Found \n" ;
165  ++countSiElFit ;
166  foundElectron = true;
167  alreadySeen[&(*strippyIter)] = true;
168 
169  reco::Electron electron((trackPtr->charge() > 0 ? 1 : -1),
170  math::XYZTLorentzVector(trackPtr->px(),
171  trackPtr->py(),
172  trackPtr->pz(),
173  trackPtr->p()),
174  math::XYZPoint(trackPtr->vx(),
175  trackPtr->vy(),
176  trackPtr->vz()));
177  electron.setSuperCluster(strippyIter->superCluster());
178  electron.setTrack(reco::TrackRef(tracks, i));
179 
180  output->push_back(electron);
181  break ; // no need to check other SiStripElectrons if this one is good
182 
183  } else {
184  LogDebug("SiStripElectronAssociator") << "Hit in Common NOT found \n" ;
185  }
186  // endif this electron belongs to this track
187 
188 
189  } // endif we haven't seen this electron before
190  LogDebug("SiStripElectronAssociator") << "Done with this electron " << "\n\n\n";
191  } // end loop over electrons
192 
193  LogDebug("SiStripElectronAssociator") << "Testing if foundElectron " << foundElectron << std::endl;
194 
195  if (!foundElectron) {
196  throw cms::Exception("Configuration")
197  << " It is possible that the trackcollection used '"
198  << trackCollection_ << "' from producer '" << trackProducer_
199  << "' is not consistent with '"<< siStripElectronCollection_
200  << "' from the producer '"<< siStripElectronProducer_
201  << "' --- Please check your cfg file " << "\n";
202  }
203 
204  LogDebug("SiStripElectronAssociator") << "At end of track loop \n" << std::endl;
205 
206  } // end loop over tracks
207 
208 
209  LogDebug("SiStripElectronAssociator") << " Number of SiStripElectrons returned with a good fit "
210  << countSiElFit << "\n"<< std::endl ;
211  iEvent.put(output, electronsLabel_.label());
212 }
#define LogDebug(id)
double p() const
momentum vector magnitude
Definition: TrackBase.h:129
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
double px() const
x coordinate of momentum vector
Definition: TrackBase.h:133
virtual void produce(edm::Event &, const edm::EventSetup &)
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:30
int iEvent
Definition: GenABIO.cc:243
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:84
trackingRecHit_iterator recHitsBegin() const
Iterator to first hit on the track.
Definition: Track.h:63
std::vector< Electron > ElectronCollection
collectin of Electron objects
Definition: ElectronFwd.h:9
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:355
double pz() const
z coordinate of momentum vector
Definition: TrackBase.h:137
double vz() const
z coordinate of the reference point on track
Definition: TrackBase.h:147
tuple tracks
Definition: testEve_cfg.py:39
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:13
edm::Ref< TrackCollection > TrackRef
persistent reference to a Track
Definition: TrackFwd.h:14
T const * product() const
Definition: Handle.h:74
std::string const & label() const
Definition: InputTag.h:25
double vy() const
y coordinate of the reference point on track
Definition: TrackBase.h:145
int charge() const
track electric charge
Definition: TrackBase.h:113
SiStripElectronAssociator(const edm::ParameterSet &)
void setSuperCluster(const reco::SuperClusterRef &r)
set refrence to Photon component
Definition: Electron.h:33
void setTrack(const reco::TrackRef &r)
set refrence to Track component
Definition: Electron.h:35
double py() const
y coordinate of momentum vector
Definition: TrackBase.h:135
double vx() const
x coordinate of the reference point on track
Definition: TrackBase.h:143
trackingRecHit_iterator recHitsEnd() const
Iterator to last hit on the track.
Definition: Track.h:65