CMS 3D CMS Logo

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