CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_5/src/RecoEgamma/EgammaElectronProducers/plugins/SiStripElectronAssociator.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    SiStripElectronAssociator
00004 // Class:      SiStripElectronAssociator
00005 // 
00013 //
00014 // Original Author:  Jim Pivarski
00015 //         Created:  Tue Aug  1 15:24:02 EDT 2006
00016 // $Id: SiStripElectronAssociator.cc,v 1.6 2010/03/06 12:45:47 sani Exp $
00017 //
00018 //
00019 
00020 
00021 #include <map>
00022 #include <sstream>
00023 
00024 #include "SiStripElectronAssociator.h"
00025 
00026 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00027 #include "DataFormats/EgammaCandidates/interface/SiStripElectronFwd.h"
00028 #include "DataFormats/EgammaCandidates/interface/SiStripElectron.h"
00029 #include "DataFormats/TrackReco/interface/Track.h"
00030 #include "DataFormats/EgammaCandidates/interface/Electron.h"
00031 #include "DataFormats/EgammaCandidates/interface/ElectronFwd.h"
00032 
00033 #include "DataFormats/DetId/interface/DetId.h"
00034 #include "FWCore/Utilities/interface/Exception.h"
00035 #include "DataFormats/GeometryVector/interface/LocalPoint.h"
00036 #include "DataFormats/TrackingRecHit/interface/TrackingRecHit.h"
00037 #include "DataFormats/Math/interface/LorentzVector.h"
00038 #include "DataFormats/Math/interface/Point3D.h"
00039 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit2DCollection.h"
00040 
00041 //
00042 // constants, enums and typedefs
00043 //
00044 
00045 //
00046 // static data member definitions
00047 //
00048 
00049 //
00050 // constructors and destructor
00051 //
00052 SiStripElectronAssociator::SiStripElectronAssociator(const edm::ParameterSet& iConfig)
00053 {
00054    //register your products
00055   electronsLabel_ = iConfig.getParameter<edm::InputTag>("electronsLabel");
00056   produces<reco::ElectronCollection>(electronsLabel_.label());
00057 
00058    //now do what ever other initialization is needed
00059   //siStripElectronProducer_ = iConfig.getParameter<edm::InputTag>("siStripElectronProducer");
00060   siStripElectronCollection_ = iConfig.getParameter<edm::InputTag>("siStripElectronCollection");
00061   //trackProducer_ = iConfig.getParameter<edm::InputTag>("trackProducer");
00062   trackCollection_ = iConfig.getParameter<edm::InputTag>("trackCollection");
00063 }
00064 
00065 
00066 SiStripElectronAssociator::~SiStripElectronAssociator()
00067 {}
00068 
00069 
00070 void
00071 SiStripElectronAssociator::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00072 {
00073 
00074   // position tolerance for equality of 2 hits set to 10 microns
00075   static const double positionTol = 1e-3 ; 
00076 
00077    edm::Handle<reco::SiStripElectronCollection> siStripElectrons;
00078    iEvent.getByLabel(siStripElectronCollection_, siStripElectrons);
00079 
00080    edm::Handle<reco::TrackCollection> tracks;
00081    iEvent.getByLabel(trackCollection_, tracks);
00082 
00083    std::map<const reco::SiStripElectron*, bool> alreadySeen;
00084    for (reco::SiStripElectronCollection::const_iterator strippyIter = siStripElectrons->begin();  strippyIter != siStripElectrons->end();  ++strippyIter) {
00085       alreadySeen[&(*strippyIter)] = false;
00086    }
00087 
00088    // Output the high-level Electrons
00089    std::auto_ptr<reco::ElectronCollection> output(new reco::ElectronCollection);
00090 
00091    LogDebug("SiStripElectronAssociator") << " About to loop over tracks " << std::endl ;
00092    LogDebug("SiStripElectronAssociator") << " Number of tracks in Associator " << tracks.product()->size() ;
00093    LogDebug("SiStripElectronAssociator") << " Number of SiStripElectron Candidates " << siStripElectrons->size();
00094 
00095    // The reco::Track's hits are a (improper?) subset of the reco::SiStripElectron's
00096    // countSiElFit counts electron candidates that return w/ a good track fit
00097    int countSiElFit = 0 ;
00098    for (unsigned int i = 0;  i < tracks.product()->size();  i++) {
00099       const reco::Track* trackPtr = &(*reco::TrackRef(tracks, i));
00100       
00101       // If the reco::Track and the reco::SiStripElectron share even
00102       // one hit in common, they belong to each other.  (Disjoint sets
00103       // of hits are assigned to electrons.)  So let's look at one hit.
00104 
00105       // But first, make sure the track's hit list is not empty.
00106       if (trackPtr->recHitsBegin() == trackPtr->recHitsEnd()) { continue; }
00107 
00108       // Detector id is not enough to completely specify a hit
00109       uint32_t id = (*trackPtr->recHitsBegin())->geographicalId().rawId();
00110       LocalPoint pos = (*trackPtr->recHitsBegin())->localPosition();
00111       
00112       LogDebug("SiStripElectronAssociator") << " New Track Candidate " << i
00113                                             << " DetId " << id
00114                                             << " pos " << pos << "\n";
00115 
00116       // Find the electron with that hit!
00117       bool foundElectron = false;
00118       for (reco::SiStripElectronCollection::const_iterator strippyIter = siStripElectrons->begin();  strippyIter != siStripElectrons->end();  ++strippyIter) {
00119         if (!alreadySeen[&(*strippyIter)]) {
00120           
00121           bool hitInCommon = false;
00122           LogDebug("SiStripElectronAssociator") << " Looping over Mono hits " << "\n" ;
00123           
00124           for (std::vector<SiStripRecHit2D>::const_iterator hitIter = strippyIter->rphiRecHits().begin();  hitIter != strippyIter->rphiRecHits().end();  ++hitIter) {
00125             
00126             LogDebug("SiStripElectronAssociator") << " SiStripCand " 
00127                                                   << " DetId " << hitIter->geographicalId().rawId()
00128                                                   << " localPos " << hitIter->localPosition()
00129                                                   << " deltasPos " << (hitIter->localPosition() - pos).mag() ;
00130             
00131             if (hitIter->geographicalId().rawId() == id   &&
00132                 (hitIter->localPosition() - pos).mag() < positionTol ) {
00133               hitInCommon = true;
00134               LogDebug("SiStripElectronAssociator") << " hitInCommon True " << "\n" ;
00135               break;
00136             } else {
00137               LogDebug("SiStripElectronAssociator") << " hitInCommon False " << "\n" ;
00138             }
00139           } // end loop over rphi hits
00140           
00141           if(!hitInCommon) {
00142             LogDebug("SiStripElectronAssociator") << " Looping over Stereo hits " << "\n" ;
00143 
00144             for (std::vector<SiStripRecHit2D>::const_iterator hitIter = strippyIter->stereoRecHits().begin();  hitIter != strippyIter->stereoRecHits().end();  ++hitIter) {
00145               
00146                LogDebug("SiStripElectronAssociator") << " SiStripCand " 
00147                                                      << " DetId " << hitIter->geographicalId().rawId()
00148                                                      << " localPos " << hitIter->localPosition()
00149                                                      << " deltasPos " << (hitIter->localPosition() - pos).mag() ;
00150                
00151                if (hitIter->geographicalId().rawId() == id   &&
00152                    (hitIter->localPosition() - pos).mag() < positionTol) {
00153                  hitInCommon = true;
00154                   LogDebug("SiStripElectronAssociator") << " hitInCommon True " << "\n"  ;
00155                  break;
00156                } else {
00157                  LogDebug("SiStripElectronAssociator") << " hitInCommon False " << "\n" ;
00158                }
00159                
00160             } // end loop over stereo hits
00161           } // end of hitInCommon check for loop over stereo hits
00162           
00163           if (hitInCommon) {
00164              LogDebug("SiStripElectronAssociator") << " Hit in Common Found \n" ;
00165              ++countSiElFit ;
00166              foundElectron = true;
00167              alreadySeen[&(*strippyIter)] = true;
00168              
00169              reco::Electron electron((trackPtr->charge() > 0 ? 1 : -1),
00170                                      math::XYZTLorentzVector(trackPtr->px(),
00171                                                              trackPtr->py(),
00172                                                              trackPtr->pz(),
00173                                                              trackPtr->p()),
00174                                      math::XYZPoint(trackPtr->vx(),
00175                                                     trackPtr->vy(),
00176                                                     trackPtr->vz()));
00177              electron.setSuperCluster(strippyIter->superCluster());
00178              electron.setTrack(reco::TrackRef(tracks, i));
00179              
00180              output->push_back(electron);
00181              break ; // no need to check other SiStripElectrons if this one is good 
00182              
00183           } else {
00184              LogDebug("SiStripElectronAssociator") << "Hit in Common NOT found \n"  ;
00185           }
00186           // endif this electron belongs to this track
00187           
00188           
00189         } // endif we haven't seen this electron before
00190         LogDebug("SiStripElectronAssociator") << "Done with this electron " << "\n\n\n";
00191       } // end loop over electrons
00192       
00193       LogDebug("SiStripElectronAssociator") << "Testing if foundElectron " << foundElectron << std::endl;
00194       
00195       if (!foundElectron) {
00196         throw cms::Exception("Configuration")
00197           << " It is possible that the trackcollection used '"
00198           << trackCollection_ << "' from producer '" << trackProducer_
00199           << "' is not consistent with '"<< siStripElectronCollection_ 
00200           << "' from the producer '"<< siStripElectronProducer_
00201           << "' --- Please check your cfg file " << "\n";
00202       }
00203       
00204       LogDebug("SiStripElectronAssociator") << "At end of track loop \n" << std::endl; 
00205       
00206    } // end loop over tracks
00207    
00208    
00209    LogDebug("SiStripElectronAssociator") << " Number of SiStripElectrons returned with a good fit " 
00210                                          << countSiElFit << "\n"<<  std::endl ;
00211    iEvent.put(output, electronsLabel_.label());
00212 }