CMS 3D CMS Logo

RecHitProcessor.cc
Go to the documentation of this file.
2 
4 
6 
8  const edm::EventSetup &iSetup,
9  const edm::EDGetToken &recHitToken,
10  std::vector<RecHitProcessor::CppfItem> &CppfVec1,
11  l1t::CPPFDigiCollection &cppfDigis,
12  const int MaxClusterSize) const {
14  iEvent.getByToken(recHitToken, recHits);
15 
17  iSetup.get<MuonGeometryRecord>().get(rpcGeom);
18 
19  // The loop is over the detector container in the rpc geometry collection. We
20  // are interested in the RPDdetID (inside of RPCChamber vectors),
21  // specifically, the RPCrechits. to assignment the CPPFDigis.
22  for (TrackingGeometry::DetContainer::const_iterator iDet = rpcGeom->dets().begin(); iDet < rpcGeom->dets().end();
23  iDet++) {
24  // we do a cast over the class RPCChamber to obtain the RPCroll vectors,
25  // inside of them, the RPCRechits are found. in other words, the method
26  // ->rolls() does not exist for other kind of vector within DetContainer
27  // and we can not obtain the rpcrechits in a suitable way.
28  if (dynamic_cast<const RPCChamber *>(*iDet) == nullptr)
29  continue;
30 
31  auto chamb = dynamic_cast<const RPCChamber *>(*iDet);
32 
33  std::vector<const RPCRoll *> rolls = (chamb->rolls());
34 
35  // Loop over rolls in the chamber
36  for (auto &iRoll : rolls) {
37  RPCDetId rpcId = (*iRoll).id();
38 
39  typedef std::pair<RPCRecHitCollection::const_iterator, RPCRecHitCollection::const_iterator> rangeRecHits;
40  rangeRecHits recHitCollection = recHits->get(rpcId);
41 
42  // Loop over the RPC digis
43  for (RPCRecHitCollection::const_iterator rechit_it = recHitCollection.first; rechit_it != recHitCollection.second;
44  rechit_it++) {
45  // const RPCDetId& rpcId = rechit_it->rpcId();
46  int rawId = rpcId.rawId();
47  // int station = rpcId.station();
48  int Bx = rechit_it->BunchX();
49  int isValid = rechit_it->isValid();
50  int firststrip = rechit_it->firstClusterStrip();
51  int clustersize = rechit_it->clusterSize();
52  LocalPoint lPos = rechit_it->localPosition();
53  const RPCRoll *roll = rpcGeom->roll(rpcId);
54  const BoundPlane &rollSurface = roll->surface();
55  GlobalPoint gPos = rollSurface.toGlobal(lPos);
56  float global_theta = emtf::rad_to_deg(gPos.theta().value());
57  float global_phi = emtf::rad_to_deg(gPos.phi().value());
58  //::::::::::::::::::::::::::::
59  // Establish the average position of the rechit
60  int rechitstrip = firststrip;
61 
62  if (clustersize > 2) {
63  int medium = 0;
64  if (clustersize % 2 == 0)
65  medium = 0.5 * (clustersize);
66  else
67  medium = 0.5 * (clustersize - 1);
68  rechitstrip += medium;
69  }
70 
71  if (clustersize > MaxClusterSize)
72  continue;
73  // This is just for test CPPFDigis with the RPC Geometry, It must be
74  // "true" in the normal runs
75  bool Geo = true;
77  // Set the EMTF Sector
78  int EMTFsector1 = 0;
79  int EMTFsector2 = 0;
80 
81  // sector 1
82  if ((global_phi > 15.) && (global_phi <= 16.3)) {
83  EMTFsector1 = 1;
84  EMTFsector2 = 6;
85  } else if ((global_phi > 16.3) && (global_phi <= 53.)) {
86  EMTFsector1 = 1;
87  EMTFsector2 = 0;
88  } else if ((global_phi > 53.) && (global_phi <= 75.)) {
89  EMTFsector1 = 1;
90  EMTFsector2 = 2;
91  }
92  // sector 2
93  else if ((global_phi > 75.) && (global_phi <= 76.3)) {
94  EMTFsector1 = 1;
95  EMTFsector2 = 2;
96  } else if ((global_phi > 76.3) && (global_phi <= 113.)) {
97  EMTFsector1 = 2;
98  EMTFsector2 = 0;
99  } else if ((global_phi > 113.) && (global_phi <= 135.)) {
100  EMTFsector1 = 2;
101  EMTFsector2 = 3;
102  }
103  // sector 3
104  // less than 180
105  else if ((global_phi > 135.) && (global_phi <= 136.3)) {
106  EMTFsector1 = 2;
107  EMTFsector2 = 3;
108  } else if ((global_phi > 136.3) && (global_phi <= 173.)) {
109  EMTFsector1 = 3;
110  EMTFsector2 = 0;
111  } else if ((global_phi > 173.) && (global_phi <= 180.)) {
112  EMTFsector1 = 3;
113  EMTFsector2 = 4;
114  }
115  // Greater than -180
116  else if ((global_phi < -165.) && (global_phi >= -180.)) {
117  EMTFsector1 = 3;
118  EMTFsector2 = 4;
119  }
120  // Fourth sector
121  else if ((global_phi > -165.) && (global_phi <= -163.7)) {
122  EMTFsector1 = 3;
123  EMTFsector2 = 4;
124  } else if ((global_phi > -163.7) && (global_phi <= -127.)) {
125  EMTFsector1 = 4;
126  EMTFsector2 = 0;
127  } else if ((global_phi > -127.) && (global_phi <= -105.)) {
128  EMTFsector1 = 4;
129  EMTFsector2 = 5;
130  }
131  // fifth sector
132  else if ((global_phi > -105.) && (global_phi <= -103.7)) {
133  EMTFsector1 = 4;
134  EMTFsector2 = 5;
135  } else if ((global_phi > -103.7) && (global_phi <= -67.)) {
136  EMTFsector1 = 5;
137  EMTFsector2 = 0;
138  } else if ((global_phi > -67.) && (global_phi <= -45.)) {
139  EMTFsector1 = 5;
140  EMTFsector2 = 6;
141  }
142  // sixth sector
143  else if ((global_phi > -45.) && (global_phi <= -43.7)) {
144  EMTFsector1 = 5;
145  EMTFsector2 = 6;
146  } else if ((global_phi > -43.7) && (global_phi <= -7.)) {
147  EMTFsector1 = 6;
148  EMTFsector2 = 0;
149  } else if ((global_phi > -7.) && (global_phi <= 15.)) {
150  EMTFsector1 = 6;
151  EMTFsector2 = 1;
152  }
153 
154  // std::vector<RecHitProcessor::CppfItem>::iterator it;
155  // for(it = CppfVec1.begin(); it != CppfVec1.end(); it++){
156  // if( (*it).rawId == rawId) if(Geo_true) std::cout << (*it).rawId
157  //<< "rawid" << rawId << std::endl;
158  // }
159  // Loop over the look up table
160  double EMTFLink1 = 0.;
161  double EMTFLink2 = 0.;
162 
163  std::vector<RecHitProcessor::CppfItem>::iterator cppf1;
164  std::vector<RecHitProcessor::CppfItem>::iterator cppf;
165  for (cppf1 = CppfVec1.begin(); cppf1 != CppfVec1.end(); cppf1++) {
166  // Condition to save the CPPFDigi
167  if (((*cppf1).rawId == rawId) && ((*cppf1).strip == rechitstrip)) {
168  int old_strip = (*cppf1).strip;
169  int before = 0;
170  int after = 0;
171 
172  if (cppf1 != CppfVec1.begin())
173  before = (*(cppf1 - 2)).strip;
174 
175  else if (cppf1 == CppfVec1.begin())
176  before = (*cppf1).strip;
177 
178  if (cppf1 != CppfVec1.end())
179  after = (*(cppf1 + 2)).strip;
180 
181  else if (cppf1 == CppfVec1.end())
182  after = (*cppf1).strip;
183 
184  cppf = cppf1;
185 
186  if (clustersize == 2) {
187  if (firststrip == 1) {
188  if (before < after)
189  cppf = (cppf1 - 1);
190  else if (before > after)
191  cppf = (cppf1 + 1);
192  } else if (firststrip > 1) {
193  if (before < after)
194  cppf = (cppf1 + 1);
195  else if (before > after)
196  cppf = (cppf1 - 1);
197  }
198  }
199  // Using the RPCGeometry
200  if (Geo) {
201  std::shared_ptr<l1t::CPPFDigi> MainVariables1(new l1t::CPPFDigi(rpcId,
202  Bx,
203  (*cppf).int_phi,
204  (*cppf).int_theta,
205  isValid,
206  (*cppf).lb,
207  (*cppf).halfchannel,
208  EMTFsector1,
209  EMTFLink1,
210  old_strip,
211  clustersize,
212  global_phi,
213  global_theta));
214  std::shared_ptr<l1t::CPPFDigi> MainVariables2(new l1t::CPPFDigi(rpcId,
215  Bx,
216  (*cppf).int_phi,
217  (*cppf).int_theta,
218  isValid,
219  (*cppf).lb,
220  (*cppf).halfchannel,
221  EMTFsector2,
222  EMTFLink2,
223  old_strip,
224  clustersize,
225  global_phi,
226  global_theta));
227 
228  if ((EMTFsector1 > 0) && (EMTFsector2 == 0)) {
229  cppfDigis.push_back(*MainVariables1.get());
230  } else if ((EMTFsector1 > 0) && (EMTFsector2 > 0)) {
231  cppfDigis.push_back(*MainVariables1.get());
232  cppfDigis.push_back(*MainVariables2.get());
233  } else if ((EMTFsector1 == 0) && (EMTFsector2 == 0)) {
234  continue;
235  }
236  } // Geo is true
237  else {
238  global_phi = 0.;
239  global_theta = 0.;
240  std::shared_ptr<l1t::CPPFDigi> MainVariables1(new l1t::CPPFDigi(rpcId,
241  Bx,
242  (*cppf).int_phi,
243  (*cppf).int_theta,
244  isValid,
245  (*cppf).lb,
246  (*cppf).halfchannel,
247  EMTFsector1,
248  EMTFLink1,
249  old_strip,
250  clustersize,
251  global_phi,
252  global_theta));
253  std::shared_ptr<l1t::CPPFDigi> MainVariables2(new l1t::CPPFDigi(rpcId,
254  Bx,
255  (*cppf).int_phi,
256  (*cppf).int_theta,
257  isValid,
258  (*cppf).lb,
259  (*cppf).halfchannel,
260  EMTFsector2,
261  EMTFLink2,
262  old_strip,
263  clustersize,
264  global_phi,
265  global_theta));
266  if ((EMTFsector1 > 0) && (EMTFsector2 == 0)) {
267  cppfDigis.push_back(*MainVariables1.get());
268  } else if ((EMTFsector1 > 0) && (EMTFsector2 > 0)) {
269  cppfDigis.push_back(*MainVariables1.get());
270  cppfDigis.push_back(*MainVariables2.get());
271  } else if ((EMTFsector1 == 0) && (EMTFsector2 == 0)) {
272  continue;
273  }
274  }
275  } // Condition to save the CPPFDigi
276  } // Loop over the LUTVector
277  } // Loop over the recHits
278  } // End loop: for (std::vector<const RPCRoll*>::const_iterator r =
279  // rolls.begin(); r != rolls.end(); ++r)
280  } // End loop: for (TrackingGeometry::DetContainer::const_iterator iDet =
281  // rpcGeom->dets().begin(); iDet < rpcGeom->dets().end(); iDet++)
282 }
283 
285  const edm::EventSetup &iSetup,
286  const edm::EDGetToken &recHitToken,
287  l1t::CPPFDigiCollection &cppfDigis) const {
288  // Get the RPC Geometry
290  iSetup.get<MuonGeometryRecord>().get(rpcGeom);
291 
292  // Get the RecHits from the event
294  iEvent.getByToken(recHitToken, recHits);
295 
296  // The loop is over the detector container in the rpc geometry collection. We
297  // are interested in the RPDdetID (inside of RPCChamber vectors),
298  // specifically, the RPCrechits. to assignment the CPPFDigis.
299  for (TrackingGeometry::DetContainer::const_iterator iDet = rpcGeom->dets().begin(); iDet < rpcGeom->dets().end();
300  iDet++) {
301  // we do a cast over the class RPCChamber to obtain the RPCroll vectors,
302  // inside of them, the RPCRechits are found. in other words, the method
303  // ->rolls() does not exist for other kind of vector within DetContainer
304  // and we can not obtain the rpcrechits in a suitable way.
305  if (dynamic_cast<const RPCChamber *>(*iDet) == nullptr)
306  continue;
307 
308  auto chamb = dynamic_cast<const RPCChamber *>(*iDet);
309  std::vector<const RPCRoll *> rolls = (chamb->rolls());
310 
311  // Loop over rolls in the chamber
312  for (auto &iRoll : rolls) {
313  RPCDetId rpcId = (*iRoll).id();
314 
315  typedef std::pair<RPCRecHitCollection::const_iterator, RPCRecHitCollection::const_iterator> rangeRecHits;
316  rangeRecHits recHitCollection = recHits->get(rpcId);
317 
318  for (RPCRecHitCollection::const_iterator rechit_it = recHitCollection.first; rechit_it != recHitCollection.second;
319  rechit_it++) {
320  // const RPCDetId& rpcId = rechit_it->rpcId();
321  // int rawId = rpcId.rawId();
322  int region = rpcId.region();
323  // int station = rpcId.station();
324  int Bx = rechit_it->BunchX();
325  int isValid = rechit_it->isValid();
326  int firststrip = rechit_it->firstClusterStrip();
327  int clustersize = rechit_it->clusterSize();
328  LocalPoint lPos = rechit_it->localPosition();
329  const RPCRoll *roll = rpcGeom->roll(rpcId);
330  const BoundPlane &rollSurface = roll->surface();
331  GlobalPoint gPos = rollSurface.toGlobal(lPos);
332  float global_theta = emtf::rad_to_deg(gPos.theta().value());
333  float global_phi = emtf::rad_to_deg(gPos.phi().value());
334  // Endcap region only
335 
336  if (region != 0) {
337  int int_theta =
338  (region == -1 ? 180. * 32. / 36.5 : 0.) + (float)region * global_theta * 32. / 36.5 - 8.5 * 32 / 36.5;
339 
340  if (region == 1) {
341  if (global_theta < 8.5)
342  int_theta = 0;
343  if (global_theta > 45.)
344  int_theta = 31;
345  } else if (region == -1) {
346  if (global_theta < 135.)
347  int_theta = 31;
348  if (global_theta > 171.5)
349  int_theta = 0;
350  }
351 
352  // Local EMTF
353  double local_phi = 0.;
354  int EMTFsector1 = 0;
355  int EMTFsector2 = 0;
356 
357  // sector 1
358  if ((global_phi > 15.) && (global_phi <= 16.3)) {
359  local_phi = global_phi - 15.;
360  EMTFsector1 = 1;
361  EMTFsector2 = 6;
362  } else if ((global_phi > 16.3) && (global_phi <= 53.)) {
363  local_phi = global_phi - 15.;
364  EMTFsector1 = 1;
365  EMTFsector2 = 0;
366  } else if ((global_phi > 53.) && (global_phi <= 75.)) {
367  local_phi = global_phi - 15.;
368  EMTFsector1 = 1;
369  EMTFsector2 = 2;
370  }
371  // sector 2
372  else if ((global_phi > 75.) && (global_phi <= 76.3)) {
373  local_phi = global_phi - 15.;
374  EMTFsector1 = 1;
375  EMTFsector2 = 2;
376  } else if ((global_phi > 76.3) && (global_phi <= 113.)) {
377  local_phi = global_phi - 75.;
378  EMTFsector1 = 2;
379  EMTFsector2 = 0;
380  } else if ((global_phi > 113.) && (global_phi <= 135.)) {
381  local_phi = global_phi - 75.;
382  EMTFsector1 = 2;
383  EMTFsector2 = 3;
384  }
385  // sector 3
386  // less than 180
387  else if ((global_phi > 135.) && (global_phi <= 136.3)) {
388  local_phi = global_phi - 75.;
389  EMTFsector1 = 2;
390  EMTFsector2 = 3;
391  } else if ((global_phi > 136.3) && (global_phi <= 173.)) {
392  local_phi = global_phi - 135.;
393  EMTFsector1 = 3;
394  EMTFsector2 = 0;
395  } else if ((global_phi > 173.) && (global_phi <= 180.)) {
396  local_phi = global_phi - 135.;
397  EMTFsector1 = 3;
398  EMTFsector2 = 4;
399  }
400  // Greater than -180
401  else if ((global_phi < -165.) && (global_phi >= -180.)) {
402  local_phi = global_phi + 225.;
403  EMTFsector1 = 3;
404  EMTFsector2 = 4;
405  }
406  // Fourth sector
407  else if ((global_phi > -165.) && (global_phi <= -163.7)) {
408  local_phi = global_phi + 225.;
409  EMTFsector1 = 3;
410  EMTFsector2 = 4;
411  } else if ((global_phi > -163.7) && (global_phi <= -127.)) {
412  local_phi = global_phi + 165.;
413  EMTFsector1 = 4;
414  EMTFsector2 = 0;
415  } else if ((global_phi > -127.) && (global_phi <= -105.)) {
416  local_phi = global_phi + 165.;
417  EMTFsector1 = 4;
418  EMTFsector2 = 5;
419  }
420  // fifth sector
421  else if ((global_phi > -105.) && (global_phi <= -103.7)) {
422  local_phi = global_phi + 165.;
423  EMTFsector1 = 4;
424  EMTFsector2 = 5;
425  } else if ((global_phi > -103.7) && (global_phi <= -67.)) {
426  local_phi = global_phi + 105.;
427  EMTFsector1 = 5;
428  EMTFsector2 = 0;
429  } else if ((global_phi > -67.) && (global_phi <= -45.)) {
430  local_phi = global_phi + 105.;
431  EMTFsector1 = 5;
432  EMTFsector2 = 6;
433  }
434  // sixth sector
435  else if ((global_phi > -45.) && (global_phi <= -43.7)) {
436  local_phi = global_phi + 105.;
437  EMTFsector1 = 5;
438  EMTFsector2 = 6;
439  } else if ((global_phi > -43.7) && (global_phi <= -7.)) {
440  local_phi = global_phi + 45.;
441  EMTFsector1 = 6;
442  EMTFsector2 = 0;
443  } else if ((global_phi > -7.) && (global_phi <= 15.)) {
444  local_phi = global_phi + 45.;
445  EMTFsector1 = 6;
446  EMTFsector2 = 1;
447  }
448 
449  int int_phi = int((local_phi + 22.0) * 15. + .5);
450 
451  double EMTFLink1 = 0.;
452  double EMTFLink2 = 0.;
453  double lb = 0.;
454  double halfchannel = 0.;
455 
456  // Invalid hit
457  if (isValid == 0)
458  int_phi = 2047;
459  // Right integers range
460  assert(0 <= int_phi && int_phi < 1250);
461  assert(0 <= int_theta && int_theta < 32);
462 
463  std::shared_ptr<l1t::CPPFDigi> MainVariables1(new l1t::CPPFDigi(rpcId,
464  Bx,
465  int_phi,
466  int_theta,
467  isValid,
468  lb,
469  halfchannel,
470  EMTFsector1,
471  EMTFLink1,
472  firststrip,
473  clustersize,
474  global_phi,
475  global_theta));
476  std::shared_ptr<l1t::CPPFDigi> MainVariables2(new l1t::CPPFDigi(rpcId,
477  Bx,
478  int_phi,
479  int_theta,
480  isValid,
481  lb,
482  halfchannel,
483  EMTFsector2,
484  EMTFLink2,
485  firststrip,
486  clustersize,
487  global_phi,
488  global_theta));
489  if (int_theta == 31)
490  continue;
491  if ((EMTFsector1 > 0) && (EMTFsector2 == 0)) {
492  cppfDigis.push_back(*MainVariables1.get());
493  }
494  if ((EMTFsector1 > 0) && (EMTFsector2 > 0)) {
495  cppfDigis.push_back(*MainVariables1.get());
496  cppfDigis.push_back(*MainVariables2.get());
497  }
498  if ((EMTFsector1 == 0) && (EMTFsector2 == 0)) {
499  continue;
500  }
501  } // No barrel rechits
502 
503  } // End loop: for (RPCRecHitCollection::const_iterator recHit =
504  // recHitCollection.first; recHit != recHitCollection.second; recHit++)
505 
506  } // End loop: for (std::vector<const RPCRoll*>::const_iterator r =
507  // rolls.begin(); r != rolls.end(); ++r)
508  } // End loop: for (TrackingGeometry::DetContainer::const_iterator iDet =
509  // rpcGeom->dets().begin(); iDet < rpcGeom->dets().end(); iDet++)
510 } // End function: void RecHitProcessor::process()
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
T1 value() const
Explicit access to value in case implicit conversion not OK.
Definition: Phi.h:75
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
Geom::Theta< T > theta() const
Definition: PV3DBase.h:72
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:43
int iEvent
Definition: GenABIO.cc:224
void processLook(const edm::Event &iEvent, const edm::EventSetup &iSetup, const edm::EDGetToken &recHitToken, std::vector< RecHitProcessor::CppfItem > &CppfVec1, l1t::CPPFDigiCollection &cppfDigis, const int MaxClusterSize) const
const DetContainer & dets() const override
Returm a vector of all GeomDet (including all GeomDetUnits)
Definition: RPCGeometry.cc:24
void process(const edm::Event &iEvent, const edm::EventSetup &iSetup, const edm::EDGetToken &recHitToken, l1t::CPPFDigiCollection &cppfDigis) const
T get() const
Definition: EventSetup.h:73
double rad_to_deg(double rad)
Definition: TrackTools.h:62
std::vector< CPPFDigi > CPPFDigiCollection
Definition: CPPFDigi.h:83
const RPCRoll * roll(RPCDetId id) const
Return a roll given its id.
Definition: RPCGeometry.cc:50
int region() const
Region id: 0 for Barrel, +/-1 For +/- Endcap.
Definition: RPCDetId.h:53