File indexing completed on 2024-04-06 12:02:58
0001 #include "FWCore/ServiceRegistry/interface/Service.h"
0002 #include "CondCore/DBOutputService/interface/PoolDBOutputService.h"
0003 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0004 #include "FWCore/Framework/interface/ESHandle.h"
0005 #include "FWCore/Framework/interface/EventSetup.h"
0006 #include "FWCore/Framework/interface/Event.h"
0007 #include "CalibCalorimetry/EcalLaserCorrection/interface/EcalLaserDbService.h"
0008 #include "CalibCalorimetry/EcalLaserCorrection/interface/EcalLaserDbRecord.h"
0009 #include "CondTools/Ecal/interface/EcalPFRecHitThresholdsMaker.h"
0010 #include "CondFormats/EcalObjects/interface/EcalPFRecHitThresholds.h"
0011 #include "CondFormats/DataRecord/interface/EcalPFRecHitThresholdsRcd.h"
0012 #include "CondFormats/DataRecord/interface/EcalPedestalsRcd.h"
0013 #include "CondFormats/DataRecord/interface/EcalIntercalibConstantsRcd.h"
0014 #include "CondFormats/EcalObjects/interface/EcalADCToGeVConstant.h"
0015 #include "CondFormats/DataRecord/interface/EcalADCToGeVConstantRcd.h"
0016
0017 #include <vector>
0018
0019 EcalPFRecHitThresholdsMaker::EcalPFRecHitThresholdsMaker(const edm::ParameterSet& iConfig)
0020 : m_timetype(iConfig.getParameter<std::string>("timetype")),
0021 ecalPedestalsToken_(esConsumes()),
0022 ecalADCToGeVConstantToken_(esConsumes()),
0023 ecalIntercalibConstantsToken_(esConsumes()),
0024 ecalLaserDbServiceToken_(esConsumes()) {
0025 m_nsigma = iConfig.getParameter<double>("NSigma");
0026 }
0027
0028 EcalPFRecHitThresholdsMaker::~EcalPFRecHitThresholdsMaker() {}
0029
0030 void EcalPFRecHitThresholdsMaker::analyze(const edm::Event& evt, const edm::EventSetup& evtSetup) {
0031 edm::Service<cond::service::PoolDBOutputService> dbOutput;
0032 if (!dbOutput.isAvailable()) {
0033 throw cms::Exception("PoolDBOutputService is not available");
0034 }
0035
0036 const EcalPedestals* ped_db = &evtSetup.getData(ecalPedestalsToken_);
0037 edm::LogInfo("EcalPFRecHitThresholdsMaker") << "ped pointer is: " << ped_db << std::endl;
0038
0039 const EcalADCToGeVConstant* adc_db = &evtSetup.getData(ecalADCToGeVConstantToken_);
0040 edm::LogInfo("EcalPFRecHitThresholdsMaker") << "adc pointer is: " << adc_db << std::endl;
0041
0042 const EcalIntercalibConstants* ical_db = &evtSetup.getData(ecalIntercalibConstantsToken_);
0043 edm::LogInfo("EcalPFRecHitThresholdsMaker") << "inter pointer is: " << ical_db << std::endl;
0044
0045 const auto laser = evtSetup.getHandle(ecalLaserDbServiceToken_);
0046
0047 EcalPFRecHitThresholds pfthresh;
0048
0049
0050
0051 float adc_EB = float(adc_db->getEEValue());
0052 float adc_EE = float(adc_db->getEBValue());
0053
0054
0055
0056 for (int iEta = -EBDetId::MAX_IETA; iEta <= EBDetId::MAX_IETA; ++iEta) {
0057 if (iEta == 0)
0058 continue;
0059 for (int iPhi = EBDetId::MIN_IPHI; iPhi <= EBDetId::MAX_IPHI; ++iPhi) {
0060
0061 if (EBDetId::validDetId(iEta, iPhi)) {
0062 EBDetId ebdetid(iEta, iPhi, EBDetId::ETAPHIMODE);
0063 EcalPedestals::const_iterator it = ped_db->find(ebdetid.rawId());
0064 EcalPedestals::Item aped = (*it);
0065
0066 EcalIntercalibConstants::const_iterator itc = ical_db->find(ebdetid.rawId());
0067 float calib = (*itc);
0068
0069
0070 float lasercalib = 1.;
0071 lasercalib = laser->getLaserCorrection(ebdetid, evt.time());
0072
0073 EcalPFRecHitThreshold thresh = aped.rms_x12 * calib * adc_EB * lasercalib * m_nsigma;
0074
0075 if (iPhi == 100)
0076 edm::LogInfo("EcalPFRecHitThresholdsMaker") << "Thresh(GeV)=" << thresh << std::endl;
0077
0078 pfthresh.insert(std::make_pair(ebdetid.rawId(), thresh));
0079 }
0080 }
0081 }
0082
0083 for (int iX = EEDetId::IX_MIN; iX <= EEDetId::IX_MAX; ++iX) {
0084 for (int iY = EEDetId::IY_MIN; iY <= EEDetId::IY_MAX; ++iY) {
0085
0086 if (EEDetId::validDetId(iX, iY, 1)) {
0087 EEDetId eedetid(iX, iY, 1);
0088
0089 EcalPedestals::const_iterator it = ped_db->find(eedetid.rawId());
0090 EcalPedestals::Item aped = (*it);
0091
0092 EcalIntercalibConstants::const_iterator itc = ical_db->find(eedetid.rawId());
0093 float calib = (*itc);
0094
0095
0096 float lasercalib = 1.;
0097 lasercalib = laser->getLaserCorrection(eedetid, evt.time());
0098
0099 EcalPFRecHitThreshold thresh = aped.rms_x12 * calib * adc_EE * lasercalib * m_nsigma;
0100 pfthresh.insert(std::make_pair(eedetid.rawId(), thresh));
0101 }
0102 if (EEDetId::validDetId(iX, iY, -1)) {
0103 EEDetId eedetid(iX, iY, -1);
0104
0105 EcalPedestals::const_iterator it = ped_db->find(eedetid.rawId());
0106 EcalPedestals::Item aped = (*it);
0107
0108 EcalIntercalibConstants::const_iterator itc = ical_db->find(eedetid.rawId());
0109 float calib = (*itc);
0110
0111
0112 float lasercalib = 1.;
0113 lasercalib = laser->getLaserCorrection(eedetid, evt.time());
0114
0115 EcalPFRecHitThreshold thresh = aped.rms_x12 * calib * adc_EE * lasercalib * m_nsigma;
0116 pfthresh.insert(std::make_pair(eedetid.rawId(), thresh));
0117
0118 if (iX == 50)
0119 edm::LogInfo("EcalPFRecHitThresholdsMaker") << "Thresh(GeV)=" << thresh << std::endl;
0120 }
0121 }
0122 }
0123
0124 dbOutput->createOneIOV<const EcalPFRecHitThresholds>(pfthresh, dbOutput->beginOfTime(), "EcalPFRecHitThresholdsRcd");
0125
0126 edm::LogInfo("EcalPFRecHitThresholdsMaker") << "EcalPFRecHitThresholdsMaker wrote it " << std::endl;
0127 }