Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 11:58:03

0001 //emacs settings:-*- mode: c++; c-basic-offset: 2; indent-tabs-mode: nil -*-
0002 /*
0003  * $Id: EcalLaserCondTools.cc,v 1.2 2010/06/14 10:45:17 pgras Exp $
0004  *
0005  * author: Ph Gras. June, 2010
0006  */
0007 
0008 #include "CalibCalorimetry/EcalTrivialCondModules/interface/EcalLaserCondTools.h"
0009 
0010 #include "FWCore/Framework/interface/Event.h"
0011 #include "FWCore/Framework/interface/EventSetup.h"
0012 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0013 
0014 #include "DataFormats/EcalDetId/interface/EBDetId.h"
0015 #include "DataFormats/EcalDetId/interface/EEDetId.h"
0016 
0017 #include "hdf5.h"
0018 
0019 #include <string>
0020 #include <fstream>
0021 #include <algorithm>
0022 #include <memory>
0023 #include <cmath>
0024 
0025 EcalLaserCondTools::EcalLaserCondTools(const edm::ParameterSet& ps)
0026     : fout_(nullptr),
0027       eventList_(nullptr),
0028       eventListFileName_(ps.getParameter<std::string>("eventListFile")),
0029       verb_(ps.getParameter<int>("verbosity")),
0030       mode_(ps.getParameter<std::string>("mode")),
0031       fnames_(ps.getParameter<std::vector<std::string> >("inputFiles")),
0032       skipIov_(ps.getParameter<int>("skipIov")),
0033       nIovs_(ps.getParameter<int>("nIovs")),
0034       fromTime_(ps.getParameter<int>("fromTime")),
0035       toTime_(ps.getParameter<int>("toTime")),
0036       minP_(ps.getParameter<double>("transparencyMin")),
0037       maxP_(ps.getParameter<double>("transparencyMax")) {
0038   if (mode_ == "db_to_ascii_file") {
0039     laserAPDPNRatiosToken_ = esConsumes();
0040   }
0041 
0042   ferr_ = fopen("corr_errors.txt", "w");
0043   fprintf(ferr_, "#t1\tdetid\tp1\tp2\tp3");
0044 
0045   if (!eventListFileName_.empty()) {
0046     eventList_ = fopen(eventListFileName_.c_str(), "r");
0047     if (eventList_ == nullptr)
0048       throw cms::Exception("User") << "Failed to open file " << eventListFileName_ << "\n";
0049   }
0050 }
0051 
0052 EcalLaserCondTools::~EcalLaserCondTools() {
0053   if (ferr_)
0054     fclose(ferr_);
0055   if (fout_)
0056     fclose(fout_);
0057 }
0058 
0059 void EcalLaserCondTools::analyze(const edm::Event& event, const edm::EventSetup& es) {
0060   if (mode_ == "ascii_file_to_db") {
0061     if (verb_ > 2)
0062       edm::LogPrint("EcalLaserCondTools") << "ascii_file_to_db mode\n";
0063 
0064     if (!db_.isAvailable()) {
0065       throw cms::Exception("CondDBAccess") << "Failed to connect to PoolDBOutputService\n";
0066     }
0067     FileReader corrReader(fnames_);
0068     corrReader.setVerbosity(verb_);
0069     fillDb(corrReader);
0070   } else if (mode_ == "hdf_file_to_db") {
0071     from_hdf_to_db();
0072   } else if (mode_ == "db_to_ascii_file") {
0073     dbToAscii(es);
0074   } else {
0075     cms::Exception("InvalidParam") << "Value of parameter mode is not valid. Expecting ascii_file_to_db or read";
0076   }
0077 }
0078 
0079 void EcalLaserCondTools::from_hdf_to_db() {
0080   cond::Time_t iovStart = 0;
0081 
0082   hid_t file, space, memspace;
0083   hid_t dset_rawid, dset_t2, dset;
0084 
0085   hsize_t dims[2] = {};
0086 
0087   for (unsigned int ifile = 0; ifile < fnames_.size(); ++ifile) {
0088     if (verb_) {
0089       edm::LogPrint("EcalLaserCondTools") << " - converting file: " << fnames_[ifile] << "\n";
0090     }
0091 
0092     file = H5Fopen(fnames_[ifile].c_str(), H5F_ACC_RDONLY, H5P_DEFAULT);
0093 
0094     dset_rawid = H5Dopen(file, "cmssw_id", H5P_DEFAULT);
0095     space = H5Dget_space(dset_rawid);
0096     H5Sget_simple_extent_dims(space, dims, nullptr);
0097 
0098     unsigned int nCrystals = dims[0];
0099     int rawid[nCrystals];
0100     herr_t status;
0101 
0102     status = H5Dread(dset_rawid, H5T_NATIVE_INT, H5S_ALL, H5S_ALL, H5P_DEFAULT, rawid);
0103     if (status < 0)
0104       throw cms::Exception("EcalLaserCondTool:HDF") << "Error while reading HD file.";
0105 
0106     H5Dclose(dset_rawid);
0107     H5Sclose(space);
0108 
0109     dset_t2 = H5Dopen(file, "t2", H5P_DEFAULT);
0110     space = H5Dget_space(dset_t2);
0111     H5Sget_simple_extent_dims(space, dims, nullptr);
0112 
0113     unsigned int nIovs = dims[0];
0114     unsigned int nLME = dims[1];
0115 
0116     if (verb_) {
0117       edm::LogPrint("EcalLaserCondTools") << "Number of crystals: " << nCrystals << "\n";
0118       edm::LogPrint("EcalLaserCondTools") << "Number of IOVs: " << nIovs << "\n";
0119       edm::LogPrint("EcalLaserCondTools") << "Number of Monitoring regions: " << nLME << "\n";
0120     }
0121 
0122     int t1[nIovs], t3[nIovs], t2[nIovs][nLME];
0123 
0124     // -- reading data (cmsswid, t2, t1, t3, p2, p1, p3
0125     if (verb_ > 1)
0126       edm::LogPrint("EcalLaserCondTools") << " * reading t2 table "
0127                                           << "\n";
0128     status = H5Dread(dset_t2, H5T_NATIVE_INT, H5S_ALL, H5S_ALL, H5P_DEFAULT, t2[0]);
0129     if (status < 0)
0130       throw cms::Exception("EcalLaserCondTool:HDF") << "Error while reading HD file.";
0131 
0132     H5Dclose(dset_t2);
0133     //H5Sclose(space);
0134 
0135     if (verb_ > 1)
0136       edm::LogPrint("EcalLaserCondTools") << " * reading t1 table "
0137                                           << "\n";
0138     dset = H5Dopen(file, "t1", H5P_DEFAULT);
0139     status = H5Dread(dset, H5T_NATIVE_INT, H5S_ALL, H5S_ALL, H5P_DEFAULT, t1);
0140     if (status < 0)
0141       throw cms::Exception("EcalLaserCondTool:HDF") << "Error while reading HD file.";
0142 
0143     H5Dclose(dset);
0144 
0145     if (verb_ > 1)
0146       edm::LogPrint("EcalLaserCondTools") << " * reading t3 table "
0147                                           << "\n";
0148     dset = H5Dopen(file, "t3", H5P_DEFAULT);
0149     status = H5Dread(dset, H5T_NATIVE_INT, H5S_ALL, H5S_ALL, H5P_DEFAULT, t3);
0150     if (status < 0)
0151       throw cms::Exception("EcalLaserCondTool:HDF") << "Error while reading HD file.";
0152 
0153     H5Dclose(dset);
0154 
0155     assert(EcalLaserCondTools::nLmes == nLME);
0156 
0157     // read crystal info IOV by IOV (otherwise too large)
0158     float p1[nCrystals], p2[nCrystals], p3[nCrystals];
0159     hsize_t iov_dim[1] = {nCrystals};
0160     memspace = H5Screate_simple(1, iov_dim, nullptr);
0161 
0162     EcalLaserAPDPNRatios corrSet;
0163     for (unsigned int iIov = skipIov_; iIov < nIovs && iIov < unsigned(nIovs_); ++iIov) {
0164       EcalLaserAPDPNRatios::EcalLaserTimeStamp t;
0165       iovStart = uint64_t(t1[iIov]) << 32;
0166       for (size_t iLme = 0; iLme < EcalLaserCondTools::nLmes; ++iLme) {
0167         t.t1 = edm::Timestamp(uint64_t(t1[iIov]) << 32);
0168         t.t2 = edm::Timestamp(uint64_t(t2[iIov][iLme]) << 32);
0169         t.t3 = edm::Timestamp(uint64_t(t3[iIov]) << 32);
0170         corrSet.setTime(iLme, t);
0171       }
0172 
0173       hsize_t offset[2] = {iIov, 0};      // shift rows: iIov, columns: 0
0174       hsize_t count[2] = {1, nCrystals};  // 1 row, nXtal columns
0175 
0176       dset = H5Dopen(file, "p1", H5P_DEFAULT);
0177       space = H5Dget_space(dset);
0178       status = H5Sselect_hyperslab(space, H5S_SELECT_SET, offset, nullptr, count, nullptr);
0179       if (status < 0)
0180         throw cms::Exception("EcalLaserCondTool:HDF") << "Error while reading HD file.";
0181 
0182       status = H5Dread(dset, H5T_NATIVE_FLOAT, memspace, space, H5P_DEFAULT, p1);
0183       if (status < 0)
0184         throw cms::Exception("EcalLaserCondTool:HDF") << "Error while reading HD file.";
0185 
0186       H5Dclose(dset);
0187       //H5Sclose(space);
0188 
0189       dset = H5Dopen(file, "p2", H5P_DEFAULT);
0190       space = H5Dget_space(dset);
0191       status = H5Sselect_hyperslab(space, H5S_SELECT_SET, offset, nullptr, count, nullptr);
0192       if (status < 0)
0193         throw cms::Exception("EcalLaserCondTool:HDF") << "Error while reading HD file.";
0194 
0195       status = H5Dread(dset, H5T_NATIVE_FLOAT, memspace, space, H5P_DEFAULT, p2);
0196       if (status < 0)
0197         throw cms::Exception("EcalLaserCondTool:HDF") << "Error while reading HD file.";
0198 
0199       H5Dclose(dset);
0200       //      H5Sclose(space);
0201 
0202       dset = H5Dopen(file, "p3", H5P_DEFAULT);
0203       space = H5Dget_space(dset);
0204       status = H5Sselect_hyperslab(space, H5S_SELECT_SET, offset, nullptr, count, nullptr);
0205       if (status < 0)
0206         throw cms::Exception("EcalLaserCondTool:HDF") << "Error while reading HD file.";
0207 
0208       status = H5Dread(dset, H5T_NATIVE_FLOAT, memspace, space, H5P_DEFAULT, p3);
0209       if (status < 0)
0210         throw cms::Exception("EcalLaserCondTool:HDF") << "Error while reading HD file.";
0211       H5Dclose(dset);
0212       H5Sclose(space);
0213 
0214       for (size_t iXtal = 0; iXtal < nCrystals; ++iXtal) {
0215         DetId detid = rawid[iXtal];
0216 
0217         EcalLaserAPDPNRatios::EcalLaserAPDPNpair corr = EcalLaserAPDPNRatios::EcalLaserAPDPNpair();
0218         corr.p1 = p1[iXtal];
0219         corr.p2 = p2[iXtal];
0220         corr.p3 = p3[iXtal];
0221 
0222         if (!std::isfinite(corr.p1) || !std::isfinite(corr.p2) || !std::isfinite(corr.p3) || corr.p1 < minP_ ||
0223             corr.p1 > maxP_ || corr.p2 < minP_ || corr.p2 > maxP_ || corr.p3 < minP_ || corr.p3 > maxP_) {
0224           fprintf(ferr_, "%d %d %f %f %f\n", t1[iIov], (int)detid, corr.p1, corr.p2, corr.p3);
0225           corr.p1 = corr.p2 = corr.p3 = 1;
0226         }
0227         corrSet.setValue((int)detid, corr);
0228       }
0229 
0230       try {
0231         //Write correction set in DB (one IOV):
0232         //if (db_->isNewTagRequest("EcalLaserAPDPNRatiosRcd")) {
0233         //  if (verb_)
0234         //      edm::LogPrint("EcalLaserCondTools") << "First IOV, extending starting time.\n";
0235         //  iovStart = db_->beginOfTime();
0236         //}
0237         timeval t;
0238         gettimeofday(&t, nullptr);
0239         if (verb_ > 1)
0240           edm::LogPrint("EcalLaserCondTools")
0241               << "[" << timeToString(t.tv_sec) << "] "
0242               << "Write IOV " << iIov << " starting from " << timeToString(iovStart >> 32) << "... ";
0243         db_->writeOneIOV(corrSet, iovStart, "EcalLaserAPDPNRatiosRcd");
0244       } catch (const cms::Exception& e) {
0245         if (verb_ > 1)
0246           edm::LogPrint("EcalLaserCondTools") << "Failed. ";
0247         edm::LogPrint("EcalLaserCondTools") << "Exception catched while writting to cond DB" << e.what() << "\n";
0248       }
0249       if (verb_ > 1)
0250         edm::LogPrint("EcalLaserCondTools") << "Suceeded.\n";
0251 
0252     }  // loop over IOVs
0253 
0254     H5Sclose(memspace);
0255     H5Fclose(file);
0256   }  // loop over input files
0257 }
0258 
0259 void EcalLaserCondTools::fillDb(CorrReader& r) {
0260   int iIov = 0;
0261   int processedIovs = 0;
0262   if (verb_ > 2)
0263     edm::LogPrint("EcalLaserCondTools") << "Starting filling DB...\n";
0264   int t1 = 0;
0265   int t3 = 0;
0266   int t2[nLmes];
0267 
0268   int prevT1 = 0;
0269   int prevT3 = 0;
0270   int prevT = 0;
0271   int t = 0;
0272   if (eventList_) {
0273     int iline = 0;
0274     while (!feof(eventList_)) {
0275       //skips comment lines:
0276       char c[2];
0277       while (fscanf(eventList_, " %1[#]%*[^\n]\n", &c[0]) == 1)
0278         ++iline;
0279 
0280       int n = fscanf(eventList_, "%*d %*d %*d %d%*[^\n]\n", &t);
0281       if (verb_ > 1)
0282         edm::LogPrint("EcalLaserCondTools") << "Event time: t = " << t << ", " << timeToString(t) << "\n";
0283       ++iline;
0284       if (n != 1)
0285         throw cms::Exception("User") << "Syntax error in event list file " << eventListFileName_ << " at line " << iline
0286                                      << " " << n << " "
0287                                      << ".\n";
0288       if (t < prevT)
0289         throw cms::Exception("User") << "Events in the event list file " << eventListFileName_
0290                                      << " are not ordered in increased time as required! See line " << iline << "\n";
0291       if (t == 0)
0292         throw cms::Exception("User") << "Found an unexpected t = 0 time value "
0293                                         "in the event list file"
0294                                      << eventListFileName_ << " at line " << iline << "!\n";
0295       //Look for IOV:
0296       bool iovFound = true;
0297       if (t <= t3) {  //IOV already inserted for previous event.
0298         if (verb_ > 1)
0299           edm::LogPrint("EcalLaserCondTools") << "Event in same IOV than previous one.\n";
0300         continue;
0301       }
0302 
0303       while ((iovFound = r.readTime(t1, t2, t3)) && t3 < t) /*NOP*/
0304         ;
0305 
0306       if (iovFound) {
0307         if (t1 < prevT1 && t3 < prevT3)
0308           throw cms::Exception("User")
0309               << "IOVs in the correction ascii file are not ordered in increased time as required!\n";
0310         else if (t1 < prevT1 || t3 < prevT3)
0311           throw cms::Exception("User") << "Found interleaved IOVs in the correction ascii file!\n";
0312         processIov(r, t1, t2, t3);
0313       } else {
0314         edm::LogPrint("EcalLaserCondTools") << "Warning: event beyond last IOV t3. Event time: " << timeToString(t)
0315                                             << ". Last IOV t3: " << timeToString(t3) << "\n";
0316       }
0317     }
0318   } else
0319     while (r.readTime(t1, t2, t3)) {
0320       ++iIov;
0321       if (iIov <= skipIov_) {
0322         edm::LogPrint("EcalLaserCondTools") << "Skipping IOV " << iIov << "\n";
0323         continue;
0324       } else if (processedIovs >= nIovs_ && nIovs_ >= 0) {
0325         edm::LogPrint("EcalLaserCondTools") << "Requested number of IOVs, " << nIovs_ << ", processed.\n";
0326         return;
0327       } else {
0328         processIov(r, t1, t2, t3);
0329         ++processedIovs;
0330       }
0331     }
0332 }
0333 
0334 void EcalLaserCondTools::processIov(CorrReader& r, int t1, int t2[EcalLaserCondTools::nLmes], int t3) {
0335   static int iIov = 0;
0336   ++iIov;
0337 
0338   //   FILE* fdebug = 0;
0339   //   if(i==) fdebug = fopen("debug.txt", "w");
0340 
0341   //  if(iIov <= skipIov_) { std::cout << "Skipping IOV " << iIov << "\n"; return; }
0342 
0343   cond::Time_t iovStart = 0;
0344 
0345   if (verb_ > 1) {
0346     edm::LogPrint("EcalLaserCondTools") << "t1:" << t1 << "(" << timeToString(t1) << ") \n"
0347                                         << "t3: " << t3 << "(" << timeToString(t3) << ")\nt2-t1: ";
0348     for (int i = 0; i < EcalLaserCondTools::nLmes; ++i) {
0349       edm::LogPrint("EcalLaserCondTools") << t2[i] - t1 << "\t";
0350     }
0351 
0352     edm::LogPrint("EcalLaserCondTools") << "\n";
0353   }
0354   if (t1 < fromTime_) {
0355     edm::LogPrint("EcalLaserCondTools") << "Skipping IOV " << iIov << ", "
0356                                         << ", which is before 'fromTime'," << timeToString(fromTime_) << "("
0357                                         << fromTime_ << ").\n";
0358     return;
0359   }
0360 
0361   if (toTime_ != -1 && t3 < toTime_) {
0362     edm::LogPrint("EcalLaserCondTools") << "Skipping IOV " << iIov << ", "
0363                                         << ", which is beyond 'toTime'," << timeToString(toTime_) << "(" << toTime_
0364                                         << ").\n";
0365     return;
0366   }
0367 
0368   if (t1 == 0) {
0369     edm::LogPrint("EcalLaserCondTools") << "Skipping IOV with t1 = 0"
0370                                         << "\n";
0371     return;
0372   }
0373 
0374   EcalLaserAPDPNRatios corrSet;
0375 
0376   EcalLaserAPDPNRatios::EcalLaserTimeStamp t;
0377   iovStart = uint64_t(t1) << 32;
0378   for (size_t i = 0; i < EcalLaserCondTools::nLmes; ++i) {
0379     t.t1 = edm::Timestamp(uint64_t(t1) << 32);
0380     t.t2 = edm::Timestamp(uint64_t(t2[i]) << 32);
0381     t.t3 = edm::Timestamp(uint64_t(t3) << 32);
0382     corrSet.setTime(i, t);
0383   }
0384 
0385   constexpr int ncrystals = 75848;
0386   std::set<int> detidList;
0387   for (int i = 0; i < ncrystals; ++i) {
0388     DetId detid;
0389     //EcalLaserAPDPNRatios::EcalLaserAPDPNpair corr = {0, 0, 0};
0390     EcalLaserAPDPNRatios::EcalLaserAPDPNpair corr = EcalLaserAPDPNRatios::EcalLaserAPDPNpair();
0391     if (verb_ > 2)
0392       edm::LogPrint("EcalLaserCondTools") << "Reading " << toNth(i + 1) << " crystal\n";
0393     if (!r.readPs(detid, corr)) {
0394       throw cms::Exception("LasCor") << "Failed to read " << toNth(i + 1) << " crystal correction.\n";
0395     }
0396 
0397     std::pair<std::set<int>::iterator, bool> res = detidList.insert(int(detid));
0398 
0399     if (!res.second) {  //detid already processed
0400       edm::LogPrint("EcalLaserCondTools")
0401           << "Duplicate det id, for IOV " << iIov << " t1 = " << t1 << " detid = " << int(detid) << "\n";
0402     }
0403 
0404     if (!std::isfinite(corr.p1) || !std::isfinite(corr.p2) || !std::isfinite(corr.p3) || corr.p1 < minP_ ||
0405         corr.p1 > maxP_ || corr.p2 < minP_ || corr.p2 > maxP_ || corr.p3 < minP_ || corr.p3 > maxP_) {
0406       fprintf(ferr_, "%d %d %f %f %f\n", t1, (int)detid, corr.p1, corr.p2, corr.p3);
0407       corr.p1 = corr.p2 = corr.p3 = 1;
0408     }
0409 
0410     if (verb_ > 2) {
0411       if (detid.subdetId() == EcalBarrel) {
0412         edm::LogPrint("EcalLaserCondTools") << EBDetId(detid);
0413       } else if (detid.subdetId() == EcalEndcap) {
0414         edm::LogPrint("EcalLaserCondTools") << EEDetId(detid);
0415       } else {
0416         edm::LogPrint("EcalLaserCondTools") << (int)detid;
0417       }
0418       edm::LogPrint("EcalLaserCondTools") << ": "
0419                                           << "p1 = " << corr.p1 << "\t"
0420                                           << "p2 = " << corr.p2 << "\t"
0421                                           << "p3 = " << corr.p3 << "\n";
0422     }
0423 
0424     corrSet.setValue((int)detid, corr);
0425   }
0426 
0427   try {
0428     //Write correction set in DB (one IOV):
0429     if (db_->isNewTagRequest("EcalLaserAPDPNRatiosRcd")) {
0430       if (verb_)
0431         edm::LogPrint("EcalLaserCondTools") << "First IOV, extending starting time.\n";
0432 
0433       iovStart = db_->beginOfTime();
0434     }
0435     timeval t;
0436     gettimeofday(&t, nullptr);
0437     if (verb_ > 1)
0438       edm::LogPrint("EcalLaserCondTools")
0439           << "[" << timeToString(t.tv_sec) << "] "
0440           << "Write IOV " << iIov << " starting from " << timeToString(iovStart >> 32) << "... ";
0441     db_->writeOneIOV(corrSet, iovStart, "EcalLaserAPDPNRatiosRcd");
0442   } catch (const cms::Exception& e) {
0443     edm::LogPrint("EcalLaserCondTools") << "Failed.\nException cathed while writting to cond DB" << e.what() << "\n";
0444   }
0445   edm::LogPrint("EcalLaserCondTools") << "Suceeded.\n";
0446 }
0447 
0448 bool EcalLaserCondTools::FileReader::nextFile() {
0449   for (;;) {
0450     ++ifile_;
0451     if (ifile_ >= fnames_.size()) {
0452       if (verb_ > 1)
0453         edm::LogPrint("EcalLaserCondTools") << "No more correction files.\n";
0454 
0455       return false;
0456     }
0457     if (verb_ > 1)
0458       edm::LogPrint("EcalLaserCondTools") << "Opening file " << fnames_[ifile_] << "\n";
0459 
0460     f_ = fopen(fnames_[ifile_].c_str(), "r");
0461     iline_ = 0;
0462     if (f_ == nullptr) {
0463       std::cerr << "Failed to open file " << fnames_[ifile_] << ". File skipped!\n";
0464     } else {
0465       return true;
0466     }
0467   }
0468 }
0469 
0470 bool EcalLaserCondTools::FileReader::readTime(int& t1, int t2[EcalLaserCondTools::nLmes], int& t3) {
0471   trim();
0472   if ((f_ == nullptr || feof(f_)) && !nextFile()) {
0473     if (verb_ > 1)
0474       edm::LogPrint("EcalLaserCondTools") << "No more record\n";
0475 
0476     return false;
0477   }
0478   int i;
0479   char* buf = nullptr;
0480   size_t s = 0;
0481   while ((i = fgetc(f_)) != 'T' && i != 'L' && i >= 0)
0482     getline(&buf, &s, f_);
0483   if (buf)
0484     free(buf);
0485   buf = nullptr;
0486 
0487   if (i == 'L') {  //last record put 3 consecutive times starting from end of prev. IOV
0488     t1 = t3;
0489     for (int i = 0; i < EcalLaserCondTools::nLmes; ++i)
0490       t2[i] = t1 + 1;
0491     t3 = t1 + 2;
0492     return true;
0493   }
0494 
0495   if (i != 'T') {
0496     if (verb_ > 1)
0497       edm::LogPrint("EcalLaserCondTools") << "No more record or bad line type/marker (getc returned " << i << ")\n";
0498 
0499     return false;
0500   }
0501 
0502   EcalLaserAPDPNRatios::EcalLaserTimeStamp t;
0503   int n = fscanf(f_, "%d %d", &t1, &t3);
0504   for (int i = 0; i < EcalLaserCondTools::nLmes; ++i) {
0505     int nn = fscanf(f_, "%d", &t2[i]);
0506     if (nn != 1)
0507       break;
0508     n += nn;
0509   }
0510 
0511   int nnn = fscanf(f_, " ");
0512 
0513   if (n != (2 + EcalLaserCondTools::nLmes) || nnn != 0)
0514     throw cms::Exception("LasCorFile") << "File " << fnames_[ifile_] << " line " << iline_
0515                                        << ": syntax error. Expecting 'T' marker followed by 94 values: "
0516                                        << "t1 t2 t3(lme 1) t3(lme 2) ... t3(lme " << EcalLaserCondTools::nLmes << ")\n";
0517 
0518   return true;
0519 }
0520 
0521 bool EcalLaserCondTools::FileReader::readPs(DetId& detid, EcalLaserAPDPNRatios::EcalLaserAPDPNpair& corr) {
0522   if (f_ == nullptr) {
0523     if (verb_)
0524       edm::LogPrint("EcalLaserCondTools") << "Requested to read p1..p3 parameter line while no file is closed.\n";
0525 
0526     return false;
0527   }
0528 
0529   trim();
0530   int i = fgetc(f_);
0531 
0532   if (i != 'P') {
0533     if (verb_ && i >= 0)
0534       edm::LogPrint("EcalLaserCondTools") << "File " << fnames_[ifile_] << " line " << iline_
0535                                           << ": unexpected line type, '" << (char)i << "' while expecting 'P'\n";
0536 
0537     if (verb_ && i < 0)
0538       edm::LogPrint("EcalLaserCondTools") << "Failed to read p1..p3 parameter line\n";
0539 
0540     return false;
0541   }
0542 
0543   int rawdetid;
0544   int n = fscanf(f_, "%d %f %f %f\n", &rawdetid, &corr.p1, &corr.p2, &corr.p3);
0545   ++iline_;
0546 
0547   if (n != 4) {
0548     //    corr.p2=corr.p1;
0549     //    corr.p3=corr.p1;
0550     throw cms::Exception("I/O") << "Syntax error at line " << iline_ << "of file " << fnames_[ifile_] << " read " << n
0551                                 << " values,"
0552                                 << " raw id" << rawdetid << ": " << corr.p1 << ", " << corr.p2;
0553   }
0554   detid = rawdetid;
0555   constexpr int ECALID = 3;
0556   if (detid.det() != ECALID)
0557     throw cms::Exception("InvalidValue") << "Line " << iline_ << "of file " << fnames_[ifile_]
0558                                          << " contains an invalid det ID (detector code is not ECAL!)\n";
0559   if (detid.subdetId() == EcalBarrel) {
0560     EBDetId ebDetId(detid);
0561     if (!EBDetId::validDetId(ebDetId.ietaAbs(), ebDetId.iphi()))
0562       throw cms::Exception("InvalidValue") << "Line " << iline_ << "of file " << fnames_[ifile_]
0563                                            << " contains an invalid det ID (detector code is not ECAL!)\n";
0564   }
0565   if (detid.subdetId() == EcalEndcap) {
0566     EEDetId eeDetId(detid);
0567     if (!EEDetId::validDetId(eeDetId.ix(), eeDetId.iy(), eeDetId.zside()))
0568       throw cms::Exception("InvalidValue") << "Line " << iline_ << "of file " << fnames_[ifile_]
0569                                            << " contains an invalid det ID (detector code is not ECAL!)\n";
0570   }
0571   ++iline_;
0572   return true;
0573 }
0574 
0575 void EcalLaserCondTools::FileReader::trim() {
0576   if (f_ == nullptr)
0577     return;
0578   bool skipLine = false;
0579   int c;
0580   while ((c = fgetc(f_)) >= 0 && (c == ' ' || c == '\t' || c == '\n' || c == '#' || skipLine)) {
0581     if (c == '#')
0582       skipLine = true;
0583     if (c == '\n') {
0584       ++iline_;
0585       skipLine = false;
0586     }
0587   }
0588   ungetc(c, f_);
0589 }
0590 
0591 std::string EcalLaserCondTools::toNth(int n) {
0592   std::stringstream s;
0593   s << n;
0594   if (n % 100 < 10 || n % 100 > 20) {
0595     switch (n % 10) {
0596       case 1:
0597         s << "st";
0598         break;
0599       case 2:
0600         s << "nd";
0601         break;
0602       case 3:
0603         s << "rd";
0604         break;
0605       default:
0606         s << "th";
0607     }
0608   } else {
0609     s << "th";
0610   }
0611   return s.str();
0612 }
0613 
0614 std::string EcalLaserCondTools::timeToString(time_t t) {
0615   char buf[256];
0616   struct tm lt;
0617   localtime_r(&t, &lt);
0618   strftime(buf, sizeof(buf), "%F %R:%S", &lt);
0619   buf[sizeof(buf) - 1] = 0;
0620   return std::string(buf);
0621 }
0622 
0623 void EcalLaserCondTools::dbToAscii(const edm::EventSetup& es) {
0624   const auto& laserAPDPNRatios = es.getData(laserAPDPNRatiosToken_);
0625 
0626   const EcalLaserAPDPNRatios::EcalLaserAPDPNRatiosMap& p = laserAPDPNRatios.getLaserMap();
0627   const EcalLaserAPDPNRatios::EcalLaserTimeStampMap& t = laserAPDPNRatios.getTimeMap();
0628 
0629   if (t.size() != EcalLaserCondTools::nLmes)
0630     throw cms::Exception("LasCor") << "Unexpected number time parameter triplets\n";
0631 
0632   if (fout_ == nullptr) {
0633     fout_ = fopen("corr_dump.txt", "w");
0634     if (fout_ == nullptr)
0635       throw cms::Exception("LasCor") << "Failed to create file corr_dump.txt\n";
0636   }
0637 
0638   unsigned t1 = t[0].t1.unixTime();
0639   unsigned t3 = t[0].t3.unixTime();
0640   fprintf(fout_, "T %d\t%d", t1, t3);
0641 
0642   if (verb_)
0643     edm::LogPrint("EcalLaserCondTools") << "Processing IOV " << t1 << " - " << t3 << "(" << timeToString(t1) << " - "
0644                                         << timeToString(t3) << "\n";
0645 
0646   for (unsigned i = 0; i < t.size(); ++i) {
0647     if (t[i].t1.unixTime() != t1 || t[i].t3.unixTime() != t3) {
0648       throw cms::Exception("LasCor") << "Inconsitency in t1, t3: "
0649                                      << "t1(lme 1) =" << t1 << ", t1(lme " << (i + 1) << ") = " << t[i].t1.unixTime()
0650                                      << ", t3(lme 1) =" << t3 << ", t3(lme " << (i + 1) << ") = " << t[i].t3.unixTime()
0651                                      << "\n";
0652     }
0653     fprintf(fout_, "\t%d", t[i].t2.unixTime());
0654   }
0655   fputc('\n', fout_);
0656   fflush(fout_);
0657 
0658   for (int ieta = -EBDetId::MAX_IETA; ieta <= EBDetId::MAX_IETA; ++ieta) {
0659     if (ieta == 0)
0660       continue;
0661     for (int iphi = EBDetId::MIN_IPHI; iphi <= EBDetId::MAX_IPHI; ++iphi) {
0662       if (EBDetId::validDetId(ieta, iphi)) {
0663         EBDetId detId(ieta, iphi);
0664         EcalLaserAPDPNRatios::EcalLaserAPDPNpair corr = p.barrel(detId.hashedIndex());
0665         fprintf(fout_, "P %d\t%f\t%f\t%f\n", (int)detId, corr.p1, corr.p2, corr.p3);
0666       }
0667     }
0668   }
0669 
0670   for (int iZ = 1; iZ >= -1; --iZ) {
0671     for (int iX = EEDetId::IX_MIN; iX <= EEDetId::IX_MAX; ++iX) {
0672       for (int iY = EEDetId::IY_MIN; iY <= EEDetId::IY_MAX; ++iY) {
0673         if (EEDetId::validDetId(iX, iY, iZ)) {
0674           EEDetId detId(iX, iY, iZ);
0675           EcalLaserAPDPNRatios::EcalLaserAPDPNpair corr = p.endcap(detId.hashedIndex());
0676           fprintf(fout_, "P %d\t%f\t%f\t%f\n", (int)detId, corr.p1, corr.p2, corr.p3);
0677         }
0678       }
0679     }
0680   }
0681 }
0682 //DEFINE_FWK_MODULE(EcalLaserCondTools);