File indexing completed on 2024-04-06 12:02:49
0001
0002
0003
0004
0005
0006
0007
0008 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0009 #include "FWCore/Framework/interface/ESHandle.h"
0010 #include "FWCore/Framework/interface/Event.h"
0011 #include "FWCore/Framework/interface/MakerMacros.h"
0012 #include <cmath>
0013 #include <fstream>
0014 #include <iostream>
0015 #include <map>
0016 #include <sstream>
0017 #include <stdexcept>
0018 #include <string>
0019
0020 #include "FWCore/Framework/interface/EventSetup.h"
0021 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0022
0023 #include "CondTools/DT/test/validate/DTROMapValidateDBRead.h"
0024 #include "CondFormats/DTObjects/interface/DTReadOutMapping.h"
0025 #include "CondFormats/DataRecord/interface/DTReadOutMappingRcd.h"
0026
0027 DTROMapValidateDBRead::DTROMapValidateDBRead(edm::ParameterSet const& p)
0028 : dataFileName(p.getParameter<std::string>("chkFile")),
0029 elogFileName(p.getParameter<std::string>("logFile")),
0030 dtreadoutmappingToken_(esConsumes()) {}
0031
0032 DTROMapValidateDBRead::DTROMapValidateDBRead(int i) : dtreadoutmappingToken_(esConsumes()) {}
0033
0034 void DTROMapValidateDBRead::analyze(const edm::Event& e, const edm::EventSetup& context) {
0035 using namespace edm::eventsetup;
0036
0037 std::cout << " I AM IN RUN NUMBER " << e.id().run() << std::endl;
0038 std::cout << " ---EVENT NUMBER " << e.id().event() << std::endl;
0039 std::stringstream run_fn;
0040
0041
0042 std::ifstream chkFile(dataFileName.c_str());
0043 std::ofstream logFile(elogFileName.c_str(), std::ios_base::app);
0044 auto ro = context.getHandle(dtreadoutmappingToken_);
0045 std::cout << ro->mapRobRos() << " " << ro->mapCellTdc() << std::endl;
0046 std::cout << std::distance(ro->begin(), ro->end()) << " data in the container" << std::endl;
0047 int whe;
0048 int sta;
0049 int sec;
0050 int qua;
0051 int lay;
0052 int cel;
0053 int ddu;
0054 int ros;
0055 int rob;
0056 int tdc;
0057 int cha;
0058 int ckwhe;
0059 int cksta;
0060 int cksec;
0061 int ckqua;
0062 int cklay;
0063 int ckcel;
0064 int ckddu;
0065 int ckros;
0066 int ckrob;
0067 int cktdc;
0068 int ckcha;
0069 int status;
0070 DTReadOutMapping::const_iterator iter = ro->begin();
0071 DTReadOutMapping::const_iterator iend = ro->end();
0072 while (iter != iend) {
0073 const DTReadOutGeometryLink& link = *iter;
0074 status = ro->readOutToGeometry(
0075 link.dduId, link.rosId, link.robId, link.tdcId, link.channelId, whe, sta, sec, qua, lay, cel);
0076 if (status)
0077 logFile << "ERROR while getting chan->cell map " << link.dduId << " " << link.rosId << " " << link.robId << " "
0078 << link.tdcId << " " << link.channelId << " , status = " << status << std::endl;
0079 if ((link.wheelId != whe) || (link.stationId != sta) || (link.sectorId != sec) || (link.slId != qua) ||
0080 (link.layerId != lay) || (link.cellId != cel))
0081 logFile << "MISMATCH WHEN READING chan->cell " << link.dduId << " " << link.rosId << " " << link.robId << " "
0082 << link.tdcId << " " << link.channelId << " : " << link.wheelId << " " << link.stationId << " "
0083 << link.sectorId << " " << link.slId << " " << link.layerId << " " << link.cellId << " -> " << whe << " "
0084 << sta << " " << sec << " " << qua << " " << lay << " " << cel << std::endl;
0085 status = ro->geometryToReadOut(
0086 link.wheelId, link.stationId, link.sectorId, link.slId, link.layerId, link.cellId, ddu, ros, rob, tdc, cha);
0087 if (status)
0088 logFile << "ERROR while getting cell->chan map " << link.wheelId << " " << link.stationId << " " << link.sectorId
0089 << " " << link.slId << " " << link.layerId << " " << link.cellId << " , status = " << status << std::endl;
0090 if ((link.dduId != ddu) || (link.rosId != ros) || (link.robId != rob) || (link.tdcId != tdc) ||
0091 (link.channelId != cha))
0092 logFile << "MISMATCH WHEN READING cell->chan " << link.wheelId << " " << link.stationId << " " << link.sectorId
0093 << " " << link.slId << " " << link.layerId << " " << link.cellId << " : " << link.dduId << " "
0094 << link.rosId << " " << link.robId << " " << link.tdcId << " " << link.channelId << " -> " << ddu << " "
0095 << ros << " " << rob << " " << tdc << " " << cha << std::endl;
0096 iter++;
0097 }
0098
0099 while (chkFile >> ckddu >> ckros >> ckrob >> cktdc >> ckcha >> ckwhe >> cksta >> cksec >> ckqua >> cklay >> ckcel) {
0100 status = ro->readOutToGeometry(ckddu, ckros, ckrob, cktdc, ckcha, whe, sta, sec, qua, lay, cel);
0101 if ((ckwhe != whe) || (cksta != sta) || (cksec != sec) || (ckqua != qua) || (cklay != lay) || (ckcel != cel))
0102 logFile << "MISMATCH IN WRITING AND READING chan->cell map " << ckddu << " " << ckros << " " << ckrob << " "
0103 << cktdc << " " << ckcha << " : " << ckwhe << " " << cksta << " " << cksec << " " << ckqua << " " << cklay
0104 << " " << ckcel << " -> " << whe << " " << sta << " " << sec << " " << qua << " " << lay << " " << cel
0105 << std::endl;
0106 status = ro->geometryToReadOut(ckwhe, cksta, cksec, ckqua, cklay, ckcel, ddu, ros, rob, tdc, cha);
0107 if ((ckddu != ddu) || (ckros != ros) || (ckrob != rob) || (cktdc != tdc) || (ckcha != cha))
0108 logFile << "MISMATCH IN WRITING AND READING cell->chan map " << ckwhe << " " << cksta << " " << cksec << " "
0109 << ckqua << " " << cklay << " " << ckcel << " : " << ckddu << " " << ckros << " " << ckrob << " " << cktdc
0110 << " " << ckcha << " -> " << ddu << " " << ros << " " << rob << " " << tdc << " " << cha << std::endl;
0111 }
0112 }
0113
0114 void DTROMapValidateDBRead::endJob() {
0115 std::ifstream logFile(elogFileName.c_str());
0116 char* line = new char[1000];
0117 int errors = 0;
0118 std::cout << "ReadOut Map validation result:" << std::endl;
0119 while (logFile.getline(line, 1000)) {
0120 std::cout << line << std::endl;
0121 errors++;
0122 }
0123 if (!errors) {
0124 std::cout << " ********************************* " << std::endl;
0125 std::cout << " *** *** " << std::endl;
0126 std::cout << " *** NO ERRORS FOUND *** " << std::endl;
0127 std::cout << " *** *** " << std::endl;
0128 std::cout << " ********************************* " << std::endl;
0129 }
0130 return;
0131 }
0132
0133 DEFINE_FWK_MODULE(DTROMapValidateDBRead);