CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
RPCCalibSetUp.cc
Go to the documentation of this file.
10 
11 #include <cmath>
12 #include <cstdlib>
13 #include <cstring>
14 #include <fstream>
15 #include <iostream>
16 #include <map>
17 #include <sstream>
18 #include <string>
19 #include <utility>
20 #include <vector>
21 
22 using namespace std;
23 
25  _mapDetIdNoise.clear();
26  _mapDetIdEff.clear();
27  _bxmap.clear();
28  _mapDetClsMap.clear();
29 
30  //------------------------ Noise Reading ----------------------------
31 
32  edm::FileInPath fp1 = ps.getParameter<edm::FileInPath>("noisemapfile");
33  std::ifstream _infile1(fp1.fullPath().c_str(), std::ios::in);
34 
35  std::vector<float> vnoise;
36 
37  int rpcdetid = 0;
38  std::string buff;
39 
40  std::vector<std::string> words;
41 
42  int count = 0;
43  while (getline(_infile1, buff, '\n')) {
44  words.clear();
45  vnoise.clear();
46 
47  stringstream ss;
48  std::string chname;
49  ss << buff;
50  ss >> chname >> rpcdetid;
51 
52  std::string::size_type pos = 0, prev_pos = 0;
53 
54  while ((pos = buff.find(" ", pos)) != string::npos) {
55  words.push_back(buff.substr(prev_pos, pos - prev_pos));
56  prev_pos = ++pos;
57  }
58  words.push_back(buff.substr(prev_pos, pos - prev_pos));
59 
60  for (unsigned int i = 2; i < words.size(); ++i) {
61  float value = atof(((words)[i]).c_str());
62  vnoise.push_back(value);
63  }
64 
65  _mapDetIdNoise.insert(make_pair(static_cast<uint32_t>(rpcdetid), vnoise));
66 
67  count++;
68  }
69  _infile1.close();
70 
71  //------------------------ Eff Reading ----------------------------
72 
73  edm::FileInPath fp2 = ps.getParameter<edm::FileInPath>("effmapfile");
74  std::ifstream _infile2(fp2.fullPath().c_str(), std::ios::in);
75 
76  std::vector<float> veff;
77  rpcdetid = 0;
78 
79  while (getline(_infile2, buff, '\n')) {
80  words.clear();
81  veff.clear();
82 
83  stringstream ss;
84  std::string chname;
85  ss << buff;
86  ss >> chname >> rpcdetid;
87 
88  std::string::size_type pos = 0, prev_pos = 0;
89  while ((pos = buff.find(" ", pos)) != string::npos) {
90  words.push_back(buff.substr(prev_pos, pos - prev_pos));
91  prev_pos = ++pos;
92  }
93  words.push_back(buff.substr(prev_pos, pos - prev_pos));
94 
95  for (unsigned int i = 2; i < words.size(); ++i) {
96  float value = atof(((words)[i]).c_str());
97  veff.push_back(value);
98  }
99  _mapDetIdEff.insert(make_pair(static_cast<uint32_t>(rpcdetid), veff));
100  }
101  _infile2.close();
102 
103  //---------------------- Timing reading ------------------------------------
104 
105  edm::FileInPath fp3 = ps.getParameter<edm::FileInPath>("timingMap");
106  std::ifstream _infile3(fp3.fullPath().c_str(), std::ios::in);
107 
108  uint32_t detUnit = 0;
109  float timing = 0.;
110  while (!_infile3.eof()) {
111  _infile3 >> detUnit >> timing;
112  _bxmap[RPCDetId(detUnit)] = timing;
113  }
114  _infile3.close();
115 
116  //---------------------- Cluster size --------------------------------------
117 
118  edm::FileInPath fp4 = ps.getParameter<edm::FileInPath>("clsmapfile");
119  std::ifstream _infile4(fp4.fullPath().c_str(), ios::in);
120 
121  string buffer;
122  double sum = 0;
123  unsigned int counter = 1;
124  unsigned int row = 1;
125  std::vector<double> sum_clsize;
126 
127  while (_infile4 >> buffer) {
128  const char *buffer1 = buffer.c_str();
129  double dato = atof(buffer1);
130  sum += dato;
131  sum_clsize.push_back(sum);
132 
133  if (counter == row * 20) {
134  _clsMap[row] = sum_clsize;
135  row++;
136  sum = 0;
137  sum_clsize.clear();
138  }
139  counter++;
140  }
141  _infile4.close();
142 
143  //---------------------- Cluster size Chamber by Chamber -------------------
144 
145  edm::FileInPath fp5 = ps.getParameter<edm::FileInPath>("clsidmapfile");
146  std::ifstream _infile5(fp5.fullPath().c_str(), ios::in);
147 
148  std::vector<double> vClsDistrib;
149  rpcdetid = 0;
150 
151  while (getline(_infile5, buff, '\n')) {
152  words.clear();
153  vClsDistrib.clear();
154 
155  stringstream ss1;
156  ss1 << buff;
157  ss1 >> rpcdetid;
158 
159  std::string::size_type pos = 0, prev_pos = 0;
160  while ((pos = buff.find(" ", pos)) != string::npos) {
161  words.push_back(buff.substr(prev_pos, pos - prev_pos));
162  prev_pos = ++pos;
163  }
164  words.push_back(buff.substr(prev_pos, pos - prev_pos));
165 
166  float clusterSizeSumData(0.);
167 
168  for (unsigned int i = 1; i < words.size(); ++i) {
169  float value = atof(((words)[i]).c_str());
170 
171  clusterSizeSumData += value;
172  vClsDistrib.push_back(clusterSizeSumData);
173  if (!(i % 20)) {
174  clusterSizeSumData = 0.;
175  }
176  }
177  if (vClsDistrib.size() != 100) {
178  throw cms::Exception("DataCorrupt") << "Exception comming from RPCCalibSetUp - cluster size - a wrong "
179  "format "
180  << std::endl;
181  }
182  _mapDetClsMap.insert(make_pair(static_cast<uint32_t>(rpcdetid), vClsDistrib));
183  std::cout << "_mapDetClsMap.size()\t" << _mapDetClsMap.size() << std::endl;
184  }
185 
186  _infile5.close();
187 }
188 
189 std::vector<float> RPCCalibSetUp::getNoise(uint32_t id) {
190  map<uint32_t, std::vector<float>>::iterator iter = _mapDetIdNoise.find(id);
191  if (iter == _mapDetIdNoise.end()) {
192  throw cms::Exception("DataCorrupt") << "Exception comming from RPCCalibSetUp - no noise information for "
193  "DetId\t"
194  << id << std::endl;
195  }
196  return (iter->second);
197 }
198 
199 std::vector<float> RPCCalibSetUp::getEff(uint32_t id) {
200  map<uint32_t, std::vector<float>>::iterator iter = _mapDetIdEff.find(id);
201  if (iter == _mapDetIdEff.end()) {
202  throw cms::Exception("DataCorrupt") << "Exception comming from RPCCalibSetUp - no efficiency information "
203  "for DetId\t"
204  << id << std::endl;
205  }
206  if ((iter->second).size() != 96) {
207  throw cms::Exception("DataCorrupt") << "Exception comming from RPCCalibSetUp - efficiency information in a "
208  "wrong format for DetId\t"
209  << id << std::endl;
210  }
211  return iter->second;
212 }
213 
214 float RPCCalibSetUp::getTime(uint32_t id) {
215  RPCDetId rpcid(id);
216 
217  std::map<RPCDetId, float>::iterator iter = _bxmap.find(rpcid);
218  if (iter == _bxmap.end()) {
219  throw cms::Exception("DataCorrupt") << "Exception comming from RPCCalibSetUp - no timing information for "
220  "rpcid.rawId()\t"
221  << rpcid.rawId() << std::endl;
222  }
223  return iter->second;
224 }
225 
226 std::map<int, std::vector<double>> RPCCalibSetUp::getClsMap() {
227  if (_clsMap.size() != 5) {
228  throw cms::Exception("DataCorrupt") << "Exception comming from RPCCalibSetUp - cluster size - a wrong "
229  "format "
230  << std::endl;
231  }
232  return _clsMap;
233 }
234 
235 std::vector<double> RPCCalibSetUp::getCls(uint32_t id) {
236  std::map<uint32_t, std::vector<double>>::iterator iter = _mapDetClsMap.find(id);
237  if (iter == _mapDetClsMap.end()) {
238  throw cms::Exception("DataCorrupt") << "Exception comming from RPCCalibSetUp - no cluster size information "
239  "for DetId\t"
240  << id << std::endl;
241  }
242  if ((iter->second).size() != 100) {
243  throw cms::Exception("DataCorrupt") << "Exception comming from RPCCalibSetUp - cluster size information in "
244  "a wrong format for DetId\t"
245  << id << std::endl;
246  }
247  return iter->second;
248 }
249 
std::vector< double > getCls(uint32_t id)
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
std::vector< float > getEff(uint32_t id)
uint16_t size_type
std::map< int, std::vector< double > > getClsMap()
virtual ~RPCCalibSetUp()
float getTime(uint32_t id)
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
static std::atomic< unsigned int > counter
std::string fullPath() const
Definition: FileInPath.cc:161
RPCCalibSetUp(const edm::ParameterSet &ps)
tuple cout
Definition: gather_cfg.py:144
std::vector< float > getNoise(uint32_t id)