Go to the documentation of this file.00001 #include "RecoHI/HiEgammaAlgos/interface/RxCalculator.h"
00002
00003 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00004
00005 #include "Geometry/Records/interface/IdealGeometryRecord.h"
00006 #include "Geometry/Records/interface/CaloGeometryRecord.h"
00007 #include "Geometry/CaloGeometry/interface/CaloSubdetectorGeometry.h"
00008
00009 #include "DataFormats/Common/interface/Handle.h"
00010 #include "DataFormats/EgammaReco/interface/BasicCluster.h"
00011 #include "DataFormats/EgammaReco/interface/SuperCluster.h"
00012
00013
00014 using namespace edm;
00015 using namespace reco;
00016
00017 #define PI 3.141592653589793238462643383279502884197169399375105820974945
00018
00019
00020 RxCalculator::RxCalculator (const edm::Event &iEvent, const edm::EventSetup &iSetup, edm::InputTag hbheLabel, edm::InputTag hfLabel, edm::InputTag hoLabel)
00021 {
00022 Handle<HFRecHitCollection> hfhandle;
00023 iEvent.getByLabel(hfLabel, hfhandle);
00024 fHFRecHits_ = hfhandle.product();
00025
00026 Handle<HORecHitCollection> hohandle;
00027 iEvent.getByLabel(hoLabel, hohandle);
00028 fHORecHits_ = hohandle.product();
00029
00030 Handle<HBHERecHitCollection> hehbhandle;
00031 iEvent.getByLabel(hbheLabel, hehbhandle);
00032 fHBHERecHits_ = hehbhandle.product();
00033
00034 ESHandle<CaloGeometry> geometryHandle;
00035 iSetup.get<CaloGeometryRecord>().get(geometryHandle);
00036 geometry_ = geometryHandle.product();
00037
00038 }
00039
00040
00041 double RxCalculator::getRx(const reco::SuperClusterRef cluster, double x, double threshold )
00042 {
00043 using namespace edm;
00044 using namespace reco;
00045
00046
00047 if(!fHBHERecHits_) {
00048 LogError("RxCalculator") << "Error! Can't get HBHERecHits for event.";
00049 return -100;
00050 }
00051
00052 double SClusterEta = cluster->eta();
00053 double SClusterPhi = cluster->phi();
00054 double TotalEt = 0;
00055
00056 for(size_t index = 0; index < fHBHERecHits_->size(); index++) {
00057 const HBHERecHit &rechit = (*fHBHERecHits_)[index];
00058 const DetId &detid = rechit.id();
00059 const GlobalPoint& hitpoint = geometry_->getPosition(detid);
00060 double eta = hitpoint.eta();
00061 double phi = hitpoint.phi();
00062 double dEta = fabs(eta-SClusterEta);
00063 double dPhi = fabs(phi-SClusterPhi);
00064 while (dPhi>2*PI) dPhi-=2*PI;
00065 if (dPhi>PI) dPhi=2*PI-dPhi;
00066
00067 if (dPhi>PI) dPhi=2*PI-dPhi;
00068
00069 double dR = sqrt(dEta * dEta + dPhi * dPhi);
00070
00071 if (dR<x*0.1) {
00072 double et = rechit.energy()/cosh(eta);
00073 if (et<threshold) et=0;
00074 TotalEt += et;
00075 }
00076 }
00077
00078 return TotalEt;
00079 }
00080
00081 double RxCalculator::getROx(const reco::SuperClusterRef cluster, double x,double threshold)
00082 {
00083 using namespace edm;
00084 using namespace reco;
00085
00086 if(!fHORecHits_) {
00087 LogError("RxCalculator") << "Error! Can't get HORecHits for event.";
00088 return -100;
00089 }
00090
00091 double SClusterEta = cluster->eta();
00092 double SClusterPhi = cluster->phi();
00093 double TotalEt = 0;
00094
00095 for(size_t index = 0; index < fHORecHits_->size(); index++) {
00096 const HORecHit &rechit = (*fHORecHits_)[index];
00097 const DetId &detid = rechit.id();
00098 const GlobalPoint& hitpoint = geometry_->getPosition(detid);
00099 double eta = hitpoint.eta();
00100 double phi = hitpoint.phi();
00101 double dEta = fabs(eta-SClusterEta);
00102 double dPhi = fabs(phi-SClusterPhi);
00103 while (dPhi>2*PI) dPhi-=2*PI;
00104 if (dPhi>PI) dPhi=2*PI-dPhi;
00105
00106 double dR = sqrt(dEta * dEta + dPhi * dPhi);
00107 if (dR<x*0.1) {
00108 double et = rechit.energy()/cosh(eta);
00109 if (et<threshold) et=0;
00110 TotalEt += et;
00111 }
00112 }
00113 return TotalEt;
00114 }
00115
00116 double RxCalculator::getRFx(const reco::SuperClusterRef cluster, double x, double threshold)
00117 {
00118 using namespace edm;
00119 using namespace reco;
00120
00121 if(!fHFRecHits_) {
00122 LogError("RxCalculator") << "Error! Can't get HFRecHits for event.";
00123 return -100;
00124 }
00125
00126 double SClusterEta = cluster->eta();
00127 double SClusterPhi = cluster->phi();
00128 double TotalEt = 0;
00129
00130 for(size_t index = 0; index < fHFRecHits_->size(); index++) {
00131 const HFRecHit &rechit = (*fHFRecHits_)[index];
00132 const DetId &detid = rechit.id();
00133 const GlobalPoint& hitpoint = geometry_->getPosition(detid);
00134 double eta = hitpoint.eta();
00135 double phi = hitpoint.phi();
00136 double dEta = fabs(eta-SClusterEta);
00137 double dPhi = fabs(phi-SClusterPhi);
00138 while (dPhi>2*PI) dPhi-=2*PI;
00139 if (dPhi>PI) dPhi=2*PI-dPhi;
00140
00141
00142 double dR = sqrt(dEta * dEta + dPhi * dPhi);
00143 if (dR<x*0.1) {
00144 double et = rechit.energy()/cosh(eta);
00145 if (et<threshold) et=0;
00146 TotalEt += et;
00147 }
00148 }
00149
00150
00151
00152 return TotalEt;
00153 }
00154
00155
00156 double RxCalculator::getCRx(const reco::SuperClusterRef cluster, double x, double threshold )
00157 {
00158 using namespace edm;
00159 using namespace reco;
00160
00161
00162 if(!fHBHERecHits_) {
00163 LogError("RxCalculator") << "Error! Can't get HBHERecHits for event.";
00164 return -100;
00165 }
00166
00167 double SClusterEta = cluster->eta();
00168 double SClusterPhi = cluster->phi();
00169 double TotalEt = 0;
00170
00171 for(size_t index = 0; index < fHBHERecHits_->size(); index++) {
00172 const HBHERecHit &rechit = (*fHBHERecHits_)[index];
00173 const DetId &detid = rechit.id();
00174 const GlobalPoint& hitpoint = geometry_->getPosition(detid);
00175 double eta = hitpoint.eta();
00176 double phi = hitpoint.phi();
00177 double dEta = fabs(eta-SClusterEta);
00178 double dPhi = fabs(phi-SClusterPhi);
00179 while (dPhi>2*PI) dPhi-=2*PI;
00180 if (dPhi>PI) dPhi=2*PI-dPhi;
00181
00182 if (dEta<x*0.1) {
00183 double et = rechit.energy()/cosh(eta);
00184 if (et<threshold) et=0;
00185 TotalEt += et;
00186 }
00187 }
00188
00189 double Rx = getRx(cluster,x,threshold);
00190 double CRx = Rx - TotalEt / 40.0 * x;
00191
00192 return CRx;
00193 }