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