CMS 3D CMS Logo

/data/git/CMSSW_5_3_11_patch5/src/SimMuon/RPCDigitizer/src/RPCSimSetUp.cc

Go to the documentation of this file.
00001 #include "FWCore/Framework/interface/EDProducer.h"
00002 #include "DataFormats/Common/interface/Handle.h"
00003 #include "FWCore/Framework/interface/ESHandle.h"
00004 #include "SimMuon/RPCDigitizer/src/RPCDigiProducer.h"
00005 #include "SimMuon/RPCDigitizer/src/RPCDigitizer.h"
00006 #include "Geometry/Records/interface/MuonGeometryRecord.h"
00007 #include "SimDataFormats/CrossingFrame/interface/CrossingFrame.h"
00008 #include "SimDataFormats/CrossingFrame/interface/MixCollection.h"
00009 #include "DataFormats/Common/interface/Handle.h"
00010 #include "FWCore/Framework/interface/ESHandle.h"
00011 #include "Geometry/Records/interface/MuonGeometryRecord.h"
00012 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
00013 #include "SimDataFormats/CrossingFrame/interface/MixCollection.h"
00014 #include "FWCore/ServiceRegistry/interface/Service.h"
00015 #include "FWCore/Utilities/interface/RandomNumberGenerator.h"
00016 #include "FWCore/Framework/interface/Event.h"
00017 #include "FWCore/Framework/interface/EventSetup.h"
00018 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00019 #include "SimDataFormats/TrackingHit/interface/PSimHitContainer.h"
00020 #include "SimMuon/RPCDigitizer/src/RPCSimSetUp.h"
00021 #include "DataFormats/MuonDetId/interface/RPCDetId.h"
00022 
00023 #include <cmath>
00024 #include <math.h>
00025 #include <fstream>
00026 #include <sstream>
00027 #include <iostream>
00028 #include<cstring>
00029 #include<string>
00030 #include<vector>
00031 #include<stdlib.h>
00032 #include <utility>
00033 #include <map>
00034 
00035 using namespace std;
00036 
00037 RPCSimSetUp::RPCSimSetUp(const edm::ParameterSet& ps) {
00038 
00039   _mapDetIdNoise.clear();
00040   _mapDetIdEff.clear();
00041   _bxmap.clear();
00042   _clsMap.clear();
00043 
00044 }
00045 
00046 void RPCSimSetUp::setRPCSetUp(std::vector<RPCStripNoises::NoiseItem> vnoise, std::vector<float> vcls){
00047 
00048   unsigned int counter = 1;
00049   unsigned int row = 1;
00050   std::vector<double> sum_clsize;
00051 
00052   for(unsigned int n = 0; n < vcls.size(); ++n){
00053 
00054     sum_clsize.push_back(vcls[n]);
00055 
00056     if(counter == row*20) {
00057 
00058       _clsMap[row] = sum_clsize;
00059       row++;
00060       sum_clsize.clear();
00061     }
00062     counter++;
00063   }
00064 
00065   unsigned int n = 0; 
00066   uint32_t temp = 0; 
00067   std::vector<float> veff, vvnoise;
00068   veff.clear();
00069   vvnoise.clear();
00070 
00071   for(std::vector<RPCStripNoises::NoiseItem>::iterator it = vnoise.begin(); it != vnoise.end(); ++it){
00072     if(n%96 == 0) {
00073       if(n > 0 ){
00074         _mapDetIdNoise[temp]= vvnoise;
00075         _mapDetIdEff[temp] = veff;
00076         _bxmap[RPCDetId(it->dpid)] = it->time;
00077         
00078         veff.clear();
00079         vvnoise.clear();
00080         vvnoise.push_back((it->noise));
00081         veff.push_back((it->eff));
00082       }
00083       else if(n == 0 ){
00084         vvnoise.push_back((it->noise));
00085         veff.push_back((it->eff));
00086         _bxmap[RPCDetId(it->dpid)] = it->time;
00087       }
00088     } else if (n == vnoise.size()-1 ){
00089       temp = it->dpid;
00090       vvnoise.push_back((it->noise));
00091       veff.push_back((it->eff));
00092       _mapDetIdNoise[temp]= vvnoise;
00093       _mapDetIdEff[temp] = veff;
00094     } else {
00095       temp = it->dpid;
00096       vvnoise.push_back((it->noise));
00097       veff.push_back((it->eff));
00098     }
00099     n++;
00100   }
00101 }
00102 
00103 void RPCSimSetUp::setRPCSetUp(std::vector<RPCStripNoises::NoiseItem> vnoise, std::vector<RPCClusterSize::ClusterSizeItem> vClusterSize){
00104 
00105   std::vector<RPCClusterSize::ClusterSizeItem>::iterator itCls;
00106   uint32_t  detId;
00107   int clsCounter(1);
00108   std::vector<double> clsVect;
00109 
00110   for(itCls = vClusterSize.begin(); itCls != vClusterSize.end(); ++itCls){
00111     clsVect.push_back(((double)(itCls->clusterSize)));
00112     if((!(clsCounter%100)) && (clsCounter!=0)){
00113       detId=itCls->dpid;
00114       _mapDetClsMap[detId]=clsVect;
00115       clsVect.clear();
00116       clsCounter=0;
00117     }
00118     ++clsCounter;
00119   }
00120 
00121   unsigned int n = 0; 
00122   uint32_t temp = 0; 
00123   std::vector<float> veff, vvnoise;
00124   veff.clear();
00125   vvnoise.clear();
00126 
00127   for(std::vector<RPCStripNoises::NoiseItem>::iterator it = vnoise.begin(); it != vnoise.end(); ++it){
00128     if(n%96 == 0) {
00129       if(n > 0 ){
00130         _mapDetIdNoise[temp]= vvnoise;
00131         _mapDetIdEff[temp] = veff;
00132         _bxmap[RPCDetId(it->dpid)] = it->time;
00133         
00134         veff.clear();
00135         vvnoise.clear();
00136         vvnoise.push_back((it->noise));
00137         veff.push_back((it->eff));
00138       }
00139       else if(n == 0 ){
00140         vvnoise.push_back((it->noise));
00141         veff.push_back((it->eff));
00142         _bxmap[RPCDetId(it->dpid)] = it->time;
00143       }
00144     } else if (n == vnoise.size()-1 ){
00145       temp = it->dpid;
00146       vvnoise.push_back((it->noise));
00147       veff.push_back((it->eff));
00148       _mapDetIdNoise[temp]= vvnoise;
00149       _mapDetIdEff[temp] = veff;
00150     } else {
00151       temp = it->dpid;
00152       vvnoise.push_back((it->noise));
00153       veff.push_back((it->eff));
00154     }
00155     n++;
00156   }
00157 }
00158 
00159 
00160 const std::vector<float>& RPCSimSetUp::getNoise(uint32_t id)
00161 {
00162   map<uint32_t,std::vector<float> >::iterator iter = _mapDetIdNoise.find(id);
00163   if(iter == _mapDetIdNoise.end()){
00164     throw cms::Exception("DataCorrupt") 
00165       << "Exception comming from RPCSimSetUp - no noise information for DetId\t"<<id<< std::endl;
00166   }
00167   return iter->second;
00168 }
00169 
00170 const std::vector<float>& RPCSimSetUp::getEff(uint32_t id)
00171 {
00172   map<uint32_t,std::vector<float> >::iterator iter = _mapDetIdEff.find(id);
00173   if(iter == _mapDetIdEff.end()){
00174     throw cms::Exception("DataCorrupt") 
00175       << "Exception comming from RPCSimSetUp - no efficiency information for DetId\t"<<id<< std::endl;
00176   }
00177   if((iter->second).size() != 96){
00178     throw cms::Exception("DataCorrupt") 
00179       << "Exception comming from RPCSimSetUp - efficiency information in a wrong format for DetId\t"<<id<< std::endl;
00180   }
00181   return iter->second;
00182 }
00183 
00184 float RPCSimSetUp::getTime(uint32_t id)
00185 {
00186   RPCDetId rpcid(id);
00187   std::map<RPCDetId, float>::iterator iter = _bxmap.find(rpcid);
00188   if(iter == _bxmap.end()){
00189     throw cms::Exception("DataCorrupt") 
00190       << "Exception comming from RPCSimSetUp - no timing information for rpcid.rawId()\t"<<rpcid.rawId()<< std::endl;
00191   }
00192   return iter->second;
00193 }
00194 
00195 const std::map< int, std::vector<double> >& RPCSimSetUp::getClsMap()
00196 {
00197   if(_clsMap.size()!=5){
00198     throw cms::Exception("DataCorrupt") 
00199       << "Exception comming from RPCSimSetUp - cluster size - a wrong format "<< std::endl;
00200   }
00201   return _clsMap;
00202 }
00203 
00204 
00205 //const std::map<int, std::vector<double> >& RPCSimSetUp::getClsMap(uint32_t id)
00206 const std::vector<double>& RPCSimSetUp::getCls(uint32_t id)
00207 {
00208 
00209   map<uint32_t,std::vector<double> >::iterator iter = _mapDetClsMap.find(id);
00210   if(iter == _mapDetClsMap.end()){
00211     throw cms::Exception("DataCorrupt") 
00212       << "Exception comming from RPCSimSetUp - no cluster size information for DetId\t"<<id<< std::endl;
00213   }
00214   if((iter->second).size() != 100){
00215     throw cms::Exception("DataCorrupt") 
00216       << "Exception comming from RPCSimSetUp - cluster size information in a wrong format for DetId\t"<<id<< std::endl;
00217   }
00218   return iter->second;
00219 }
00220 
00221 RPCSimSetUp::~RPCSimSetUp(){}