Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2022-02-18 08:23:02

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