CMS 3D CMS Logo

EgHLTOffEle.h

Go to the documentation of this file.
00001 #ifndef DQMOFFLINE_TRIGGER_EGHLTOFFELE
00002 #define DQMOFFLINE_TRIGGER_EGHLTOFFELE
00003 
00004 //class: EgHLTOffEle
00005 //
00006 //author: Sam Harper (July 2008)
00007 //
00008 //
00009 //aim: to allow easy access to electron ID variables
00010 //     currently the CMSSW electron classes are a mess with key electron selection variables not being accessable from GsfElectron
00011 //     this a stop gap to produce a simple electron class with all variables easily accessable via methods 
00012 //     note as this is meant for HLT Offline DQM, I do not want the overhead of converting to pat
00013 //
00014 //implimentation: aims to be a wrapper for GsfElectron methods, it is hoped that in time these methods will be directly added to GsfElectron and so
00015 //                make this class obsolute
00016 //                unfortunately can not be a pure wrapper as needs to store isol and cluster shape
00017 //
00018 
00019 
00020 #include "DataFormats/EgammaCandidates/interface/GsfElectronFwd.h"
00021 #include "DataFormats/EgammaReco/interface/ClusterShapeFwd.h"
00022 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
00023 #include "DataFormats/EgammaReco/interface/ClusterShape.h"
00024 
00025 #include "DQMOffline/Trigger/interface/CutCodes.h"
00026 #include "DQMOffline/Trigger/interface/TrigCodes.h"
00027 
00028 class EgHLTOffEle { 
00029   
00030  public:
00031   //helper struct to store the isolations
00032   struct IsolData {
00033     float nrTrks;
00034     float ptTrks;
00035     float em;
00036     float had;
00037   };
00038 
00039 public:
00040   //helper struct to store the cluster shapes
00041   struct ClusShapeData {
00042     float sigmaEtaEta;
00043     //float sigmaIEtaIEta;
00044     //float e2x5MaxOver5x5;
00045     float sigmaPhiPhi;
00046     //float sigmaIPhiIPhi;
00047   };
00048 
00049 
00050  private:
00051   const reco::GsfElectron* gsfEle_; //pointers to the underlying electron (we do not own this)
00052   // const reco::ClusterShape* clusShape_; //pointers to the underlying cluster shape (we do not own this and it may be null)
00053   ClusShapeData clusShapeData_;
00054   IsolData isolData_;
00055 
00056   //these are bit-packed words telling me which cuts the electron fail (ie 0x0 is passed all cuts)
00057   int tagCutCode_;
00058   int probeCutCode_;
00059   int cutCode_;
00060   
00061   //and these are the trigger bits stored
00062   //note that the trigger bits are defined at the begining of each job
00063   //and do not necessaryly map between jobs
00064   TrigCodes::TrigBitSet trigBits_;
00065 
00066  public:
00067   
00068  EgHLTOffEle(const reco::GsfElectron& ele,const ClusShapeData& shapeData,const IsolData& isolData):
00069    gsfEle_(&ele),clusShapeData_(shapeData),isolData_(isolData),
00070    tagCutCode_(CutCodes::INVALID),probeCutCode_(CutCodes::INVALID),cutCode_(CutCodes::INVALID){}
00071   ~EgHLTOffEle(){}
00072   
00073   //modifiers  
00074   void setTagCutCode(int code){tagCutCode_=code;}
00075   void setProbeCutCode(int code){probeCutCode_=code;} 
00076   void setCutCode(int code){cutCode_=code;}
00077   void setTrigBits(TrigCodes::TrigBitSet bits){trigBits_=bits;}
00078 
00079   //kinematic and geometric methods
00080   float et()const{return gsfEle_->et();}
00081   float energy()const{return gsfEle_->energy();}
00082   float eta()const{return gsfEle_->eta();}
00083   float phi()const{return gsfEle_->phi();}
00084   float etSC()const{return gsfEle_->superCluster()->position().rho()/gsfEle_->superCluster()->position().r()*caloEnergy();}
00085   float caloEnergy()const{return gsfEle_->caloEnergy();}
00086   float etaSC()const{return gsfEle_->superCluster()->eta();}
00087   float detEta()const{return etaSC();}
00088   float phiSC()const{return gsfEle_->superCluster()->phi();}
00089   float zVtx()const{return gsfEle_->TrackPositionAtVtx().z();}
00090   const math::XYZTLorentzVector& p4()const{return gsfEle_->p4();}
00091 
00092   //classification (couldnt they have just named it 'type')
00093   int classification()const{return gsfEle_->classification();}
00094 
00095   //track methods
00096   int charge()const{return gsfEle_->charge();}
00097   float pVtx()const{return gsfEle_->trackMomentumAtVtx().R();}
00098   float pCalo()const{return gsfEle_->trackMomentumAtCalo().R();}
00099   float ptVtx()const{return gsfEle_->trackMomentumAtVtx().rho();}
00100   float ptCalo()const{return gsfEle_->trackMomentumAtCalo().rho();}
00101  
00102   
00103   //abreviations of overly long GsfElectron methods, I'm sorry but if you cant figure out what hOverE() means, you shouldnt be using this class
00104   float hOverE()const{return gsfEle_->hadronicOverEm();}
00105   float dEtaIn()const{return gsfEle_->deltaEtaSuperClusterTrackAtVtx();}
00106   float dPhiIn()const{return gsfEle_->deltaPhiSuperClusterTrackAtVtx();}
00107   float dPhiOut()const{return gsfEle_->deltaPhiSeedClusterTrackAtCalo();}
00108   float epIn()const{return gsfEle_->eSuperClusterOverP();}
00109   float epOut()const{return gsfEle_->eSeedClusterOverPout();}
00110   
00111   //variables with no direct method
00112   float sigmaEtaEta()const;
00113   float sigmaEtaEtaUnCorr()const{return clusShapeData_.sigmaEtaEta;}
00114   //float sigmaIEtaIEta()const;                                         
00115   float sigmaPhiPhi()const{return clusShapeData_.sigmaPhiPhi;}
00116   //float sigmaIPhiIPhi()const{return clusShapeData_.sigmaIPhiIPhi;}
00117   //float e2x5MaxOver5x5()const{return clusShapeData_.e2x5MaxOver5x5;}
00118   //float sigmaPhiPhi()const{return clusShape_!=NULL ? sqrt(clusShape_->covPhiPhi()) : 999;}
00119   float bremFrac()const{return (pVtx()-pCalo())/pVtx();}
00120   float invEOverInvP()const{return 1./gsfEle_->caloEnergy() - 1./gsfEle_->trackMomentumAtVtx().R();}
00121   //float e9OverE25()const{return clusShape_!=NULL ? clusShape_->e3x3()/clusShape_->e5x5() : -999;}
00122   
00123   //isolation
00124   float isolEm()const{return isolData_.em;}
00125   float isolHad()const{return isolData_.had;}
00126   float isolNrTrks()const{return isolData_.nrTrks;}
00127   float isolPtTrks()const{return isolData_.ptTrks;}
00128 
00129   //selection cuts
00130   int tagCutCode()const{return tagCutCode_;}
00131   int probeCutCode()const{return probeCutCode_;}
00132   int cutCode()const{return cutCode_;}
00133 
00134   //trigger
00135   TrigCodes::TrigBitSet trigBits()const{return trigBits_;}
00136   
00137 
00138 };
00139 
00140 #endif

Generated on Tue Jun 9 17:34:09 2009 for CMSSW by  doxygen 1.5.4