CMS 3D CMS Logo

RPCSimSetUp.cc
Go to the documentation of this file.
19 
20 #include <cmath>
21 #include <cmath>
22 #include <fstream>
23 #include <sstream>
24 #include <iostream>
25 #include <cstring>
26 #include <string>
27 #include <vector>
28 #include <cstdlib>
29 #include <utility>
30 #include <map>
31 
32 using namespace std;
33 
35  _mapDetIdNoise.clear();
36  _mapDetIdEff.clear();
37  _bxmap.clear();
38  _clsMap.clear();
39 }
40 
41 void RPCSimSetUp::setRPCSetUp(const std::vector<RPCStripNoises::NoiseItem>& vnoise, const std::vector<float>& vcls) {
42  unsigned int counter = 1;
43  unsigned int row = 1;
44  std::vector<double> sum_clsize;
45 
46  for (unsigned int n = 0; n < vcls.size(); ++n) {
47  sum_clsize.push_back(vcls[n]);
48 
49  if (counter == row * 20) {
50  _clsMap[row] = sum_clsize;
51  row++;
52  sum_clsize.clear();
53  }
54  counter++;
55  }
56 
57  unsigned int n = 0;
58  uint32_t temp = 0;
59  std::vector<float> veff, vvnoise;
60  veff.clear();
61  vvnoise.clear();
62 
63  for (std::vector<RPCStripNoises::NoiseItem>::const_iterator it = vnoise.begin(); it != vnoise.end(); ++it) {
64  if (n % 96 == 0) {
65  if (n > 0) {
66  _mapDetIdNoise[temp] = vvnoise;
67  _mapDetIdEff[temp] = veff;
68  _bxmap[RPCDetId(it->dpid)] = it->time;
69 
70  veff.clear();
71  vvnoise.clear();
72  vvnoise.push_back((it->noise));
73  veff.push_back((it->eff));
74  } else if (n == 0) {
75  vvnoise.push_back((it->noise));
76  veff.push_back((it->eff));
77  _bxmap[RPCDetId(it->dpid)] = it->time;
78  }
79  } else if (n == vnoise.size() - 1) {
80  temp = it->dpid;
81  vvnoise.push_back((it->noise));
82  veff.push_back((it->eff));
83  _mapDetIdNoise[temp] = vvnoise;
84  _mapDetIdEff[temp] = veff;
85  } else {
86  temp = it->dpid;
87  vvnoise.push_back((it->noise));
88  veff.push_back((it->eff));
89  }
90  n++;
91  }
92 }
93 
94 void RPCSimSetUp::setRPCSetUp(const std::vector<RPCStripNoises::NoiseItem>& vnoise,
95  const std::vector<RPCClusterSize::ClusterSizeItem>& vClusterSize) {
96  LogDebug("RPCSimSetup") << "RPCSimSetUp::setRPCSetUp(vector<NoiseItem>, vector<ClusterSizeItem>)" << std::endl;
97 
98  uint32_t detId = 0, current_detId, this_detId;
99  RPCDetId rpcId, current_rpcId, this_rpcId;
100  const RPCRoll* current_roll = nullptr;
101  const RPCRoll* this_roll = nullptr;
102  unsigned int current_nStrips;
103 
104  LogDebug("RPCSimSetup") << "RPCSimSetUp::setRPCSetUp :: ClusterSizeItem :: begin" << std::endl;
105 #ifdef EDM_ML_DEBUG
106  std::stringstream sslogclsitem;
107 #endif
108  // ### ClusterSizeItem #######################################################
109  std::vector<RPCClusterSize::ClusterSizeItem>::const_iterator itCls;
110  int clsCounter(1);
111  std::vector<double> clsVect;
112  // ### loop for New Format (120 entries)
113  for (itCls = vClusterSize.begin(); itCls != vClusterSize.end(); ++itCls) {
114  clsVect.push_back(((double)(itCls->clusterSize)));
115 #ifdef EDM_ML_DEBUG
116  sslogclsitem << " Push back clustersize = " << itCls->clusterSize << std::endl;
117  sslogclsitem << "Filling cls in _mapDetCls[detId,clsVect] :: detId = " << detId;
118  sslogclsitem << " --> will it be accepted? clsCounter = " << clsCounter << " accepted?";
119  sslogclsitem << " New Format ::" << ((!(clsCounter % 120)) && (clsCounter != 0)); // <<std::endl;
120  sslogclsitem << " Old Format ::" << ((!(clsCounter % 100)) && (clsCounter != 0)); // <<std::endl;
121  sslogclsitem << std::endl;
122 #endif
123  // New Format :: loop until 120
124  if ((!(clsCounter % 120)) && (clsCounter != 0)) {
125  detId = itCls->dpid;
126  _mapDetClsMap[detId] = clsVect;
127 #ifdef EDM_ML_DEBUG
128  std::stringstream LogDebugClsVectString;
129  LogDebugClsVectString << "[";
130  for (std::vector<double>::iterator itClsVect = clsVect.begin(); itClsVect != clsVect.end(); ++itClsVect) {
131  LogDebugClsVectString << *itClsVect << ",";
132  }
133  LogDebugClsVectString << "]";
134  std::string LogDebugClsVectStr = LogDebugClsVectString.str();
135  LogDebug("RPCSimSetup") << "Filling clsVect in _mapDetCls[detId,clsVect] :: detId = " << RPCDetId(detId) << " = "
136  << detId << " clsVec = " << LogDebugClsVectStr;
137 
138  sslogclsitem << " --> New Method ";
139  sslogclsitem << " --> saved in map " << std::endl;
140  sslogclsitem << "Filling cls in _mapDetClsMap[detId,clsVect] :: detId = " << detId;
141  sslogclsitem << " --> will it be accepted? clsCounter = " << clsCounter << " accepted? "
142  << ((!(clsCounter % 120)) && (clsCounter != 0)) << std::endl;
143 #endif
144  clsVect.clear();
145  clsCounter = 0;
146  } else {
147 #ifdef EDM_ML_DEBUG
148  sslogclsitem << " --> not saved in map " << std::endl;
149 #endif
150  }
151  ++clsCounter;
152  }
153  // ### loop for Old Format (100 entries)
154  for (itCls = vClusterSize.begin(); itCls != vClusterSize.end(); ++itCls) {
155  clsVect.push_back(((double)(itCls->clusterSize)));
156 #ifdef EDM_ML_DEBUG
157  sslogclsitem << " Push back clustersize = " << itCls->clusterSize << std::endl;
158  sslogclsitem << "Filling cls in _mapDetClsMapLegacy[detId,clsVect] :: detId = " << detId;
159  sslogclsitem << " --> will it be accepted? clsCounter = " << clsCounter << " accepted?";
160  sslogclsitem << " New Format ::" << ((!(clsCounter % 120)) && (clsCounter != 0)); // <<std::endl;
161  sslogclsitem << " Old Format ::" << ((!(clsCounter % 100)) && (clsCounter != 0)); // <<std::endl;
162  sslogclsitem << std::endl;
163 #endif
164  // Old Format :: same until 100
165  if ((!(clsCounter % 100)) && (clsCounter != 0)) {
166  detId = itCls->dpid;
167  _mapDetClsMapLegacy[detId] = clsVect;
168 #ifdef EDM_ML_DEBUG
169  std::stringstream LogDebugClsVectString;
170  LogDebugClsVectString << "[";
171  for (std::vector<double>::iterator itClsVect = clsVect.begin(); itClsVect != clsVect.end(); ++itClsVect) {
172  LogDebugClsVectString << *itClsVect << ",";
173  }
174  LogDebugClsVectString << "]";
175  std::string LogDebugClsVectStr = LogDebugClsVectString.str();
176  LogDebug("RPCSimSetup") << "Filling clsVect in _mapDetClsLegacy[detId,clsVect] :: detId = " << RPCDetId(detId)
177  << " = " << detId << " clsVec = " << LogDebugClsVectStr;
178 
179  sslogclsitem << " --> Old Method ";
180  sslogclsitem << " --> saved in map " << std::endl;
181  sslogclsitem << "Filling cls in _mapDetClsMapLegacy[detId,clsVect] :: detId = " << detId;
182  sslogclsitem << " --> will it be accepted? clsCounter = " << clsCounter << " accepted? "
183  << ((!(clsCounter % 120)) && (clsCounter != 0)) << std::endl;
184 #endif
185  clsVect.clear();
186  clsCounter = 0;
187  } else {
188 #ifdef EDM_ML_DEBUG
189  sslogclsitem << " --> not saved in map " << std::endl;
190 #endif
191  }
192  ++clsCounter;
193  }
194  // ###########################################################################
195 #ifdef EDM_ML_DEBUG
196  std::string logclsitem = sslogclsitem.str();
197  sslogclsitem.clear();
198  LogDebug("RPCSimSetupClsLoopDetails") << logclsitem << std::endl;
199  LogDebug("RPCSimSetup") << "RPCSimSetUp::setRPCSetUp :: ClusterSizeItem :: end" << std::endl;
200 
201  LogDebug("RPCSimSetup") << "RPCSimSetUp::setRPCSetUp :: NoiseItem :: begin" << std::endl;
202  std::stringstream sslognoiseitem;
203 #endif
204  // ### NoiseItem #############################################################
205  unsigned int count_strips = 1;
206 #ifdef EDM_ML_DEBUG
207  unsigned int count_all = 1;
208 #endif
209  std::vector<float> vveff, vvnoise;
210 
211  // DetId to start with needs to be a DetId inside the Geometry used
212  // Therefore loop on the NoiseItems and search for the first valid roll in the Geometry
213  // Assign this as the DetId to start with (so called current_roll) and quit the loop
214  bool quitLoop = false;
215  current_detId = 0;
216  current_nStrips = 0; // current_rpcId = 0; current_roll = 0;
217  for (std::vector<RPCStripNoises::NoiseItem>::const_iterator it = vnoise.begin(); it != vnoise.end() && !quitLoop;
218  ++it) {
219  // roll associated to the conditions of this strip (iterator)
220  current_detId = it->dpid;
221  current_rpcId = RPCDetId(current_detId);
222  // Test whether this roll (picked up from the conditions) is inside the RPC Geometry
223  const RPCRoll* roll = theGeometry->roll(current_rpcId);
224  if (roll == nullptr) {
225 #ifdef EDM_ML_DEBUG
226  sslognoiseitem << "Searching for first valid detid :: current_detId = " << current_detId;
227  sslognoiseitem << " aka " << current_rpcId << " is not in current Geometry --> Skip " << std::endl;
228 #endif
229  continue;
230  } else {
231 #ifdef EDM_ML_DEBUG
232  sslognoiseitem << "Searching for first valid detid :: current_detId = " << current_detId;
233  sslognoiseitem << " aka " << current_rpcId
234  << " is the first (valid) roll in the current Geometry --> Accept, Assign & Quit Loop"
235  << std::endl;
236 #endif
237  current_roll = theGeometry->roll(current_rpcId);
238  current_nStrips = current_roll->nstrips();
239  quitLoop = true;
240  }
241  }
242 
243 #ifdef EDM_ML_DEBUG
244  sslognoiseitem << "Start Position :: current_detId = " << current_detId << " aka " << current_rpcId;
245  sslognoiseitem << " is a valid roll with pointer " << current_roll << " and has "
246  << (current_roll ? current_roll->nstrips() : 0) << " strips" << std::endl;
247  sslognoiseitem << " -------------------------------------------------------------------------------------------------"
248  "------------------------------------ "
249  << std::endl;
250 #endif
251  for (std::vector<RPCStripNoises::NoiseItem>::const_iterator it = vnoise.begin(); it != vnoise.end(); ++it) {
252  // roll associated to the conditions of this strip (iterator)
253  this_detId = it->dpid;
254  this_rpcId = RPCDetId(this_detId);
255  // Test whether this roll (picked up from the conditions) is inside the RPC Geometry
256  const RPCRoll* roll = theGeometry->roll(this_rpcId);
257  if (roll == nullptr) {
258 #ifdef EDM_ML_DEBUG
259  sslognoiseitem << "Inside Loop :: [" << std::setw(6) << count_all << "][" << std::setw(3) << count_strips
260  << "] :: this_detId = " << this_detId << " aka " << this_rpcId
261  << " which is not in current Geometry --> Skip " << std::endl;
262 #endif
263  continue;
264  }
265 
266  // Case 1 :: FIRST ENTRY
267  // ---------------------
268  if (this_detId == current_detId && count_strips == 1) {
269  // fill bx in map
270  _bxmap[current_detId] = it->time;
271  // clear vectors
272  vveff.clear();
273  vvnoise.clear();
274  // fill the vectors
275  vvnoise.push_back((it->noise));
276  vveff.push_back((it->eff));
277 #ifdef EDM_ML_DEBUG
278  sslognoiseitem << "RPCSimSetUp::setRPCSetUp :: NoiseItem :: case 1" << std::endl;
279  sslognoiseitem << this_detId << " = " << this_rpcId << " with " << roll->nstrips() << " strips" << std::endl;
280  sslognoiseitem << "[NoiseItem :: n = " << count_all
281  << "] Filling time in _bxmap[detId] :: detId = " << RPCDetId(it->dpid) << " time = " << it->time
282  << std::endl;
283  sslognoiseitem << "First Value :: [" << std::setw(6) << count_all << "][" << std::setw(3) << count_strips
284  << "] :: this_detId = " << this_detId << " aka " << this_rpcId;
285  sslognoiseitem << " Strip " << std::setw(3) << count_strips << " Noise = " << it->noise << " Hz/cm2" << std::endl;
286  // update one counter
287  ++count_all;
288 #endif
289  // update the other counter
290  ++count_strips;
291  }
292  // Case 2 :: 2ND ENTRY --> LAST-1 ENTRY
293  // ------------------------------------
294  else if (this_detId == current_detId && count_strips > 1 && count_strips < current_nStrips) {
295 #ifdef EDM_ML_DEBUG
296  sslognoiseitem << "RPCSimSetUp::setRPCSetUp :: NoiseItem :: case 2" << std::endl;
297  sslognoiseitem << "Inside Loop :: [" << std::setw(6) << count_all << "][" << std::setw(3) << count_strips
298  << "] :: this_detId = " << this_detId << " aka " << this_rpcId;
299  sslognoiseitem << " Strip " << std::setw(3) << count_strips << " Noise = " << it->noise << " Hz/cm2" << std::endl;
300  // update one counter
301  ++count_all;
302 #endif
303  // fill the vectors
304  vvnoise.push_back((it->noise));
305  vveff.push_back((it->eff));
306  // update the other counter
307  ++count_strips;
308  }
309 
310  // Case 3 :: LAST ENTRY
311  // --------------------
312  else if (this_detId == current_detId && count_strips == current_nStrips) {
313 #ifdef EDM_ML_DEBUG
314  sslognoiseitem << "RPCSimSetUp::setRPCSetUp :: NoiseItem :: case 3" << std::endl;
315  sslognoiseitem << "Last Value :: [" << std::setw(6) << count_all << "][" << std::setw(3) << count_strips
316  << "] :: this_detId = " << this_detId << " aka " << this_rpcId;
317  sslognoiseitem << " Strip " << std::setw(3) << count_strips << " Noise = " << it->noise << " Hz/cm2" << std::endl;
318  // update one counter
319  ++count_all;
320 #endif
321  // fill last value in the vector
322  vvnoise.push_back((it->noise));
323  vveff.push_back((it->eff));
324  // update the other counter
325  ++count_strips;
326  // fill vectors into map
327  _mapDetIdNoise[current_detId] = vvnoise;
328  _mapDetIdEff[current_detId] = vveff;
329 
330 #ifdef EDM_ML_DEBUG
331  sslognoiseitem << " fill vectors into map" << std::endl;
332  std::stringstream LogDebugNoiVectString, LogDebugEffVectString;
333  LogDebugNoiVectString << "[";
334  for (std::vector<float>::iterator itNoiVect = vvnoise.begin(); itNoiVect != vvnoise.end(); ++itNoiVect) {
335  LogDebugNoiVectString << (*itNoiVect) << ",";
336  }
337  LogDebugNoiVectString << "]";
338  std::string LogDebugNoiVectStr = LogDebugNoiVectString.str();
339  LogDebugEffVectString << "[";
340  for (std::vector<float>::iterator itEffVect = vveff.begin(); itEffVect != vveff.end(); ++itEffVect) {
341  LogDebugEffVectString << (*itEffVect) << ",";
342  }
343  LogDebugEffVectString << "]";
344  std::string LogDebugEffVectStr = LogDebugEffVectString.str();
345  LogDebug("RPCSimSetup") << "Filling vvnoise in _mapDetIdNoise[detId] :: detId = " << RPCDetId(it->dpid) << " = "
346  << (RPCDetId(it->dpid)).rawId() << " vvnoise = " << LogDebugNoiVectStr;
347  LogDebug("RPCSimSetup") << "Filling veff in _mapDetIdEff[detId] :: detId = " << RPCDetId(it->dpid) << " = "
348  << (RPCDetId(it->dpid)).rawId() << " veff = " << LogDebugEffVectStr;
349 #endif
350  // look for next different detId and rename it to the current_detId
351  // at this point we skip all the conditions for the strips that are not in this roll
352  // and we will go to the conditions for the first strip of the next roll
353  bool next_detId_found = false;
354 #ifdef EDM_ML_DEBUG
355  sslognoiseitem << "look for next different detId" << std::endl;
356 #endif
357  while (next_detId_found == 0 && it != vnoise.end() - 1) {
358  ++it;
359  this_detId = it->dpid;
360  this_rpcId = RPCDetId(this_detId);
361  this_roll = theGeometry->roll(this_rpcId);
362  if (!this_roll)
363  continue;
364 #ifdef EDM_ML_DEBUG
365  sslognoiseitem << "Inside While:: [" << std::setw(6) << count_all << "][" << std::setw(3) << count_strips
366  << "] :: this_detId = " << this_detId << " aka " << this_rpcId << " Noise = " << it->noise
367  << " Hz/cm2" << std::endl;
368 // ++count_all;
369 #endif
370  ++count_strips;
371  if (this_detId != current_detId) {
372 #ifdef EDM_ML_DEBUG
373  sslognoiseitem << "Different detId is found :: " << this_detId << " aka " << this_rpcId
374  << " Noise = " << it->noise << " Hz/cm2";
375 #endif
376  // next roll is found. update current_detId to this newly found detId
377  // and update also the number of strips
378  current_detId = this_detId;
379  current_rpcId = RPCDetId(current_detId);
380  next_detId_found = true;
381  current_nStrips = (theGeometry->roll(current_rpcId))->nstrips();
382 #ifdef EDM_ML_DEBUG
383  sslognoiseitem << " with " << current_nStrips << " strips" << std::endl;
384 #endif
385  --it; // subtract one, because at the end of the loop the iterator will be increased with one
386  // in fact the treatment for roll N stops when we find the first occurence of roll N+1
387  // however we want to start the treatment for roll N+1 with the first occurence of roll N+1
388  // so the first entry of each new roll N+1 is manipulated twice in the loop (once as a stop, once as a start)
389  // therefore we have to manipulate the iterator here, subtracting one, to treat again this entry
390  }
391  }
392  // reset count_strips
393  count_strips = 1;
394  }
395  // There should be no Case 4
396  // -------------------------
397  else {
398  }
399  }
400  // ###########################################################################
401 #ifdef EDM_ML_DEBUG
402  std::string lognoiseitem = sslognoiseitem.str();
403  sslognoiseitem.clear();
404  LogDebug("RPCSimSetupNoiseLoopDetails") << lognoiseitem << std::endl;
405  LogDebug("RPCSimSetup") << "RPCSimSetUp::setRPCSetUp :: NoiseItem :: end" << std::endl;
406 
407  LogDebug("RPCSimSetup") << "RPCSimSetUp::setRPCSetUp :: end" << std::endl;
408 #endif
409 }
410 
411 const std::vector<float>& RPCSimSetUp::getNoise(uint32_t id) {
412  std::map<uint32_t, std::vector<float> >::iterator iter = _mapDetIdNoise.find(id);
413  if (iter == _mapDetIdNoise.end()) {
414  throw cms::Exception("DataCorrupt") << "Exception from RPCSimSetUp - no noise information for DetId\t" << id
415  << std::endl;
416  }
417  LogDebug("RPCSimSetupChecks") << "All OK from RPCSimSetUp - noise information for DetId\t" << id << std::endl;
418  return iter->second;
419 }
420 
421 const std::vector<float>& RPCSimSetUp::getEff(uint32_t id) {
422  std::map<uint32_t, std::vector<float> >::iterator iter = _mapDetIdEff.find(id);
423  if (iter == _mapDetIdEff.end()) {
424  throw cms::Exception("DataCorrupt") << "Exception from RPCSimSetUp - no efficiency information for DetId\t" << id
425  << std::endl;
426  }
427 
428  RPCDetId rpcId = RPCDetId(id);
429  const RPCRoll* roll = theGeometry->roll(rpcId);
430  unsigned int numbStrips = roll->nstrips();
431 
432  if ((iter->second).size() < numbStrips) {
433  LogDebug("RPCSimSetup") << "Exception from RPCSimSetUp - efficiency information in a wrong format for DetId\t" << id
434  << " aka " << RPCDetId(id) << std::endl;
435  LogDebug("RPCSimSetup") << " number of strips in Conditions\t" << (iter->second).size()
436  << " number of strips in Geometry\t" << numbStrips << std::endl;
437  throw cms::Exception("DataCorrupt")
438  << "Exception from RPCSimSetUp - efficiency information in a wrong format for DetId\t" << id << std::endl;
439  }
440 
441  return iter->second;
442 }
443 
444 float RPCSimSetUp::getTime(uint32_t id) {
445  RPCDetId rpcid(id);
446  std::map<RPCDetId, float>::iterator iter = _bxmap.find(rpcid);
447  if (iter == _bxmap.end()) {
448  throw cms::Exception("DataCorrupt") << "Exception from RPCSimSetUp - no timing information for rpcid.rawId()\t"
449  << rpcid.rawId() << std::endl;
450  }
451  return iter->second;
452 }
453 
454 const std::map<int, std::vector<double> >& RPCSimSetUp::getClsMap() {
455  if (_clsMap.size() != 5) {
456  throw cms::Exception("DataCorrupt") << "Exception from RPCSimSetUp - cluster size - a wrong format " << std::endl;
457  }
458  return _clsMap;
459 }
460 
461 //const std::map<int, std::vector<double> >& RPCSimSetUp::getClsMap(uint32_t id)
462 const std::vector<double>& RPCSimSetUp::getCls(uint32_t id) //legacy member function
463 {
464  LogDebug("RPCSimSetupChecks") << "RPCSimSetUp::getCls" << std::endl;
465 
466  std::map<uint32_t, std::vector<double> >::iterator iter = _mapDetClsMapLegacy.find(id);
467  if (iter == _mapDetClsMapLegacy.end()) {
468  throw cms::Exception("DataCorrupt") << "Exception from RPCSimSetUp - no cluster size information for DetId\t" << id
469  << std::endl;
470  }
471  if ((iter->second).size() != 100) {
472  throw cms::Exception("DataCorrupt")
473  << "Exception from RPCSimSetUp - _mapDetClsMapLegacy - cluster size information in a wrong format for DetId\t"
474  << id << std::endl;
475  }
476  LogDebug("RPCSimSetupChecks")
477  << "All OK from RPCSimSetUp - _mapDetClsMapLegacy - cluster size information for DetId\t" << id << std::endl;
478  return iter->second;
479 }
480 
481 const std::vector<double>& RPCSimSetUp::getAsymmetricClsDistribution(uint32_t id, uint32_t slice) {
482  LogDebug("RPCSimSetupChecks") << "RPCSimSetUp::getAsymmetricClsDistribution" << std::endl;
483 
484  std::map<uint32_t, std::vector<double> >::const_iterator iter = _mapDetClsMap.find(id);
485  if (iter == _mapDetClsMap.end()) {
486  throw cms::Exception("DataCorrupt")
487  << "Exception from RPCSimSetUp - _mapDetClsMap - no cluster size information for DetId\t" << id << std::endl;
488  }
489  if ((iter->second).size() != 120) {
490  throw cms::Exception("DataCorrupt")
491  << "Exception from RPCSimSetUp - _mapDetClsMap - cluster size information in a wrong format for DetId\t" << id
492  << std::endl;
493  }
494  // return iter->second;
495 
496  std::vector<double> dataForAsymmCls = iter->second;
497  if (slice > 4) {
498  throw cms::Exception("DataCorrupt") << "Exception from RPCSimSetUp - slice variable not in the range" << std::endl;
499  }
500 
501  _DetClsAsymmetric.clear();
502 
503  std::vector<double> clsFewStripsDistribution;
504  std::vector<double> clsDistribution;
505  std::vector<double> clsAccumulativeDistribution;
506 
507  std::map<int, std::vector<double> > mapSliceVsDistribution;
508 
509  const int slices = 5;
510  const int distributionFewStrips = 24;
511 
512  double sliceVsFewStripsDistribution[slices][distributionFewStrips];
513 
514  for (int j = 0; j < distributionFewStrips; j++) {
515  for (int i = 0; i < slices; i++) {
516  sliceVsFewStripsDistribution[i][j] = dataForAsymmCls[j * slices + i];
517  }
518  }
519 
520  int i = slice;
521  double sum = 0;
522  int counter = 0;
523  for (int j = 0; j < distributionFewStrips; j++) {
524  counter++;
525  sum += sliceVsFewStripsDistribution[i][j];
526  if (counter % 4 == 0) {
527  _DetClsAsymmetric.push_back(sum);
528  }
529  }
530  return _DetClsAsymmetric;
531 }
532 
533 const std::vector<double>& RPCSimSetUp::getAsymmetryForCls(uint32_t id, uint32_t slice, uint32_t cls) {
534  LogDebug("RPCSimSetupChecks") << "RPCSimSetUp::getAsymmetryForCls" << std::endl;
535 
536  std::map<uint32_t, std::vector<double> >::const_iterator iter = _mapDetClsMap.find(id);
537  if (iter == _mapDetClsMap.end()) {
538  throw cms::Exception("DataCorrupt")
539  << "Exception from RPCSimSetUp - _mapDetClsMap - no cluster size information for DetId\t" << id << std::endl;
540  }
541  if ((iter->second).size() != 120) {
542  throw cms::Exception("DataCorrupt")
543  << "Exception from RPCSimSetUp - _mapDetClsMap - cluster size information in a wrong format for DetId\t" << id
544  << '\t' << (iter->second).size() << std::endl;
545  }
546 
547  std::vector<double> dataForAsymmCls = iter->second;
548 
549  if (slice > 4) {
550  throw cms::Exception("DataCorrupt") << "Exception from RPCSimSetUp - slice variable not in the range" << std::endl;
551  }
552 
553  _DetAsymmetryForCls.clear();
554 
555  std::vector<double> clsFewStripsDistribution;
556  std::vector<double> clsDistribution;
557  std::vector<double> clsAccumulativeDistribution;
558  std::vector<double> clsDetAsymmetryForCls;
559  clsDetAsymmetryForCls.clear();
560 
561  std::map<int, std::vector<double> > mapSliceVsDistribution;
562 
563  const int slices = 5;
564  const int distributionFewStrips = 24;
565 
566  double sliceVsFewStripsDistribution[slices][distributionFewStrips];
567 
568  for (int j = 0; j < distributionFewStrips; j++) {
569  for (int i = 0; i < slices; i++) {
570  sliceVsFewStripsDistribution[i][j] = dataForAsymmCls[j * slices + i];
571  }
572  }
573 
574  int vector_lenght;
575  switch (cls) {
576  case 1:
577  case 3:
578  case 5:
579  vector_lenght = 3;
580  break;
581  case 2:
582  case 4:
583  vector_lenght = 4;
584  break;
585  case 6:
586  default:
587  vector_lenght = 1;
588  break;
589  }
590 
591  float sum = 0;
592  float value;
593  for (int i = 0; i < vector_lenght; i++) {
594  value = sliceVsFewStripsDistribution[slice][(cls - 1) * 4 + i];
595  clsDetAsymmetryForCls.push_back(value);
596  sum += value;
597  // LogDebug ("RPCSimSetup")<<"value\t"<<value<<std::endl;
598  // LogDebug ("RPCSimSetup")<<"sum\t"<<sum<<std::endl;
599  }
600 
601  float accum = 0;
602  for (int i = clsDetAsymmetryForCls.size() - 1; i > -1; i--) {
603  accum += clsDetAsymmetryForCls[i];
604  _DetAsymmetryForCls.push_back(accum / sum);
605  }
606  return _DetAsymmetryForCls;
607 }
608 
float getTime(uint32_t id)
Definition: RPCSimSetUp.cc:444
RPCSimSetUp(const edm::ParameterSet &ps)
Definition: RPCSimSetUp.cc:34
const std::vector< double > & getAsymmetricClsDistribution(uint32_t id, uint32_t slice)
Definition: RPCSimSetUp.cc:481
const std::vector< float > & getEff(uint32_t id)
Definition: RPCSimSetUp.cc:421
int nstrips() const
Definition: RPCRoll.cc:24
const std::vector< float > & getNoise(uint32_t id)
Definition: RPCSimSetUp.cc:411
virtual ~RPCSimSetUp()
Definition: RPCSimSetUp.cc:609
const std::vector< double > & getCls(uint32_t id)
Definition: RPCSimSetUp.cc:462
const std::map< int, std::vector< double > > & getClsMap()
Definition: RPCSimSetUp.cc:454
Definition: value.py:1
const std::vector< double > & getAsymmetryForCls(uint32_t id, uint32_t slice, uint32_t cls)
Definition: RPCSimSetUp.cc:533
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
void setRPCSetUp(const std::vector< RPCStripNoises::NoiseItem > &vnoise, const std::vector< float > &vcls)
Definition: RPCSimSetUp.cc:41
#define LogDebug(id)