CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/RecoEcal/EgammaCoreTools/interface/EcalClusterLazyTools.h

Go to the documentation of this file.
00001 #ifndef RecoEcal_EgammaCoreTools_EcalClusterLazyTools_h
00002 #define RecoEcal_EgammaCoreTools_EcalClusterLazyTools_h
00003 
00014 #include "FWCore/Framework/interface/Event.h"
00015 #include "DataFormats/EgammaReco/interface/BasicCluster.h"
00016 #include "DataFormats/EgammaReco/interface/SuperCluster.h"
00017 #include "DataFormats/EcalRecHit/interface/EcalRecHitCollections.h"
00018 
00019 #include "FWCore/Framework/interface/ESHandle.h"
00020 
00021 #include "CondFormats/EcalObjects/interface/EcalIntercalibConstants.h"
00022 #include "CondFormats/EcalObjects/interface/EcalADCToGeVConstant.h"
00023 #include "CalibCalorimetry/EcalLaserCorrection/interface/EcalLaserDbService.h"
00024 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00025 #include "RecoLocalCalo/EcalRecAlgos/interface/EcalSeverityLevelAlgoRcd.h"
00026 #include "RecoLocalCalo/EcalRecAlgos/interface/EcalSeverityLevelAlgo.h"
00027 
00028 
00029 class CaloTopology;
00030 class CaloGeometry;
00031 class CaloSubdetectorTopology;
00032 
00033 class EcalClusterLazyTools {
00034         public:
00035        EcalClusterLazyTools( const edm::Event &ev, const edm::EventSetup &es, edm::InputTag redEBRecHits, edm::InputTag redEERecHits,const edm::ParameterSet& config );
00036        EcalClusterLazyTools( const edm::Event &ev, const edm::EventSetup &es, edm::InputTag redEBRecHits, edm::InputTag redEERecHits);
00037        EcalClusterLazyTools( const edm::Event &ev, const edm::EventSetup &es, edm::InputTag redEBRecHits, edm::InputTag redEERecHits, edm::InputTag redESRecHits);
00038 
00039                 ~EcalClusterLazyTools();
00040 
00041                 // various energies in the matrix nxn surrounding the maximum energy crystal of the input cluster  
00042                 //NOTE (29/10/08): we now use an eta/phi coordinate system rather than phi/eta
00043                 //to minmise possible screwups, for now e5x1 isnt defined all the majority of people who call it actually want e1x5 and 
00044                 //it is thought it is better that their code doesnt compile rather than pick up the wrong function
00045                 //therefore in this version and later e1x5 = e5x1 in the old version 
00046                 //so 1x5 is 1 crystal in eta and 5 crystals in phi
00047                 //note e3x2 does not have a definate eta/phi geometry, it takes the maximum 3x2 block containing the 
00048                 //seed regardless of whether that 3 in eta or phi
00049                 float e1x3( const reco::BasicCluster &cluster );
00050                 float e1x3( const reco::BasicCluster &cluster, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv );
00051 
00052                 float e3x1( const reco::BasicCluster &cluster );
00053                 float e3x1( const reco::BasicCluster &cluster, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv );
00054 
00055                 float e1x5( const reco::BasicCluster &cluster );
00056                 float e1x5( const reco::BasicCluster &cluster, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv );
00057 
00058                 float e5x1( const reco::BasicCluster &cluster );
00059                 float e5x1( const reco::BasicCluster &cluster, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv );
00060 
00061                 float e2x2( const reco::BasicCluster &cluster );
00062                 float e2x2( const reco::BasicCluster &cluster, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv );
00063 
00064                 float e3x2( const reco::BasicCluster &cluster );
00065                 float e3x2( const reco::BasicCluster &cluster, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv );
00066 
00067                 float e3x3( const reco::BasicCluster &cluster );
00068                 float e3x3( const reco::BasicCluster &cluster, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv );
00069 
00070                 float e4x4( const reco::BasicCluster &cluster );
00071                 float e4x4( const reco::BasicCluster &cluster, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv );
00072 
00073                 float e5x5( const reco::BasicCluster &cluster );
00074                 float e5x5( const reco::BasicCluster &cluster, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv );
00075                 // energy in the 2x5 strip right of the max crystal (does not contain max crystal)
00076                 // 2 crystals wide in eta, 5 wide in phi.
00077                 float e2x5Right( const reco::BasicCluster &cluster );
00078                 float e2x5Right( const reco::BasicCluster &cluster, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv );
00079                 // energy in the 2x5 strip left of the max crystal (does not contain max crystal)
00080                 float e2x5Left( const reco::BasicCluster &cluster );
00081                 float e2x5Left( const reco::BasicCluster &cluster, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv );
00082                 // energy in the 5x2 strip above the max crystal (does not contain max crystal)
00083                 // 5 crystals wide in eta, 2 wide in phi.
00084                 float e2x5Top( const reco::BasicCluster &cluster );
00085                 float e2x5Top( const reco::BasicCluster &cluster, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv );
00086 
00087                 // energy in the 5x2 strip below the max crystal (does not contain max crystal)
00088                 float e2x5Bottom( const reco::BasicCluster &cluster );
00089                 float e2x5Bottom( const reco::BasicCluster &cluster, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv );
00090                 // energy in a 2x5 strip containing the seed (max) crystal.
00091                 // 2 crystals wide in eta, 5 wide in phi.
00092                 // it is the maximum of either (1x5left + 1x5center) or (1x5right + 1x5center)
00093                 float e2x5Max( const reco::BasicCluster &cluster );
00094                 float e2x5Max( const reco::BasicCluster &cluster, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv );
00095                 // energies in the crystal left, right, top, bottom w.r.t. to the most energetic crystal
00096                 float eLeft( const reco::BasicCluster &cluster );
00097                 float eLeft( const reco::BasicCluster &cluster, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv );
00098 
00099                 float eRight( const reco::BasicCluster &cluster );
00100                 float eRight( const reco::BasicCluster &cluster, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv );
00101 
00102                 float eTop( const reco::BasicCluster &cluster );
00103                 float eTop( const reco::BasicCluster &cluster, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv );
00104 
00105                 float eBottom( const reco::BasicCluster &cluster );
00106                 float eBottom( const reco::BasicCluster &cluster, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv );
00107                 // the energy of the most energetic crystal in the cluster
00108                 float eMax( const reco::BasicCluster &cluster );
00109                 float eMax( const reco::BasicCluster &cluster, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv );
00110                 // the energy of the second most energetic crystal in the cluster
00111                 float e2nd( const reco::BasicCluster &cluster );
00112                 float e2nd( const reco::BasicCluster &cluster, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv );
00113 
00114                 // get the DetId and the energy of the maximum energy crystal of the input cluster
00115                 std::pair<DetId, float> getMaximum( const reco::BasicCluster &cluster );
00116                 std::pair<DetId, float> getMaximum( const reco::BasicCluster &cluster, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv );
00117 
00118                 std::vector<float> energyBasketFractionEta( const reco::BasicCluster &cluster );
00119                 std::vector<float> energyBasketFractionEta( const reco::BasicCluster &cluster,std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv  );
00120 
00121                 std::vector<float> energyBasketFractionPhi( const reco::BasicCluster &cluster );
00122                 std::vector<float> energyBasketFractionPhi( const reco::BasicCluster &cluster, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv );
00123 
00124                 // return a vector v with v[0] = etaLat, v[1] = phiLat, v[2] = lat
00125                 std::vector<float> lat( const reco::BasicCluster &cluster, bool logW = true, float w0 = 4.7 );
00126                 std::vector<float> lat( const reco::BasicCluster &cluster, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv, bool logW = true, float w0 = 4.7 );
00127 
00128                 // return a vector v with v[0] = covEtaEta, v[1] = covEtaPhi, v[2] = covPhiPhi
00129                 std::vector<float> covariances(const reco::BasicCluster &cluster, float w0 = 4.7 );
00130                 std::vector<float> covariances(const reco::BasicCluster &cluster, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv, float w0 = 4.7 );
00131 
00132                 // return a vector v with v[0] = covIEtaIEta, v[1] = covIEtaIPhi, v[2] = covIPhiIPhi
00133                 //this function calculates differences in eta/phi in units of crystals not global eta/phi
00134                 //this is gives better performance in the crack regions of the calorimeter but gives otherwise identical results to covariances function
00135                 //this is only defined for the barrel, it returns covariances when the cluster is in the endcap
00136                 //Warning: covIEtaIEta has been studied by egamma, but so far covIPhiIPhi hasnt been studied extensively so there could be a bug in 
00137                 //         the covIPhiIEta or covIPhiIPhi calculations. I dont think there is but as it hasnt been heavily used, there might be one
00138                 std::vector<float> localCovariances(const reco::BasicCluster &cluster, float w0 = 4.7);
00139                 std::vector<float> localCovariances(const reco::BasicCluster &cluster, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv, float w0 = 4.7);
00140 
00141                 std::vector<float> scLocalCovariances(const reco::SuperCluster &cluster, float w0 = 4.7);
00142                 std::vector<float> scLocalCovariances(const reco::SuperCluster &cluster, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv, float w0 = 4.7);
00143 
00144                 double zernike20( const reco::BasicCluster &cluster, double R0 = 6.6, bool logW = true, float w0 = 4.7 );
00145                 double zernike42( const reco::BasicCluster &cluster, double R0 = 6.6, bool logW = true, float w0 = 4.7 );
00146 
00147                 // get the detId's of a matrix centered in the maximum energy crystal = (0,0)
00148                 // the size is specified by ixMin, ixMax, iyMin, iyMax in unit of crystals
00149                 std::vector<DetId> matrixDetId( DetId id, int ixMin, int ixMax, int iyMin, int iyMax );
00150                 // get the energy deposited in a matrix centered in the maximum energy crystal = (0,0)
00151                 // the size is specified by ixMin, ixMax, iyMin, iyMax in unit of crystals
00152                 float matrixEnergy( const reco::BasicCluster &cluster, DetId id, int ixMin, int ixMax, int iyMin, int iyMax );
00153                 float matrixEnergy( const reco::BasicCluster &cluster, DetId id, int ixMin, int ixMax, int iyMin, int iyMax, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv );
00154 
00155                 // get time of basic cluster seed crystal 
00156                 float BasicClusterSeedTime(const reco::BasicCluster &cluster);
00157                 // error-weighted average of time from constituents of basic cluster 
00158                 float BasicClusterTime(const reco::BasicCluster &cluster, const edm::Event &ev);
00159                 // get BasicClusterSeedTime of the seed basic cluser of the supercluster
00160                 float SuperClusterSeedTime(const reco::SuperCluster &cluster);
00161                 // get BasicClusterTime of the seed basic cluser of the supercluster
00162                 float SuperClusterTime(const reco::SuperCluster &cluster, const edm::Event &ev);
00163 
00164                 // mapping for preshower rechits
00165                 std::map<DetId, EcalRecHit> rechits_map_;
00166                 // get Preshower hit array
00167                 std::vector<float> getESHits(double X, double Y, double Z, std::map<DetId, EcalRecHit> rechits_map, const CaloGeometry* geometry, CaloSubdetectorTopology *topology_p, int row=0, int plane=1);
00168                 // get Preshower hit shape
00169                 float getESShape(std::vector<float> ESHits0);
00170                 // get Preshower effective sigmaRR
00171                 float eseffsirir( const reco::SuperCluster &cluster );
00172                 float eseffsixix( const reco::SuperCluster &cluster );
00173                 float eseffsiyiy( const reco::SuperCluster &cluster );
00174 
00175 //  std::vector<int> flagsexcl_;
00176   //std::vector<int> severitiesexcl_;
00177  // const EcalSeverityLevelAlgo *sevLv;
00178 
00179         private:
00180                 void getGeometry( const edm::EventSetup &es );
00181                 void getTopology( const edm::EventSetup &es );
00182                 void getEBRecHits( const edm::Event &ev, edm::InputTag redEBRecHits );
00183                 void getEERecHits( const edm::Event &ev, edm::InputTag redEERecHits );
00184                 void getESRecHits( const edm::Event &ev, edm::InputTag redESRecHits );
00185                 const EcalRecHitCollection * getEcalRecHitCollection( const reco::BasicCluster &cluster );
00186 
00187                 const CaloGeometry *geometry_;
00188                 const CaloTopology *topology_;
00189                 const EcalRecHitCollection *ebRecHits_;
00190                 const EcalRecHitCollection *eeRecHits_;
00191                 const EcalRecHitCollection *esRecHits_;
00192                 
00193                 //const EcalIntercalibConstantMap& icalMap;
00194                 edm::ESHandle<EcalIntercalibConstants> ical;
00195                 EcalIntercalibConstantMap        icalMap;
00196                 edm::ESHandle<EcalADCToGeVConstant>    agc;
00197                 edm::ESHandle<EcalLaserDbService>      laser;
00198                 void getIntercalibConstants( const edm::EventSetup &es );
00199                 void getADCToGeV           ( const edm::EventSetup &es );
00200                 void getLaserDbService     ( const edm::EventSetup &es );
00201 
00202 //  std::vector<int> flagsexcl_;
00203 //  std::vector<int> severitiesexcl_;
00204                 
00205 };
00206 
00207 #endif