Electron

IPTYPE

Macros

Line Code
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381
//
//

#ifndef DataFormats_PatCandidates_Electron_h
#define DataFormats_PatCandidates_Electron_h

/**
  \class    pat::Electron Electron.h "DataFormats/PatCandidates/interface/Electron.h"
  \brief    Analysis-level electron class

   pat::Electron implements the analysis-level electron class within the
   'pat' namespace.

   Please post comments and questions to the Physics Tools hypernews:
   https://hypernews.cern.ch/HyperNews/CMS/get/physTools.html

  \author   Steven Lowette, Giovanni Petrucciani, Frederic Ronga
*/

#include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
#include "DataFormats/EgammaCandidates/interface/GsfElectronFwd.h"
#include "DataFormats/EgammaCandidates/interface/GsfElectronCore.h"
#include "DataFormats/EgammaCandidates/interface/GsfElectronCoreFwd.h"
#include "DataFormats/EgammaReco/interface/BasicCluster.h"
#include "DataFormats/GsfTrackReco/interface/GsfTrack.h"
#include "DataFormats/PatCandidates/interface/Lepton.h"

#include "DataFormats/EcalRecHit/interface/EcalRecHitCollections.h"
#include "DataFormats/ParticleFlowCandidate/interface/PFCandidateFwd.h"
#include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
#include "DataFormats/PatCandidates/interface/PackedCandidate.h"
#include "DataFormats/Common/interface/AtomicPtrCache.h"

// Define typedefs for convenience
namespace pat {
  class Electron;
  typedef std::vector<Electron> ElectronCollection;
  typedef edm::Ref<ElectronCollection> ElectronRef;
  typedef edm::RefVector<ElectronCollection> ElectronRefVector;
}  // namespace pat

namespace reco {
  /// pipe operator (introduced to use pat::Electron with PFTopProjectors)
  std::ostream& operator<<(std::ostream& out, const pat::Electron& obj);
}  // namespace reco

// Class definition
namespace pat {
  class PATElectronSlimmer;

  class Electron : public Lepton<reco::GsfElectron> {
  public:
    typedef std::pair<std::string, float> IdPair;

    /// default constructor
    Electron();
    /// constructor from reco::GsfElectron
    Electron(const reco::GsfElectron& anElectron);
    /// constructor from a RefToBase to a reco::GsfElectron (to be superseded by Ptr counterpart)
    Electron(const edm::RefToBase<reco::GsfElectron>& anElectronRef);
    /// constructor from a Ptr to a reco::GsfElectron
    Electron(const edm::Ptr<reco::GsfElectron>& anElectronRef);
    /// destructor
    ~Electron() override;

    /// required reimplementation of the Candidate's clone method
    Electron* clone() const override { return new Electron(*this); }

    // ---- methods for content embedding ----
    /// override the virtual reco::GsfElectron::core method, so that the embedded core can be used by GsfElectron client methods
    reco::GsfElectronCoreRef core() const override;
    /// override the reco::GsfElectron::gsfTrack method, to access the internal storage of the supercluster
    reco::GsfTrackRef gsfTrack() const override;
    /// override the reco::GsfElectron::superCluster method, to access the internal storage of the supercluster
    reco::SuperClusterRef superCluster() const override;
    /// override the reco::GsfElectron::pflowSuperCluster method, to access the internal storage of the pflowSuperCluster
    reco::SuperClusterRef parentSuperCluster() const override;
    /// returns nothing. Use either gsfTrack or closestCtfTrack
    reco::TrackRef track() const override;
    /// override the reco::GsfElectron::closestCtfTrackRef method, to access the internal storage of the track
    reco::TrackRef closestCtfTrackRef() const override;
    /// direct access to the seed cluster
    reco::CaloClusterPtr seed() const;

    //method to access the basic clusters
    const std::vector<reco::CaloCluster>& basicClusters() const { return basicClusters_; }
    //method to access the preshower clusters
    const std::vector<reco::CaloCluster>& preshowerClusters() const { return preshowerClusters_; }
    //method to access the pflow basic clusters
    const std::vector<reco::CaloCluster>& pflowBasicClusters() const { return pflowBasicClusters_; }
    //method to access the pflow preshower clusters
    const std::vector<reco::CaloCluster>& pflowPreshowerClusters() const { return pflowPreshowerClusters_; }

    using reco::RecoCandidate::track;  // avoid hiding the base implementation
    /// method to store the electron's core internally
    void embedGsfElectronCore();
    /// method to store the electron's GsfTrack internally
    void embedGsfTrack();
    /// method to store the electron's SuperCluster internally
    void embedSuperCluster();
    /// method to store the electron's PflowSuperCluster internally
    void embedPflowSuperCluster();
    /// method to store the electron's seedcluster internally
    void embedSeedCluster();
    /// method to store the electron's basic clusters
    void embedBasicClusters();
    /// method to store the electron's preshower clusters
    void embedPreshowerClusters();
    /// method to store the electron's pflow basic clusters
    void embedPflowBasicClusters();
    /// method to store the electron's pflow preshower clusters
    void embedPflowPreshowerClusters();
    /// method to store the electron's Track internally
    void embedTrack();
    /// method to store the RecHits internally - can be called from the PATElectronProducer
    void embedRecHits(const EcalRecHitCollection* rechits);

    // ---- methods for electron ID ----
    /// Returns a specific electron ID associated to the pat::Electron given its name
    // For cut-based IDs, the value map has the following meaning:
    // 0: fails,
    // 1: passes electron ID only,
    // 2: passes electron Isolation only,
    // 3: passes electron ID and Isolation only,
    // 4: passes conversion rejection,
    // 5: passes conversion rejection and ID,
    // 6: passes conversion rejection and Isolation,
    // 7: passes the whole selection.
    // For more details have a look at:
    // https://twiki.cern.ch/twiki/bin/view/CMS/SimpleCutBasedEleID
    // https://twiki.cern.ch/twiki/bin/view/CMS/SWGuideCategoryBasedElectronID
    // Note: an exception is thrown if the specified ID is not available
    float electronID(const std::string& name) const;
    float electronID(const char* name) const { return electronID(std::string(name)); }
    /// Returns true if a specific ID is available in this pat::Electron
    bool isElectronIDAvailable(const std::string& name) const;
    bool isElectronIDAvailable(const char* name) const { return isElectronIDAvailable(std::string(name)); }
    /// Returns all the electron IDs in the form of <name,value> pairs. The 'default' ID is the first in the list
    const std::vector<IdPair>& electronIDs() const { return electronIDs_; }
    /// Store multiple electron ID values, discarding existing ones. The first one in the list becomes the 'default' electron id
    void setElectronIDs(const std::vector<IdPair>& ids) { electronIDs_ = ids; }

    // ---- overload of isolation functions ----
    /// Overload of pat::Lepton::trackIso(); returns the value of the summed track pt in a cone of deltaR<0.4
    float trackIso() const { return dr04TkSumPt(); }
    /// Overload of pat::Lepton::ecalIso(); returns the value of the summed Et of all recHits in the ecal in a cone of deltaR<0.4
    float ecalIso() const { return dr04EcalRecHitSumEt(); }
    /// Overload of pat::Lepton::hcalIso(); returns the value of the summed Et of all caloTowers in the hcal in a cone of deltaR<0.4
    float hcalIso() const { return dr04HcalTowerSumEt(); }
    /// Overload of pat::Lepton::caloIso(); returns the sum of ecalIso() and hcalIso
    float caloIso() const { return ecalIso() + hcalIso(); }
    /// returns PUPPI isolations
    float puppiChargedHadronIso() const { return puppiChargedHadronIso_; }
    float puppiNeutralHadronIso() const { return puppiNeutralHadronIso_; }
    float puppiPhotonIso() const { return puppiPhotonIso_; }
    /// returns PUPPINoLeptons isolations
    float puppiNoLeptonsChargedHadronIso() const { return puppiNoLeptonsChargedHadronIso_; }
    float puppiNoLeptonsNeutralHadronIso() const { return puppiNoLeptonsNeutralHadronIso_; }
    float puppiNoLeptonsPhotonIso() const { return puppiNoLeptonsPhotonIso_; }
    /// sets PUPPI isolations
    void setIsolationPUPPI(float chargedhadrons_, float neutralhadrons_, float photons_) {
      puppiChargedHadronIso_ = chargedhadrons_;
      puppiNeutralHadronIso_ = neutralhadrons_;
      puppiPhotonIso_ = photons_;
    }
    /// sets PUPPINoLeptons isolations
    void setIsolationPUPPINoLeptons(float chargedhadrons_, float neutralhadrons_, float photons_) {
      puppiNoLeptonsChargedHadronIso_ = chargedhadrons_;
      puppiNoLeptonsNeutralHadronIso_ = neutralhadrons_;
      puppiNoLeptonsPhotonIso_ = photons_;
    }
    // ---- PF specific methods ----
    bool isPF() const { return isPF_; }
    void setIsPF(bool hasPFCandidate) { isPF_ = hasPFCandidate; }

    /// reference to the source PFCandidates; null if this has been built from a standard electron
    reco::PFCandidateRef pfCandidateRef() const;
    /// add a reference to the source IsolatedPFCandidate
    void setPFCandidateRef(const reco::PFCandidateRef& ref) { pfCandidateRef_ = ref; }
    /// embed the PFCandidate pointed to by pfCandidateRef_
    void embedPFCandidate();
    /// get the number of non-null PFCandidates
    size_t numberOfSourceCandidatePtrs() const override {
      return (pfCandidateRef_.isNonnull() ? 1 : 0) + associatedPackedFCandidateIndices_.size();
    }
    /// get the source candidate pointer with index i
    reco::CandidatePtr sourceCandidatePtr(size_type i) const override;

    // ---- embed various impact parameters with errors ----
    typedef enum IPTYPE { PV2D = 0, PV3D = 1, BS2D = 2, BS3D = 3, PVDZ = 4, IpTypeSize = 5 } IpType;
    /// Impact parameter wrt primary vertex or beamspot
    double dB(IPTYPE type) const;
    /// Uncertainty on the corresponding impact parameter
    double edB(IPTYPE type) const;
    /// the version without arguments returns PD2D, but with an absolute value (for backwards compatibility)
    double dB() const { return std::abs(dB(PV2D)); }
    /// the version without arguments returns PD2D, but with an absolute value (for backwards compatibility)
    double edB() const { return std::abs(edB(PV2D)); }
    /// Set impact parameter of a certain type and its uncertainty
    void setDB(double dB, double edB, IPTYPE type);

    // ---- Momentum estimate specific methods ----
    const LorentzVector& ecalDrivenMomentum() const { return ecalDrivenMomentum_; }
    void setEcalDrivenMomentum(const Candidate::LorentzVector& mom) { ecalDrivenMomentum_ = mom; }

    /// pipe operator (introduced to use pat::Electron with PFTopProjectors)
    friend std::ostream& reco::operator<<(std::ostream& out, const pat::Electron& obj);

    /// additional mva input variables
    /// sigmaIEtaIPhi
    float sigmaIetaIphi() const { return sigmaIetaIphi_; }
    /// sigmaIEtaIPhi (from full 5x5 non-ZS clusters without fractions, a la 5.3.X)
    float full5x5_sigmaIetaIphi() const { return full5x5_sigmaIetaIphi_; }
    /// ip3d
    double ip3d() const { return ip3d_; }
    /// set missing mva input variables
    void setMvaVariables(double sigmaIetaIphi, double ip3d);
    void full5x5_setSigmaIetaIphi(float sigmaIetaIphi) { full5x5_sigmaIetaIphi_ = sigmaIetaIphi; }

    const EcalRecHitCollection* recHits() const { return &recHits_; }

    /// additional regression variables
    /// regression1
    double ecalRegressionEnergy() const { return ecalRegressionEnergy_; }
    double ecalRegressionError() const { return ecalRegressionError_; }
    /// regression2
    double ecalTrackRegressionEnergy() const { return ecalTrackRegressionEnergy_; }
    double ecalTrackRegressionError() const { return ecalTrackRegressionError_; }
    /// set regression1
    void setEcalRegressionEnergy(double val, double err) {
      ecalRegressionEnergy_ = val;
      ecalRegressionError_ = err;
    }
    /// set regression2
    void setEcalTrackRegressionEnergy(double val, double err) {
      ecalTrackRegressionEnergy_ = val;
      ecalTrackRegressionError_ = err;
    }

    /// set scale corrections / smearings
    void setEcalScale(double val) { ecalScale_ = val; }
    void setEcalSmear(double val) { ecalSmear_ = val; }
    void setEcalRegressionScale(double val) { ecalRegressionScale_ = val; }
    void setEcalRegressionSmear(double val) { ecalRegressionSmear_ = val; }
    void setEcalTrackRegressionScale(double val) { ecalTrackRegressionScale_ = val; }
    void setEcalTrackRegressionSmear(double val) { ecalTrackRegressionSmear_ = val; }

    /// get scale corrections /smearings
    double ecalScale() const { return ecalScale_; }
    double ecalSmear() const { return ecalSmear_; }
    double ecalRegressionScale() const { return ecalRegressionScale_; }
    double ecalRegressionSmear() const { return ecalRegressionSmear_; }
    double ecalTrackRegressionScale() const { return ecalTrackRegressionScale_; }
    double ecalTrackRegressionSmear() const { return ecalTrackRegressionSmear_; }
    /// vertex fit combined with missing number of hits method
    bool passConversionVeto() const { return passConversionVeto_; }
    void setPassConversionVeto(bool flag) { passConversionVeto_ = flag; }

    /// References to PFCandidates linked to this object (e.g. for isolation vetos or masking before jet reclustering)
    edm::RefVector<pat::PackedCandidateCollection> associatedPackedPFCandidates() const;
    /// References to PFCandidates linked to this object (e.g. for isolation vetos or masking before jet reclustering)
    template <typename T>
    void setAssociatedPackedPFCandidates(const edm::RefProd<pat::PackedCandidateCollection>& refprod,
                                         T beginIndexItr,
                                         T endIndexItr) {
      packedPFCandidates_ = refprod;
      associatedPackedFCandidateIndices_.clear();
      associatedPackedFCandidateIndices_.insert(associatedPackedFCandidateIndices_.end(), beginIndexItr, endIndexItr);
    }

    friend class PATElectronSlimmer;

  protected:
    /// init impact parameter defaults (for use in a constructor)
    void initImpactParameters();

    // ---- for content embedding ----
    /// True if electron's gsfElectronCore is stored internally
    bool embeddedGsfElectronCore_;
    /// Place to store electron's gsfElectronCore internally
    std::vector<reco::GsfElectronCore> gsfElectronCore_;
    /// True if electron's gsfTrack is stored internally
    bool embeddedGsfTrack_;
    /// Place to store electron's gsfTrack internally
    std::vector<reco::GsfTrack> gsfTrack_;
    /// True if electron's supercluster is stored internally
    bool embeddedSuperCluster_;
    /// True if electron's pflowsupercluster is stored internally
    bool embeddedPflowSuperCluster_;
    /// Place to store electron's supercluster internally
    std::vector<reco::SuperCluster> superCluster_;
    /// Place to temporarily store the electron's supercluster after relinking the seed to it
    edm::AtomicPtrCache<std::vector<reco::SuperCluster> > superClusterRelinked_;
    /// Place to store electron's basic clusters internally
    std::vector<reco::CaloCluster> basicClusters_;
    /// Place to store electron's preshower clusters internally
    std::vector<reco::CaloCluster> preshowerClusters_;
    /// Place to store electron's pflow basic clusters internally
    std::vector<reco::CaloCluster> pflowBasicClusters_;
    /// Place to store electron's pflow preshower clusters internally
    std::vector<reco::CaloCluster> pflowPreshowerClusters_;
    /// Place to store electron's pflow supercluster internally
    std::vector<reco::SuperCluster> pflowSuperCluster_;
    /// True if electron's track is stored internally
    bool embeddedTrack_;
    /// Place to store electron's track internally
    std::vector<reco::Track> track_;
    /// True if seed cluster is stored internally
    bool embeddedSeedCluster_;
    /// Place to store electron's seed cluster internally
    std::vector<reco::CaloCluster> seedCluster_;
    /// True if RecHits stored internally
    bool embeddedRecHits_;
    /// Place to store electron's RecHits internally (5x5 around seed+ all RecHits)
    EcalRecHitCollection recHits_;

    // ---- electron ID's holder ----
    /// Electron IDs
    std::vector<IdPair> electronIDs_;

    // ---- PF specific members ----
    bool isPF_;
    /// true if the IsolatedPFCandidate is embedded
    bool embeddedPFCandidate_;
    /// A copy of the source IsolatedPFCandidate is stored in this vector if embeddedPFCandidate_ if True
    reco::PFCandidateCollection pfCandidate_;
    /// reference to the IsolatedPFCandidate this has been built from; null if this has been built from a standard electron
    reco::PFCandidateRef pfCandidateRef_;

    // ---- specific members : Momentum estimates ----
    /// ECAL-driven momentum
    LorentzVector ecalDrivenMomentum_;

    /// additional missing mva variables : 14/04/2012
    float sigmaIetaIphi_, full5x5_sigmaIetaIphi_;
    double ip3d_;

    /// output of regression
    double ecalRegressionEnergy_;
    double ecalTrackRegressionEnergy_;
    double ecalRegressionError_;
    double ecalTrackRegressionError_;

    /// scale corrections and smearing applied or to be be applied. Initialized to -99999.
    double ecalScale_;
    double ecalSmear_;

    double ecalRegressionScale_;
    double ecalRegressionSmear_;

    double ecalTrackRegressionScale_;
    double ecalTrackRegressionSmear_;

    /// PUPPI isolations
    float puppiChargedHadronIso_;
    float puppiNeutralHadronIso_;
    float puppiPhotonIso_;

    /// PUPPINoLeptons isolations
    float puppiNoLeptonsChargedHadronIso_;
    float puppiNoLeptonsNeutralHadronIso_;
    float puppiNoLeptonsPhotonIso_;

    /// conversion veto
    bool passConversionVeto_;

    // ---- cached impact parameters ----
    /// True if the IP (former dB) has been cached
    uint8_t cachedIP_;
    /// Impact parameter at the primary vertex,
    float ip_[IpTypeSize];
    /// Impact parameter uncertainty as recommended by the tracking group
    float eip_[IpTypeSize];

    // ---- link to PackedPFCandidates
    edm::RefProd<pat::PackedCandidateCollection> packedPFCandidates_;
    std::vector<uint16_t> associatedPackedFCandidateIndices_;
  };
}  // namespace pat

#endif