CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
RPCSimSetUp.cc
Go to the documentation of this file.
22 
23 #include <cmath>
24 #include <math.h>
25 #include <fstream>
26 #include <sstream>
27 #include <iostream>
28 #include<cstring>
29 #include<string>
30 #include<vector>
31 #include<stdlib.h>
32 #include <utility>
33 #include <map>
34 
35 using namespace std;
36 
38 
39  _mapDetIdNoise.clear();
40  _mapDetIdEff.clear();
41  _bxmap.clear();
42  _clsMap.clear();
43 
44 }
45 
46 void RPCSimSetUp::setRPCSetUp(std::vector<RPCStripNoises::NoiseItem> vnoise, std::vector<float> vcls){
47 
48  double sum = 0;
49  unsigned int counter = 1;
50  unsigned int row = 1;
51  std::vector<double> sum_clsize;
52 
53  for(unsigned int n = 0; n < vcls.size(); ++n){
54 
55  sum_clsize.push_back(vcls[n]);
56 
57  if(counter == row*20) {
58 
59  _clsMap[row] = sum_clsize;
60  row++;
61  sum = 0;
62  sum_clsize.clear();
63  }
64  counter++;
65  }
66 
67  unsigned int n = 0;
68  uint32_t temp = 0;
69  std::vector<float> veff, vvnoise;
70  veff.clear();
71  vvnoise.clear();
72 
73  for(std::vector<RPCStripNoises::NoiseItem>::iterator it = vnoise.begin(); it != vnoise.end(); ++it){
74  if(n%96 == 0) {
75  if(n > 0 ){
76  _mapDetIdNoise[temp]= vvnoise;
77  _mapDetIdEff[temp] = veff;
78  _bxmap[RPCDetId(it->dpid)] = it->time;
79 
80  veff.clear();
81  vvnoise.clear();
82  vvnoise.push_back((it->noise));
83  veff.push_back((it->eff));
84  }
85  else if(n == 0 ){
86  vvnoise.push_back((it->noise));
87  veff.push_back((it->eff));
88  _bxmap[RPCDetId(it->dpid)] = it->time;
89  }
90  } else if (n == vnoise.size()-1 ){
91  temp = it->dpid;
92  vvnoise.push_back((it->noise));
93  veff.push_back((it->eff));
94  _mapDetIdNoise[temp]= vvnoise;
95  _mapDetIdEff[temp] = veff;
96  } else {
97  temp = it->dpid;
98  vvnoise.push_back((it->noise));
99  veff.push_back((it->eff));
100  }
101  n++;
102  }
103 }
104 
105 void RPCSimSetUp::setRPCSetUp(std::vector<RPCStripNoises::NoiseItem> vnoise, std::vector<RPCClusterSize::ClusterSizeItem> vClusterSize){
106 
107  std::vector<RPCClusterSize::ClusterSizeItem>::iterator itCls;
108  uint32_t detId;
109  int clsCounter(1);
110  std::vector<double> clsVect;
111 
112  for(itCls = vClusterSize.begin(); itCls != vClusterSize.end(); ++itCls){
113  clsVect.push_back(((double)(itCls->clusterSize)));
114  if((!(clsCounter%100)) && (clsCounter!=0)){
115  detId=itCls->dpid;
116  _mapDetClsMap[detId]=clsVect;
117  clsVect.clear();
118  clsCounter=0;
119  }
120  ++clsCounter;
121  }
122 
123  unsigned int n = 0;
124  uint32_t temp = 0;
125  std::vector<float> veff, vvnoise;
126  veff.clear();
127  vvnoise.clear();
128 
129  for(std::vector<RPCStripNoises::NoiseItem>::iterator it = vnoise.begin(); it != vnoise.end(); ++it){
130  if(n%96 == 0) {
131  if(n > 0 ){
132  _mapDetIdNoise[temp]= vvnoise;
133  _mapDetIdEff[temp] = veff;
134  _bxmap[RPCDetId(it->dpid)] = it->time;
135 
136  veff.clear();
137  vvnoise.clear();
138  vvnoise.push_back((it->noise));
139  veff.push_back((it->eff));
140  }
141  else if(n == 0 ){
142  vvnoise.push_back((it->noise));
143  veff.push_back((it->eff));
144  _bxmap[RPCDetId(it->dpid)] = it->time;
145  }
146  } else if (n == vnoise.size()-1 ){
147  temp = it->dpid;
148  vvnoise.push_back((it->noise));
149  veff.push_back((it->eff));
150  _mapDetIdNoise[temp]= vvnoise;
151  _mapDetIdEff[temp] = veff;
152  } else {
153  temp = it->dpid;
154  vvnoise.push_back((it->noise));
155  veff.push_back((it->eff));
156  }
157  n++;
158  }
159 }
160 
161 
162 const std::vector<float>& RPCSimSetUp::getNoise(uint32_t id)
163 {
164  map<uint32_t,std::vector<float> >::iterator iter = _mapDetIdNoise.find(id);
165  if(iter == _mapDetIdNoise.end()){
166  throw cms::Exception("DataCorrupt")
167  << "Exception comming from RPCSimSetUp - no noise information for DetId\t"<<id<< std::endl;
168  }
169  return iter->second;
170 }
171 
172 const std::vector<float>& RPCSimSetUp::getEff(uint32_t id)
173 {
174  map<uint32_t,std::vector<float> >::iterator iter = _mapDetIdEff.find(id);
175  if(iter == _mapDetIdEff.end()){
176  throw cms::Exception("DataCorrupt")
177  << "Exception comming from RPCSimSetUp - no efficiency information for DetId\t"<<id<< std::endl;
178  }
179  if((iter->second).size() != 96){
180  throw cms::Exception("DataCorrupt")
181  << "Exception comming from RPCSimSetUp - efficiency information in a wrong format for DetId\t"<<id<< std::endl;
182  }
183  return iter->second;
184 }
185 
186 float RPCSimSetUp::getTime(uint32_t id)
187 {
188  RPCDetId rpcid(id);
189  std::map<RPCDetId, float>::iterator iter = _bxmap.find(rpcid);
190  if(iter == _bxmap.end()){
191  throw cms::Exception("DataCorrupt")
192  << "Exception comming from RPCSimSetUp - no timing information for rpcid.rawId()\t"<<rpcid.rawId()<< std::endl;
193  }
194  return iter->second;
195 }
196 
197 const std::map< int, std::vector<double> >& RPCSimSetUp::getClsMap()
198 {
199  if(_clsMap.size()!=5){
200  throw cms::Exception("DataCorrupt")
201  << "Exception comming from RPCSimSetUp - cluster size - a wrong format "<< std::endl;
202  }
203  return _clsMap;
204 }
205 
206 
207 //const std::map<int, std::vector<double> >& RPCSimSetUp::getClsMap(uint32_t id)
208 const std::vector<double>& RPCSimSetUp::getCls(uint32_t id)
209 {
210 
211  map<uint32_t,std::vector<double> >::iterator iter = _mapDetClsMap.find(id);
212  if(iter == _mapDetClsMap.end()){
213  throw cms::Exception("DataCorrupt")
214  << "Exception comming from RPCSimSetUp - no cluster size information for DetId\t"<<id<< std::endl;
215  }
216  if((iter->second).size() != 100){
217  throw cms::Exception("DataCorrupt")
218  << "Exception comming from RPCSimSetUp - cluster size information in a wrong format for DetId\t"<<id<< std::endl;
219  }
220  return iter->second;
221 }
222 
float getTime(uint32_t id)
Definition: RPCSimSetUp.cc:186
RPCSimSetUp(const edm::ParameterSet &ps)
Definition: RPCSimSetUp.cc:37
const std::vector< float > & getEff(uint32_t id)
Definition: RPCSimSetUp.cc:172
uint32_t rawId() const
get the raw id
Definition: DetId.h:45
const std::vector< float > & getNoise(uint32_t id)
Definition: RPCSimSetUp.cc:162
virtual ~RPCSimSetUp()
Definition: RPCSimSetUp.cc:223
const std::vector< double > & getCls(uint32_t id)
Definition: RPCSimSetUp.cc:208
const std::map< int, std::vector< double > > & getClsMap()
Definition: RPCSimSetUp.cc:197
void setRPCSetUp(std::vector< RPCStripNoises::NoiseItem > vnoise, std::vector< float > vcls)
Definition: RPCSimSetUp.cc:46