CMS 3D CMS Logo

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 //
17 //
18 
19 
20 #include <map>
21 #include <sstream>
22 
24 
30 
38 
39 //
40 // constants, enums and typedefs
41 //
42 
43 //
44 // static data member definitions
45 //
46 
47 //
48 // constructors and destructor
49 //
51 {
52  //register your products
53  electronsLabel_ = iConfig.getParameter<edm::InputTag>("electronsLabel");
54  produces<reco::ElectronCollection>(electronsLabel_.label());
55 
56  //now do what ever other initialization is needed
57  siStripElectronCollection_ = consumes<reco::SiStripElectronCollection>(iConfig.getParameter<edm::InputTag>("siStripElectronCollection"));
58  trackCollection_ = consumes<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("trackCollection"));
59 }
60 
61 
63 {}
64 
65 
66 void
68 {
69 
70  // position tolerance for equality of 2 hits set to 10 microns
71  static const double positionTol = 1e-3 ;
72 
74  iEvent.getByToken(siStripElectronCollection_, siStripElectrons);
75 
77  iEvent.getByToken(trackCollection_, tracks);
78 
79  std::map<const reco::SiStripElectron*, bool> alreadySeen;
80  for (reco::SiStripElectronCollection::const_iterator strippyIter = siStripElectrons->begin(); strippyIter != siStripElectrons->end(); ++strippyIter) {
81  alreadySeen[&(*strippyIter)] = false;
82  }
83 
84  // Output the high-level Electrons
85  auto output = std::make_unique<reco::ElectronCollection>();
86 
87  LogDebug("SiStripElectronAssociator") << " About to loop over tracks " << std::endl ;
88  LogDebug("SiStripElectronAssociator") << " Number of tracks in Associator " << tracks.product()->size() ;
89  LogDebug("SiStripElectronAssociator") << " Number of SiStripElectron Candidates " << siStripElectrons->size();
90 
91  // The reco::Track's hits are a (improper?) subset of the reco::SiStripElectron's
92  // countSiElFit counts electron candidates that return w/ a good track fit
93  int countSiElFit = 0 ;
94  for (unsigned int i = 0; i < tracks.product()->size(); i++) {
95  const reco::Track* trackPtr = &(*reco::TrackRef(tracks, i));
96 
97  // If the reco::Track and the reco::SiStripElectron share even
98  // one hit in common, they belong to each other. (Disjoint sets
99  // of hits are assigned to electrons.) So let's look at one hit.
100 
101  // But first, make sure the track's hit list is not empty.
102  if (trackPtr->recHitsBegin() == trackPtr->recHitsEnd()) { continue; }
103 
104  // Detector id is not enough to completely specify a hit
105  uint32_t id = (*trackPtr->recHitsBegin())->geographicalId().rawId();
106  LocalPoint pos = (*trackPtr->recHitsBegin())->localPosition();
107 
108  LogDebug("SiStripElectronAssociator") << " New Track Candidate " << i
109  << " DetId " << id
110  << " pos " << pos << "\n";
111 
112  // Find the electron with that hit!
113  bool foundElectron = false;
114  for (reco::SiStripElectronCollection::const_iterator strippyIter = siStripElectrons->begin(); strippyIter != siStripElectrons->end(); ++strippyIter) {
115  if (!alreadySeen[&(*strippyIter)]) {
116 
117  bool hitInCommon = false;
118  LogDebug("SiStripElectronAssociator") << " Looping over Mono hits " << "\n" ;
119 
120  for (std::vector<SiStripRecHit2D>::const_iterator hitIter = strippyIter->rphiRecHits().begin(); hitIter != strippyIter->rphiRecHits().end(); ++hitIter) {
121 
122  LogDebug("SiStripElectronAssociator") << " SiStripCand "
123  << " DetId " << hitIter->geographicalId().rawId()
124  << " localPos " << hitIter->localPosition()
125  << " deltasPos " << (hitIter->localPosition() - pos).mag() ;
126 
127  if (hitIter->geographicalId().rawId() == id &&
128  (hitIter->localPosition() - pos).mag() < positionTol ) {
129  hitInCommon = true;
130  LogDebug("SiStripElectronAssociator") << " hitInCommon True " << "\n" ;
131  break;
132  } else {
133  LogDebug("SiStripElectronAssociator") << " hitInCommon False " << "\n" ;
134  }
135  } // end loop over rphi hits
136 
137  if(!hitInCommon) {
138  LogDebug("SiStripElectronAssociator") << " Looping over Stereo hits " << "\n" ;
139 
140  for (std::vector<SiStripRecHit2D>::const_iterator hitIter = strippyIter->stereoRecHits().begin(); hitIter != strippyIter->stereoRecHits().end(); ++hitIter) {
141 
142  LogDebug("SiStripElectronAssociator") << " SiStripCand "
143  << " DetId " << hitIter->geographicalId().rawId()
144  << " localPos " << hitIter->localPosition()
145  << " deltasPos " << (hitIter->localPosition() - pos).mag() ;
146 
147  if (hitIter->geographicalId().rawId() == id &&
148  (hitIter->localPosition() - pos).mag() < positionTol) {
149  hitInCommon = true;
150  LogDebug("SiStripElectronAssociator") << " hitInCommon True " << "\n" ;
151  break;
152  } else {
153  LogDebug("SiStripElectronAssociator") << " hitInCommon False " << "\n" ;
154  }
155 
156  } // end loop over stereo hits
157  } // end of hitInCommon check for loop over stereo hits
158 
159  if (hitInCommon) {
160  LogDebug("SiStripElectronAssociator") << " Hit in Common Found \n" ;
161  ++countSiElFit ;
162  foundElectron = true;
163  alreadySeen[&(*strippyIter)] = true;
164 
165  reco::Electron electron((trackPtr->charge() > 0 ? 1 : -1),
166  math::XYZTLorentzVector(trackPtr->px(),
167  trackPtr->py(),
168  trackPtr->pz(),
169  trackPtr->p()),
170  math::XYZPoint(trackPtr->vx(),
171  trackPtr->vy(),
172  trackPtr->vz()));
173  electron.setSuperCluster(strippyIter->superCluster());
174  electron.setTrack(reco::TrackRef(tracks, i));
175 
176  output->push_back(electron);
177  break ; // no need to check other SiStripElectrons if this one is good
178 
179  } else {
180  LogDebug("SiStripElectronAssociator") << "Hit in Common NOT found \n" ;
181  }
182  // endif this electron belongs to this track
183 
184 
185  } // endif we haven't seen this electron before
186  LogDebug("SiStripElectronAssociator") << "Done with this electron " << "\n\n\n";
187  } // end loop over electrons
188 
189  LogDebug("SiStripElectronAssociator") << "Testing if foundElectron " << foundElectron << std::endl;
190 
191  if (!foundElectron) {
192  throw cms::Exception("Configuration")
193  << " Inconsistent track collection!";
194  }
195 
196  LogDebug("SiStripElectronAssociator") << "At end of track loop \n" << std::endl;
197 
198  } // end loop over tracks
199 
200 
201  LogDebug("SiStripElectronAssociator") << " Number of SiStripElectrons returned with a good fit "
202  << countSiElFit << "\n"<< std::endl ;
204 }
#define LogDebug(id)
double p() const
momentum vector magnitude
Definition: TrackBase.h:610
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:122
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
double px() const
x coordinate of momentum vector
Definition: TrackBase.h:622
virtual void produce(edm::Event &, const edm::EventSetup &)
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
int iEvent
Definition: GenABIO.cc:230
edm::EDGetTokenT< reco::TrackCollection > trackCollection_
trackingRecHit_iterator recHitsBegin() const
Iterator to first hit on the track.
Definition: Track.h:104
double pz() const
z coordinate of momentum vector
Definition: TrackBase.h:634
double vz() const
z coordinate of the reference point on track
Definition: TrackBase.h:664
T const * product() const
Definition: Handle.h:81
hitCont::const_iterator hitIter
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
edm::Ref< TrackCollection > TrackRef
persistent reference to a Track
Definition: TrackFwd.h:20
std::string const & label() const
Definition: InputTag.h:36
double vy() const
y coordinate of the reference point on track
Definition: TrackBase.h:658
int charge() const
track electric charge
Definition: TrackBase.h:562
SiStripElectronAssociator(const edm::ParameterSet &)
void setSuperCluster(const reco::SuperClusterRef &r)
set refrence to Photon component
Definition: Electron.h:35
void setTrack(const reco::TrackRef &r)
set refrence to Track component
Definition: Electron.h:37
def move(src, dest)
Definition: eostools.py:510
edm::EDGetTokenT< reco::SiStripElectronCollection > siStripElectronCollection_
double py() const
y coordinate of momentum vector
Definition: TrackBase.h:628
double vx() const
x coordinate of the reference point on track
Definition: TrackBase.h:652
trackingRecHit_iterator recHitsEnd() const
Iterator to last hit on the track.
Definition: Track.h:109