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  unsigned int count_all = 1;
207  std::vector<float> vveff, vvnoise;
208 
209  // DetId to start with needs to be a DetId inside the Geometry used
210  // Therefore loop on the NoiseItems and search for the first valid roll in the Geometry
211  // Assign this as the DetId to start with (so called current_roll) and quit the loop
212  bool quitLoop = false;
213  current_detId = 0;
214  current_nStrips = 0; // current_rpcId = 0; current_roll = 0;
215  for (std::vector<RPCStripNoises::NoiseItem>::const_iterator it = vnoise.begin(); it != vnoise.end() && !quitLoop;
216  ++it) {
217  // roll associated to the conditions of this strip (iterator)
218  current_detId = it->dpid;
219  current_rpcId = RPCDetId(current_detId);
220  // Test whether this roll (picked up from the conditions) is inside the RPC Geometry
221  const RPCRoll* roll = theGeometry->roll(current_rpcId);
222  if (roll == nullptr) {
223 #ifdef EDM_ML_DEBUG
224  sslognoiseitem << "Searching for first valid detid :: current_detId = " << current_detId;
225  sslognoiseitem << " aka " << current_rpcId << " is not in current Geometry --> Skip " << std::endl;
226 #endif
227  continue;
228  } else {
229 #ifdef EDM_ML_DEBUG
230  sslognoiseitem << "Searching for first valid detid :: current_detId = " << current_detId;
231  sslognoiseitem << " aka " << current_rpcId
232  << " is the first (valid) roll in the current Geometry --> Accept, Assign & Quit Loop"
233  << std::endl;
234 #endif
235  current_roll = theGeometry->roll(current_rpcId);
236  current_nStrips = current_roll->nstrips();
237  quitLoop = true;
238  }
239  }
240 
241 #ifdef EDM_ML_DEBUG
242  sslognoiseitem << "Start Position :: current_detId = " << current_detId << " aka " << current_rpcId;
243  sslognoiseitem << " is a valid roll with pointer " << current_roll << " and has "
244  << (current_roll ? current_roll->nstrips() : 0) << " strips" << std::endl;
245  sslognoiseitem << " -------------------------------------------------------------------------------------------------"
246  "------------------------------------ "
247  << std::endl;
248 #endif
249  for (std::vector<RPCStripNoises::NoiseItem>::const_iterator it = vnoise.begin(); it != vnoise.end(); ++it) {
250  // roll associated to the conditions of this strip (iterator)
251  this_detId = it->dpid;
252  this_rpcId = RPCDetId(this_detId);
253  // Test whether this roll (picked up from the conditions) is inside the RPC Geometry
254  const RPCRoll* roll = theGeometry->roll(this_rpcId);
255  if (roll == nullptr) {
256 #ifdef EDM_ML_DEBUG
257  sslognoiseitem << "Inside Loop :: [" << std::setw(6) << count_all << "][" << std::setw(3) << count_strips
258  << "] :: this_detId = " << this_detId << " aka " << this_rpcId
259  << " which is not in current Geometry --> Skip " << std::endl;
260 #endif
261  continue;
262  }
263 
264  // Case 1 :: FIRST ENTRY
265  // ---------------------
266  if (this_detId == current_detId && count_strips == 1) {
267  // fill bx in map
268  _bxmap[current_detId] = it->time;
269  // clear vectors
270  vveff.clear();
271  vvnoise.clear();
272  // fill the vectors
273  vvnoise.push_back((it->noise));
274  vveff.push_back((it->eff));
275 #ifdef EDM_ML_DEBUG
276  sslognoiseitem << "RPCSimSetUp::setRPCSetUp :: NoiseItem :: case 1" << std::endl;
277  sslognoiseitem << this_detId << " = " << this_rpcId << " with " << roll->nstrips() << " strips" << std::endl;
278  sslognoiseitem << "[NoiseItem :: n = " << count_all
279  << "] Filling time in _bxmap[detId] :: detId = " << RPCDetId(it->dpid) << " time = " << it->time
280  << std::endl;
281  sslognoiseitem << "First Value :: [" << std::setw(6) << count_all << "][" << std::setw(3) << count_strips
282  << "] :: this_detId = " << this_detId << " aka " << this_rpcId;
283  sslognoiseitem << " Strip " << std::setw(3) << count_strips << " Noise = " << it->noise << " Hz/cm2" << std::endl;
284 #endif
285  // update counter
286  ++count_strips;
287  ++count_all;
288  }
289  // Case 2 :: 2ND ENTRY --> LAST-1 ENTRY
290  // ------------------------------------
291  else if (this_detId == current_detId && count_strips > 1 && count_strips < current_nStrips) {
292 #ifdef EDM_ML_DEBUG
293  sslognoiseitem << "RPCSimSetUp::setRPCSetUp :: NoiseItem :: case 2" << std::endl;
294  sslognoiseitem << "Inside Loop :: [" << std::setw(6) << count_all << "][" << std::setw(3) << count_strips
295  << "] :: this_detId = " << this_detId << " aka " << this_rpcId;
296  sslognoiseitem << " Strip " << std::setw(3) << count_strips << " Noise = " << it->noise << " Hz/cm2" << std::endl;
297 #endif
298  // fill the vectors
299  vvnoise.push_back((it->noise));
300  vveff.push_back((it->eff));
301  // update counter
302  ++count_strips;
303  ++count_all;
304  }
305 
306  // Case 3 :: LAST ENTRY
307  // --------------------
308  else if (this_detId == current_detId && count_strips == current_nStrips) {
309 #ifdef EDM_ML_DEBUG
310  sslognoiseitem << "RPCSimSetUp::setRPCSetUp :: NoiseItem :: case 3" << std::endl;
311  sslognoiseitem << "Last Value :: [" << std::setw(6) << count_all << "][" << std::setw(3) << count_strips
312  << "] :: this_detId = " << this_detId << " aka " << this_rpcId;
313  sslognoiseitem << " Strip " << std::setw(3) << count_strips << " Noise = " << it->noise << " Hz/cm2" << std::endl;
314 #endif
315  // fill last value in the vector
316  vvnoise.push_back((it->noise));
317  vveff.push_back((it->eff));
318  // update counter
319  ++count_strips;
320  ++count_all;
321  // fill vectors into map
322  _mapDetIdNoise[current_detId] = vvnoise;
323  _mapDetIdEff[current_detId] = vveff;
324 
325 #ifdef EDM_ML_DEBUG
326  sslognoiseitem << " fill vectors into map" << std::endl;
327  std::stringstream LogDebugNoiVectString, LogDebugEffVectString;
328  LogDebugNoiVectString << "[";
329  for (std::vector<float>::iterator itNoiVect = vvnoise.begin(); itNoiVect != vvnoise.end(); ++itNoiVect) {
330  LogDebugNoiVectString << (*itNoiVect) << ",";
331  }
332  LogDebugNoiVectString << "]";
333  std::string LogDebugNoiVectStr = LogDebugNoiVectString.str();
334  LogDebugEffVectString << "[";
335  for (std::vector<float>::iterator itEffVect = vveff.begin(); itEffVect != vveff.end(); ++itEffVect) {
336  LogDebugEffVectString << (*itEffVect) << ",";
337  }
338  LogDebugEffVectString << "]";
339  std::string LogDebugEffVectStr = LogDebugEffVectString.str();
340  LogDebug("RPCSimSetup") << "Filling vvnoise in _mapDetIdNoise[detId] :: detId = " << RPCDetId(it->dpid) << " = "
341  << (RPCDetId(it->dpid)).rawId() << " vvnoise = " << LogDebugNoiVectStr;
342  LogDebug("RPCSimSetup") << "Filling veff in _mapDetIdEff[detId] :: detId = " << RPCDetId(it->dpid) << " = "
343  << (RPCDetId(it->dpid)).rawId() << " veff = " << LogDebugEffVectStr;
344 #endif
345  // look for next different detId and rename it to the current_detId
346  // at this point we skip all the conditions for the strips that are not in this roll
347  // and we will go to the conditions for the first strip of the next roll
348  bool next_detId_found = false;
349 #ifdef EDM_ML_DEBUG
350  sslognoiseitem << "look for next different detId" << std::endl;
351 #endif
352  while (next_detId_found == 0 && it != vnoise.end() - 1) {
353  ++it;
354  this_detId = it->dpid;
355  this_rpcId = RPCDetId(this_detId);
356  this_roll = theGeometry->roll(this_rpcId);
357  if (!this_roll)
358  continue;
359 #ifdef EDM_ML_DEBUG
360  sslognoiseitem << "Inside While:: [" << std::setw(6) << count_all << "][" << std::setw(3) << count_strips
361  << "] :: this_detId = " << this_detId << " aka " << this_rpcId << " Noise = " << it->noise
362  << " Hz/cm2" << std::endl;
363 #endif
364  ++count_strips;
365  // ++count_all;
366  if (this_detId != current_detId) {
367 #ifdef EDM_ML_DEBUG
368  sslognoiseitem << "Different detId is found :: " << this_detId << " aka " << this_rpcId
369  << " Noise = " << it->noise << " Hz/cm2";
370 #endif
371  // next roll is found. update current_detId to this newly found detId
372  // and update also the number of strips
373  current_detId = this_detId;
374  current_rpcId = RPCDetId(current_detId);
375  next_detId_found = true;
376  current_nStrips = (theGeometry->roll(current_rpcId))->nstrips();
377 #ifdef EDM_ML_DEBUG
378  sslognoiseitem << " with " << current_nStrips << " strips" << std::endl;
379 #endif
380  --it; // subtract one, because at the end of the loop the iterator will be increased with one
381  // in fact the treatment for roll N stops when we find the first occurence of roll N+1
382  // however we want to start the treatment for roll N+1 with the first occurence of roll N+1
383  // so the first entry of each new roll N+1 is manipulated twice in the loop (once as a stop, once as a start)
384  // therefore we have to manipulate the iterator here, subtracting one, to treat again this entry
385  }
386  }
387  // reset count_strips
388  count_strips = 1;
389  }
390  // There should be no Case 4
391  // -------------------------
392  else {
393  }
394  }
395  // ###########################################################################
396 #ifdef EDM_ML_DEBUG
397  std::string lognoiseitem = sslognoiseitem.str();
398  sslognoiseitem.clear();
399  LogDebug("RPCSimSetupNoiseLoopDetails") << lognoiseitem << std::endl;
400  LogDebug("RPCSimSetup") << "RPCSimSetUp::setRPCSetUp :: NoiseItem :: end" << std::endl;
401 
402  LogDebug("RPCSimSetup") << "RPCSimSetUp::setRPCSetUp :: end" << std::endl;
403 #endif
404 }
405 
406 const std::vector<float>& RPCSimSetUp::getNoise(uint32_t id) {
407  map<uint32_t, std::vector<float> >::iterator iter = _mapDetIdNoise.find(id);
408  if (iter == _mapDetIdNoise.end()) {
409  throw cms::Exception("DataCorrupt") << "Exception from RPCSimSetUp - no noise information for DetId\t" << id
410  << std::endl;
411  }
412  LogDebug("RPCSimSetupChecks") << "All OK from RPCSimSetUp - noise information for DetId\t" << id << std::endl;
413  return iter->second;
414 }
415 
416 const std::vector<float>& RPCSimSetUp::getEff(uint32_t id) {
417  map<uint32_t, std::vector<float> >::iterator iter = _mapDetIdEff.find(id);
418  if (iter == _mapDetIdEff.end()) {
419  throw cms::Exception("DataCorrupt") << "Exception from RPCSimSetUp - no efficiency information for DetId\t" << id
420  << std::endl;
421  }
422 
423  RPCDetId rpcId = RPCDetId(id);
424  const RPCRoll* roll = theGeometry->roll(rpcId);
425  unsigned int numbStrips = roll->nstrips();
426 
427  if ((iter->second).size() < numbStrips) {
428  LogDebug("RPCSimSetup") << "Exception from RPCSimSetUp - efficiency information in a wrong format for DetId\t" << id
429  << " aka " << RPCDetId(id) << std::endl;
430  LogDebug("RPCSimSetup") << " number of strips in Conditions\t" << (iter->second).size()
431  << " number of strips in Geometry\t" << numbStrips << std::endl;
432  throw cms::Exception("DataCorrupt")
433  << "Exception from RPCSimSetUp - efficiency information in a wrong format for DetId\t" << id << std::endl;
434  }
435 
436  return iter->second;
437 }
438 
439 float RPCSimSetUp::getTime(uint32_t id) {
440  RPCDetId rpcid(id);
441  std::map<RPCDetId, float>::iterator iter = _bxmap.find(rpcid);
442  if (iter == _bxmap.end()) {
443  throw cms::Exception("DataCorrupt") << "Exception from RPCSimSetUp - no timing information for rpcid.rawId()\t"
444  << rpcid.rawId() << std::endl;
445  }
446  return iter->second;
447 }
448 
449 const std::map<int, std::vector<double> >& RPCSimSetUp::getClsMap() {
450  if (_clsMap.size() != 5) {
451  throw cms::Exception("DataCorrupt") << "Exception from RPCSimSetUp - cluster size - a wrong format " << std::endl;
452  }
453  return _clsMap;
454 }
455 
456 //const std::map<int, std::vector<double> >& RPCSimSetUp::getClsMap(uint32_t id)
457 const std::vector<double>& RPCSimSetUp::getCls(uint32_t id) //legacy member function
458 {
459  LogDebug("RPCSimSetupChecks") << "RPCSimSetUp::getCls" << std::endl;
460 
461  map<uint32_t, std::vector<double> >::iterator iter = _mapDetClsMapLegacy.find(id);
462  if (iter == _mapDetClsMapLegacy.end()) {
463  throw cms::Exception("DataCorrupt") << "Exception from RPCSimSetUp - no cluster size information for DetId\t" << id
464  << std::endl;
465  }
466  if ((iter->second).size() != 100) {
467  throw cms::Exception("DataCorrupt")
468  << "Exception from RPCSimSetUp - _mapDetClsMapLegacy - cluster size information in a wrong format for DetId\t"
469  << id << std::endl;
470  }
471  LogDebug("RPCSimSetupChecks")
472  << "All OK from RPCSimSetUp - _mapDetClsMapLegacy - cluster size information for DetId\t" << id << std::endl;
473  return iter->second;
474 }
475 
476 const std::vector<double>& RPCSimSetUp::getAsymmetricClsDistribution(uint32_t id, uint32_t slice) {
477  LogDebug("RPCSimSetupChecks") << "RPCSimSetUp::getAsymmetricClsDistribution" << std::endl;
478 
479  map<uint32_t, std::vector<double> >::const_iterator iter = _mapDetClsMap.find(id);
480  if (iter == _mapDetClsMap.end()) {
481  throw cms::Exception("DataCorrupt")
482  << "Exception from RPCSimSetUp - _mapDetClsMap - no cluster size information for DetId\t" << id << std::endl;
483  }
484  if ((iter->second).size() != 120) {
485  throw cms::Exception("DataCorrupt")
486  << "Exception from RPCSimSetUp - _mapDetClsMap - cluster size information in a wrong format for DetId\t" << id
487  << std::endl;
488  }
489  // return iter->second;
490 
491  std::vector<double> dataForAsymmCls = iter->second;
492  if (slice > 4) {
493  throw cms::Exception("DataCorrupt") << "Exception from RPCSimSetUp - slice variable not in the range" << std::endl;
494  }
495 
496  _DetClsAsymmetric.clear();
497 
498  vector<double> clsFewStripsDistribution;
499  vector<double> clsDistribution;
500  vector<double> clsAccumulativeDistribution;
501 
502  std::map<int, std::vector<double> > mapSliceVsDistribution;
503 
504  const int slices = 5;
505  const int distributionFewStrips = 24;
506 
507  double sliceVsFewStripsDistribution[slices][distributionFewStrips];
508 
509  for (int j = 0; j < distributionFewStrips; j++) {
510  for (int i = 0; i < slices; i++) {
511  sliceVsFewStripsDistribution[i][j] = dataForAsymmCls[j * slices + i];
512  }
513  }
514 
515  int i = slice;
516  double sum = 0;
517  int counter = 0;
518  for (int j = 0; j < distributionFewStrips; j++) {
519  counter++;
520  sum += sliceVsFewStripsDistribution[i][j];
521  if (counter % 4 == 0) {
522  _DetClsAsymmetric.push_back(sum);
523  }
524  }
525  return _DetClsAsymmetric;
526 }
527 
528 const std::vector<double>& RPCSimSetUp::getAsymmetryForCls(uint32_t id, uint32_t slice, uint32_t cls) {
529  LogDebug("RPCSimSetupChecks") << "RPCSimSetUp::getAsymmetryForCls" << std::endl;
530 
531  map<uint32_t, std::vector<double> >::const_iterator iter = _mapDetClsMap.find(id);
532  if (iter == _mapDetClsMap.end()) {
533  throw cms::Exception("DataCorrupt")
534  << "Exception from RPCSimSetUp - _mapDetClsMap - no cluster size information for DetId\t" << id << std::endl;
535  }
536  if ((iter->second).size() != 120) {
537  throw cms::Exception("DataCorrupt")
538  << "Exception from RPCSimSetUp - _mapDetClsMap - cluster size information in a wrong format for DetId\t" << id
539  << '\t' << (iter->second).size() << std::endl;
540  }
541 
542  std::vector<double> dataForAsymmCls = iter->second;
543 
544  if (slice > 4) {
545  throw cms::Exception("DataCorrupt") << "Exception from RPCSimSetUp - slice variable not in the range" << std::endl;
546  }
547 
548  _DetAsymmetryForCls.clear();
549 
550  vector<double> clsFewStripsDistribution;
551  vector<double> clsDistribution;
552  vector<double> clsAccumulativeDistribution;
553  vector<double> clsDetAsymmetryForCls;
554  clsDetAsymmetryForCls.clear();
555 
556  std::map<int, std::vector<double> > mapSliceVsDistribution;
557 
558  const int slices = 5;
559  const int distributionFewStrips = 24;
560 
561  double sliceVsFewStripsDistribution[slices][distributionFewStrips];
562 
563  for (int j = 0; j < distributionFewStrips; j++) {
564  for (int i = 0; i < slices; i++) {
565  sliceVsFewStripsDistribution[i][j] = dataForAsymmCls[j * slices + i];
566  }
567  }
568 
569  int vector_lenght;
570  switch (cls) {
571  case 1:
572  case 3:
573  case 5:
574  vector_lenght = 3;
575  break;
576  case 2:
577  case 4:
578  vector_lenght = 4;
579  break;
580  case 6:
581  default:
582  vector_lenght = 1;
583  break;
584  }
585 
586  float sum = 0;
587  float value;
588  for (int i = 0; i < vector_lenght; i++) {
589  value = sliceVsFewStripsDistribution[slice][(cls - 1) * 4 + i];
590  clsDetAsymmetryForCls.push_back(value);
591  sum += value;
592  // LogDebug ("RPCSimSetup")<<"value\t"<<value<<std::endl;
593  // LogDebug ("RPCSimSetup")<<"sum\t"<<sum<<std::endl;
594  }
595 
596  float accum = 0;
597  for (int i = clsDetAsymmetryForCls.size() - 1; i > -1; i--) {
598  accum += clsDetAsymmetryForCls[i];
599  _DetAsymmetryForCls.push_back(accum / sum);
600  }
601  return _DetAsymmetryForCls;
602 }
603 
size
Write out results.
float getTime(uint32_t id)
Definition: RPCSimSetUp.cc:439
RPCSimSetUp(const edm::ParameterSet &ps)
Definition: RPCSimSetUp.cc:34
const std::vector< double > & getAsymmetricClsDistribution(uint32_t id, uint32_t slice)
Definition: RPCSimSetUp.cc:476
const std::vector< float > & getEff(uint32_t id)
Definition: RPCSimSetUp.cc:416
int nstrips() const
Definition: RPCRoll.cc:24
const std::vector< float > & getNoise(uint32_t id)
Definition: RPCSimSetUp.cc:406
virtual ~RPCSimSetUp()
Definition: RPCSimSetUp.cc:604
const std::vector< double > & getCls(uint32_t id)
Definition: RPCSimSetUp.cc:457
const std::map< int, std::vector< double > > & getClsMap()
Definition: RPCSimSetUp.cc:449
Definition: value.py:1
const std::vector< double > & getAsymmetryForCls(uint32_t id, uint32_t slice, uint32_t cls)
Definition: RPCSimSetUp.cc:528
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)