CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
L1AnalysisCSCTF.cc
Go to the documentation of this file.
2 
4 
6 
8  const L1MuTriggerScales* ts,
9  const L1MuTriggerPtScale* tpts,
10  CSCSectorReceiverLUT* srLUTs_[5][2],
11  CSCTFPtLUT* ptLUTs_) {
12  //for (int i =0; i<MAXCSCTFTR; i++) csctf_.trSector[i]=0;
13 
14  csctf_.trSize = csctfTrks->size();
15  //cout << " csctf_.trSize: " << csctf_.trSize << endl;
16 
17  int nTrk = 0;
18  for (L1CSCTrackCollection::const_iterator trk = csctfTrks->begin(); trk < csctfTrks->end(); trk++) {
19  nTrk++;
20 
21  // trk->first.endcap() = 2 for - endcap
22  // = 1 for + endcap
23  csctf_.trEndcap.push_back(trk->first.endcap() == 2 ? trk->first.endcap() - 3 : trk->first.endcap());
24  //sectors: 1->6 (plus endcap), 7->12 (minus endcap)
25  csctf_.trSector.push_back(6 * (trk->first.endcap() - 1) + trk->first.sector());
26  csctf_.trBx.push_back(trk->first.BX());
27 
28  csctf_.trME1ID.push_back(trk->first.me1ID());
29  csctf_.trME2ID.push_back(trk->first.me2ID());
30  csctf_.trME3ID.push_back(trk->first.me3ID());
31  csctf_.trME4ID.push_back(trk->first.me4ID());
32  csctf_.trMB1ID.push_back(trk->first.mb1ID());
33 
34  csctf_.trME1TBin.push_back(trk->first.me1Tbin());
35  csctf_.trME2TBin.push_back(trk->first.me2Tbin());
36  csctf_.trME3TBin.push_back(trk->first.me3Tbin());
37  csctf_.trME4TBin.push_back(trk->first.me4Tbin());
38  csctf_.trMB1TBin.push_back(trk->first.mb1Tbin());
39 
40  // output link of the track
41  csctf_.trOutputLink.push_back(trk->first.outputLink());
42 
43  // PtAddress gives an handle on other parameters
44  ptadd thePtAddress(trk->first.ptLUTAddress());
45 
46  csctf_.trPhiSign.push_back(thePtAddress.delta_phi_sign);
47  csctf_.trPhi12.push_back(thePtAddress.delta_phi_12);
48  csctf_.trPhi23.push_back(thePtAddress.delta_phi_23);
49  csctf_.trMode.push_back(thePtAddress.track_mode);
50  csctf_.trForR.push_back(thePtAddress.track_fr);
51  csctf_.trCharge.push_back(trk->first.chargeValue());
52 
53  //Pt needs some more workaround since it is not in the unpacked data
55  ptdat thePtData = ptLUTs_->Pt(thePtAddress);
56 
57  // front or rear bit?
58  if (thePtAddress.track_fr) {
59  csctf_.trPtBit.push_back(thePtData.front_rank & 0x1f);
60  csctf_.trQuality.push_back((thePtData.front_rank >> 5) & 0x3);
61  csctf_.trChargeValid.push_back(thePtData.charge_valid_front);
62  } else {
63  csctf_.trPtBit.push_back(thePtData.rear_rank & 0x1f);
64  csctf_.trQuality.push_back((thePtData.rear_rank >> 5) & 0x3);
65  csctf_.trChargeValid.push_back(thePtData.charge_valid_rear);
66  }
67 
68  // convert the Pt in human readable values (GeV/c)
69  csctf_.trPt.push_back(tpts->getPtScale()->getLowEdge(csctf_.trPtBit.back()));
70 
71  // track's Eta and Phi in bit...
72  csctf_.trEtaBit.push_back(trk->first.eta_packed());
73  csctf_.trPhiBit.push_back(trk->first.localPhi());
74 
75  //... in radians
76  // Type 2 is CSC
77  //csctf_.trEta.push_back(theTriggerScales->getRegionalEtaScale(2)->getCenter(trk->first.eta_packed()));
78  //csctf_.trPhi.push_back(theTriggerScales->getPhiScale()->getLowEdge(trk->first.localPhi()));
79  csctf_.trEta.push_back(ts->getRegionalEtaScale(2)->getCenter(trk->first.eta_packed()));
80  csctf_.trPhi.push_back(ts->getPhiScale()->getLowEdge(trk->first.localPhi()));
81  //Phi in one sector varies from [0,62] degrees -> Rescale manually to global coords.
82  csctf_.trPhi_02PI.push_back(
83  fmod(csctf_.trPhi[nTrk - 1] + ((trk->first.sector() - 1) * TMath::Pi() / 3) + //sector 1 starts at 15 degrees
84  (TMath::Pi() / 12),
85  2 * TMath::Pi()));
86 
87  // //csctf lcts of tracks
88  // if( csctfLCTSource_.label() != "none" ){
89 
90  // For each trk, get the list of its LCTs
91  CSCCorrelatedLCTDigiCollection lctsOfTracks = trk->second;
92 
93  int LctTrkId_ = 0;
94 
95  for (CSCCorrelatedLCTDigiCollection::DigiRangeIterator lctOfTrks = lctsOfTracks.begin();
96  lctOfTrks != lctsOfTracks.end();
97  lctOfTrks++) {
98  int lctTrkId = 0;
99 
100  CSCCorrelatedLCTDigiCollection::Range lctRange = lctsOfTracks.get((*lctOfTrks).first);
101 
102  for (CSCCorrelatedLCTDigiCollection::const_iterator lctTrk = lctRange.first; lctTrk != lctRange.second;
103  lctTrk++, lctTrkId++) {
104  csctf_.trLctEndcap[nTrk - 1][LctTrkId_] = (*lctOfTrks).first.zendcap();
105  if ((*lctOfTrks).first.zendcap() > 0)
106  csctf_.trLctSector[nTrk - 1][LctTrkId_] = (*lctOfTrks).first.triggerSector();
107  else
108  csctf_.trLctSector[nTrk - 1][LctTrkId_] = 6 + (*lctOfTrks).first.triggerSector();
109  csctf_.trLctSubSector[nTrk - 1][LctTrkId_] =
111  ;
112  csctf_.trLctBx[nTrk - 1][LctTrkId_] = lctTrk->getBX();
113  csctf_.trLctBx0[nTrk - 1][LctTrkId_] = lctTrk->getBX0();
114 
115  csctf_.trLctStation[nTrk - 1][LctTrkId_] = (*lctOfTrks).first.station();
116  csctf_.trLctRing[nTrk - 1][LctTrkId_] = (*lctOfTrks).first.ring();
117  csctf_.trLctChamber[nTrk - 1][LctTrkId_] = (*lctOfTrks).first.chamber();
118  csctf_.trLctTriggerCSCID[nTrk - 1][LctTrkId_] = (*lctOfTrks).first.triggerCscId();
119  csctf_.trLctFpga[nTrk - 1][LctTrkId_] =
120  (csctf_.trLctSubSector[nTrk - 1][LctTrkId_] ? csctf_.trLctSubSector[nTrk - 1][LctTrkId_]
121  : (*lctOfTrks).first.station() + 1);
122 
123  // Check if DetId is within range
124  if (csctf_.trLctSector[nTrk - 1][LctTrkId_] < 1 || csctf_.trLctSector[nTrk - 1][LctTrkId_] > 12 ||
125  csctf_.trLctStation[nTrk - 1][LctTrkId_] < 1 || csctf_.trLctStation[nTrk - 1][LctTrkId_] > 4 ||
126  csctf_.trLctTriggerCSCID[nTrk - 1][LctTrkId_] < 1 || csctf_.trLctTriggerCSCID[nTrk - 1][LctTrkId_] > 9 ||
127  lctTrkId < 0 || lctTrkId > 1) {
128  edm::LogInfo("L1NtupleProducer") << " TRACK ERROR: CSC digi are out of range: ";
129 
130  continue;
131  }
132 
133  // handles not to overload the method: mostly for readability
134  int endcap = (*lctOfTrks).first.zendcap();
135  if (endcap < 0)
136  endcap = 0;
137 
138  int StationLctTrk = (*lctOfTrks).first.station();
139  int CscIdLctTrk = (*lctOfTrks).first.triggerCscId();
140  int SubSectorLctTrk = CSCTriggerNumbering::triggerSubSectorFromLabels((*lctOfTrks).first);
141 
142  int FPGALctTrk = (SubSectorLctTrk ? SubSectorLctTrk - 1 : StationLctTrk);
143 
144  // local Phi
145  lclphidat lclPhi;
146 
147  try {
148  csctf_.trLctstripNum[nTrk - 1][LctTrkId_] = lctTrk->getStrip();
149  lclPhi = srLUTs_[FPGALctTrk][endcap]->localPhi(
150  lctTrk->getStrip(), lctTrk->getPattern(), lctTrk->getQuality(), lctTrk->getBend());
151 
152  csctf_.trLctlocalPhi[nTrk - 1][LctTrkId_] = lclPhi.phi_local;
153  //csctf_.trLctlocalPhi_bend[nTrk-1][LctTrkId_] = lclPhi.phi_bend_local;
154  //csctf_.trLctCLCT_pattern[nTrk-1][LctTrkId_] = lctTrk->getPattern();
155  csctf_.trLctQuality[nTrk - 1][LctTrkId_] = lctTrk->getQuality();
156 
157  //std::cout <<"lctTrk->getPattern() = " << lctTrk->getPattern() << std::endl;
158  } catch (...) {
159  bzero(&lclPhi, sizeof(lclPhi));
160  csctf_.trLctlocalPhi[nTrk - 1][LctTrkId_] = -999;
161  //csctf_.trLctlocalPhi_bend[nTrk-1][LctTrkId_] = -999;
162  //csctf_.trLctCLCT_pattern[nTrk-1][LctTrkId_] = -999;
163  csctf_.trLctQuality[nTrk - 1][LctTrkId_] = -999;
164  }
165  // clct pattern
166  lclphiadd lclPattern;
167  try {
168  //std::cout <<"lclPattern.clct_pattern = " << lclPattern.clct_pattern << std::endl;
169  //std::cout <<"lclPattern.pattern_type = " << lclPattern.pattern_type << std::endl;
170 
171  } catch (...) {
172  }
173 
174  // Global Phi
175  gblphidat gblPhi;
176 
177  try {
178  csctf_.trLctwireGroup[nTrk - 1][LctTrkId_] = lctTrk->getKeyWG();
179  gblPhi = srLUTs_[FPGALctTrk][endcap]->globalPhiME(lclPhi.phi_local, lctTrk->getKeyWG(), CscIdLctTrk);
180 
181  csctf_.trLctglobalPhi[nTrk - 1][LctTrkId_] = gblPhi.global_phi;
182 
183  } catch (...) {
184  bzero(&gblPhi, sizeof(gblPhi));
185  csctf_.trLctglobalPhi[nTrk - 1][LctTrkId_] = -999;
186  }
187 
188  // Global Eta
189  gbletadat gblEta;
190 
191  try {
192  gblEta = srLUTs_[FPGALctTrk][endcap]->globalEtaME(
193  lclPhi.phi_bend_local, lclPhi.phi_local, lctTrk->getKeyWG(), CscIdLctTrk);
194  csctf_.trLctglobalEta[nTrk - 1][LctTrkId_] = gblEta.global_eta;
195  csctf_.trLctCLCT_pattern[nTrk - 1][LctTrkId_] = gblEta.global_bend;
196  } catch (...) {
197  bzero(&gblEta, sizeof(gblEta));
198  csctf_.trLctglobalEta[nTrk - 1][LctTrkId_] = -999;
199  csctf_.trLctCLCT_pattern[nTrk - 1][LctTrkId_] = -999;
200  }
201 
202  ++LctTrkId_;
203 
204  } // for(CSCCorrelatedLCTDigiCollection::const_iterator lctTrk
205  } // for(CSCCorrelatedLCTDigiCollection::DigiRangeIterator lctOfTrks
206 
207  csctf_.trNumLCTs.push_back(LctTrkId_);
208  //}
209  //else
210  // edm::LogInfo("L1NtupleProducer")<<" No valid CSCCorrelatedLCTDigiCollection products found";
211 
213 
214  } //for(L1CSCTrackCollection::const_iterator trk=csctfTrks->begin(); trk<csctfTrks->end(); trk++,nTrk++){
215 }
216 
217 //ALL csctf lcts
219  CSCSectorReceiverLUT* srLUTs_[5][2]) {
220  int nLCT = 0;
221  for (CSCCorrelatedLCTDigiCollection::DigiRangeIterator corrLct = corrlcts.product()->begin();
222  corrLct != corrlcts.product()->end();
223  corrLct++) {
224  nLCT++;
225 
226  int lctId = 0;
227 
228  CSCCorrelatedLCTDigiCollection::Range lctRange = corrlcts.product()->get((*corrLct).first);
229 
230  for (CSCCorrelatedLCTDigiCollection::const_iterator lct = lctRange.first; lct != lctRange.second; lct++, lctId++) {
231  csctf_.lctEndcap.push_back((*corrLct).first.zendcap());
232  if ((*corrLct).first.zendcap() > 0)
233  csctf_.lctSector.push_back((*corrLct).first.triggerSector());
234  else
235  csctf_.lctSector.push_back(6 + (*corrLct).first.triggerSector());
236 
237  csctf_.lctSubSector.push_back(CSCTriggerNumbering::triggerSubSectorFromLabels((*corrLct).first));
238  csctf_.lctBx.push_back(lct->getBX());
239  csctf_.lctBx0.push_back(lct->getBX0());
240 
241  csctf_.lctStation.push_back((*corrLct).first.station());
242  csctf_.lctRing.push_back((*corrLct).first.ring());
243  csctf_.lctChamber.push_back((*corrLct).first.chamber());
244  csctf_.lctTriggerCSCID.push_back((*corrLct).first.triggerCscId());
245  csctf_.lctFpga.push_back(
246  (csctf_.lctSubSector.back() ? csctf_.lctSubSector.back() : (*corrLct).first.station() + 1));
247 
248  // Check if DetId is within range
249  if (csctf_.lctSector.back() < 1 || csctf_.lctSector.back() > 12 || csctf_.lctStation.back() < 1 ||
250  csctf_.lctStation.back() > 4 || csctf_.lctTriggerCSCID.back() < 1 || csctf_.lctTriggerCSCID.back() > 9 ||
251  lctId < 0 || lctId > 1) {
252  edm::LogInfo("L1NtupleProducer") << " LCT ERROR: CSC digi are out of range: ";
253 
254  continue;
255  }
256 
257  // handles not to overload the method: mostly for readability
258  int endcap = (*corrLct).first.zendcap();
259  if (endcap < 0)
260  endcap = 0;
261 
262  int StationLct = (*corrLct).first.station();
263  int CscIdLct = (*corrLct).first.triggerCscId();
264  int SubSectorLct = CSCTriggerNumbering::triggerSubSectorFromLabels((*corrLct).first);
265 
266  int FPGALct = (SubSectorLct ? SubSectorLct - 1 : StationLct);
267 
268  // local Phi
269  lclphidat lclPhi;
270 
271  /*
272  try {
273 
274  csctf_.lctstripNum.push_back(lct->getStrip());
275 
276 
277  csctf_.lctlocalPhi.push_back(lclPhi.phi_local);
278  }
279  catch(...) {
280  bzero(&lclPhi,sizeof(lclPhi));
281  csctf_.lctlocalPhi.push_back(-999);
282  }
283 
284 */
285  try {
286  csctf_.lctstripNum.push_back(lct->getStrip());
287  lclPhi =
288  srLUTs_[FPGALct][endcap]->localPhi(lct->getStrip(), lct->getPattern(), lct->getQuality(), lct->getBend());
289 
290  csctf_.lctlocalPhi.push_back(lclPhi.phi_local);
291  //csctf_.lctlocalPhi_bend.push_back(lclPhi.phi_bend_local);
292  //csctf_.lctCLCT_pattern.push_back(lct->getPattern());
293  csctf_.lctQuality.push_back(lct->getQuality());
294  //std::cout <<"localPhi: lclPhi.phi_bend_local = " << lclPhi.phi_bend_local << std::endl;
295  //std::cout <<"localPhi: lct->getBend() = " << lct->getBend() << std::endl;
296 
297  } catch (...) {
298  bzero(&lclPhi, sizeof(lclPhi));
299  csctf_.lctlocalPhi.push_back(-999);
300  //csctf_.lctlocalPhi_bend.push_back(-999);
301  //csctf_.lctCLCT_pattern.push_back(-999);
302  csctf_.lctQuality.push_back(-999);
303  }
304 
305  // Global Phi
306  gblphidat gblPhi;
307 
308  try {
309  csctf_.lctwireGroup.push_back(lct->getKeyWG());
310 
311  //std::cout << "lclPhi.phi_local: " << lclPhi.phi_local << std::endl;
312  //std::cout << "lct->getKeyWG(): " << lct->getKeyWG() << std::endl;
313  //std::cout << "CscIdLct: " << CscIdLct << std::endl;
314 
315  gblPhi = srLUTs_[FPGALct][endcap]->globalPhiME(lclPhi.phi_local, lct->getKeyWG(), CscIdLct);
316  csctf_.lctglobalPhi.push_back(gblPhi.global_phi);
317 
318  } catch (...) {
319  bzero(&gblPhi, sizeof(gblPhi));
320  csctf_.lctglobalPhi.push_back(-999);
321  }
322 
323  // Global Eta
324  gbletadat gblEta;
325 
326  try {
327  gblEta =
328  srLUTs_[FPGALct][endcap]->globalEtaME(lclPhi.phi_bend_local, lclPhi.phi_local, lct->getKeyWG(), CscIdLct);
329  //std::cout <<"gblEta: lclPhi.phi_bend_local = " << lclPhi.phi_bend_local << std::endl;
330  csctf_.lctglobalEta.push_back(gblEta.global_eta);
331  csctf_.lctCLCT_pattern.push_back(gblEta.global_bend);
332  } catch (...) {
333  bzero(&gblEta, sizeof(gblEta));
334  csctf_.lctglobalEta.push_back(-999);
335  csctf_.lctCLCT_pattern.push_back(-999);
336  }
337 
338  } // for(CSCCorrelatedLCTDigiCollection::const_iterator lct
339  } // for(CSCCorrelatedLCTDigiCollection::DigiRangeIterator lct
340 
341  csctf_.lctSize = nLCT;
342 }
343 
345  int nStat = 0;
346  for (std::vector<L1CSCSPStatusDigi>::const_iterator stat = status->second.begin(); stat != status->second.end();
347  stat++) {
348  // fill the Ntuple
349  if (stat->VPs() != 0) {
350  csctf_.stSPslot.push_back(stat->slot());
351  csctf_.stL1A_BXN.push_back(stat->BXN());
352  csctf_.stTrkCounter.push_back((const_cast<L1CSCSPStatusDigi*>(&(*stat)))->track_counter());
353  csctf_.stOrbCounter.push_back((const_cast<L1CSCSPStatusDigi*>(&(*stat)))->orbit_counter());
354 
355  nStat++;
356  }
357  }
358 
359  csctf_.nsp = nStat;
360 }
361 
362 //DT Stubs added by Alex Ji
363 //Code modeled from DQM/L1TMonitor/src/L1TCSCTF.cc, v1.36
365  //get the vector of DT Stubs
366  std::vector<csctf::TrackStub> vstubs = dtStubs->get();
367  //iterate through DT Stubs
368  for (std::vector<csctf::TrackStub>::const_iterator stub = vstubs.begin(); stub != vstubs.end(); stub++) {
369  csctf_.dtBXN.push_back(stub->BX());
370  csctf_.dtFLAG.push_back(stub->getStrip()); //getStrip() is actually the "FLAG" bit
371  csctf_.dtCAL.push_back(stub->getKeyWG()); //getKeyWG() is actually the "CAL" bit
372 
373  csctf_.dtSector.push_back(6 * (stub->endcap() - 1) + stub->sector());
374  csctf_.dtSubSector.push_back(stub->subsector());
375 
376  csctf_.dtBX0.push_back(stub->getBX0()); //it is unclear what this variable is...
377  csctf_.dtPhiBend.push_back(stub->getBend());
378  csctf_.dtPhiPacked.push_back(stub->phiPacked());
379  csctf_.dtQuality.push_back(stub->getQuality());
380  }
381 
382  csctf_.dtSize = vstubs.size();
383 }
const double Pi
ptdat Pt(const ptadd &) const
Definition: CSCTFPtLUT.cc:178
lclphidat localPhi(int strip, int pattern, int quality, int lr, const bool gangedME1a=false) const
Geometry Lookup Tables.
const L1MuScale * getPtScale() const
get the Pt scale
void SetStatus(const edm::Handle< L1CSCStatusDigiCollection > status)
virtual float getLowEdge(unsigned packed) const =0
get the low edge of bin represented by packed
list status
Definition: mps_update.py:107
virtual float getCenter(unsigned packed) const =0
get the center of bin represented by packed
gblphidat globalPhiME(int phi_local, int wire_group, int cscid, const bool gangedME1a=false) const
void SetTracks(const edm::Handle< L1CSCTrackCollection > csctfTrks, const L1MuTriggerScales *ts, const L1MuTriggerPtScale *tpts, CSCSectorReceiverLUT *srLUTs_[5][2], CSCTFPtLUT *ptLUTs_)
const L1MuScale * getPhiScale() const
get the phi scale
Log< level::Info, false > LogInfo
T const * product() const
Definition: Handle.h:70
class pt_address ptadd
class pt_data ptdat
class global_phi_data gblphidat
std::pair< const_iterator, const_iterator > Range
class local_phi_address lclphiadd
std::vector< DigiType >::const_iterator const_iterator
class local_phi_data lclphidat
Data Types.
static int triggerSubSectorFromLabels(int station, int chamber)
const L1MuScale * getRegionalEtaScale(int isys) const
get the regioanl muon trigger eta scale, isys = 0(DT), 1(bRPC), 2(CSC), 3(fwdRPC) ...
gbletadat globalEtaME(int phi_bend, int phi_local, int wire_group, int cscid, const bool gangedME1a=false) const
void SetLCTs(const edm::Handle< CSCCorrelatedLCTDigiCollection > corrlcts, CSCSectorReceiverLUT *srLUTs_[5][2])
void SetDTStubs(const edm::Handle< CSCTriggerContainer< csctf::TrackStub > > dtStubs)
class global_eta_data gbletadat