CMS 3D CMS Logo

RPCSimModelTiming.cc
Go to the documentation of this file.
5 
10 
11 #include <cmath>
12 
19 
26 
27 #include <cstring>
28 #include <iostream>
29 #include <fstream>
30 #include <string>
31 #include <vector>
32 #include <cstdlib>
33 #include <utility>
34 #include <map>
35 
36 #include "CLHEP/Random/RandFlat.h"
37 #include "CLHEP/Random/RandPoissonQ.h"
38 #include "CLHEP/Random/RandGaussQ.h"
39 
41 
43  aveEff = config.getParameter<double>("averageEfficiency");
44  aveCls = config.getParameter<double>("averageClusterSize");
45  resRPC = config.getParameter<double>("timeResolution");
46  timOff = config.getParameter<double>("timingRPCOffset");
47  dtimCs = config.getParameter<double>("deltatimeAdjacentStrip");
48  resEle = config.getParameter<double>("timeJitter");
49  sspeed = config.getParameter<double>("signalPropagationSpeed");
50  lbGate = config.getParameter<double>("linkGateWidth");
51  rpcdigiprint = config.getParameter<bool>("printOutDigitizer");
52 
53  rate = config.getParameter<double>("Rate");
54  nbxing = config.getParameter<int>("Nbxing");
55  gate = config.getParameter<double>("Gate");
56  frate = config.getParameter<double>("Frate");
57  do_Y = config.getParameter<bool>("do_Y_coordinate");
58  sigmaY = config.getParameter<double>("sigmaY");
59  eledig = config.getParameter<bool>("digitizeElectrons");
60 
61  if (rpcdigiprint) {
62  edm::LogInfo("RPC digitizer parameters") << "Average Efficiency = " << aveEff;
63  edm::LogInfo("RPC digitizer parameters") << "Average Cluster Size = " << aveCls << " strips";
64  edm::LogInfo("RPC digitizer parameters") << "RPC Time Resolution = " << resRPC << " ns";
65  edm::LogInfo("RPC digitizer parameters") << "RPC Signal formation time = " << timOff << " ns";
66  edm::LogInfo("RPC digitizer parameters") << "RPC adjacent strip delay = " << dtimCs << " ns";
67  edm::LogInfo("RPC digitizer parameters") << "Electronic Jitter = " << resEle << " ns";
68  edm::LogInfo("RPC digitizer parameters") << "Signal propagation time = " << sspeed << " x c";
69  edm::LogInfo("RPC digitizer parameters") << "Link Board Gate Width = " << lbGate << " ns";
70  }
71 
72  _rpcSync = new RPCSynchronizer(config);
73 }
74 
76 
78  const edm::PSimHitContainer& rpcHits,
79  CLHEP::HepRandomEngine* engine) {
82  theDetectorHitMap.clear();
84 
85  RPCDetId rpcId = roll->id();
86  RPCGeomServ RPCname(rpcId);
87 
88  const Topology& topology = roll->specs()->topology();
89 
90  for (edm::PSimHitContainer::const_iterator _hit = rpcHits.begin(); _hit != rpcHits.end(); ++_hit) {
91  if (!eledig && _hit->particleType() == 11)
92  continue;
93  // Here I hould check if the RPC are up side down;
94  const LocalPoint& entr = _hit->entryPoint();
95 
96  int time_hit = _rpcSync->getSimHitBxAndTimingForIRPC(&(*_hit), engine);
97  double precise_time = _rpcSync->getSmearedTime();
98 
99  float posX = roll->strip(_hit->localPosition()) - static_cast<int>(roll->strip(_hit->localPosition()));
100 
101  std::vector<float> veff = (getRPCSimSetUp())->getEff(rpcId.rawId());
102 
103  // Effinciecy
104  int centralStrip = topology.channel(entr) + 1;
105  ;
106  float fire = CLHEP::RandFlat::shoot(engine);
107 
108  float smearedPositionY = CLHEP::RandGaussQ::shoot(engine, _hit->localPosition().y(), sigmaY);
109 
110  if (fire < veff[centralStrip - 1]) {
111  int fstrip = centralStrip;
112  int lstrip = centralStrip;
113 
114  // Compute the cluster size
115  int clsize = this->getClSize(rpcId.rawId(), posX, engine); // This is for cluster size chamber by chamber
116  std::vector<int> cls;
117  cls.push_back(centralStrip);
118  if (clsize > 1) {
119  for (int cl = 0; cl < (clsize - 1) / 2; cl++) {
120  if (centralStrip - cl - 1 >= 1) {
121  fstrip = centralStrip - cl - 1;
122  cls.push_back(fstrip);
123  }
124  if (centralStrip + cl + 1 <= roll->nstrips()) {
125  lstrip = centralStrip + cl + 1;
126  cls.push_back(lstrip);
127  }
128  }
129  if (clsize % 2 == 0) {
130  // insert the last strip according to the
131  // simhit position in the central strip
132  int lr = LeftRightNeighbour(*roll, entr, centralStrip);
133  if (lr == 1) {
134  if (lstrip < roll->nstrips()) {
135  lstrip++;
136  cls.push_back(lstrip);
137  }
138  } else {
139  if (fstrip > 1) {
140  fstrip--;
141  cls.push_back(fstrip);
142  }
143  }
144  }
145  }
146 
147  //digitize all the strips in the cluster
148  //in the previuos version some strips were dropped
149  //leading to un-physical "shift" of the cluster
150  for (std::vector<int>::iterator i = cls.begin(); i != cls.end(); i++) {
151  std::pair<int, int> digi(*i, time_hit);
152  RPCDigi adigi(*i, time_hit);
153  adigi.hasTime(true);
154  adigi.setTime(precise_time);
155  if (do_Y) {
156  adigi.hasY(true);
157  adigi.setY(smearedPositionY);
158  adigi.setDeltaY(sigmaY);
159  }
160  irpc_digis.insert(adigi);
161  theDetectorHitMap.insert(DetectorHitMap::value_type(digi, &(*_hit)));
162  }
163  }
164  }
165 }
166 
167 void RPCSimModelTiming::simulateNoise(const RPCRoll* roll, CLHEP::HepRandomEngine* engine) {
168  RPCDetId rpcId = roll->id();
169  RPCGeomServ RPCname(rpcId);
170  std::vector<float> vnoise = (getRPCSimSetUp())->getNoise(rpcId.rawId());
171  std::vector<float> veff = (getRPCSimSetUp())->getEff(rpcId.rawId());
172  unsigned int nstrips = roll->nstrips();
173  double area = 0.0;
174  float striplength, xmin, xmax;
175  if (rpcId.region() == 0) {
176  const RectangularStripTopology* top_ = dynamic_cast<const RectangularStripTopology*>(&(roll->topology()));
177  xmin = (top_->localPosition(0.)).x();
178  xmax = (top_->localPosition((float)roll->nstrips())).x();
179  striplength = (top_->stripLength());
180  area = striplength * (xmax - xmin);
181  } else {
182  const TrapezoidalStripTopology* top_ = dynamic_cast<const TrapezoidalStripTopology*>(&(roll->topology()));
183  xmin = (top_->localPosition(0.)).x();
184  xmax = (top_->localPosition((float)roll->nstrips())).x();
185  striplength = (top_->stripLength());
186  area = striplength * (xmax - xmin);
187  }
188 
189  for (unsigned int j = 0; j < vnoise.size(); ++j) {
190  if (j >= nstrips)
191  break;
192 
193  double ave = vnoise[j] * nbxing * gate * area * 1.0e-9 * frate / ((float)roll->nstrips());
194 
195  CLHEP::RandPoissonQ randPoissonQ(*engine, ave);
196  N_hits = randPoissonQ.fire();
197  for (int i = 0; i < N_hits; i++) {
198  double precise_time = CLHEP::RandFlat::shoot(engine, (nbxing * gate) / gate);
199  int time_hit = (static_cast<int>(precise_time)) - nbxing / 2;
200  RPCDigi adigi(j + 1, time_hit);
201  adigi.hasTime(true);
202  adigi.setTime(precise_time);
203  if (do_Y) {
204  double positionY = CLHEP::RandFlat::shoot(engine, striplength);
205  positionY -= striplength / 2;
206  adigi.hasY(true);
207  adigi.setY(positionY);
208  adigi.setDeltaY(sigmaY);
209  }
210  irpc_digis.insert(adigi);
211  }
212  }
213 }
214 
215 int RPCSimModelTiming::getClSize(uint32_t id, float posX, CLHEP::HepRandomEngine* engine) {
216  std::vector<double> clsForDetId = getRPCSimSetUp()->getCls(id);
217 
218  int cnt = 1;
219  int min = 1;
220  double func = 0.0;
221  std::vector<double> sum_clsize;
222 
223  sum_clsize.clear();
224  sum_clsize = clsForDetId;
225  int vectOffset(0);
226 
227  double rr_cl = CLHEP::RandFlat::shoot(engine);
228 
229  if (0.0 <= posX && posX < 0.2) {
230  func = clsForDetId[19] * (rr_cl);
231  vectOffset = 0;
232  }
233  if (0.2 <= posX && posX < 0.4) {
234  func = clsForDetId[39] * (rr_cl);
235  vectOffset = 20;
236  }
237  if (0.4 <= posX && posX < 0.6) {
238  func = clsForDetId[59] * (rr_cl);
239  vectOffset = 40;
240  }
241  if (0.6 <= posX && posX < 0.8) {
242  func = clsForDetId[79] * (rr_cl);
243  vectOffset = 60;
244  }
245  if (0.8 <= posX && posX < 1.0) {
246  func = clsForDetId[89] * (rr_cl);
247  vectOffset = 80;
248  }
249 
250  for (int i = vectOffset; i < (vectOffset + 20); i++) {
251  cnt++;
252  if (func > clsForDetId[i]) {
253  min = cnt;
254  } else if (func < clsForDetId[i]) {
255  break;
256  }
257  }
258  return min;
259 }
260 
261 int RPCSimModelTiming::LeftRightNeighbour(const RPCRoll& roll, const LocalPoint& hit_pos, int strip) {
262  //if left return -1
263  //if right return +1
264 
265  int leftStrip = strip - 1;
266  int rightStrip = strip + 1;
267 
268  if (leftStrip < 0)
269  return +1;
270  if (rightStrip > roll.nstrips())
271  return -1;
272 
273  double deltawL = fabs((roll.centreOfStrip(leftStrip)).x() - hit_pos.x());
274  double deltawR = fabs((roll.centreOfStrip(rightStrip)).x() - hit_pos.x());
275 
276  if (deltawL >= deltawR) {
277  return +1;
278  } else {
279  return -1;
280  }
281 }
T getParameter(std::string const &) const
float strip(const LocalPoint &lp) const
Definition: RPCRoll.cc:35
const Topology & topology() const override
Definition: RPCRollSpecs.cc:36
int LeftRightNeighbour(const RPCRoll &roll, const LocalPoint &hit_pos, int strip)
LocalPoint centreOfStrip(int strip) const
Definition: RPCRoll.cc:26
LocalPoint localPosition(float strip) const override
bool hasY() const
Definition: RPCDigi.h:34
CaloTopology const * topology(0)
DetectorHitMap theDetectorHitMap
Definition: RPCSim.h:68
void setY(double y)
Definition: RPCDigi.h:44
void setRPCSimSetUp(RPCSimSetUp *simsetup)
int nstrips() const
Definition: RPCRoll.cc:24
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
LocalPoint localPosition(float strip) const override
RPCSimSetUp * getRPCSimSetUp()
Definition: RPCSim.h:45
Definition: config.py:1
const RPCRollSpecs * specs() const
Definition: RPCRoll.cc:14
std::vector< double > clsForDetId
const Topology & topology() const override
Definition: RPCRoll.cc:18
~RPCSimModelTiming() override
void simulateNoise(const RPCRoll *, CLHEP::HepRandomEngine *) override
void simulate(const RPCRoll *roll, const edm::PSimHitContainer &rpcHits, CLHEP::HepRandomEngine *) override
edm::DetSet< RPCDigiSimLink > RPCDigiSimLinks
Definition: RPCSim.h:33
RPCDetId id() const
Definition: RPCRoll.cc:16
double getSmearedTime() const
std::vector< double > sum_clsize
const std::vector< double > & getCls(uint32_t id)
Definition: RPCSimSetUp.cc:421
float stripLength() const override
det heigth (strip length in the middle)
virtual int channel(const LocalPoint &p) const =0
Definition: RPCSim.h:30
T min(T a, T b)
Definition: MathUtil.h:58
int getClSize(uint32_t id, float posX, CLHEP::HepRandomEngine *)
RPCSimModelTiming(const edm::ParameterSet &config)
std::set< RPCDigi > irpc_digis
Definition: RPCSim.h:56
bool hasTime() const
Definition: RPCDigi.h:32
void setDeltaY(double dy)
Definition: RPCDigi.h:46
RPCSynchronizer * _rpcSync
void clear()
Definition: DetSet.h:72
std::vector< PSimHit > PSimHitContainer
T x() const
Definition: PV3DBase.h:59
RPCDigiSimLinks theRpcDigiSimLinks
Definition: RPCSim.h:70
void setTime(double time)
Definition: RPCDigi.h:41
int getSimHitBxAndTimingForIRPC(const PSimHit *, CLHEP::HepRandomEngine *)
int region() const
Region id: 0 for Barrel, +/-1 For +/- Endcap.
Definition: RPCDetId.h:53
float stripLength() const override