CMS 3D CMS Logo

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