** Warning **

Issuing rollback() due to DESTROY without explicit disconnect() of DBD::mysql::db handle dbname=lxr at /lxr/lib/LXR/Common.pm line 1113.

Last-Modified: Fri, 12 Dec 2024 03:12:41 GMT Content-Type: text/html; charset=utf-8 /CMSSW_15_0_X_2024-12-11-2300/DataFormats/PatCandidates/src/Electron.cc
Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:04:56

0001 //
0002 //
0003 
0004 #include "DataFormats/PatCandidates/interface/Electron.h"
0005 #include "FWCore/Utilities/interface/Exception.h"
0006 #include "DataFormats/Common/interface/RefToPtr.h"
0007 
0008 #include <limits>
0009 
0010 using namespace pat;
0011 
0012 /// default constructor
0013 Electron::Electron()
0014     : Lepton<reco::GsfElectron>(),
0015       embeddedGsfElectronCore_(false),
0016       embeddedGsfTrack_(false),
0017       embeddedSuperCluster_(false),
0018       embeddedPflowSuperCluster_(false),
0019       embeddedTrack_(false),
0020       embeddedSeedCluster_(false),
0021       embeddedRecHits_(false),
0022       embeddedPFCandidate_(false),
0023       ecalDrivenMomentum_(Candidate::LorentzVector(0., 0., 0., 0.)),
0024       ecalRegressionEnergy_(0.0),
0025       ecalTrackRegressionEnergy_(0.0),
0026       ecalRegressionError_(0.0),
0027       ecalTrackRegressionError_(0.0),
0028       ecalScale_(-99999.),
0029       ecalSmear_(-99999.),
0030       ecalRegressionScale_(-99999.),
0031       ecalRegressionSmear_(-99999.),
0032       ecalTrackRegressionScale_(-99999.),
0033       ecalTrackRegressionSmear_(-99999.),
0034       packedPFCandidates_(),
0035       associatedPackedFCandidateIndices_() {
0036   initImpactParameters();
0037 }
0038 
0039 /// constructor from reco::GsfElectron
0040 Electron::Electron(const reco::GsfElectron& anElectron)
0041     : Lepton<reco::GsfElectron>(anElectron),
0042       embeddedGsfElectronCore_(false),
0043       embeddedGsfTrack_(false),
0044       embeddedSuperCluster_(false),
0045       embeddedPflowSuperCluster_(false),
0046       embeddedTrack_(false),
0047       embeddedSeedCluster_(false),
0048       embeddedRecHits_(false),
0049       embeddedPFCandidate_(false),
0050       ecalDrivenMomentum_(anElectron.p4()) {
0051   initImpactParameters();
0052 }
0053 
0054 /// constructor from a RefToBase to a reco::GsfElectron (to be superseded by Ptr counterpart)
0055 Electron::Electron(const edm::RefToBase<reco::GsfElectron>& anElectronRef)
0056     : Lepton<reco::GsfElectron>(anElectronRef),
0057       embeddedGsfElectronCore_(false),
0058       embeddedGsfTrack_(false),
0059       embeddedSuperCluster_(false),
0060       embeddedPflowSuperCluster_(false),
0061       embeddedTrack_(false),
0062       embeddedSeedCluster_(false),
0063       embeddedRecHits_(false),
0064       embeddedPFCandidate_(false),
0065       ecalDrivenMomentum_(anElectronRef->p4()) {
0066   initImpactParameters();
0067 }
0068 
0069 /// constructor from a Ptr to a reco::GsfElectron
0070 Electron::Electron(const edm::Ptr<reco::GsfElectron>& anElectronRef)
0071     : Lepton<reco::GsfElectron>(anElectronRef),
0072       embeddedGsfElectronCore_(false),
0073       embeddedGsfTrack_(false),
0074       embeddedSuperCluster_(false),
0075       embeddedPflowSuperCluster_(false),
0076       embeddedTrack_(false),
0077       embeddedSeedCluster_(false),
0078       embeddedRecHits_(false),
0079       embeddedPFCandidate_(false),
0080       ecalDrivenMomentum_(anElectronRef->p4()) {
0081   initImpactParameters();
0082 }
0083 
0084 /// destructor
0085 Electron::~Electron() {}
0086 
0087 /// pipe operator (introduced to use pat::Electron with PFTopProjectors)
0088 std::ostream& reco::operator<<(std::ostream& out, const pat::Electron& obj) {
0089   if (!out)
0090     return out;
0091 
0092   out << "\tpat::Electron: ";
0093   out << std::setiosflags(std::ios::right);
0094   out << std::setiosflags(std::ios::fixed);
0095   out << std::setprecision(3);
0096   out << " E/pT/eta/phi " << obj.energy() << "/" << obj.pt() << "/" << obj.eta() << "/" << obj.phi();
0097   return out;
0098 }
0099 
0100 /// initializes the impact parameter container vars
0101 void Electron::initImpactParameters() {
0102   std::fill(ip_, ip_ + IpTypeSize, 0.0f);
0103   std::fill(eip_, eip_ + IpTypeSize, 0.0f);
0104   cachedIP_ = 0;
0105 }
0106 
0107 /// override the reco::GsfElectron::gsfTrack method, to access the internal storage of the supercluster
0108 reco::GsfTrackRef Electron::gsfTrack() const {
0109   if (embeddedGsfTrack_) {
0110     return reco::GsfTrackRef(&gsfTrack_, 0);
0111   } else {
0112     return reco::GsfElectron::gsfTrack();
0113   }
0114 }
0115 
0116 /// override the virtual reco::GsfElectron::core method, so that the embedded core can be used by GsfElectron client methods
0117 reco::GsfElectronCoreRef Electron::core() const {
0118   if (embeddedGsfElectronCore_) {
0119     return reco::GsfElectronCoreRef(&gsfElectronCore_, 0);
0120   } else {
0121     return reco::GsfElectron::core();
0122   }
0123 }
0124 
0125 /// override the reco::GsfElectron::superCluster method, to access the internal storage of the supercluster
0126 reco::SuperClusterRef Electron::superCluster() const {
0127   if (embeddedSuperCluster_) {
0128     if (embeddedSeedCluster_ || !basicClusters_.empty() || !preshowerClusters_.empty()) {
0129       if (!superClusterRelinked_.isSet()) {
0130         std::unique_ptr<std::vector<reco::SuperCluster> > sc(new std::vector<reco::SuperCluster>(superCluster_));
0131         if (embeddedSeedCluster_ && !(*sc)[0].seed().isAvailable()) {
0132           (*sc)[0].setSeed(seed());
0133         }
0134         if (!basicClusters_.empty() && !(*sc)[0].clusters().isAvailable()) {
0135           reco::CaloClusterPtrVector clusters;
0136           for (unsigned int iclus = 0; iclus < basicClusters_.size(); ++iclus) {
0137             clusters.push_back(reco::CaloClusterPtr(&basicClusters_, iclus));
0138           }
0139           (*sc)[0].setClusters(clusters);
0140         }
0141         if (!preshowerClusters_.empty() && !(*sc)[0].preshowerClusters().isAvailable()) {
0142           reco::CaloClusterPtrVector clusters;
0143           for (unsigned int iclus = 0; iclus < preshowerClusters_.size(); ++iclus) {
0144             clusters.push_back(reco::CaloClusterPtr(&preshowerClusters_, iclus));
0145           }
0146           (*sc)[0].setPreshowerClusters(clusters);
0147         }
0148         superClusterRelinked_.set(std::move(sc));
0149       }
0150       return reco::SuperClusterRef(&*superClusterRelinked_, 0);
0151     } else {
0152       return reco::SuperClusterRef(&superCluster_, 0);
0153     }
0154     //relink caloclusters if needed
0155     return reco::SuperClusterRef(&superCluster_, 0);
0156   } else {
0157     return reco::GsfElectron::superCluster();
0158   }
0159 }
0160 
0161 /// override the reco::GsfElectron::pflowSuperCluster method, to access the internal storage of the supercluster
0162 reco::SuperClusterRef Electron::parentSuperCluster() const {
0163   if (embeddedPflowSuperCluster_) {
0164     return reco::SuperClusterRef(&pflowSuperCluster_, 0);
0165   } else {
0166     return reco::GsfElectron::parentSuperCluster();
0167   }
0168 }
0169 
0170 /// direct access to the seed cluster
0171 reco::CaloClusterPtr Electron::seed() const {
0172   if (embeddedSeedCluster_) {
0173     return reco::CaloClusterPtr(&seedCluster_, 0);
0174   } else {
0175     return reco::GsfElectron::superCluster()->seed();
0176   }
0177 }
0178 
0179 /// override the reco::GsfElectron::closestCtfTrack method, to access the internal storage of the track
0180 reco::TrackRef Electron::closestCtfTrackRef() const {
0181   if (embeddedTrack_) {
0182     return reco::TrackRef(&track_, 0);
0183   } else {
0184     return reco::GsfElectron::closestCtfTrackRef();
0185   }
0186 }
0187 
0188 // the name of the method is misleading, users should use gsfTrack of closestCtfTrack
0189 reco::TrackRef Electron::track() const { return reco::TrackRef(); }
0190 
0191 /// Stores the electron's core (reco::GsfElectronCoreRef) internally
0192 void Electron::embedGsfElectronCore() {
0193   gsfElectronCore_.clear();
0194   if (reco::GsfElectron::core().isNonnull()) {
0195     gsfElectronCore_.push_back(*reco::GsfElectron::core());
0196     embeddedGsfElectronCore_ = true;
0197   }
0198 }
0199 
0200 /// Stores the electron's gsfTrack (reco::GsfTrackRef) internally
0201 void Electron::embedGsfTrack() {
0202   gsfTrack_.clear();
0203   if (reco::GsfElectron::gsfTrack().isNonnull()) {
0204     gsfTrack_.push_back(*reco::GsfElectron::gsfTrack());
0205     embeddedGsfTrack_ = true;
0206   }
0207 }
0208 
0209 /// Stores the electron's SuperCluster (reco::SuperClusterRef) internally
0210 void Electron::embedSuperCluster() {
0211   superCluster_.clear();
0212   if (reco::GsfElectron::superCluster().isNonnull()) {
0213     superCluster_.push_back(*reco::GsfElectron::superCluster());
0214     embeddedSuperCluster_ = true;
0215   }
0216 }
0217 
0218 /// Stores the electron's SuperCluster (reco::SuperClusterRef) internally
0219 void Electron::embedPflowSuperCluster() {
0220   pflowSuperCluster_.clear();
0221   if (reco::GsfElectron::parentSuperCluster().isNonnull()) {
0222     pflowSuperCluster_.push_back(*reco::GsfElectron::parentSuperCluster());
0223     embeddedPflowSuperCluster_ = true;
0224   }
0225 }
0226 
0227 /// Stores the electron's SeedCluster (reco::BasicClusterPtr) internally
0228 void Electron::embedSeedCluster() {
0229   seedCluster_.clear();
0230   if (reco::GsfElectron::superCluster().isNonnull() && reco::GsfElectron::superCluster()->seed().isNonnull()) {
0231     seedCluster_.push_back(*reco::GsfElectron::superCluster()->seed());
0232     embeddedSeedCluster_ = true;
0233   }
0234 }
0235 
0236 /// Stores the electron's BasicCluster (reco::CaloCluster) internally
0237 void Electron::embedBasicClusters() {
0238   basicClusters_.clear();
0239   if (reco::GsfElectron::superCluster().isNonnull()) {
0240     reco::CaloCluster_iterator itscl = reco::GsfElectron::superCluster()->clustersBegin();
0241     reco::CaloCluster_iterator itsclE = reco::GsfElectron::superCluster()->clustersEnd();
0242     for (; itscl != itsclE; ++itscl) {
0243       basicClusters_.push_back(**itscl);
0244     }
0245   }
0246 }
0247 
0248 /// Stores the electron's PreshowerCluster (reco::CaloCluster) internally
0249 void Electron::embedPreshowerClusters() {
0250   preshowerClusters_.clear();
0251   if (reco::GsfElectron::superCluster().isNonnull()) {
0252     reco::CaloCluster_iterator itscl = reco::GsfElectron::superCluster()->preshowerClustersBegin();
0253     reco::CaloCluster_iterator itsclE = reco::GsfElectron::superCluster()->preshowerClustersEnd();
0254     for (; itscl != itsclE; ++itscl) {
0255       preshowerClusters_.push_back(**itscl);
0256     }
0257   }
0258 }
0259 
0260 /// Stores the electron's PflowBasicCluster (reco::CaloCluster) internally
0261 void Electron::embedPflowBasicClusters() {
0262   pflowBasicClusters_.clear();
0263   if (reco::GsfElectron::parentSuperCluster().isNonnull()) {
0264     reco::CaloCluster_iterator itscl = reco::GsfElectron::parentSuperCluster()->clustersBegin();
0265     reco::CaloCluster_iterator itsclE = reco::GsfElectron::parentSuperCluster()->clustersEnd();
0266     for (; itscl != itsclE; ++itscl) {
0267       pflowBasicClusters_.push_back(**itscl);
0268     }
0269   }
0270 }
0271 
0272 /// Stores the electron's PflowPreshowerCluster (reco::CaloCluster) internally
0273 void Electron::embedPflowPreshowerClusters() {
0274   pflowPreshowerClusters_.clear();
0275   if (reco::GsfElectron::parentSuperCluster().isNonnull()) {
0276     reco::CaloCluster_iterator itscl = reco::GsfElectron::parentSuperCluster()->preshowerClustersBegin();
0277     reco::CaloCluster_iterator itsclE = reco::GsfElectron::parentSuperCluster()->preshowerClustersEnd();
0278     for (; itscl != itsclE; ++itscl) {
0279       pflowPreshowerClusters_.push_back(**itscl);
0280     }
0281   }
0282 }
0283 
0284 /// method to store the electron's track internally
0285 void Electron::embedTrack() {
0286   track_.clear();
0287   if (reco::GsfElectron::closestCtfTrackRef().isNonnull()) {
0288     track_.push_back(*reco::GsfElectron::closestCtfTrackRef());
0289     embeddedTrack_ = true;
0290   }
0291 }
0292 
0293 // method to store the RecHits internally
0294 void Electron::embedRecHits(const EcalRecHitCollection* rechits) {
0295   if (rechits != nullptr) {
0296     recHits_ = *rechits;
0297     embeddedRecHits_ = true;
0298   }
0299 }
0300 
0301 /// Returns a specific electron ID associated to the pat::Electron given its name
0302 /// For cut-based IDs, the value map has the following meaning:
0303 /// 0: fails,
0304 /// 1: passes electron ID only,
0305 /// 2: passes electron Isolation only,
0306 /// 3: passes electron ID and Isolation only,
0307 /// 4: passes conversion rejection,
0308 /// 5: passes conversion rejection and ID,
0309 /// 6: passes conversion rejection and Isolation,
0310 /// 7: passes the whole selection.
0311 /// For more details have a look at:
0312 /// https://twiki.cern.ch/twiki/bin/view/CMS/SimpleCutBasedEleID
0313 /// https://twiki.cern.ch/twiki/bin/view/CMS/SWGuideCategoryBasedElectronID
0314 /// Note: an exception is thrown if the specified ID is not available
0315 float Electron::electronID(const std::string& name) const {
0316   for (std::vector<IdPair>::const_iterator it = electronIDs_.begin(), ed = electronIDs_.end(); it != ed; ++it) {
0317     if (it->first == name)
0318       return it->second;
0319   }
0320   cms::Exception ex("Key not found");
0321   ex << "pat::Electron: the ID " << name << " can't be found in this pat::Electron.\n";
0322   ex << "The available IDs are: ";
0323   for (std::vector<IdPair>::const_iterator it = electronIDs_.begin(), ed = electronIDs_.end(); it != ed; ++it) {
0324     ex << "'" << it->first << "' ";
0325   }
0326   ex << ".\n";
0327   throw ex;
0328 }
0329 
0330 /// Checks if a specific electron ID is associated to the pat::Electron.
0331 bool Electron::isElectronIDAvailable(const std::string& name) const {
0332   for (std::vector<IdPair>::const_iterator it = electronIDs_.begin(), ed = electronIDs_.end(); it != ed; ++it) {
0333     if (it->first == name)
0334       return true;
0335   }
0336   return false;
0337 }
0338 
0339 /// reference to the source PFCandidates
0340 reco::PFCandidateRef Electron::pfCandidateRef() const {
0341   if (embeddedPFCandidate_) {
0342     return reco::PFCandidateRef(&pfCandidate_, 0);
0343   } else {
0344     return pfCandidateRef_;
0345   }
0346 }
0347 
0348 /// Stores the PFCandidate pointed to by pfCandidateRef_ internally
0349 void Electron::embedPFCandidate() {
0350   pfCandidate_.clear();
0351   if (pfCandidateRef_.isAvailable() && pfCandidateRef_.isNonnull()) {
0352     pfCandidate_.push_back(*pfCandidateRef_);
0353     embeddedPFCandidate_ = true;
0354   }
0355 }
0356 
0357 /// Returns the reference to the parent PF candidate with index i.
0358 /// For use in TopProjector.
0359 reco::CandidatePtr Electron::sourceCandidatePtr(size_type i) const {
0360   if (pfCandidateRef_.isNonnull()) {
0361     if (i == 0) {
0362       return reco::CandidatePtr(edm::refToPtr(pfCandidateRef_));
0363     } else {
0364       i--;
0365     }
0366   }
0367   if (i >= associatedPackedFCandidateIndices_.size()) {
0368     return reco::CandidatePtr();
0369   } else {
0370     return reco::CandidatePtr(edm::refToPtr(
0371         edm::Ref<pat::PackedCandidateCollection>(packedPFCandidates_, associatedPackedFCandidateIndices_[i])));
0372   }
0373 }
0374 
0375 /// dB gives the impact parameter wrt the beamline.
0376 /// If this is not cached it is not meaningful, since
0377 /// it relies on the distance to the beamline.
0378 ///
0379 /// IpType defines the type of the impact parameter
0380 /// None is default and reverts to the old functionality.
0381 ///
0382 /// Example: electron->dB(pat::Electron::PV2D)
0383 /// will return the electron transverse impact parameter
0384 /// relative to the primary vertex.
0385 double Electron::dB(IpType type_) const {
0386   // more IP types (new)
0387   if (cachedIP_ & (1 << int(type_))) {
0388     return ip_[type_];
0389   } else {
0390     return std::numeric_limits<double>::max();
0391   }
0392 }
0393 
0394 /// edB gives the uncertainty on the impact parameter wrt the beamline.
0395 /// If this is not cached it is not meaningful, since
0396 /// it relies on the distance to the beamline.
0397 ///
0398 /// IpType defines the type of the impact parameter
0399 /// None is default and reverts to the old functionality.
0400 ///
0401 /// Example: electron->edB(pat::Electron::PV2D)
0402 /// will return the electron transverse impact parameter uncertainty
0403 /// relative to the primary vertex.
0404 double Electron::edB(IpType type_) const {
0405   // more IP types (new)
0406   if (cachedIP_ & (1 << int(type_))) {
0407     return eip_[type_];
0408   } else {
0409     return std::numeric_limits<double>::max();
0410   }
0411 }
0412 
0413 /// Sets the impact parameter and its error wrt the beamline and caches it.
0414 void Electron::setDB(double dB, double edB, IpType type) {
0415   ip_[type] = dB;
0416   eip_[type] = edB;
0417   cachedIP_ |= (1 << int(type));
0418 }
0419 
0420 /// Set additional missing mva input variables for new mva ID (71X update)
0421 void Electron::setMvaVariables(double sigmaIetaIphi, double ip3d) {
0422   sigmaIetaIphi_ = sigmaIetaIphi;
0423   ip3d_ = ip3d;
0424 }
0425 
0426 edm::RefVector<pat::PackedCandidateCollection> Electron::associatedPackedPFCandidates() const {
0427   edm::RefVector<pat::PackedCandidateCollection> ret(packedPFCandidates_.id());
0428   for (uint16_t idx : associatedPackedFCandidateIndices_) {
0429     ret.push_back(edm::Ref<pat::PackedCandidateCollection>(packedPFCandidates_, idx));
0430   }
0431   return ret;
0432 }