CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_1/src/DataFormats/PatCandidates/src/Electron.cc

Go to the documentation of this file.
00001 //
00002 // $Id: Electron.cc,v 1.33 2012/10/04 11:00:10 beaudett Exp $
00003 //
00004 
00005 #include "DataFormats/PatCandidates/interface/Electron.h"
00006 #include "FWCore/Utilities/interface/Exception.h"
00007 
00008 #include <limits>
00009 
00010 using namespace pat;
00011 
00013 Electron::Electron() :
00014     Lepton<reco::GsfElectron>(),
00015     embeddedGsfElectronCore_(false),
00016     embeddedGsfTrack_(false),
00017     embeddedSuperCluster_(false),
00018     embeddedTrack_(false),
00019     embeddedSeedCluster_(false),
00020     embeddedRecHits_(false),
00021     embeddedPFCandidate_(false),
00022     ecalDrivenMomentum_(Candidate::LorentzVector(0.,0.,0.,0.)),
00023     cachedDB_(false),
00024     dB_(0.0),
00025     edB_(0.0),   
00026     ecalRegressionEnergy_(0.0),
00027     ecalTrackRegressionEnergy_(0.0),
00028     ecalRegressionError_(0.0),
00029     ecalTrackRegressionError_(0.0),
00030     ecalScale_(-99999.),
00031     ecalSmear_(-99999.),
00032     ecalRegressionScale_(-99999.),
00033     ecalRegressionSmear_(-99999.),
00034     ecalTrackRegressionScale_(-99999.),
00035     ecalTrackRegressionSmear_(-99999.)
00036 {
00037   initImpactParameters();
00038 }
00039 
00041 Electron::Electron(const reco::GsfElectron & anElectron) :
00042     Lepton<reco::GsfElectron>(anElectron),
00043     embeddedGsfElectronCore_(false),
00044     embeddedGsfTrack_(false),
00045     embeddedSuperCluster_(false),
00046     embeddedTrack_(false),
00047     embeddedSeedCluster_(false),
00048     embeddedRecHits_(false),
00049     embeddedPFCandidate_(false),
00050     ecalDrivenMomentum_(anElectron.p4()),
00051     cachedDB_(false),
00052     dB_(0.0),
00053     edB_(0.0)
00054 {
00055   initImpactParameters();
00056 }
00057 
00059 Electron::Electron(const edm::RefToBase<reco::GsfElectron> & anElectronRef) :
00060     Lepton<reco::GsfElectron>(anElectronRef),
00061     embeddedGsfElectronCore_(false),
00062     embeddedGsfTrack_(false),
00063     embeddedSuperCluster_(false),
00064     embeddedTrack_(false), 
00065     embeddedSeedCluster_(false),
00066     embeddedRecHits_(false),
00067     embeddedPFCandidate_(false),
00068     ecalDrivenMomentum_(anElectronRef->p4()),
00069     cachedDB_(false),
00070     dB_(0.0),
00071     edB_(0.0)
00072 {
00073   initImpactParameters();
00074 }
00075 
00077 Electron::Electron(const edm::Ptr<reco::GsfElectron> & anElectronRef) :
00078     Lepton<reco::GsfElectron>(anElectronRef),
00079     embeddedGsfElectronCore_(false),
00080     embeddedGsfTrack_(false),
00081     embeddedSuperCluster_(false),
00082     embeddedTrack_(false),
00083     embeddedSeedCluster_(false),
00084     embeddedRecHits_(false),
00085     embeddedPFCandidate_(false),
00086     ecalDrivenMomentum_(anElectronRef->p4()),
00087     cachedDB_(false),
00088     dB_(0.0),
00089     edB_(0.0)
00090 {
00091   initImpactParameters();
00092 }
00093 
00095 Electron::~Electron() {
00096 }
00097 
00099 std::ostream& 
00100 reco::operator<<(std::ostream& out, const pat::Electron& obj) 
00101 {
00102   if(!out) return out;
00103   
00104   out << "\tpat::Electron: ";
00105   out << std::setiosflags(std::ios::right);
00106   out << std::setiosflags(std::ios::fixed);
00107   out << std::setprecision(3);
00108   out << " E/pT/eta/phi " 
00109       << obj.energy()<<"/"
00110       << obj.pt()<<"/"
00111       << obj.eta()<<"/"
00112       << obj.phi();
00113   return out; 
00114 }
00115 
00117 void Electron::initImpactParameters() {
00118   for (int i_ = 0; i_<5; ++i_){
00119     ip_.push_back(0.0);
00120     eip_.push_back(0.0);
00121     cachedIP_.push_back(false);
00122   }
00123 }
00124 
00125 
00127 reco::GsfTrackRef Electron::gsfTrack() const {
00128   if (embeddedGsfTrack_) {
00129     return reco::GsfTrackRef(&gsfTrack_, 0);
00130   } else {
00131     return reco::GsfElectron::gsfTrack();
00132   }
00133 }
00134 
00136 reco::GsfElectronCoreRef Electron::core() const {
00137   if (embeddedGsfElectronCore_) {
00138     return reco::GsfElectronCoreRef(&gsfElectronCore_, 0);
00139   } else {
00140     return reco::GsfElectron::core();
00141   }
00142 }
00143 
00144 
00146 reco::SuperClusterRef Electron::superCluster() const {
00147   if (embeddedSuperCluster_) {
00148     return reco::SuperClusterRef(&superCluster_, 0);
00149   } else {
00150     return reco::GsfElectron::superCluster();
00151   }
00152 }
00153 
00155 reco::CaloClusterPtr Electron::seed() const {
00156   if(embeddedSeedCluster_){
00157     return reco::CaloClusterPtr(&seedCluster_,0);
00158   } else {
00159     return reco::GsfElectron::superCluster()->seed();
00160   }
00161 }
00162 
00164 reco::TrackRef Electron::closestCtfTrackRef() const {
00165   if (embeddedTrack_) {
00166     return reco::TrackRef(&track_, 0);
00167   } else {
00168     return reco::GsfElectron::closestCtfTrackRef();
00169   }
00170 }
00171 
00172 // the name of the method is misleading, users should use gsfTrack of closestCtfTrack
00173 reco::TrackRef Electron::track() const {
00174   return reco::TrackRef();
00175 }
00176 
00178 void Electron::embedGsfElectronCore() {
00179   gsfElectronCore_.clear();
00180   if (reco::GsfElectron::core().isNonnull()) {
00181       gsfElectronCore_.push_back(*reco::GsfElectron::core());
00182       embeddedGsfElectronCore_ = true;
00183   }
00184 }
00185 
00187 void Electron::embedGsfTrack() {
00188   gsfTrack_.clear();
00189   if (reco::GsfElectron::gsfTrack().isNonnull()) {
00190       gsfTrack_.push_back(*reco::GsfElectron::gsfTrack());
00191       embeddedGsfTrack_ = true;
00192   }
00193 }
00194 
00195 
00197 void Electron::embedSuperCluster() {
00198   superCluster_.clear();
00199   if (reco::GsfElectron::superCluster().isNonnull()) {
00200       superCluster_.push_back(*reco::GsfElectron::superCluster());
00201       embeddedSuperCluster_ = true;
00202   }
00203 }
00204 
00206 void Electron::embedSeedCluster() {
00207   seedCluster_.clear();
00208   if (reco::GsfElectron::superCluster().isNonnull() && reco::GsfElectron::superCluster()->seed().isNonnull()) {
00209     seedCluster_.push_back(*reco::GsfElectron::superCluster()->seed());
00210     embeddedSeedCluster_ = true;
00211   }
00212 }
00213 
00214 
00216 void Electron::embedTrack() {
00217   track_.clear();
00218   if (reco::GsfElectron::closestCtfTrackRef().isNonnull()) {
00219       track_.push_back(*reco::GsfElectron::closestCtfTrackRef());
00220       embeddedTrack_ = true;
00221   }
00222 }
00223 
00224 // method to store the RecHits internally
00225 void Electron::embedRecHits(const EcalRecHitCollection * rechits) {
00226   if (rechits!=0) {
00227     recHits_ = *rechits;
00228     embeddedRecHits_ = true;
00229   }
00230 }
00231 
00246 float Electron::electronID(const std::string& name) const {
00247     for (std::vector<IdPair>::const_iterator it = electronIDs_.begin(), ed = electronIDs_.end(); it != ed; ++it) {
00248         if (it->first == name) return it->second;
00249     }
00250     cms::Exception ex("Key not found");
00251     ex << "pat::Electron: the ID " << name << " can't be found in this pat::Electron.\n";
00252     ex << "The available IDs are: ";
00253     for (std::vector<IdPair>::const_iterator it = electronIDs_.begin(), ed = electronIDs_.end(); it != ed; ++it) {
00254         ex << "'" << it->first << "' ";
00255     }
00256     ex << ".\n";
00257     throw ex;
00258 }
00259 
00261 bool Electron::isElectronIDAvailable(const std::string& name) const {
00262     for (std::vector<IdPair>::const_iterator it = electronIDs_.begin(), ed = electronIDs_.end(); it != ed; ++it) {
00263         if (it->first == name) return true;
00264     }
00265     return false;
00266 }
00267 
00268 
00270 reco::PFCandidateRef Electron::pfCandidateRef() const {
00271   if (embeddedPFCandidate_) {
00272     return reco::PFCandidateRef(&pfCandidate_, 0);
00273   } else {
00274     return pfCandidateRef_;
00275   }
00276 }
00277 
00279 void Electron::embedPFCandidate() {
00280   pfCandidate_.clear();
00281   if ( pfCandidateRef_.isAvailable() && pfCandidateRef_.isNonnull()) {
00282     pfCandidate_.push_back( *pfCandidateRef_ );
00283     embeddedPFCandidate_ = true;
00284   }
00285 }
00286 
00289 reco::CandidatePtr Electron::sourceCandidatePtr( size_type i ) const {
00290   if (embeddedPFCandidate_) {
00291     return reco::CandidatePtr( pfCandidateRef_.id(), pfCandidateRef_.get(), pfCandidateRef_.key() );
00292   } else {
00293     return reco::CandidatePtr();
00294   }
00295 }
00296 
00297 
00298 
00309 double Electron::dB(IpType type_) const {
00310   // preserve old functionality exactly
00311   if (type_ == None){
00312     if ( cachedDB_ ) {
00313       return dB_;
00314     } else {
00315       return std::numeric_limits<double>::max();
00316     }
00317   }
00318   // more IP types (new)
00319   else if ( cachedIP_[type_] ) {
00320     return ip_[type_];
00321   } else {
00322     return std::numeric_limits<double>::max();
00323   }
00324 }
00325 
00336 double Electron::edB(IpType type_) const {
00337   // preserve old functionality exactly
00338   if (type_ == None) {
00339     if ( cachedDB_ ) {
00340       return edB_;
00341     } else {
00342       return std::numeric_limits<double>::max();
00343     }
00344   }
00345   // more IP types (new)
00346   else if ( cachedIP_[type_] ) {
00347     return eip_[type_];
00348   } else {
00349     return std::numeric_limits<double>::max();
00350   }
00351 
00352 }
00353 
00355 void Electron::setDB(double dB, double edB, IpType type){
00356   if (type == None) { // Preserve  old functionality exactly
00357     dB_ = dB; edB_ = edB;
00358     cachedDB_ = true;
00359   } else {
00360     ip_[type] = dB; 
00361     eip_[type] = edB; 
00362     cachedIP_[type] = true;
00363   }
00364 }
00365 
00367 void Electron::setMvaVariables( double r9, double sigmaIphiIphi, double sigmaIetaIphi, double ip3d){
00368   r9_ = r9;
00369   sigmaIphiIphi_ = sigmaIphiIphi;
00370   sigmaIetaIphi_ = sigmaIetaIphi;
00371   ip3d_ = ip3d;
00372 }