CMS 3D CMS Logo

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