#include <RxCalculator.h>
Public Member Functions | |
double | getCRx (const reco::SuperClusterRef clus, double i, double threshold, double innerR=0.0) |
double | getRFx (const reco::SuperClusterRef clus, double i, double threshold) |
double | getROx (const reco::SuperClusterRef clus, double i, double threshold) |
double | getRx (const reco::SuperClusterRef clus, double i, double threshold, double innerR=0.0) |
RxCalculator (const edm::Event &iEvent, const edm::EventSetup &iSetup, edm::InputTag hbheLabel, edm::InputTag hfLabel, edm::InputTag hoLabel) | |
Private Attributes | |
const HBHERecHitCollection * | fHBHERecHits_ |
const HFRecHitCollection * | fHFRecHits_ |
const HORecHitCollection * | fHORecHits_ |
const CaloGeometry * | geometry_ |
Definition at line 19 of file RxCalculator.h.
RxCalculator::RxCalculator | ( | const edm::Event & | iEvent, |
const edm::EventSetup & | iSetup, | ||
edm::InputTag | hbheLabel, | ||
edm::InputTag | hfLabel, | ||
edm::InputTag | hoLabel | ||
) |
Definition at line 20 of file RxCalculator.cc.
References edm::EventSetup::get(), edm::Event::getByLabel(), edm::HandleBase::isValid(), edm::ESHandleBase::isValid(), NULL, edm::ESHandle< T >::product(), and edm::Handle< T >::product().
{ Handle<HFRecHitCollection> hfhandle; iEvent.getByLabel(hfLabel, hfhandle); if(hfhandle.isValid()) fHFRecHits_ = hfhandle.product(); else fHFRecHits_ = NULL; Handle<HORecHitCollection> hohandle; iEvent.getByLabel(hoLabel, hohandle); if(hohandle.isValid()) fHORecHits_ = hohandle.product(); else fHORecHits_ = NULL; Handle<HBHERecHitCollection> hehbhandle; iEvent.getByLabel(hbheLabel, hehbhandle); if(hehbhandle.isValid()) fHBHERecHits_ = hehbhandle.product(); else fHBHERecHits_ = NULL; ESHandle<CaloGeometry> geometryHandle; iSetup.get<CaloGeometryRecord>().get(geometryHandle); if(geometryHandle.isValid()) geometry_ = geometryHandle.product(); else geometry_ = NULL; }
double RxCalculator::getCRx | ( | const reco::SuperClusterRef | clus, |
double | i, | ||
double | threshold, | ||
double | innerR = 0.0 |
||
) |
Definition at line 169 of file RxCalculator.cc.
References cond::rpcobgas::detid, dPhi(), CaloRecHit::energy(), eta(), PV3DBase< T, PVType, FrameType >::eta(), HBHERecHit::id(), getHLTprescales::index, PV3DBase< T, PVType, FrameType >::phi(), phi, PI, and dt_dqm_sourceclient_common_cff::reco.
Referenced by HiEgammaIsolationProducer::produce().
{ using namespace edm; using namespace reco; if(!fHBHERecHits_) { // LogError("RxCalculator") << "Error! Can't get HBHERecHits for event."; return -100; } double SClusterEta = cluster->eta(); double SClusterPhi = cluster->phi(); double TotalEt = 0; for(size_t index = 0; index < fHBHERecHits_->size(); index++) { const HBHERecHit &rechit = (*fHBHERecHits_)[index]; const DetId &detid = rechit.id(); const GlobalPoint& hitpoint = geometry_->getPosition(detid); double eta = hitpoint.eta(); double phi = hitpoint.phi(); double dEta = fabs(eta-SClusterEta); double dPhi = fabs(phi-SClusterPhi); while (dPhi>2*PI) dPhi-=2*PI; if (dPhi>PI) dPhi=2*PI-dPhi; if (dEta<x*0.1) { double et = rechit.energy()/cosh(eta); if (et<threshold) et=0; TotalEt += et; } } double Rx = getRx(cluster,x,threshold,innerR); double CRx = Rx - TotalEt * (0.01*x*x - innerR*innerR) / (2 * 2 * 0.1 * x) ; return CRx; }
double RxCalculator::getRFx | ( | const reco::SuperClusterRef | clus, |
double | i, | ||
double | threshold | ||
) |
Definition at line 129 of file RxCalculator.cc.
References cond::rpcobgas::detid, dPhi(), CaloRecHit::energy(), eta(), PV3DBase< T, PVType, FrameType >::eta(), HFRecHit::id(), getHLTprescales::index, PV3DBase< T, PVType, FrameType >::phi(), phi, PI, dt_dqm_sourceclient_common_cff::reco, and mathSSE::sqrt().
{ using namespace edm; using namespace reco; if(!fHFRecHits_) { LogError("RxCalculator") << "Error! Can't get HFRecHits for event."; return -100; } double SClusterEta = cluster->eta(); double SClusterPhi = cluster->phi(); double TotalEt = 0; for(size_t index = 0; index < fHFRecHits_->size(); index++) { const HFRecHit &rechit = (*fHFRecHits_)[index]; const DetId &detid = rechit.id(); const GlobalPoint& hitpoint = geometry_->getPosition(detid); double eta = hitpoint.eta(); double phi = hitpoint.phi(); double dEta = fabs(eta-SClusterEta); double dPhi = fabs(phi-SClusterPhi); while (dPhi>2*PI) dPhi-=2*PI; if (dPhi>PI) dPhi=2*PI-dPhi; double dR = sqrt(dEta * dEta + dPhi * dPhi); if (dR<x*0.1) { double et = rechit.energy()/cosh(eta); if (et<threshold) et=0; TotalEt += et; } } return TotalEt; }
double RxCalculator::getROx | ( | const reco::SuperClusterRef | clus, |
double | i, | ||
double | threshold | ||
) |
Definition at line 94 of file RxCalculator.cc.
References cond::rpcobgas::detid, dPhi(), CaloRecHit::energy(), eta(), PV3DBase< T, PVType, FrameType >::eta(), HORecHit::id(), getHLTprescales::index, PV3DBase< T, PVType, FrameType >::phi(), phi, PI, dt_dqm_sourceclient_common_cff::reco, and mathSSE::sqrt().
{ using namespace edm; using namespace reco; if(!fHORecHits_) { // LogError("RxCalculator") << "Error! Can't get HORecHits for event."; return -100; } double SClusterEta = cluster->eta(); double SClusterPhi = cluster->phi(); double TotalEt = 0; for(size_t index = 0; index < fHORecHits_->size(); index++) { const HORecHit &rechit = (*fHORecHits_)[index]; const DetId &detid = rechit.id(); const GlobalPoint& hitpoint = geometry_->getPosition(detid); double eta = hitpoint.eta(); double phi = hitpoint.phi(); double dEta = fabs(eta-SClusterEta); double dPhi = fabs(phi-SClusterPhi); while (dPhi>2*PI) dPhi-=2*PI; if (dPhi>PI) dPhi=2*PI-dPhi; double dR = sqrt(dEta * dEta + dPhi * dPhi); if (dR<x*0.1) { double et = rechit.energy()/cosh(eta); if (et<threshold) et=0; TotalEt += et; } } return TotalEt; }
double RxCalculator::getRx | ( | const reco::SuperClusterRef | clus, |
double | i, | ||
double | threshold, | ||
double | innerR = 0.0 |
||
) |
Definition at line 53 of file RxCalculator.cc.
References cond::rpcobgas::detid, dPhi(), CaloRecHit::energy(), eta(), PV3DBase< T, PVType, FrameType >::eta(), HBHERecHit::id(), getHLTprescales::index, PV3DBase< T, PVType, FrameType >::phi(), phi, PI, dt_dqm_sourceclient_common_cff::reco, and mathSSE::sqrt().
Referenced by HiEgammaIsolationProducer::produce().
{ using namespace edm; using namespace reco; if(!fHBHERecHits_) { // LogError("RxCalculator") << "Error! Can't get HBHERecHits for event."; return -100; } double SClusterEta = cluster->eta(); double SClusterPhi = cluster->phi(); double TotalEt = 0; for(size_t index = 0; index < fHBHERecHits_->size(); index++) { const HBHERecHit &rechit = (*fHBHERecHits_)[index]; const DetId &detid = rechit.id(); const GlobalPoint& hitpoint = geometry_->getPosition(detid); double eta = hitpoint.eta(); double phi = hitpoint.phi(); double dEta = fabs(eta-SClusterEta); double dPhi = fabs(phi-SClusterPhi); while (dPhi>2*PI) dPhi-=2*PI; if (dPhi>PI) dPhi=2*PI-dPhi; if (dPhi>PI) dPhi=2*PI-dPhi; double dR = sqrt(dEta * dEta + dPhi * dPhi); // veto inner cone/////////////// if ( dR < innerR ) continue; if (dR<x*0.1) { double et = rechit.energy()/cosh(eta); if (et<threshold) et=0; TotalEt += et; } } return TotalEt; }
const HBHERecHitCollection* RxCalculator::fHBHERecHits_ [private] |
Definition at line 32 of file RxCalculator.h.
const HFRecHitCollection* RxCalculator::fHFRecHits_ [private] |
Definition at line 34 of file RxCalculator.h.
const HORecHitCollection* RxCalculator::fHORecHits_ [private] |
Definition at line 33 of file RxCalculator.h.
const CaloGeometry* RxCalculator::geometry_ [private] |
Definition at line 35 of file RxCalculator.h.