Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:30:50

0001 #include "Geometry/RPCGeometry/interface/RPCRoll.h"
0002 #include "Geometry/RPCGeometry/interface/RPCRollSpecs.h"
0003 #include "SimMuon/RPCDigitizer/src/RPCSimParam.h"
0004 #include "Geometry/CommonTopologies/interface/RectangularStripTopology.h"
0005 #include "Geometry/CommonTopologies/interface/TrapezoidalStripTopology.h"
0006 
0007 #include <cmath>
0008 
0009 #include "CLHEP/Random/RandFlat.h"
0010 #include "CLHEP/Random/RandPoissonQ.h"
0011 
0012 RPCSimParam::RPCSimParam(const edm::ParameterSet& config) : RPCSim(config) {
0013   aveEff = config.getParameter<double>("averageEfficiency");
0014   aveCls = config.getParameter<double>("averageClusterSize");
0015   resRPC = config.getParameter<double>("timeResolution");
0016   timOff = config.getParameter<double>("timingRPCOffset");
0017   dtimCs = config.getParameter<double>("deltatimeAdjacentStrip");
0018   resEle = config.getParameter<double>("timeJitter");
0019   sspeed = config.getParameter<double>("signalPropagationSpeed");
0020   lbGate = config.getParameter<double>("linkGateWidth");
0021   rpcdigiprint = config.getParameter<bool>("printOutDigitizer");
0022 
0023   rate = config.getParameter<double>("Rate");
0024   nbxing = config.getParameter<int>("Nbxing");
0025   gate = config.getParameter<double>("Gate");
0026 
0027   if (rpcdigiprint) {
0028     std::cout << "Average Efficiency        = " << aveEff << std::endl;
0029     std::cout << "Average Cluster Size      = " << aveCls << " strips" << std::endl;
0030     std::cout << "RPC Time Resolution       = " << resRPC << " ns" << std::endl;
0031     std::cout << "RPC Signal formation time = " << timOff << " ns" << std::endl;
0032     std::cout << "RPC adjacent strip delay  = " << dtimCs << " ns" << std::endl;
0033     std::cout << "Electronic Jitter         = " << resEle << " ns" << std::endl;
0034     std::cout << "Signal propagation time   = " << sspeed << " x c" << std::endl;
0035     std::cout << "Link Board Gate Width     = " << lbGate << " ns" << std::endl;
0036   }
0037 
0038   _rpcSync = new RPCSynchronizer(config);
0039 }
0040 
0041 RPCSimParam::~RPCSimParam() { delete _rpcSync; }
0042 
0043 void RPCSimParam::simulate(const RPCRoll* roll, const edm::PSimHitContainer& rpcHits, CLHEP::HepRandomEngine* engine) {
0044   _rpcSync->setRPCSimSetUp(getRPCSimSetUp());
0045   theRpcDigiSimLinks.clear();
0046   theDetectorHitMap.clear();
0047   theRpcDigiSimLinks = RPCDigiSimLinks(roll->id().rawId());
0048 
0049   const Topology& topology = roll->specs()->topology();
0050 
0051   for (edm::PSimHitContainer::const_iterator _hit = rpcHits.begin(); _hit != rpcHits.end(); ++_hit) {
0052     // Here I hould check if the RPC are up side down;
0053     const LocalPoint& entr = _hit->entryPoint();
0054     int time_hit = _rpcSync->getSimHitBx(&(*_hit), engine);
0055 
0056     // Effinciecy
0057     float eff = CLHEP::RandFlat::shoot(engine);
0058     if (eff < aveEff) {
0059       int centralStrip = topology.channel(entr) + 1;
0060       int fstrip = centralStrip;
0061       int lstrip = centralStrip;
0062       // Compute the cluster size
0063       double w = CLHEP::RandFlat::shoot(engine);
0064       if (w < 1.e-10)
0065         w = 1.e-10;
0066       int clsize = static_cast<int>(-1. * aveCls * log(w) + 1.);
0067       std::vector<int> cls;
0068       cls.push_back(centralStrip);
0069       if (clsize > 1) {
0070         for (int cl = 0; cl < (clsize - 1) / 2; cl++)
0071           if (centralStrip - cl - 1 >= 1) {
0072             fstrip = centralStrip - cl - 1;
0073             cls.push_back(fstrip);
0074           }
0075         for (int cl = 0; cl < (clsize - 1) / 2; cl++)
0076           if (centralStrip + cl + 1 <= roll->nstrips()) {
0077             lstrip = centralStrip + cl + 1;
0078             cls.push_back(lstrip);
0079           }
0080         if (clsize % 2 == 0) {
0081           // insert the last strip according to the
0082           // simhit position in the central strip
0083           double deltaw = roll->centreOfStrip(centralStrip).x() - entr.x();
0084           if (deltaw < 0.) {
0085             if (lstrip < roll->nstrips()) {
0086               lstrip++;
0087               cls.push_back(lstrip);
0088             }
0089           } else {
0090             if (fstrip > 1) {
0091               fstrip--;
0092               cls.push_back(fstrip);
0093             }
0094           }
0095         }
0096       }
0097 
0098       for (std::vector<int>::iterator i = cls.begin(); i != cls.end(); i++) {
0099         // Check the timing of the adjacent strip
0100         std::pair<unsigned int, int> digi(*i, time_hit);
0101 
0102         theDetectorHitMap.insert(DetectorHitMap::value_type(digi, &(*_hit)));
0103         strips.insert(digi);
0104       }
0105     }
0106   }
0107 }
0108 
0109 void RPCSimParam::simulateNoise(const RPCRoll* roll, CLHEP::HepRandomEngine* engine) {
0110   RPCDetId rpcId = roll->id();
0111   int nstrips = roll->nstrips();
0112   double area = 0.0;
0113 
0114   if (rpcId.region() == 0) {
0115     const RectangularStripTopology* top_ = dynamic_cast<const RectangularStripTopology*>(&(roll->topology()));
0116     float xmin = (top_->localPosition(0.)).x();
0117     float xmax = (top_->localPosition((float)roll->nstrips())).x();
0118     float striplength = (top_->stripLength());
0119     area = striplength * (xmax - xmin);
0120   } else {
0121     const TrapezoidalStripTopology* top_ = dynamic_cast<const TrapezoidalStripTopology*>(&(roll->topology()));
0122     float xmin = (top_->localPosition(0.)).x();
0123     float xmax = (top_->localPosition((float)roll->nstrips())).x();
0124     float striplength = (top_->stripLength());
0125     area = striplength * (xmax - xmin);
0126   }
0127 
0128   double ave = rate * nbxing * gate * area * 1.0e-9;
0129 
0130   CLHEP::RandPoissonQ randPoissonQ(*engine, ave);
0131   N_hits = randPoissonQ.fire();
0132 
0133   for (int i = 0; i < N_hits; i++) {
0134     int strip = static_cast<int>(CLHEP::RandFlat::shoot(engine, 1, nstrips));
0135     int time_hit;
0136     time_hit = (static_cast<int>(CLHEP::RandFlat::shoot(engine, (nbxing * gate) / gate))) - nbxing / 2;
0137     std::pair<int, int> digi(strip, time_hit);
0138     strips.insert(digi);
0139   }
0140 }