CMS 3D CMS Logo

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