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/RPCSimModelTiming.h"
0004 #include "SimMuon/RPCDigitizer/src/RPCSimSetUp.h"
0005
0006 #include "SimMuon/RPCDigitizer/src/RPCSynchronizer.h"
0007 #include "Geometry/CommonTopologies/interface/RectangularStripTopology.h"
0008 #include "Geometry/CommonTopologies/interface/TrapezoidalStripTopology.h"
0009 #include "Geometry/RPCGeometry/interface/RPCGeomServ.h"
0010
0011 #include <cmath>
0012
0013 #include "FWCore/Framework/interface/Frameworkfwd.h"
0014 #include "FWCore/Framework/interface/EventSetup.h"
0015 #include "FWCore/Framework/interface/Event.h"
0016 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0017
0018 #include "SimDataFormats/TrackingHit/interface/PSimHitContainer.h"
0019 #include "SimDataFormats/TrackingHit/interface/PSimHit.h"
0020 #include "Geometry/RPCGeometry/interface/RPCGeometry.h"
0021 #include "Geometry/Records/interface/MuonGeometryRecord.h"
0022 #include "DataFormats/MuonDetId/interface/RPCDetId.h"
0023 #include "SimMuon/RPCDigitizer/src/RPCSimSetUp.h"
0024
0025 #include <cstring>
0026 #include <iostream>
0027 #include <fstream>
0028 #include <string>
0029 #include <vector>
0030 #include <cstdlib>
0031 #include <utility>
0032 #include <map>
0033
0034 #include "CLHEP/Random/RandFlat.h"
0035 #include "CLHEP/Random/RandPoissonQ.h"
0036 #include "CLHEP/Random/RandGaussQ.h"
0037
0038 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0039
0040 RPCSimModelTiming::RPCSimModelTiming(const edm::ParameterSet& config) : RPCSim(config) {
0041 aveEff = config.getParameter<double>("averageEfficiency");
0042 aveCls = config.getParameter<double>("averageClusterSize");
0043 resRPC = config.getParameter<double>("timeResolution");
0044 timOff = config.getParameter<double>("timingRPCOffset");
0045 dtimCs = config.getParameter<double>("deltatimeAdjacentStrip");
0046 resEle = config.getParameter<double>("timeJitter");
0047 sspeed = config.getParameter<double>("signalPropagationSpeed");
0048 lbGate = config.getParameter<double>("linkGateWidth");
0049 rpcdigiprint = config.getParameter<bool>("printOutDigitizer");
0050
0051 rate = config.getParameter<double>("Rate");
0052 nbxing = config.getParameter<int>("Nbxing");
0053 gate = config.getParameter<double>("Gate");
0054 frate = config.getParameter<double>("Frate");
0055 do_Y = config.getParameter<bool>("do_Y_coordinate");
0056 sigmaY = config.getParameter<double>("sigmaY");
0057 eledig = config.getParameter<bool>("digitizeElectrons");
0058
0059 if (rpcdigiprint) {
0060 edm::LogInfo("RPC digitizer parameters") << "Average Efficiency = " << aveEff;
0061 edm::LogInfo("RPC digitizer parameters") << "Average Cluster Size = " << aveCls << " strips";
0062 edm::LogInfo("RPC digitizer parameters") << "RPC Time Resolution = " << resRPC << " ns";
0063 edm::LogInfo("RPC digitizer parameters") << "RPC Signal formation time = " << timOff << " ns";
0064 edm::LogInfo("RPC digitizer parameters") << "RPC adjacent strip delay = " << dtimCs << " ns";
0065 edm::LogInfo("RPC digitizer parameters") << "Electronic Jitter = " << resEle << " ns";
0066 edm::LogInfo("RPC digitizer parameters") << "Signal propagation time = " << sspeed << " x c";
0067 edm::LogInfo("RPC digitizer parameters") << "Link Board Gate Width = " << lbGate << " ns";
0068 }
0069
0070 _rpcSync = new RPCSynchronizer(config);
0071 }
0072
0073 RPCSimModelTiming::~RPCSimModelTiming() { delete _rpcSync; }
0074
0075 void RPCSimModelTiming::simulate(const RPCRoll* roll,
0076 const edm::PSimHitContainer& rpcHits,
0077 CLHEP::HepRandomEngine* engine) {
0078 _rpcSync->setRPCSimSetUp(getRPCSimSetUp());
0079 theRpcDigiSimLinks.clear();
0080 theDetectorHitMap.clear();
0081 theRpcDigiSimLinks = RPCDigiSimLinks(roll->id().rawId());
0082
0083 RPCDetId rpcId = roll->id();
0084 RPCGeomServ RPCname(rpcId);
0085
0086 const Topology& topology = roll->specs()->topology();
0087
0088 for (edm::PSimHitContainer::const_iterator _hit = rpcHits.begin(); _hit != rpcHits.end(); ++_hit) {
0089 if (!eledig && _hit->particleType() == 11)
0090 continue;
0091
0092 const LocalPoint& entr = _hit->entryPoint();
0093
0094 int time_hit = _rpcSync->getSimHitBxAndTimingForIRPC(&(*_hit), engine);
0095 double precise_time = _rpcSync->getSmearedTime();
0096
0097 float posX = roll->strip(_hit->localPosition()) - static_cast<int>(roll->strip(_hit->localPosition()));
0098
0099 std::vector<float> veff = (getRPCSimSetUp())->getEff(rpcId.rawId());
0100
0101
0102 int centralStrip = topology.channel(entr) + 1;
0103 ;
0104 float fire = CLHEP::RandFlat::shoot(engine);
0105
0106 float smearedPositionY = CLHEP::RandGaussQ::shoot(engine, _hit->localPosition().y(), sigmaY);
0107
0108 if (fire < veff[centralStrip - 1]) {
0109 int fstrip = centralStrip;
0110 int lstrip = centralStrip;
0111
0112
0113 int clsize = this->getClSize(rpcId.rawId(), posX, engine);
0114 std::vector<int> cls;
0115 cls.push_back(centralStrip);
0116 if (clsize > 1) {
0117 for (int cl = 0; cl < (clsize - 1) / 2; cl++) {
0118 if (centralStrip - cl - 1 >= 1) {
0119 fstrip = centralStrip - cl - 1;
0120 cls.push_back(fstrip);
0121 }
0122 if (centralStrip + cl + 1 <= roll->nstrips()) {
0123 lstrip = centralStrip + cl + 1;
0124 cls.push_back(lstrip);
0125 }
0126 }
0127 if (clsize % 2 == 0) {
0128
0129
0130 int lr = LeftRightNeighbour(*roll, entr, centralStrip);
0131 if (lr == 1) {
0132 if (lstrip < roll->nstrips()) {
0133 lstrip++;
0134 cls.push_back(lstrip);
0135 }
0136 } else {
0137 if (fstrip > 1) {
0138 fstrip--;
0139 cls.push_back(fstrip);
0140 }
0141 }
0142 }
0143 }
0144
0145
0146
0147
0148 for (std::vector<int>::iterator i = cls.begin(); i != cls.end(); i++) {
0149 std::pair<int, int> digi(*i, time_hit);
0150 RPCDigi adigi(*i, time_hit);
0151 adigi.hasTime(true);
0152 adigi.setTime(precise_time);
0153 if (do_Y) {
0154 adigi.hasY(true);
0155 adigi.setY(smearedPositionY);
0156 adigi.setDeltaY(sigmaY);
0157 }
0158 irpc_digis.insert(adigi);
0159 theDetectorHitMap.insert(DetectorHitMap::value_type(digi, &(*_hit)));
0160 }
0161 }
0162 }
0163 }
0164
0165 void RPCSimModelTiming::simulateNoise(const RPCRoll* roll, CLHEP::HepRandomEngine* engine) {
0166 RPCDetId rpcId = roll->id();
0167 RPCGeomServ RPCname(rpcId);
0168 std::vector<float> vnoise = (getRPCSimSetUp())->getNoise(rpcId.rawId());
0169 std::vector<float> veff = (getRPCSimSetUp())->getEff(rpcId.rawId());
0170 unsigned int nstrips = roll->nstrips();
0171 double area = 0.0;
0172 float striplength, xmin, xmax;
0173 if (rpcId.region() == 0) {
0174 const RectangularStripTopology* top_ = dynamic_cast<const RectangularStripTopology*>(&(roll->topology()));
0175 xmin = (top_->localPosition(0.)).x();
0176 xmax = (top_->localPosition((float)roll->nstrips())).x();
0177 striplength = (top_->stripLength());
0178 area = striplength * (xmax - xmin);
0179 } else {
0180 const TrapezoidalStripTopology* top_ = dynamic_cast<const TrapezoidalStripTopology*>(&(roll->topology()));
0181 xmin = (top_->localPosition(0.)).x();
0182 xmax = (top_->localPosition((float)roll->nstrips())).x();
0183 striplength = (top_->stripLength());
0184 area = striplength * (xmax - xmin);
0185 }
0186
0187 for (unsigned int j = 0; j < vnoise.size(); ++j) {
0188 if (j >= nstrips)
0189 break;
0190
0191 double ave = vnoise[j] * nbxing * gate * area * 1.0e-9 * frate / ((float)roll->nstrips());
0192
0193 CLHEP::RandPoissonQ randPoissonQ(*engine, ave);
0194 N_hits = randPoissonQ.fire();
0195 for (int i = 0; i < N_hits; i++) {
0196 double precise_time = CLHEP::RandFlat::shoot(engine, (nbxing * gate) / gate);
0197 int time_hit = (static_cast<int>(precise_time)) - nbxing / 2;
0198 RPCDigi adigi(j + 1, time_hit);
0199 adigi.hasTime(true);
0200 adigi.setTime(precise_time);
0201 if (do_Y) {
0202 double positionY = CLHEP::RandFlat::shoot(engine, striplength);
0203 positionY -= striplength / 2;
0204 adigi.hasY(true);
0205 adigi.setY(positionY);
0206 adigi.setDeltaY(sigmaY);
0207 }
0208 irpc_digis.insert(adigi);
0209 }
0210 }
0211 }
0212
0213 int RPCSimModelTiming::getClSize(uint32_t id, float posX, CLHEP::HepRandomEngine* engine) {
0214 std::vector<double> clsForDetId = getRPCSimSetUp()->getCls(id);
0215
0216 int cnt = 1;
0217 int min = 1;
0218 double func = 0.0;
0219 std::vector<double> sum_clsize;
0220
0221 sum_clsize.clear();
0222 sum_clsize = clsForDetId;
0223 int vectOffset(0);
0224
0225 double rr_cl = CLHEP::RandFlat::shoot(engine);
0226
0227 if (0.0 <= posX && posX < 0.2) {
0228 func = clsForDetId[19] * (rr_cl);
0229 vectOffset = 0;
0230 }
0231 if (0.2 <= posX && posX < 0.4) {
0232 func = clsForDetId[39] * (rr_cl);
0233 vectOffset = 20;
0234 }
0235 if (0.4 <= posX && posX < 0.6) {
0236 func = clsForDetId[59] * (rr_cl);
0237 vectOffset = 40;
0238 }
0239 if (0.6 <= posX && posX < 0.8) {
0240 func = clsForDetId[79] * (rr_cl);
0241 vectOffset = 60;
0242 }
0243 if (0.8 <= posX && posX < 1.0) {
0244 func = clsForDetId[89] * (rr_cl);
0245 vectOffset = 80;
0246 }
0247
0248 for (int i = vectOffset; i < (vectOffset + 20); i++) {
0249 cnt++;
0250 if (func > clsForDetId[i]) {
0251 min = cnt;
0252 } else if (func < clsForDetId[i]) {
0253 break;
0254 }
0255 }
0256 return min;
0257 }
0258
0259 int RPCSimModelTiming::LeftRightNeighbour(const RPCRoll& roll, const LocalPoint& hit_pos, int strip) {
0260
0261
0262
0263 int leftStrip = strip - 1;
0264 int rightStrip = strip + 1;
0265
0266 if (leftStrip < 0)
0267 return +1;
0268 if (rightStrip > roll.nstrips())
0269 return -1;
0270
0271 double deltawL = fabs((roll.centreOfStrip(leftStrip)).x() - hit_pos.x());
0272 double deltawR = fabs((roll.centreOfStrip(rightStrip)).x() - hit_pos.x());
0273
0274 if (deltawL >= deltawR) {
0275 return +1;
0276 } else {
0277 return -1;
0278 }
0279 }