CMS 3D CMS Logo

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