CMS 3D CMS Logo

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