Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 11:25:30

0001 #include "DataFormats/Common/interface/Handle.h"
0002 #include "SimMuon/RPCDigitizer/src/RPCDigiProducer.h"
0003 #include "SimMuon/RPCDigitizer/src/RPCDigitizer.h"
0004 #include "Geometry/Records/interface/MuonGeometryRecord.h"
0005 #include "SimDataFormats/CrossingFrame/interface/CrossingFrame.h"
0006 #include "SimDataFormats/CrossingFrame/interface/MixCollection.h"
0007 #include "DataFormats/Common/interface/Handle.h"
0008 #include "Geometry/Records/interface/MuonGeometryRecord.h"
0009 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
0010 #include "SimDataFormats/CrossingFrame/interface/MixCollection.h"
0011 #include "FWCore/ServiceRegistry/interface/Service.h"
0012 #include "FWCore/Utilities/interface/RandomNumberGenerator.h"
0013 #include "FWCore/Framework/interface/Event.h"
0014 #include "FWCore/Framework/interface/EventSetup.h"
0015 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0016 #include "SimDataFormats/TrackingHit/interface/PSimHitContainer.h"
0017 #include "SimMuon/RPCDigitizer/src/RPCSimSetUp.h"
0018 #include "DataFormats/MuonDetId/interface/RPCDetId.h"
0019 
0020 #include <cmath>
0021 #include <cmath>
0022 #include <fstream>
0023 #include <sstream>
0024 #include <iostream>
0025 #include <cstring>
0026 #include <string>
0027 #include <vector>
0028 #include <cstdlib>
0029 #include <utility>
0030 #include <map>
0031 
0032 using namespace std;
0033 
0034 RPCSimSetUp::RPCSimSetUp(const edm::ParameterSet& ps) {
0035   _mapDetIdNoise.clear();
0036   _mapDetIdEff.clear();
0037   _bxmap.clear();
0038   _clsMap.clear();
0039 }
0040 
0041 void RPCSimSetUp::setRPCSetUp(const std::vector<RPCStripNoises::NoiseItem>& vnoise, const std::vector<float>& vcls) {
0042   unsigned int counter = 1;
0043   unsigned int row = 1;
0044   std::vector<double> sum_clsize;
0045 
0046   for (unsigned int n = 0; n < vcls.size(); ++n) {
0047     sum_clsize.push_back(vcls[n]);
0048 
0049     if (counter == row * 20) {
0050       _clsMap[row] = sum_clsize;
0051       row++;
0052       sum_clsize.clear();
0053     }
0054     counter++;
0055   }
0056 
0057   unsigned int n = 0;
0058   uint32_t temp = 0;
0059   std::vector<float> veff, vvnoise;
0060   veff.clear();
0061   vvnoise.clear();
0062 
0063   for (std::vector<RPCStripNoises::NoiseItem>::const_iterator it = vnoise.begin(); it != vnoise.end(); ++it) {
0064     if (n % 96 == 0) {
0065       if (n > 0) {
0066         _mapDetIdNoise[temp] = vvnoise;
0067         _mapDetIdEff[temp] = veff;
0068         _bxmap[RPCDetId(it->dpid)] = it->time;
0069 
0070         veff.clear();
0071         vvnoise.clear();
0072         vvnoise.push_back((it->noise));
0073         veff.push_back((it->eff));
0074       } else if (n == 0) {
0075         vvnoise.push_back((it->noise));
0076         veff.push_back((it->eff));
0077         _bxmap[RPCDetId(it->dpid)] = it->time;
0078       }
0079     } else if (n == vnoise.size() - 1) {
0080       temp = it->dpid;
0081       vvnoise.push_back((it->noise));
0082       veff.push_back((it->eff));
0083       _mapDetIdNoise[temp] = vvnoise;
0084       _mapDetIdEff[temp] = veff;
0085     } else {
0086       temp = it->dpid;
0087       vvnoise.push_back((it->noise));
0088       veff.push_back((it->eff));
0089     }
0090     n++;
0091   }
0092 }
0093 
0094 void RPCSimSetUp::setRPCSetUp(const std::vector<RPCStripNoises::NoiseItem>& vnoise,
0095                               const std::vector<RPCClusterSize::ClusterSizeItem>& vClusterSize) {
0096   LogDebug("RPCSimSetup") << "RPCSimSetUp::setRPCSetUp(vector<NoiseItem>, vector<ClusterSizeItem>)" << std::endl;
0097 
0098   uint32_t detId = 0, current_detId, this_detId;
0099   RPCDetId rpcId, current_rpcId, this_rpcId;
0100   const RPCRoll* current_roll = nullptr;
0101   const RPCRoll* this_roll = nullptr;
0102   unsigned int current_nStrips;
0103 
0104   LogDebug("RPCSimSetup") << "RPCSimSetUp::setRPCSetUp :: ClusterSizeItem :: begin" << std::endl;
0105 #ifdef EDM_ML_DEBUG
0106   std::stringstream sslogclsitem;
0107 #endif
0108   // ### ClusterSizeItem #######################################################
0109   std::vector<RPCClusterSize::ClusterSizeItem>::const_iterator itCls;
0110   int clsCounter(1);
0111   std::vector<double> clsVect;
0112   // ### loop for New Format (120 entries)
0113   for (itCls = vClusterSize.begin(); itCls != vClusterSize.end(); ++itCls) {
0114     clsVect.push_back(((double)(itCls->clusterSize)));
0115 #ifdef EDM_ML_DEBUG
0116     sslogclsitem << " Push back clustersize = " << itCls->clusterSize << std::endl;
0117     sslogclsitem << "Filling cls in _mapDetCls[detId,clsVect] :: detId = " << detId;
0118     sslogclsitem << " --> will it be accepted? clsCounter = " << clsCounter << " accepted?";
0119     sslogclsitem << " New Format ::" << ((!(clsCounter % 120)) && (clsCounter != 0));  // <<std::endl;
0120     sslogclsitem << " Old Format ::" << ((!(clsCounter % 100)) && (clsCounter != 0));  // <<std::endl;
0121     sslogclsitem << std::endl;
0122 #endif
0123     // New Format :: loop until 120
0124     if ((!(clsCounter % 120)) && (clsCounter != 0)) {
0125       detId = itCls->dpid;
0126       _mapDetClsMap[detId] = clsVect;
0127 #ifdef EDM_ML_DEBUG
0128       std::stringstream LogDebugClsVectString;
0129       LogDebugClsVectString << "[";
0130       for (std::vector<double>::iterator itClsVect = clsVect.begin(); itClsVect != clsVect.end(); ++itClsVect) {
0131         LogDebugClsVectString << *itClsVect << ",";
0132       }
0133       LogDebugClsVectString << "]";
0134       std::string LogDebugClsVectStr = LogDebugClsVectString.str();
0135       LogDebug("RPCSimSetup") << "Filling clsVect in _mapDetCls[detId,clsVect] :: detId = " << RPCDetId(detId) << " = "
0136                               << detId << " clsVec = " << LogDebugClsVectStr;
0137 
0138       sslogclsitem << " --> New Method ";
0139       sslogclsitem << " --> saved in map " << std::endl;
0140       sslogclsitem << "Filling cls in _mapDetClsMap[detId,clsVect] :: detId = " << detId;
0141       sslogclsitem << " --> will it be accepted? clsCounter = " << clsCounter << " accepted? "
0142                    << ((!(clsCounter % 120)) && (clsCounter != 0)) << std::endl;
0143 #endif
0144       clsVect.clear();
0145       clsCounter = 0;
0146     } else {
0147 #ifdef EDM_ML_DEBUG
0148       sslogclsitem << " --> not saved in map " << std::endl;
0149 #endif
0150     }
0151     ++clsCounter;
0152   }
0153   // ### loop for Old Format (100 entries)
0154   for (itCls = vClusterSize.begin(); itCls != vClusterSize.end(); ++itCls) {
0155     clsVect.push_back(((double)(itCls->clusterSize)));
0156 #ifdef EDM_ML_DEBUG
0157     sslogclsitem << " Push back clustersize = " << itCls->clusterSize << std::endl;
0158     sslogclsitem << "Filling cls in _mapDetClsMapLegacy[detId,clsVect] :: detId = " << detId;
0159     sslogclsitem << " --> will it be accepted? clsCounter = " << clsCounter << " accepted?";
0160     sslogclsitem << " New Format ::" << ((!(clsCounter % 120)) && (clsCounter != 0));  // <<std::endl;
0161     sslogclsitem << " Old Format ::" << ((!(clsCounter % 100)) && (clsCounter != 0));  // <<std::endl;
0162     sslogclsitem << std::endl;
0163 #endif
0164     // Old Format :: same until 100
0165     if ((!(clsCounter % 100)) && (clsCounter != 0)) {
0166       detId = itCls->dpid;
0167       _mapDetClsMapLegacy[detId] = clsVect;
0168 #ifdef EDM_ML_DEBUG
0169       std::stringstream LogDebugClsVectString;
0170       LogDebugClsVectString << "[";
0171       for (std::vector<double>::iterator itClsVect = clsVect.begin(); itClsVect != clsVect.end(); ++itClsVect) {
0172         LogDebugClsVectString << *itClsVect << ",";
0173       }
0174       LogDebugClsVectString << "]";
0175       std::string LogDebugClsVectStr = LogDebugClsVectString.str();
0176       LogDebug("RPCSimSetup") << "Filling clsVect in _mapDetClsLegacy[detId,clsVect] :: detId = " << RPCDetId(detId)
0177                               << " = " << detId << " clsVec = " << LogDebugClsVectStr;
0178 
0179       sslogclsitem << " --> Old Method ";
0180       sslogclsitem << " --> saved in map " << std::endl;
0181       sslogclsitem << "Filling cls in _mapDetClsMapLegacy[detId,clsVect] :: detId = " << detId;
0182       sslogclsitem << " --> will it be accepted? clsCounter = " << clsCounter << " accepted? "
0183                    << ((!(clsCounter % 120)) && (clsCounter != 0)) << std::endl;
0184 #endif
0185       clsVect.clear();
0186       clsCounter = 0;
0187     } else {
0188 #ifdef EDM_ML_DEBUG
0189       sslogclsitem << " --> not saved in map " << std::endl;
0190 #endif
0191     }
0192     ++clsCounter;
0193   }
0194   // ###########################################################################
0195 #ifdef EDM_ML_DEBUG
0196   std::string logclsitem = sslogclsitem.str();
0197   sslogclsitem.clear();
0198   LogDebug("RPCSimSetupClsLoopDetails") << logclsitem << std::endl;
0199   LogDebug("RPCSimSetup") << "RPCSimSetUp::setRPCSetUp :: ClusterSizeItem :: end" << std::endl;
0200 
0201   LogDebug("RPCSimSetup") << "RPCSimSetUp::setRPCSetUp :: NoiseItem :: begin" << std::endl;
0202   std::stringstream sslognoiseitem;
0203 #endif
0204   // ### NoiseItem #############################################################
0205   unsigned int count_strips = 1;
0206   unsigned int count_all = 1;
0207   std::vector<float> vveff, vvnoise;
0208 
0209   // DetId to start with needs to be a DetId inside the Geometry used
0210   // Therefore loop on the NoiseItems and search for the first valid roll in the Geometry
0211   // Assign this as the DetId to start with (so called current_roll) and quit the loop
0212   bool quitLoop = false;
0213   current_detId = 0;
0214   current_nStrips = 0;  // current_rpcId = 0; current_roll = 0;
0215   for (std::vector<RPCStripNoises::NoiseItem>::const_iterator it = vnoise.begin(); it != vnoise.end() && !quitLoop;
0216        ++it) {
0217     // roll associated to the conditions of this strip (iterator)
0218     current_detId = it->dpid;
0219     current_rpcId = RPCDetId(current_detId);
0220     // Test whether this roll (picked up from the conditions) is inside the RPC Geometry
0221     const RPCRoll* roll = theGeometry->roll(current_rpcId);
0222     if (roll == nullptr) {
0223 #ifdef EDM_ML_DEBUG
0224       sslognoiseitem << "Searching for first valid detid :: current_detId = " << current_detId;
0225       sslognoiseitem << " aka " << current_rpcId << " is not in current Geometry --> Skip " << std::endl;
0226 #endif
0227       continue;
0228     } else {
0229 #ifdef EDM_ML_DEBUG
0230       sslognoiseitem << "Searching for first valid detid :: current_detId = " << current_detId;
0231       sslognoiseitem << " aka " << current_rpcId
0232                      << " is the first (valid) roll in the current Geometry --> Accept, Assign & Quit Loop"
0233                      << std::endl;
0234 #endif
0235       current_roll = theGeometry->roll(current_rpcId);
0236       current_nStrips = current_roll->nstrips();
0237       quitLoop = true;
0238     }
0239   }
0240 
0241 #ifdef EDM_ML_DEBUG
0242   sslognoiseitem << "Start Position ::            current_detId = " << current_detId << " aka " << current_rpcId;
0243   sslognoiseitem << " is a valid roll with pointer " << current_roll << " and has "
0244                  << (current_roll ? current_roll->nstrips() : 0) << " strips" << std::endl;
0245   sslognoiseitem << " -------------------------------------------------------------------------------------------------"
0246                     "------------------------------------ "
0247                  << std::endl;
0248 #endif
0249   for (std::vector<RPCStripNoises::NoiseItem>::const_iterator it = vnoise.begin(); it != vnoise.end(); ++it) {
0250     // roll associated to the conditions of this strip (iterator)
0251     this_detId = it->dpid;
0252     this_rpcId = RPCDetId(this_detId);
0253     // Test whether this roll (picked up from the conditions) is inside the RPC Geometry
0254     const RPCRoll* roll = theGeometry->roll(this_rpcId);
0255     if (roll == nullptr) {
0256 #ifdef EDM_ML_DEBUG
0257       sslognoiseitem << "Inside Loop :: [" << std::setw(6) << count_all << "][" << std::setw(3) << count_strips
0258                      << "] :: this_detId = " << this_detId << " aka " << this_rpcId
0259                      << " which is not in current Geometry --> Skip " << std::endl;
0260 #endif
0261       continue;
0262     }
0263 
0264     // Case 1 :: FIRST ENTRY
0265     // ---------------------
0266     if (this_detId == current_detId && count_strips == 1) {
0267       // fill bx in map
0268       _bxmap[current_detId] = it->time;
0269       // clear vectors
0270       vveff.clear();
0271       vvnoise.clear();
0272       // fill the vectors
0273       vvnoise.push_back((it->noise));
0274       vveff.push_back((it->eff));
0275 #ifdef EDM_ML_DEBUG
0276       sslognoiseitem << "RPCSimSetUp::setRPCSetUp :: NoiseItem :: case 1" << std::endl;
0277       sslognoiseitem << this_detId << " = " << this_rpcId << " with " << roll->nstrips() << " strips" << std::endl;
0278       sslognoiseitem << "[NoiseItem :: n = " << count_all
0279                      << "] Filling time in _bxmap[detId] :: detId = " << RPCDetId(it->dpid) << " time = " << it->time
0280                      << std::endl;
0281       sslognoiseitem << "First Value :: [" << std::setw(6) << count_all << "][" << std::setw(3) << count_strips
0282                      << "] :: this_detId = " << this_detId << " aka " << this_rpcId;
0283       sslognoiseitem << " Strip " << std::setw(3) << count_strips << " Noise = " << it->noise << " Hz/cm2" << std::endl;
0284 #endif
0285       // update counter
0286       ++count_strips;
0287       ++count_all;
0288     }
0289     // Case 2 :: 2ND ENTRY --> LAST-1 ENTRY
0290     // ------------------------------------
0291     else if (this_detId == current_detId && count_strips > 1 && count_strips < current_nStrips) {
0292 #ifdef EDM_ML_DEBUG
0293       sslognoiseitem << "RPCSimSetUp::setRPCSetUp :: NoiseItem :: case 2" << std::endl;
0294       sslognoiseitem << "Inside Loop :: [" << std::setw(6) << count_all << "][" << std::setw(3) << count_strips
0295                      << "] :: this_detId = " << this_detId << " aka " << this_rpcId;
0296       sslognoiseitem << " Strip " << std::setw(3) << count_strips << " Noise = " << it->noise << " Hz/cm2" << std::endl;
0297 #endif
0298       // fill the vectors
0299       vvnoise.push_back((it->noise));
0300       vveff.push_back((it->eff));
0301       // update counter
0302       ++count_strips;
0303       ++count_all;
0304     }
0305 
0306     // Case 3 :: LAST ENTRY
0307     // --------------------
0308     else if (this_detId == current_detId && count_strips == current_nStrips) {
0309 #ifdef EDM_ML_DEBUG
0310       sslognoiseitem << "RPCSimSetUp::setRPCSetUp :: NoiseItem :: case 3" << std::endl;
0311       sslognoiseitem << "Last Value ::  [" << std::setw(6) << count_all << "][" << std::setw(3) << count_strips
0312                      << "] :: this_detId = " << this_detId << " aka " << this_rpcId;
0313       sslognoiseitem << " Strip " << std::setw(3) << count_strips << " Noise = " << it->noise << " Hz/cm2" << std::endl;
0314 #endif
0315       // fill last value in the vector
0316       vvnoise.push_back((it->noise));
0317       vveff.push_back((it->eff));
0318       // update counter
0319       ++count_strips;
0320       ++count_all;
0321       // fill vectors into map
0322       _mapDetIdNoise[current_detId] = vvnoise;
0323       _mapDetIdEff[current_detId] = vveff;
0324 
0325 #ifdef EDM_ML_DEBUG
0326       sslognoiseitem << " fill vectors into map" << std::endl;
0327       std::stringstream LogDebugNoiVectString, LogDebugEffVectString;
0328       LogDebugNoiVectString << "[";
0329       for (std::vector<float>::iterator itNoiVect = vvnoise.begin(); itNoiVect != vvnoise.end(); ++itNoiVect) {
0330         LogDebugNoiVectString << (*itNoiVect) << ",";
0331       }
0332       LogDebugNoiVectString << "]";
0333       std::string LogDebugNoiVectStr = LogDebugNoiVectString.str();
0334       LogDebugEffVectString << "[";
0335       for (std::vector<float>::iterator itEffVect = vveff.begin(); itEffVect != vveff.end(); ++itEffVect) {
0336         LogDebugEffVectString << (*itEffVect) << ",";
0337       }
0338       LogDebugEffVectString << "]";
0339       std::string LogDebugEffVectStr = LogDebugEffVectString.str();
0340       LogDebug("RPCSimSetup") << "Filling vvnoise in _mapDetIdNoise[detId] :: detId = " << RPCDetId(it->dpid) << " = "
0341                               << (RPCDetId(it->dpid)).rawId() << " vvnoise = " << LogDebugNoiVectStr;
0342       LogDebug("RPCSimSetup") << "Filling veff    in _mapDetIdEff[detId]   :: detId = " << RPCDetId(it->dpid) << " = "
0343                               << (RPCDetId(it->dpid)).rawId() << " veff    = " << LogDebugEffVectStr;
0344 #endif
0345       // look for next different detId and rename it to the current_detId
0346       // at this point we skip all the conditions for the strips that are not in this roll
0347       // and we will go to the conditions for the first strip of the next roll
0348       bool next_detId_found = false;
0349 #ifdef EDM_ML_DEBUG
0350       sslognoiseitem << "look for next different detId" << std::endl;
0351 #endif
0352       while (next_detId_found == 0 && it != vnoise.end() - 1) {
0353         ++it;
0354         this_detId = it->dpid;
0355         this_rpcId = RPCDetId(this_detId);
0356         this_roll = theGeometry->roll(this_rpcId);
0357         if (!this_roll)
0358           continue;
0359 #ifdef EDM_ML_DEBUG
0360         sslognoiseitem << "Inside While:: [" << std::setw(6) << count_all << "][" << std::setw(3) << count_strips
0361                        << "] :: this_detId = " << this_detId << " aka " << this_rpcId << " Noise = " << it->noise
0362                        << " Hz/cm2" << std::endl;
0363 #endif
0364         ++count_strips;
0365         // ++count_all;
0366         if (this_detId != current_detId) {
0367 #ifdef EDM_ML_DEBUG
0368           sslognoiseitem << "Different detId is found ::                  " << this_detId << " aka " << this_rpcId
0369                          << " Noise = " << it->noise << " Hz/cm2";
0370 #endif
0371           // next roll is found. update current_detId to this newly found detId
0372           // and update also the number of strips
0373           current_detId = this_detId;
0374           current_rpcId = RPCDetId(current_detId);
0375           next_detId_found = true;
0376           current_nStrips = (theGeometry->roll(current_rpcId))->nstrips();
0377 #ifdef EDM_ML_DEBUG
0378           sslognoiseitem << " with " << current_nStrips << " strips" << std::endl;
0379 #endif
0380           --it;  // subtract one, because at the end of the loop the iterator will be increased with one
0381                  // in fact the treatment for roll N stops when we find the first occurence of roll N+1
0382                  // however we want to start the treatment for roll N+1 with the first occurence of roll N+1
0383           // so the first entry of each new roll N+1 is manipulated twice in the loop (once as a stop, once as a start)
0384           // therefore we have to manipulate the iterator here, subtracting one, to treat again this entry
0385         }
0386       }
0387       // reset count_strips
0388       count_strips = 1;
0389     }
0390     // There should be no Case 4
0391     // -------------------------
0392     else {
0393     }
0394   }
0395   // ###########################################################################
0396 #ifdef EDM_ML_DEBUG
0397   std::string lognoiseitem = sslognoiseitem.str();
0398   sslognoiseitem.clear();
0399   LogDebug("RPCSimSetupNoiseLoopDetails") << lognoiseitem << std::endl;
0400   LogDebug("RPCSimSetup") << "RPCSimSetUp::setRPCSetUp :: NoiseItem :: end" << std::endl;
0401 
0402   LogDebug("RPCSimSetup") << "RPCSimSetUp::setRPCSetUp :: end" << std::endl;
0403 #endif
0404 }
0405 
0406 const std::vector<float>& RPCSimSetUp::getNoise(uint32_t id) {
0407   map<uint32_t, std::vector<float> >::iterator iter = _mapDetIdNoise.find(id);
0408   if (iter == _mapDetIdNoise.end()) {
0409     throw cms::Exception("DataCorrupt") << "Exception from RPCSimSetUp - no noise information for DetId\t" << id
0410                                         << std::endl;
0411   }
0412   LogDebug("RPCSimSetupChecks") << "All OK from RPCSimSetUp - noise information for DetId\t" << id << std::endl;
0413   return iter->second;
0414 }
0415 
0416 const std::vector<float>& RPCSimSetUp::getEff(uint32_t id) {
0417   map<uint32_t, std::vector<float> >::iterator iter = _mapDetIdEff.find(id);
0418   if (iter == _mapDetIdEff.end()) {
0419     throw cms::Exception("DataCorrupt") << "Exception from RPCSimSetUp - no efficiency information for DetId\t" << id
0420                                         << std::endl;
0421   }
0422 
0423   RPCDetId rpcId = RPCDetId(id);
0424   const RPCRoll* roll = theGeometry->roll(rpcId);
0425   unsigned int numbStrips = roll->nstrips();
0426 
0427   if ((iter->second).size() < numbStrips) {
0428     LogDebug("RPCSimSetup") << "Exception from RPCSimSetUp - efficiency information in a wrong format for DetId\t" << id
0429                             << " aka " << RPCDetId(id) << std::endl;
0430     LogDebug("RPCSimSetup") << " number of strips in Conditions\t" << (iter->second).size()
0431                             << " number of strips in Geometry\t" << numbStrips << std::endl;
0432     throw cms::Exception("DataCorrupt")
0433         << "Exception from RPCSimSetUp - efficiency information in a wrong format for DetId\t" << id << std::endl;
0434   }
0435 
0436   return iter->second;
0437 }
0438 
0439 float RPCSimSetUp::getTime(uint32_t id) {
0440   RPCDetId rpcid(id);
0441   std::map<RPCDetId, float>::iterator iter = _bxmap.find(rpcid);
0442   if (iter == _bxmap.end()) {
0443     throw cms::Exception("DataCorrupt") << "Exception from RPCSimSetUp - no timing information for rpcid.rawId()\t"
0444                                         << rpcid.rawId() << std::endl;
0445   }
0446   return iter->second;
0447 }
0448 
0449 const std::map<int, std::vector<double> >& RPCSimSetUp::getClsMap() {
0450   if (_clsMap.size() != 5) {
0451     throw cms::Exception("DataCorrupt") << "Exception from RPCSimSetUp - cluster size - a wrong format " << std::endl;
0452   }
0453   return _clsMap;
0454 }
0455 
0456 //const std::map<int, std::vector<double> >& RPCSimSetUp::getClsMap(uint32_t id)
0457 const std::vector<double>& RPCSimSetUp::getCls(uint32_t id)  //legacy member function
0458 {
0459   LogDebug("RPCSimSetupChecks") << "RPCSimSetUp::getCls" << std::endl;
0460 
0461   map<uint32_t, std::vector<double> >::iterator iter = _mapDetClsMapLegacy.find(id);
0462   if (iter == _mapDetClsMapLegacy.end()) {
0463     throw cms::Exception("DataCorrupt") << "Exception from RPCSimSetUp - no cluster size information for DetId\t" << id
0464                                         << std::endl;
0465   }
0466   if ((iter->second).size() != 100) {
0467     throw cms::Exception("DataCorrupt")
0468         << "Exception from RPCSimSetUp - _mapDetClsMapLegacy - cluster size information in a wrong format for DetId\t"
0469         << id << std::endl;
0470   }
0471   LogDebug("RPCSimSetupChecks")
0472       << "All OK from RPCSimSetUp - _mapDetClsMapLegacy - cluster size information for DetId\t" << id << std::endl;
0473   return iter->second;
0474 }
0475 
0476 const std::vector<double>& RPCSimSetUp::getAsymmetricClsDistribution(uint32_t id, uint32_t slice) {
0477   LogDebug("RPCSimSetupChecks") << "RPCSimSetUp::getAsymmetricClsDistribution" << std::endl;
0478 
0479   map<uint32_t, std::vector<double> >::const_iterator iter = _mapDetClsMap.find(id);
0480   if (iter == _mapDetClsMap.end()) {
0481     throw cms::Exception("DataCorrupt")
0482         << "Exception from RPCSimSetUp - _mapDetClsMap - no cluster size information for DetId\t" << id << std::endl;
0483   }
0484   if ((iter->second).size() != 120) {
0485     throw cms::Exception("DataCorrupt")
0486         << "Exception from RPCSimSetUp - _mapDetClsMap - cluster size information in a wrong format for DetId\t" << id
0487         << std::endl;
0488   }
0489   //  return iter->second;
0490 
0491   std::vector<double> dataForAsymmCls = iter->second;
0492   if (slice > 4) {
0493     throw cms::Exception("DataCorrupt") << "Exception from RPCSimSetUp - slice variable not in the range" << std::endl;
0494   }
0495 
0496   _DetClsAsymmetric.clear();
0497 
0498   vector<double> clsFewStripsDistribution;
0499   vector<double> clsDistribution;
0500   vector<double> clsAccumulativeDistribution;
0501 
0502   std::map<int, std::vector<double> > mapSliceVsDistribution;
0503 
0504   const int slices = 5;
0505   const int distributionFewStrips = 24;
0506 
0507   double sliceVsFewStripsDistribution[slices][distributionFewStrips];
0508 
0509   for (int j = 0; j < distributionFewStrips; j++) {
0510     for (int i = 0; i < slices; i++) {
0511       sliceVsFewStripsDistribution[i][j] = dataForAsymmCls[j * slices + i];
0512     }
0513   }
0514 
0515   int i = slice;
0516   double sum = 0;
0517   int counter = 0;
0518   for (int j = 0; j < distributionFewStrips; j++) {
0519     counter++;
0520     sum += sliceVsFewStripsDistribution[i][j];
0521     if (counter % 4 == 0) {
0522       _DetClsAsymmetric.push_back(sum);
0523     }
0524   }
0525   return _DetClsAsymmetric;
0526 }
0527 
0528 const std::vector<double>& RPCSimSetUp::getAsymmetryForCls(uint32_t id, uint32_t slice, uint32_t cls) {
0529   LogDebug("RPCSimSetupChecks") << "RPCSimSetUp::getAsymmetryForCls" << std::endl;
0530 
0531   map<uint32_t, std::vector<double> >::const_iterator iter = _mapDetClsMap.find(id);
0532   if (iter == _mapDetClsMap.end()) {
0533     throw cms::Exception("DataCorrupt")
0534         << "Exception from RPCSimSetUp - _mapDetClsMap - no cluster size information for DetId\t" << id << std::endl;
0535   }
0536   if ((iter->second).size() != 120) {
0537     throw cms::Exception("DataCorrupt")
0538         << "Exception from RPCSimSetUp - _mapDetClsMap - cluster size information in a wrong format for DetId\t" << id
0539         << '\t' << (iter->second).size() << std::endl;
0540   }
0541 
0542   std::vector<double> dataForAsymmCls = iter->second;
0543 
0544   if (slice > 4) {
0545     throw cms::Exception("DataCorrupt") << "Exception from RPCSimSetUp - slice variable not in the range" << std::endl;
0546   }
0547 
0548   _DetAsymmetryForCls.clear();
0549 
0550   vector<double> clsFewStripsDistribution;
0551   vector<double> clsDistribution;
0552   vector<double> clsAccumulativeDistribution;
0553   vector<double> clsDetAsymmetryForCls;
0554   clsDetAsymmetryForCls.clear();
0555 
0556   std::map<int, std::vector<double> > mapSliceVsDistribution;
0557 
0558   const int slices = 5;
0559   const int distributionFewStrips = 24;
0560 
0561   double sliceVsFewStripsDistribution[slices][distributionFewStrips];
0562 
0563   for (int j = 0; j < distributionFewStrips; j++) {
0564     for (int i = 0; i < slices; i++) {
0565       sliceVsFewStripsDistribution[i][j] = dataForAsymmCls[j * slices + i];
0566     }
0567   }
0568 
0569   int vector_lenght;
0570   switch (cls) {
0571     case 1:
0572     case 3:
0573     case 5:
0574       vector_lenght = 3;
0575       break;
0576     case 2:
0577     case 4:
0578       vector_lenght = 4;
0579       break;
0580     case 6:
0581     default:
0582       vector_lenght = 1;
0583       break;
0584   }
0585 
0586   float sum = 0;
0587   float value;
0588   for (int i = 0; i < vector_lenght; i++) {
0589     value = sliceVsFewStripsDistribution[slice][(cls - 1) * 4 + i];
0590     clsDetAsymmetryForCls.push_back(value);
0591     sum += value;
0592     //     LogDebug ("RPCSimSetup")<<"value\t"<<value<<std::endl;
0593     //    LogDebug ("RPCSimSetup")<<"sum\t"<<sum<<std::endl;
0594   }
0595 
0596   float accum = 0;
0597   for (int i = clsDetAsymmetryForCls.size() - 1; i > -1; i--) {
0598     accum += clsDetAsymmetryForCls[i];
0599     _DetAsymmetryForCls.push_back(accum / sum);
0600   }
0601   return _DetAsymmetryForCls;
0602 }
0603 
0604 RPCSimSetUp::~RPCSimSetUp() {}