CMS 3D CMS Logo

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