CMS 3D CMS Logo

RPCSimAverageNoise.cc

Go to the documentation of this file.
00001 #include "Geometry/RPCGeometry/interface/RPCRoll.h"
00002 #include "Geometry/RPCGeometry/interface/RPCRollSpecs.h"
00003 #include "SimMuon/RPCDigitizer/src/RPCSimAverageNoise.h"
00004 
00005 #include "SimMuon/RPCDigitizer/src/RPCSynchronizer.h"
00006 #include "Geometry/CommonTopologies/interface/RectangularStripTopology.h"
00007 #include "Geometry/CommonTopologies/interface/TrapezoidalStripTopology.h"
00008 
00009 #include <cmath>
00010 #include "FWCore/ServiceRegistry/interface/Service.h"
00011 #include "FWCore/Utilities/interface/RandomNumberGenerator.h"
00012 #include "FWCore/Utilities/interface/Exception.h"
00013 #include "CLHEP/Random/RandomEngine.h"
00014 #include "CLHEP/Random/RandFlat.h"
00015 #include <CLHEP/Random/RandGaussQ.h>
00016 #include <CLHEP/Random/RandFlat.h>
00017 
00018 #include <FWCore/Framework/interface/Frameworkfwd.h>
00019 #include <FWCore/Framework/interface/EventSetup.h>
00020 #include <FWCore/Framework/interface/EDAnalyzer.h>
00021 #include <FWCore/Framework/interface/Event.h>
00022 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00023 #include <FWCore/Framework/interface/ESHandle.h>
00024 
00025 #include "SimDataFormats/TrackingHit/interface/PSimHitContainer.h"
00026 #include "SimDataFormats/TrackingHit/interface/PSimHit.h"
00027 #include "Geometry/RPCGeometry/interface/RPCGeometry.h"
00028 #include <Geometry/Records/interface/MuonGeometryRecord.h>
00029 #include "DataFormats/MuonDetId/interface/RPCDetId.h"
00030 #include "SimMuon/RPCDigitizer/src/RPCSimSetUp.h"
00031 
00032 #include<cstring>
00033 #include<iostream>
00034 #include<fstream>
00035 #include<string>
00036 #include<vector>
00037 #include<stdlib.h>
00038 #include <utility>
00039 #include <map>
00040 
00041 #include "CLHEP/config/CLHEP.h"
00042 #include "CLHEP/Random/Random.h"
00043 #include "CLHEP/Random/RandFlat.h"
00044 #include "CLHEP/Random/RandPoissonQ.h"
00045 
00046 using namespace std;
00047 
00048 RPCSimAverageNoise::RPCSimAverageNoise(const edm::ParameterSet& config) : 
00049   RPCSim(config)
00050 {
00051 
00052   _rpcSync = new RPCSynchronizer(config);
00053 
00054   aveEff = config.getParameter<double>("averageEfficiency");
00055   aveCls = config.getParameter<double>("averageClusterSize");
00056   resRPC = config.getParameter<double>("timeResolution");
00057   timOff = config.getParameter<double>("timingRPCOffset");
00058   dtimCs = config.getParameter<double>("deltatimeAdjacentStrip");
00059   resEle = config.getParameter<double>("timeJitter");
00060   sspeed = config.getParameter<double>("signalPropagationSpeed");
00061   lbGate = config.getParameter<double>("linkGateWidth");
00062   rpcdigiprint = config.getParameter<bool>("printOutDigitizer");
00063   rate=config.getParameter<double>("Rate");
00064   nbxing=config.getParameter<int>("Nbxing");
00065   gate=config.getParameter<double>("Gate");
00066   frate=config.getParameter<double>("Frate");
00067 
00068   if (rpcdigiprint) {
00069     std::cout <<"Average Efficiency        = "<<aveEff<<std::endl;
00070     std::cout <<"Average Cluster Size      = "<<aveCls<<" strips"<<std::endl;
00071     std::cout <<"RPC Time Resolution       = "<<resRPC<<" ns"<<std::endl;
00072     std::cout <<"RPC Signal formation time = "<<timOff<<" ns"<<std::endl;
00073     std::cout <<"RPC adjacent strip delay  = "<<dtimCs<<" ns"<<std::endl;
00074     std::cout <<"Electronic Jitter         = "<<resEle<<" ns"<<std::endl;
00075     std::cout <<"Signal propagation time   = "<<sspeed<<" x c"<<std::endl;
00076     std::cout <<"Link Board Gate Width     = "<<lbGate<<" ns"<<std::endl;
00077   }
00078 
00079   edm::Service<edm::RandomNumberGenerator> rng;
00080   if ( ! rng.isAvailable()) {
00081     throw cms::Exception("Configuration")
00082       << "RPCDigitizer requires the RandomNumberGeneratorService\n"
00083       "which is not present in the configuration file.  You must add the service\n"
00084       "in the configuration file or remove the modules that require it.";
00085   }
00086   
00087   rndEngine = &(rng->getEngine());
00088   flatDistribution = new CLHEP::RandFlat(rndEngine);
00089 }
00090 
00091 RPCSimAverageNoise::~RPCSimAverageNoise(){}
00092 
00093 int RPCSimAverageNoise::getClSize(float posX)
00094 {
00095 
00096   std::map< int, std::vector<double> > clsMap = getRPCSimSetUp()->getClsMap();
00097 
00098   int cnt = 1;
00099   int min = 1;
00100   int max = 1;
00101   double func=0.0;
00102   std::vector<double> sum_clsize;
00103 
00104   double rr_cl = flatDistribution->fire();
00105   if(0.0 <= posX && posX < 0.2)  {
00106     func = (clsMap[1])[(clsMap[1]).size()-1]*(rr_cl);
00107     sum_clsize = clsMap[1];
00108   }
00109   if(0.2 <= posX && posX < 0.4) {
00110     func = (clsMap[2])[(clsMap[2]).size()-1]*(rr_cl);
00111     sum_clsize = clsMap[2];
00112   }
00113   if(0.4 <= posX && posX < 0.6) {
00114     func = (clsMap[3])[(clsMap[3]).size()-1]*(rr_cl);
00115     sum_clsize = clsMap[3];
00116   }
00117   if(0.6 <= posX && posX < 0.8) {
00118     func = (clsMap[4])[(clsMap[4]).size()-1]*(rr_cl);
00119     sum_clsize = clsMap[4];
00120   }
00121   if(0.8 <= posX && posX < 1.0)  {
00122     func = (clsMap[5])[(clsMap[5]).size()-1]*(rr_cl);
00123     sum_clsize = clsMap[5];
00124   }
00125 
00126   for(vector<double>::iterator iter = sum_clsize.begin();
00127       iter != sum_clsize.end(); ++iter){
00128     cnt++;
00129     if(func > (*iter)){
00130       min = cnt;
00131     }
00132     else if(func < (*iter)){
00133       max = cnt;
00134       break;
00135     }
00136   }
00137   return min;
00138 }
00139 
00140 
00141 void
00142 RPCSimAverageNoise::simulate(const RPCRoll* roll,
00143                         const edm::PSimHitContainer& rpcHits)
00144 {
00145   _rpcSync->setRPCSimSetUp(getRPCSimSetUp());
00146   theRpcDigiSimLinks.clear();
00147   theDetectorHitMap.clear();
00148   theRpcDigiSimLinks = RPCDigiSimLinks(roll->id().rawId());
00149 
00150   const Topology& topology=roll->specs()->topology();
00151 
00152   for (edm::PSimHitContainer::const_iterator _hit = rpcHits.begin();
00153        _hit != rpcHits.end(); ++_hit){
00154 
00155     // Here I hould check if the RPC are up side down;
00156     const LocalPoint& entr=_hit->entryPoint();
00157     int time_hit = _rpcSync->getSimHitBx(&(*_hit));
00158     float posX = roll->strip(_hit->localPosition()) - static_cast<int>(roll->strip(_hit->localPosition()));
00159 
00160     // Effinciecy
00161 
00162     if (flatDistribution->fire() < aveEff) {
00163 
00164       int centralStrip = topology.channel(entr)+1;  
00165       int fstrip=centralStrip;
00166       int lstrip=centralStrip;
00167       // Compute the cluster size
00168       double w = flatDistribution->fire(1);
00169       if (w < 1.e-10) w=1.e-10;
00170       int clsize = this->getClSize(posX);
00171 
00172       std::vector<int> cls;
00173       cls.push_back(centralStrip);
00174       if (clsize > 1){
00175         for (int cl = 0; cl < (clsize-1)/2; cl++)
00176           if (centralStrip - cl -1 >= 1  ){
00177             fstrip = centralStrip-cl-1;
00178             cls.push_back(fstrip);
00179           }
00180         for (int cl = 0; cl < (clsize-1)/2; cl++)
00181           if (centralStrip + cl + 1 <= roll->nstrips() ){
00182             lstrip = centralStrip+cl+1;
00183             cls.push_back(lstrip);
00184           }
00185         if (clsize%2 == 0 ){
00186           // insert the last strip according to the 
00187           // simhit position in the central strip 
00188           double deltaw=roll->centreOfStrip(centralStrip).x()-entr.x();
00189           if (deltaw<0.) {
00190             if (lstrip < roll->nstrips() ){
00191               lstrip++;
00192               cls.push_back(lstrip);
00193             }
00194           }else{
00195             if (fstrip > 1 ){
00196               fstrip--;
00197               cls.push_back(fstrip);
00198             }
00199           }
00200         }
00201       }
00202 
00203       for (std::vector<int>::iterator i=cls.begin(); i!=cls.end();i++){
00204         // Check the timing of the adjacent strip
00205         std::pair<int, int> digi(*i,time_hit );
00206 
00207         theDetectorHitMap.insert(DetectorHitMap::value_type(digi,&(*_hit)));
00208         strips.insert(digi);
00209       }
00210     }
00211   }
00212 }
00213 
00214 void RPCSimAverageNoise::simulateNoise(const RPCRoll* roll)
00215 {
00216   RPCDetId rpcId = roll->id();
00217   std::vector<float> vnoise = (getRPCSimSetUp())->getNoise(rpcId.rawId());
00218   unsigned int nstrips = roll->nstrips();
00219 
00220   double area = 0.0;
00221   
00222   if ( rpcId.region() == 0 )
00223     {
00224       const RectangularStripTopology* top_ = dynamic_cast<const
00225         RectangularStripTopology*>(&(roll->topology()));
00226       float xmin = (top_->localPosition(0.)).x();
00227       float xmax = (top_->localPosition((float)roll->nstrips())).x();
00228       float striplength = (top_->stripLength());
00229       area = striplength*(xmax-xmin);
00230     }
00231   else
00232     {
00233       const TrapezoidalStripTopology* top_=dynamic_cast<const TrapezoidalStripTopology*>(&(roll->topology()));
00234       float xmin = (top_->localPosition(0.)).x();
00235       float xmax = (top_->localPosition((float)roll->nstrips())).x();
00236       float striplength = (top_->stripLength());
00237       area = striplength*(xmax-xmin);
00238     }
00239 
00240   for(unsigned int j = 0; j < vnoise.size(); ++j){
00241     
00242     if(j >= nstrips) break; 
00243 
00244     double ave = frate*vnoise[j]*nbxing*gate*area*1.0e-9;
00245     poissonDistribution_ = new CLHEP::RandPoissonQ(rndEngine, ave);
00246     N_hits = poissonDistribution_->fire();
00247 
00248     for (int i = 0; i < N_hits; i++ ){
00249       flatDistribution = new CLHEP::RandFlat(rndEngine, (nbxing*gate)/gate);
00250       int time_hit = (static_cast<int>(flatDistribution->fire())) - nbxing/2;
00251       std::pair<int, int> digi(j+1,time_hit);
00252       strips.insert(digi);
00253     }
00254   }
00255 }

Generated on Tue Jun 9 17:47:46 2009 for CMSSW by  doxygen 1.5.4