File indexing completed on 2024-04-06 12:02:56
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019 #include <memory>
0020
0021
0022 #include "FWCore/Framework/interface/Frameworkfwd.h"
0023 #include "FWCore/Framework/interface/Event.h"
0024 #include "FWCore/Framework/interface/MakerMacros.h"
0025 #include "FWCore/Framework/interface/EventSetup.h"
0026 #include "FWCore/Framework/interface/EventSetupRecordKey.h"
0027 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0028 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0029 #include "FWCore/ServiceRegistry/interface/Service.h"
0030 #include "FWCore/Framework/interface/IOVSyncValue.h"
0031
0032 #include "CondCore/DBOutputService/interface/PoolDBOutputService.h"
0033
0034 #include "CondFormats/EcalObjects/interface/EcalLaserAPDPNRatios.h"
0035 #include "CondFormats/DataRecord/interface/EcalLaserAPDPNRatiosRcd.h"
0036 #include "CondFormats/DataRecord/interface/EcalLaserAPDPNRatiosRefRcd.h"
0037 #include "CondFormats/DataRecord/interface/EcalLaserAlphasRcd.h"
0038
0039 #include "OnlineDB/EcalCondDB/interface/all_monitoring_types.h"
0040 #include "OnlineDB/Oracle/interface/Oracle.h"
0041 #include "OnlineDB/EcalCondDB/interface/EcalCondDBInterface.h"
0042
0043 #include "DataFormats/EcalDetId/interface/EBDetId.h"
0044 #include "DataFormats/EcalDetId/interface/EEDetId.h"
0045
0046 #include "CondTools/Ecal/interface/EcalGetLaserData.h"
0047
0048
0049
0050
0051
0052
0053
0054
0055
0056
0057
0058
0059 EcalGetLaserData::EcalGetLaserData(const edm::ParameterSet& iConfig)
0060 :
0061 m_cacheIDs(),
0062 m_records(),
0063 ecalLaserAPDPNRatiosToken_(esConsumes()),
0064 ecalLaserAPDPNRatiosRefToken_(esConsumes()),
0065 ecalLaserAlphasToken_(esConsumes()) {
0066
0067
0068 std::string container;
0069 std::string record;
0070 typedef std::vector<edm::ParameterSet> Parameters;
0071 Parameters toGet = iConfig.getParameter<Parameters>("toGet");
0072 for (const auto& iparam : toGet) {
0073 container = iparam.getParameter<std::string>("container");
0074 record = iparam.getParameter<std::string>("record");
0075 m_cacheIDs.emplace(container, 0);
0076 m_records.emplace(container, record);
0077
0078 }
0079 }
0080
0081 EcalGetLaserData::~EcalGetLaserData() {
0082
0083
0084 }
0085
0086
0087
0088
0089
0090
0091 void EcalGetLaserData::analyze(const edm::Event& evt, const edm::EventSetup& evtSetup) {
0092 using namespace edm;
0093
0094
0095 for (const auto& irec : m_records) {
0096 const std::string& container = irec.first;
0097
0098
0099 if (container == "EcalLaserAPDPNRatios") {
0100
0101 const EcalLaserAPDPNRatios* laserapdpnrRatios = &evtSetup.getData(ecalLaserAPDPNRatiosToken_);
0102
0103 EcalLaserAPDPNRatios::EcalLaserTimeStamp timestamp;
0104 EcalLaserAPDPNRatios::EcalLaserAPDPNpair apdpnpair;
0105
0106 const EcalLaserAPDPNRatios::EcalLaserAPDPNRatiosMap& laserRatiosMap = laserapdpnrRatios->getLaserMap();
0107 const EcalLaserAPDPNRatios::EcalLaserTimeStampMap& laserTimeMap = laserapdpnrRatios->getTimeMap();
0108
0109
0110 for (int iEta = -EBDetId::MAX_IETA; iEta <= EBDetId::MAX_IETA; ++iEta) {
0111 if (iEta == 0)
0112 continue;
0113 for (int iPhi = EBDetId::MIN_IPHI; iPhi <= EBDetId::MAX_IPHI; ++iPhi) {
0114 EBDetId ebdetid(iEta, iPhi);
0115 int hi = ebdetid.hashedIndex();
0116
0117 if (hi < static_cast<int>(laserRatiosMap.size())) {
0118 apdpnpair = laserRatiosMap[hi];
0119 edm::LogInfo("EcalGetLaserData") << "A sample value of APDPN pair EB : " << hi << " : " << apdpnpair.p1
0120 << " , " << apdpnpair.p2 << std::endl;
0121 } else {
0122 edm::LogError("EcalGetLaserData") << "error with laserRatiosMap!" << std::endl;
0123 }
0124 }
0125 }
0126
0127
0128 for (int iX = EEDetId::IX_MIN; iX <= EEDetId::IX_MAX; ++iX) {
0129 for (int iY = EEDetId::IY_MIN; iY <= EEDetId::IY_MAX; ++iY) {
0130 if (!EEDetId::validDetId(iX, iY, 1))
0131 continue;
0132
0133 EEDetId eedetidpos(iX, iY, 1);
0134 int hi = eedetidpos.hashedIndex();
0135
0136 if (hi < static_cast<int>(laserRatiosMap.size())) {
0137 apdpnpair = laserRatiosMap[hi];
0138 edm::LogInfo("EcalGetLaserData") << "A sample value of APDPN pair EE+ : " << hi << " : " << apdpnpair.p1
0139 << " , " << apdpnpair.p2 << std::endl;
0140 } else {
0141 edm::LogError("EcalGetLaserData") << "error with laserRatiosMap!" << std::endl;
0142 }
0143
0144 if (!EEDetId::validDetId(iX, iY, -1))
0145 continue;
0146 EEDetId eedetidneg(iX, iY, 1);
0147 hi = eedetidneg.hashedIndex();
0148
0149 if (hi < static_cast<int>(laserRatiosMap.size())) {
0150 apdpnpair = laserRatiosMap[hi];
0151 edm::LogInfo("EcalGetLaserData") << "A sample value of APDPN pair EE- : " << hi << " : " << apdpnpair.p1
0152 << " , " << apdpnpair.p2 << std::endl;
0153 } else {
0154 edm::LogError("EcalGetLaserData") << "error with laserRatiosMap!" << std::endl;
0155 }
0156 }
0157 }
0158
0159 for (int i = 0; i < 92; i++) {
0160 timestamp = laserTimeMap[i];
0161 edm::LogInfo("EcalGetLaserData") << "A value of timestamp pair : " << i << " " << timestamp.t1.value() << " , "
0162 << timestamp.t2.value() << std::endl;
0163 }
0164
0165 edm::LogInfo("EcalGetLaserData") << ".. just retrieved the last valid record from DB " << std::endl;
0166
0167 } else if (container == "EcalLaserAPDPNRatiosRef") {
0168
0169 EcalLaserAPDPNref apdpnref;
0170 const EcalLaserAPDPNRatiosRefMap& laserRefMap = (&evtSetup.getData(ecalLaserAPDPNRatiosRefToken_))->getMap();
0171
0172
0173 for (int iEta = -EBDetId::MAX_IETA; iEta <= EBDetId::MAX_IETA; ++iEta) {
0174 if (iEta == 0)
0175 continue;
0176 for (int iPhi = EBDetId::MIN_IPHI; iPhi <= EBDetId::MAX_IPHI; ++iPhi) {
0177 EBDetId ebdetid(iEta, iPhi);
0178 int hi = ebdetid.hashedIndex();
0179
0180 if (hi < static_cast<int>(laserRefMap.size())) {
0181 apdpnref = laserRefMap[hi];
0182 edm::LogInfo("EcalGetLaserData")
0183 << "A sample value of APDPN Reference value EB : " << hi << " : " << apdpnref << std::endl;
0184 } else {
0185 edm::LogError("EcalGetLaserData") << "error with laserRefMap!" << std::endl;
0186 }
0187 }
0188 }
0189
0190
0191 for (int iX = EEDetId::IX_MIN; iX <= EEDetId::IX_MAX; ++iX) {
0192 for (int iY = EEDetId::IY_MIN; iY <= EEDetId::IY_MAX; ++iY) {
0193 if (!EEDetId::validDetId(iX, iY, 1))
0194 continue;
0195
0196 EEDetId eedetidpos(iX, iY, 1);
0197 int hi = eedetidpos.hashedIndex();
0198
0199 if (hi < static_cast<int>(laserRefMap.size())) {
0200 apdpnref = laserRefMap[hi];
0201 edm::LogInfo("EcalGetLaserData")
0202 << "A sample value of APDPN Reference value EE+ : " << hi << " : " << apdpnref << std::endl;
0203
0204 } else {
0205 edm::LogError("EcalGetLaserData") << "error with laserRefMap!" << std::endl;
0206 }
0207
0208 if (!EEDetId::validDetId(iX, iY, -1))
0209 continue;
0210 EEDetId eedetidneg(iX, iY, -1);
0211 hi = eedetidneg.hashedIndex();
0212
0213 if (hi < static_cast<int>(laserRefMap.size())) {
0214 apdpnref = laserRefMap[hi];
0215 edm::LogInfo("EcalGetLaserData")
0216 << "A sample value of APDPN Reference value EE- : " << hi << " : " << apdpnref << std::endl;
0217 } else {
0218 edm::LogError("EcalGetLaserData") << "error with laserRefMap!" << std::endl;
0219 }
0220 }
0221 }
0222
0223 edm::LogInfo("EcalGetLaserData") << "... just retrieved the last valid record from DB " << std::endl;
0224
0225 } else if (container == "EcalLaserAlphas") {
0226
0227
0228 EcalLaserAlpha alpha;
0229 const EcalLaserAlphaMap& laserAlphaMap = (&evtSetup.getData(ecalLaserAlphasToken_))->getMap();
0230
0231
0232 for (int iEta = -EBDetId::MAX_IETA; iEta <= EBDetId::MAX_IETA; ++iEta) {
0233 if (iEta == 0)
0234 continue;
0235 for (int iPhi = EBDetId::MIN_IPHI; iPhi <= EBDetId::MAX_IPHI; ++iPhi) {
0236 EBDetId ebdetid(iEta, iPhi);
0237 int hi = ebdetid.hashedIndex();
0238
0239 if (hi < static_cast<int>(laserAlphaMap.size())) {
0240 alpha = laserAlphaMap[hi];
0241 edm::LogInfo("EcalGetLaserData")
0242 << " A sample value of Alpha value EB : " << hi << " : " << alpha << std::endl;
0243 } else {
0244 edm::LogError("EcalGetLaserData") << "error with laserAlphaMap!" << std::endl;
0245 }
0246 }
0247 }
0248
0249
0250 for (int iX = EEDetId::IX_MIN; iX <= EEDetId::IX_MAX; ++iX) {
0251 for (int iY = EEDetId::IY_MIN; iY <= EEDetId::IY_MAX; ++iY) {
0252 if (!EEDetId::validDetId(iX, iY, 1))
0253 continue;
0254
0255 EEDetId eedetidpos(iX, iY, 1);
0256 int hi = eedetidpos.hashedIndex();
0257
0258 if (hi < static_cast<int>(laserAlphaMap.size())) {
0259 alpha = laserAlphaMap[hi];
0260 edm::LogInfo("EcalGetLaserData")
0261 << " A sample value of Alpha value EE+ : " << hi << " : " << alpha << std::endl;
0262 } else {
0263 edm::LogError("EcalGetLaserData") << "error with laserAlphaMap!" << std::endl;
0264 }
0265
0266 if (!EEDetId::validDetId(iX, iY, -1))
0267 continue;
0268 EEDetId eedetidneg(iX, iY, -1);
0269 hi = eedetidneg.hashedIndex();
0270
0271 if (hi < static_cast<int>(laserAlphaMap.size())) {
0272 alpha = laserAlphaMap[hi];
0273 edm::LogInfo("EcalGetLaserData")
0274 << " A sample value of Alpha value EE- : " << hi << " : " << alpha << std::endl;
0275 } else {
0276 edm::LogError("EcalGetLaserData") << "error with laserAlphaMap!" << std::endl;
0277 }
0278 }
0279 }
0280
0281 edm::LogInfo("EcalGetLaserData") << "... just retrieved the last valid record from DB " << std::endl;
0282
0283 } else {
0284 edm::LogError("EcalGetLaserData") << "Cannot retrieve for container: " << container << std::endl;
0285 }
0286 }
0287 }
0288
0289
0290 void EcalGetLaserData::beginJob() {}
0291
0292
0293 void EcalGetLaserData::endJob() {}