CMS 3D CMS Logo

RPCSynchronizer.cc
Go to the documentation of this file.
6 
11 
19 
20 #include "CLHEP/Random/RandGaussQ.h"
21 
22 #include <cstring>
23 #include <iostream>
24 #include <fstream>
25 #include <string>
26 #include <vector>
27 #include <cstdlib>
28 #include <cmath>
29 
30 using namespace std;
31 
33  resRPC = config.getParameter<double>("timeResolution");
34  timOff = config.getParameter<double>("timingRPCOffset");
35  dtimCs = config.getParameter<double>("deltatimeAdjacentStrip");
36  resEle = config.getParameter<double>("timeJitter");
37  sspeed = config.getParameter<double>("signalPropagationSpeed");
38  lbGate = config.getParameter<double>("linkGateWidth");
39  LHCGate = config.getParameter<double>("Gate");
40  cosmics = config.getParameter<bool>("cosmics");
41  irpc_timing_res = config.getParameter<double>("IRPC_time_resolution");
42  irpc_electronics_jitter = config.getParameter<double>("IRPC_electronics_jitter");
43  N_BX = config.getParameter<int>("BX_range");
44  //"magic" parameter for cosmics
45  cosmicPar = 37.62;
46 
47  double c = 299792458; // [m/s]
48  //light speed in [cm/ns]
49  cspeed = c * 1e+2 * 1e-9;
50  //signal propagation speed [cm/ns]
51  sspeed = sspeed * cspeed;
52 }
53 
55 
56 int RPCSynchronizer::getSimHitBx(const PSimHit* simhit, CLHEP::HepRandomEngine* engine) {
57  RPCSimSetUp* simsetup = this->getRPCSimSetUp();
58  const RPCGeometry* geometry = simsetup->getGeometry();
59  float timeref = simsetup->getTime(simhit->detUnitId());
60 
61  int bx = -999;
62  LocalPoint simHitPos = simhit->localPosition();
63  float tof = simhit->timeOfFlight();
64 
65  //automatic variable to prevent memory leak
66 
67  float rr_el = CLHEP::RandGaussQ::shoot(engine, 0., resEle);
68 
69  RPCDetId SimDetId(simhit->detUnitId());
70 
71  const RPCRoll* SimRoll = nullptr;
72 
73  for (TrackingGeometry::DetContainer::const_iterator it = geometry->dets().begin(); it != geometry->dets().end();
74  it++) {
75  if (dynamic_cast<const RPCChamber*>(*it) != nullptr) {
76  auto ch = dynamic_cast<const RPCChamber*>(*it);
77 
78  std::vector<const RPCRoll*> rollsRaf = (ch->rolls());
79  for (std::vector<const RPCRoll*>::iterator r = rollsRaf.begin(); r != rollsRaf.end(); ++r) {
80  if ((*r)->id() == SimDetId) {
81  SimRoll = &(*(*r));
82  break;
83  }
84  }
85  }
86  }
87 
88  if (SimRoll != nullptr) {
89  float distanceFromEdge = 0;
90  float half_stripL = 0.;
91 
92  if (SimRoll->id().region() == 0) {
93  const RectangularStripTopology* top_ = dynamic_cast<const RectangularStripTopology*>(&(SimRoll->topology()));
94  half_stripL = top_->stripLength() / 2;
95  distanceFromEdge = half_stripL + simHitPos.y();
96  } else {
97  const TrapezoidalStripTopology* top_ = dynamic_cast<const TrapezoidalStripTopology*>(&(SimRoll->topology()));
98  half_stripL = top_->stripLength() / 2;
99  distanceFromEdge = half_stripL - simHitPos.y();
100  }
101 
102  float prop_time = distanceFromEdge / sspeed;
103 
104  double rr_tim1 = CLHEP::RandGaussQ::shoot(engine, 0., resRPC);
105  double total_time = tof + prop_time + timOff + rr_tim1 + rr_el;
106 
107  // Bunch crossing assignment
108  double time_differ = 0.;
109 
110  if (cosmics) {
111  time_differ = (total_time - (timeref + ((half_stripL / sspeed) + timOff))) / cosmicPar;
112  } else if (!cosmics) {
113  time_differ = total_time - (timeref + (half_stripL / sspeed) + timOff);
114  }
115 
116  double inf_time = 0;
117  double sup_time = 0;
118 
119  for (int n = -N_BX; n <= N_BX; ++n) {
120  if (cosmics) {
121  inf_time = (-lbGate / 2 + n * LHCGate) / cosmicPar;
122  sup_time = (lbGate / 2 + n * LHCGate) / cosmicPar;
123  } else if (!cosmics) {
124  inf_time = -lbGate / 2 + n * LHCGate;
125  sup_time = lbGate / 2 + n * LHCGate;
126  }
127 
128  if (inf_time < time_differ && time_differ < sup_time) {
129  bx = n;
130  break;
131  }
132  }
133  }
134 
135  return bx;
136 }
137 
138 int RPCSynchronizer::getSimHitBxAndTimingForIRPC(const PSimHit* simhit, CLHEP::HepRandomEngine* engine) {
139  RPCSimSetUp* simsetup = this->getRPCSimSetUp();
140  const RPCGeometry* geometry = simsetup->getGeometry();
141  float timeref = simsetup->getTime(simhit->detUnitId());
142 
143  int bx = -999;
144  LocalPoint simHitPos = simhit->localPosition();
145  float tof = simhit->timeOfFlight();
146 
147  //automatic variable to prevent memory leak
148 
149  // float rr_el = CLHEP::RandGaussQ::shoot(engine, 0.,resEle);
150  float rr_el = CLHEP::RandGaussQ::shoot(engine, 0., irpc_electronics_jitter);
151 
152  RPCDetId SimDetId(simhit->detUnitId());
153 
154  const RPCRoll* SimRoll = nullptr;
155 
156  for (TrackingGeometry::DetContainer::const_iterator it = geometry->dets().begin(); it != geometry->dets().end();
157  it++) {
158  if (dynamic_cast<const RPCChamber*>(*it) != nullptr) {
159  auto ch = dynamic_cast<const RPCChamber*>(*it);
160 
161  std::vector<const RPCRoll*> rollsRaf = (ch->rolls());
162  for (std::vector<const RPCRoll*>::iterator r = rollsRaf.begin(); r != rollsRaf.end(); ++r) {
163  if ((*r)->id() == SimDetId) {
164  SimRoll = &(*(*r));
165  break;
166  }
167  }
168  }
169  }
170 
171  if (SimRoll != nullptr) {
172  float distanceFromEdge = 0;
173  float half_stripL = 0.;
174 
175  if (SimRoll->id().region() == 0) {
176  const RectangularStripTopology* top_ = dynamic_cast<const RectangularStripTopology*>(&(SimRoll->topology()));
177  half_stripL = top_->stripLength() / 2;
178  distanceFromEdge = half_stripL + simHitPos.y();
179  } else {
180  const TrapezoidalStripTopology* top_ = dynamic_cast<const TrapezoidalStripTopology*>(&(SimRoll->topology()));
181  half_stripL = top_->stripLength() / 2;
182  distanceFromEdge = half_stripL - simHitPos.y();
183  }
184 
185  float prop_time = distanceFromEdge / sspeed;
186 
187  // double rr_tim1 = CLHEP::RandGaussQ::shoot(engine, 0.,resRPC);
188  double rr_tim1 = CLHEP::RandGaussQ::shoot(engine, 0., irpc_timing_res);
189 
190  double total_time = tof + prop_time + timOff + rr_tim1 + rr_el;
191 
192  // Bunch crossing assignment
193  double time_differ = 0.;
194 
195  if (cosmics) {
196  time_differ = (total_time - (timeref + ((half_stripL / sspeed) + timOff))) / cosmicPar;
197  } else if (!cosmics) {
198  time_differ = total_time - (timeref + (half_stripL / sspeed) + timOff);
199  }
200 
201  double exact_total_time = tof + prop_time + timOff;
202  double exact_time_differ = 0.;
203 
204  if (cosmics) {
205  exact_time_differ = (exact_total_time - (timeref + ((half_stripL / sspeed) + timOff))) / cosmicPar;
206  } else if (!cosmics) {
207  exact_time_differ = exact_total_time - (timeref + (half_stripL / sspeed) + timOff);
208  }
209 
210  double inf_time = 0;
211  double sup_time = 0;
212 
213  for (int n = -N_BX; n <= N_BX; ++n) {
214  if (cosmics) {
215  inf_time = (-lbGate / 2 + n * LHCGate) / cosmicPar;
216  sup_time = (lbGate / 2 + n * LHCGate) / cosmicPar;
217  } else if (!cosmics) {
218  inf_time = -lbGate / 2 + n * LHCGate;
219  sup_time = lbGate / 2 + n * LHCGate;
220  }
221 
222  if (inf_time < time_differ && time_differ < sup_time) {
223  bx = n;
224  the_exact_time = exact_time_differ;
225  the_smeared_time = time_differ;
226 
227  //cout<<"Debug\t"<<inf_time<<'\t'<<sup_time<<endl;
229  // cout<<"Bingo\t"<<time_differ<<'\t'<<bx<<'\t'<<exact_time_differ<<'\t'<<exact_time_differ-time_differ<<'\t'<<exact_time_differ-bx*25.<<endl;
230  break;
231  }
232  }
233  }
234  return bx;
235 }
unsigned int detUnitId() const
Definition: PSimHit.h:97
float getTime(uint32_t id)
Definition: RPCSimSetUp.cc:444
Definition: config.py:1
RPCSynchronizer(const edm::ParameterSet &config)
int getSimHitBx(const PSimHit *, CLHEP::HepRandomEngine *)
float stripLength() const override
T y() const
Definition: PV3DBase.h:60
const RPCGeometry * getGeometry()
Definition: RPCSimSetUp.h:48
const unsigned int N_BX
Definition: ConstantDef.h:5
Local3DPoint localPosition() const
Definition: PSimHit.h:52
float timeOfFlight() const
Definition: PSimHit.h:73
float stripLength() const override
det heigth (strip length in the middle)
int getSimHitBxAndTimingForIRPC(const PSimHit *, CLHEP::HepRandomEngine *)