CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch13/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   double sum = 0;
00049   unsigned int counter = 1;
00050   unsigned int row = 1;
00051   std::vector<double> sum_clsize;
00052 
00053   for(unsigned int n = 0; n < vcls.size(); ++n){
00054 
00055     sum_clsize.push_back(vcls[n]);
00056 
00057     if(counter == row*20) {
00058 
00059       _clsMap[row] = sum_clsize;
00060       row++;
00061       sum = 0;
00062       sum_clsize.clear();
00063     }
00064     counter++;
00065   }
00066 
00067   unsigned int n = 0; 
00068   uint32_t temp = 0; 
00069   std::vector<float> veff, vvnoise;
00070   veff.clear();
00071   vvnoise.clear();
00072 
00073   for(std::vector<RPCStripNoises::NoiseItem>::iterator it = vnoise.begin(); it != vnoise.end(); ++it){
00074     if(n%96 == 0) {
00075       if(n > 0 ){
00076         _mapDetIdNoise[temp]= vvnoise;
00077         _mapDetIdEff[temp] = veff;
00078         _bxmap[RPCDetId(it->dpid)] = it->time;
00079         
00080         veff.clear();
00081         vvnoise.clear();
00082         vvnoise.push_back((it->noise));
00083         veff.push_back((it->eff));
00084       }
00085       else if(n == 0 ){
00086         vvnoise.push_back((it->noise));
00087         veff.push_back((it->eff));
00088         _bxmap[RPCDetId(it->dpid)] = it->time;
00089       }
00090     } else if (n == vnoise.size()-1 ){
00091       temp = it->dpid;
00092       vvnoise.push_back((it->noise));
00093       veff.push_back((it->eff));
00094       _mapDetIdNoise[temp]= vvnoise;
00095       _mapDetIdEff[temp] = veff;
00096     } else {
00097       temp = it->dpid;
00098       vvnoise.push_back((it->noise));
00099       veff.push_back((it->eff));
00100     }
00101     n++;
00102   }
00103 }
00104 
00105 void RPCSimSetUp::setRPCSetUp(std::vector<RPCStripNoises::NoiseItem> vnoise, std::vector<RPCClusterSize::ClusterSizeItem> vClusterSize){
00106 
00107   std::vector<RPCClusterSize::ClusterSizeItem>::iterator itCls;
00108   uint32_t  detId;
00109   int clsCounter(1);
00110   std::vector<double> clsVect;
00111 
00112   for(itCls = vClusterSize.begin(); itCls != vClusterSize.end(); ++itCls){
00113     clsVect.push_back(((double)(itCls->clusterSize)));
00114     if((!(clsCounter%100)) && (clsCounter!=0)){
00115       detId=itCls->dpid;
00116       _mapDetClsMap[detId]=clsVect;
00117       clsVect.clear();
00118       clsCounter=0;
00119     }
00120     ++clsCounter;
00121   }
00122 
00123   unsigned int n = 0; 
00124   uint32_t temp = 0; 
00125   std::vector<float> veff, vvnoise;
00126   veff.clear();
00127   vvnoise.clear();
00128 
00129   for(std::vector<RPCStripNoises::NoiseItem>::iterator it = vnoise.begin(); it != vnoise.end(); ++it){
00130     if(n%96 == 0) {
00131       if(n > 0 ){
00132         _mapDetIdNoise[temp]= vvnoise;
00133         _mapDetIdEff[temp] = veff;
00134         _bxmap[RPCDetId(it->dpid)] = it->time;
00135         
00136         veff.clear();
00137         vvnoise.clear();
00138         vvnoise.push_back((it->noise));
00139         veff.push_back((it->eff));
00140       }
00141       else if(n == 0 ){
00142         vvnoise.push_back((it->noise));
00143         veff.push_back((it->eff));
00144         _bxmap[RPCDetId(it->dpid)] = it->time;
00145       }
00146     } else if (n == vnoise.size()-1 ){
00147       temp = it->dpid;
00148       vvnoise.push_back((it->noise));
00149       veff.push_back((it->eff));
00150       _mapDetIdNoise[temp]= vvnoise;
00151       _mapDetIdEff[temp] = veff;
00152     } else {
00153       temp = it->dpid;
00154       vvnoise.push_back((it->noise));
00155       veff.push_back((it->eff));
00156     }
00157     n++;
00158   }
00159 }
00160 
00161 
00162 const std::vector<float>& RPCSimSetUp::getNoise(uint32_t id)
00163 {
00164   map<uint32_t,std::vector<float> >::iterator iter = _mapDetIdNoise.find(id);
00165   if(iter == _mapDetIdNoise.end()){
00166     throw cms::Exception("DataCorrupt") 
00167       << "Exception comming from RPCSimSetUp - no noise information for DetId\t"<<id<< std::endl;
00168   }
00169   return iter->second;
00170 }
00171 
00172 const std::vector<float>& RPCSimSetUp::getEff(uint32_t id)
00173 {
00174   map<uint32_t,std::vector<float> >::iterator iter = _mapDetIdEff.find(id);
00175   if(iter == _mapDetIdEff.end()){
00176     throw cms::Exception("DataCorrupt") 
00177       << "Exception comming from RPCSimSetUp - no efficiency information for DetId\t"<<id<< std::endl;
00178   }
00179   if((iter->second).size() != 96){
00180     throw cms::Exception("DataCorrupt") 
00181       << "Exception comming from RPCSimSetUp - efficiency information in a wrong format for DetId\t"<<id<< std::endl;
00182   }
00183   return iter->second;
00184 }
00185 
00186 float RPCSimSetUp::getTime(uint32_t id)
00187 {
00188   RPCDetId rpcid(id);
00189   std::map<RPCDetId, float>::iterator iter = _bxmap.find(rpcid);
00190   if(iter == _bxmap.end()){
00191     throw cms::Exception("DataCorrupt") 
00192       << "Exception comming from RPCSimSetUp - no timing information for rpcid.rawId()\t"<<rpcid.rawId()<< std::endl;
00193   }
00194   return iter->second;
00195 }
00196 
00197 const std::map< int, std::vector<double> >& RPCSimSetUp::getClsMap()
00198 {
00199   if(_clsMap.size()!=5){
00200     throw cms::Exception("DataCorrupt") 
00201       << "Exception comming from RPCSimSetUp - cluster size - a wrong format "<< std::endl;
00202   }
00203   return _clsMap;
00204 }
00205 
00206 
00207 //const std::map<int, std::vector<double> >& RPCSimSetUp::getClsMap(uint32_t id)
00208 const std::vector<double>& RPCSimSetUp::getCls(uint32_t id)
00209 {
00210 
00211   map<uint32_t,std::vector<double> >::iterator iter = _mapDetClsMap.find(id);
00212   if(iter == _mapDetClsMap.end()){
00213     throw cms::Exception("DataCorrupt") 
00214       << "Exception comming from RPCSimSetUp - no cluster size information for DetId\t"<<id<< std::endl;
00215   }
00216   if((iter->second).size() != 100){
00217     throw cms::Exception("DataCorrupt") 
00218       << "Exception comming from RPCSimSetUp - cluster size information in a wrong format for DetId\t"<<id<< std::endl;
00219   }
00220   return iter->second;
00221 }
00222 
00223 RPCSimSetUp::~RPCSimSetUp(){}