Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-05-10 03:53:10

0001 #include "CalibMuon/RPCCalibration/interface/RPCCalibSetUp.h"
0002 #include "DataFormats/Common/interface/Handle.h"
0003 #include "DataFormats/MuonDetId/interface/RPCDetId.h"
0004 #include "FWCore/Framework/interface/ESHandle.h"
0005 #include "FWCore/Framework/interface/Event.h"
0006 #include "FWCore/Framework/interface/EventSetup.h"
0007 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0008 #include "FWCore/ServiceRegistry/interface/Service.h"
0009 #include "Geometry/Records/interface/MuonGeometryRecord.h"
0010 
0011 #include <cmath>
0012 #include <cstdlib>
0013 #include <cstring>
0014 #include <fstream>
0015 #include <iostream>
0016 #include <map>
0017 #include <sstream>
0018 #include <string>
0019 #include <utility>
0020 #include <vector>
0021 
0022 using namespace std;
0023 
0024 RPCCalibSetUp::RPCCalibSetUp(const edm::ParameterSet &ps) {
0025   _mapDetIdNoise.clear();
0026   _mapDetIdEff.clear();
0027   _bxmap.clear();
0028   _mapDetClsMap.clear();
0029 
0030   //------------------------ Noise Reading ----------------------------
0031 
0032   edm::FileInPath fp1 = ps.getParameter<edm::FileInPath>("noisemapfile");
0033   std::ifstream _infile1(fp1.fullPath().c_str(), std::ios::in);
0034 
0035   std::vector<float> vnoise;
0036 
0037   int rpcdetid = 0;
0038   std::string buff;
0039 
0040   std::vector<std::string> words;
0041 
0042   while (getline(_infile1, buff, '\n')) {
0043     words.clear();
0044     vnoise.clear();
0045 
0046     stringstream ss;
0047     std::string chname;
0048     ss << buff;
0049     ss >> chname >> rpcdetid;
0050 
0051     std::string::size_type pos = 0, prev_pos = 0;
0052 
0053     while ((pos = buff.find("  ", pos)) != string::npos) {
0054       words.push_back(buff.substr(prev_pos, pos - prev_pos));
0055       prev_pos = ++pos;
0056     }
0057     words.push_back(buff.substr(prev_pos, pos - prev_pos));
0058 
0059     for (unsigned int i = 2; i < words.size(); ++i) {
0060       float value = atof(((words)[i]).c_str());
0061       vnoise.push_back(value);
0062     }
0063 
0064     _mapDetIdNoise.insert(make_pair(static_cast<uint32_t>(rpcdetid), vnoise));
0065   }
0066   _infile1.close();
0067 
0068   //------------------------ Eff Reading ----------------------------
0069 
0070   edm::FileInPath fp2 = ps.getParameter<edm::FileInPath>("effmapfile");
0071   std::ifstream _infile2(fp2.fullPath().c_str(), std::ios::in);
0072 
0073   std::vector<float> veff;
0074   rpcdetid = 0;
0075 
0076   while (getline(_infile2, buff, '\n')) {
0077     words.clear();
0078     veff.clear();
0079 
0080     stringstream ss;
0081     std::string chname;
0082     ss << buff;
0083     ss >> chname >> rpcdetid;
0084 
0085     std::string::size_type pos = 0, prev_pos = 0;
0086     while ((pos = buff.find("  ", pos)) != string::npos) {
0087       words.push_back(buff.substr(prev_pos, pos - prev_pos));
0088       prev_pos = ++pos;
0089     }
0090     words.push_back(buff.substr(prev_pos, pos - prev_pos));
0091 
0092     for (unsigned int i = 2; i < words.size(); ++i) {
0093       float value = atof(((words)[i]).c_str());
0094       veff.push_back(value);
0095     }
0096     _mapDetIdEff.insert(make_pair(static_cast<uint32_t>(rpcdetid), veff));
0097   }
0098   _infile2.close();
0099 
0100   //---------------------- Timing reading ------------------------------------
0101 
0102   edm::FileInPath fp3 = ps.getParameter<edm::FileInPath>("timingMap");
0103   std::ifstream _infile3(fp3.fullPath().c_str(), std::ios::in);
0104 
0105   uint32_t detUnit = 0;
0106   float timing = 0.;
0107   while (!_infile3.eof()) {
0108     _infile3 >> detUnit >> timing;
0109     _bxmap[RPCDetId(detUnit)] = timing;
0110   }
0111   _infile3.close();
0112 
0113   //---------------------- Cluster size --------------------------------------
0114 
0115   edm::FileInPath fp4 = ps.getParameter<edm::FileInPath>("clsmapfile");
0116   std::ifstream _infile4(fp4.fullPath().c_str(), ios::in);
0117 
0118   string buffer;
0119   double sum = 0;
0120   unsigned int counter = 1;
0121   unsigned int row = 1;
0122   std::vector<double> sum_clsize;
0123 
0124   while (_infile4 >> buffer) {
0125     const char *buffer1 = buffer.c_str();
0126     double dato = atof(buffer1);
0127     sum += dato;
0128     sum_clsize.push_back(sum);
0129 
0130     if (counter == row * 20) {
0131       _clsMap[row] = sum_clsize;
0132       row++;
0133       sum = 0;
0134       sum_clsize.clear();
0135     }
0136     counter++;
0137   }
0138   _infile4.close();
0139 
0140   //---------------------- Cluster size Chamber by Chamber -------------------
0141 
0142   edm::FileInPath fp5 = ps.getParameter<edm::FileInPath>("clsidmapfile");
0143   std::ifstream _infile5(fp5.fullPath().c_str(), ios::in);
0144 
0145   std::vector<double> vClsDistrib;
0146   rpcdetid = 0;
0147 
0148   while (getline(_infile5, buff, '\n')) {
0149     words.clear();
0150     vClsDistrib.clear();
0151 
0152     stringstream ss1;
0153     ss1 << buff;
0154     ss1 >> rpcdetid;
0155 
0156     std::string::size_type pos = 0, prev_pos = 0;
0157     while ((pos = buff.find("  ", pos)) != string::npos) {
0158       words.push_back(buff.substr(prev_pos, pos - prev_pos));
0159       prev_pos = ++pos;
0160     }
0161     words.push_back(buff.substr(prev_pos, pos - prev_pos));
0162 
0163     float clusterSizeSumData(0.);
0164 
0165     for (unsigned int i = 1; i < words.size(); ++i) {
0166       float value = atof(((words)[i]).c_str());
0167 
0168       clusterSizeSumData += value;
0169       vClsDistrib.push_back(clusterSizeSumData);
0170       if (!(i % 20)) {
0171         clusterSizeSumData = 0.;
0172       }
0173     }
0174     if (vClsDistrib.size() != 100) {
0175       throw cms::Exception("DataCorrupt") << "Exception comming from RPCCalibSetUp - cluster size - a wrong "
0176                                              "format "
0177                                           << std::endl;
0178     }
0179     _mapDetClsMap.insert(make_pair(static_cast<uint32_t>(rpcdetid), vClsDistrib));
0180     std::cout << "_mapDetClsMap.size()\t" << _mapDetClsMap.size() << std::endl;
0181   }
0182 
0183   _infile5.close();
0184 }
0185 
0186 std::vector<float> RPCCalibSetUp::getNoise(uint32_t id) {
0187   map<uint32_t, std::vector<float>>::iterator iter = _mapDetIdNoise.find(id);
0188   if (iter == _mapDetIdNoise.end()) {
0189     throw cms::Exception("DataCorrupt") << "Exception comming from RPCCalibSetUp - no noise information for "
0190                                            "DetId\t"
0191                                         << id << std::endl;
0192   }
0193   return (iter->second);
0194 }
0195 
0196 std::vector<float> RPCCalibSetUp::getEff(uint32_t id) {
0197   map<uint32_t, std::vector<float>>::iterator iter = _mapDetIdEff.find(id);
0198   if (iter == _mapDetIdEff.end()) {
0199     throw cms::Exception("DataCorrupt") << "Exception comming from RPCCalibSetUp - no efficiency information "
0200                                            "for DetId\t"
0201                                         << id << std::endl;
0202   }
0203   if ((iter->second).size() != 96) {
0204     throw cms::Exception("DataCorrupt") << "Exception comming from RPCCalibSetUp - efficiency information in a "
0205                                            "wrong format for DetId\t"
0206                                         << id << std::endl;
0207   }
0208   return iter->second;
0209 }
0210 
0211 float RPCCalibSetUp::getTime(uint32_t id) {
0212   RPCDetId rpcid(id);
0213 
0214   std::map<RPCDetId, float>::iterator iter = _bxmap.find(rpcid);
0215   if (iter == _bxmap.end()) {
0216     throw cms::Exception("DataCorrupt") << "Exception comming from RPCCalibSetUp - no timing information for "
0217                                            "rpcid.rawId()\t"
0218                                         << rpcid.rawId() << std::endl;
0219   }
0220   return iter->second;
0221 }
0222 
0223 std::map<int, std::vector<double>> RPCCalibSetUp::getClsMap() {
0224   if (_clsMap.size() != 5) {
0225     throw cms::Exception("DataCorrupt") << "Exception comming from RPCCalibSetUp - cluster size - a wrong "
0226                                            "format "
0227                                         << std::endl;
0228   }
0229   return _clsMap;
0230 }
0231 
0232 std::vector<double> RPCCalibSetUp::getCls(uint32_t id) {
0233   std::map<uint32_t, std::vector<double>>::iterator iter = _mapDetClsMap.find(id);
0234   if (iter == _mapDetClsMap.end()) {
0235     throw cms::Exception("DataCorrupt") << "Exception comming from RPCCalibSetUp - no cluster size information "
0236                                            "for DetId\t"
0237                                         << id << std::endl;
0238   }
0239   if ((iter->second).size() != 100) {
0240     throw cms::Exception("DataCorrupt") << "Exception comming from RPCCalibSetUp - cluster size information in "
0241                                            "a wrong format for DetId\t"
0242                                         << id << std::endl;
0243   }
0244   return iter->second;
0245 }
0246 
0247 RPCCalibSetUp::~RPCCalibSetUp() {}