CMS 3D CMS Logo

L1MuDTSectorReceiver.cc
Go to the documentation of this file.
1 //-------------------------------------------------
2 //
3 // Class: L1MuDTSectorReceiver
4 //
5 // Description: Sector Receiver
6 //
7 //
8 //
9 // Author :
10 // N. Neumeister CERN EP
11 // J. Troconiz UAM Madrid
12 //
13 //--------------------------------------------------
14 
15 //-----------------------
16 // This Class's Header --
17 //-----------------------
18 
20 
21 //---------------
22 // C++ Headers --
23 //---------------
24 
25 #include <iostream>
26 #include <cmath>
27 
28 //-------------------------------
29 // Collaborating Class Headers --
30 //-------------------------------
31 
46 
47 using namespace std;
48 
49 // --------------------------------
50 // class L1MuDTSectorReceiver
51 //---------------------------------
52 
53 //----------------
54 // Constructors --
55 //----------------
57  m_sp(sp),
58  m_DTDigiToken(iC.consumes<L1MuDTChambPhContainer>(m_sp.tf().config()->getDTDigiInputTag())),
59  m_CSCTrSToken(iC.mayConsume<CSCTriggerContainer<csctf::TrackStub> >(m_sp.tf().config()->getCSCTrSInputTag())) {
60 }
61 
62 
63 //--------------
64 // Destructor --
65 //--------------
67 
68 // reset();
69 
70 }
71 
72 
73 //--------------
74 // Operations --
75 //--------------
76 
77 //
78 // receive track segment data from the DTBX and CSC chamber triggers
79 //
80 void L1MuDTSectorReceiver::run(int bx, const edm::Event& e, const edm::EventSetup& c) {
81 
82  c.get< L1MuDTTFParametersRcd >().get( pars );
83  c.get< L1MuDTTFMasksRcd >().get( msks );
84 
85  // get track segments from DTBX chamber trigger
86  receiveDTBXData(bx, e, c);
87 
88  // get track segments from CSC chamber trigger
89  if ( m_sp.tf().config()->overlap() && m_sp.ovl() ) {
90  receiveCSCData(bx, e, c);
91  }
92 
93 }
94 
95 
96 //
97 // clear
98 //
100 
101 }
102 
103 
104 //
105 // receive track segment data from the DTBX chamber trigger
106 //
108 
110  e.getByToken(m_DTDigiToken,dttrig);
111 
112  L1MuDTChambPhDigi const* ts=nullptr;
113 
114  // const int bx_offset = dttrig->correctBX();
115  int bx_offset=0;
116  bx = bx + bx_offset;
117 
118  // get DTBX phi track segments
119  int address = 0;
120  for ( int station = 1; station <= 4; station++ ) {
121  int max_address = (station == 1) ? 2 : 12;
122  for (int reladr =0; reladr < max_address; reladr++) {
123  address++;
124  if ( m_sp.ovl() && (reladr/2)%2 != 0 ) continue;
125  int wheel = address2wheel(reladr);
126  int sector = address2sector(reladr);
127  if ( reladr%2 == 0 ) ts = dttrig->chPhiSegm1(wheel,station,sector,bx);
128  if ( reladr%2 == 1 ) ts = dttrig->chPhiSegm2(wheel,station,sector,bx);
129  if ( ts ) {
130  int phi = ts->phi();
131  int phib = ts->phiB();
132  int qual = ts->code();
133  bool tag = (reladr%2 == 1) ? true : false;
134 
135  int lwheel = m_sp.id().wheel();
136  lwheel = abs(lwheel)/lwheel*(abs(wheel)+1);
137 
138  if ( station == 1 ) {
139  if ( msks->get_inrec_chdis_st1(lwheel, sector) ) continue;
140  if ( qual < pars->get_inrec_qual_st1(lwheel, sector) ) continue;
141  }
142  else if ( station == 2 ) {
143  if ( msks->get_inrec_chdis_st2(lwheel, sector) ) continue;
144  if ( qual < pars->get_inrec_qual_st2(lwheel, sector) ) continue;
145  }
146  else if ( station == 3 ) {
147  if ( msks->get_inrec_chdis_st3(lwheel, sector) ) continue;
148  if ( qual < pars->get_inrec_qual_st3(lwheel, sector) ) continue;
149  }
150  else if ( station == 4 ) {
151  if ( msks->get_inrec_chdis_st4(lwheel, sector) ) continue;
152  if ( qual < pars->get_inrec_qual_st4(lwheel, sector) ) continue;
153  }
154 
155  if ( reladr/2 == 1 && qual < pars->get_soc_stdis_n(m_sp.id().wheel(), m_sp.id().sector()) ) continue;
156  if ( reladr/2 == 2 && qual < pars->get_soc_stdis_wl(m_sp.id().wheel(), m_sp.id().sector()) ) continue;
157  if ( reladr/2 == 3 && qual < pars->get_soc_stdis_zl(m_sp.id().wheel(), m_sp.id().sector()) ) continue;
158  if ( reladr/2 == 4 && qual < pars->get_soc_stdis_wr(m_sp.id().wheel(), m_sp.id().sector()) ) continue;
159  if ( reladr/2 == 5 && qual < pars->get_soc_stdis_zr(m_sp.id().wheel(), m_sp.id().sector()) ) continue;
160 
161  //
162  // out-of-time TS filter (compare TS at +-1 bx)
163  //
164  bool skipTS = false;
165 
166  bool nbx_del = pars->get_soc_nbx_del(m_sp.id().wheel(), m_sp.id().sector());
167  if ( m_sp.tf().config()->getTSOutOfTimeFilter() || nbx_del ) {
168 
169  int sh_phi = 12 - m_sp.tf().config()->getNbitsExtPhi();
171 
172  L1MuDTChambPhDigi const * tsPreviousBX_1 = dttrig->chPhiSegm1(wheel,station,sector,bx-1);
173  if ( tsPreviousBX_1 ) {
174  int phiBX = tsPreviousBX_1->phi();
175  int qualBX = tsPreviousBX_1->code();
176  if ( abs( (phi >> sh_phi) - (phiBX >> sh_phi) ) <= tolerance &&
177  qualBX > qual ) skipTS = true;
178  }
179 
180  L1MuDTChambPhDigi const * tsPreviousBX_2 = dttrig->chPhiSegm2(wheel,station,sector,bx-1);
181  if ( tsPreviousBX_2 ) {
182  int phiBX = tsPreviousBX_2->phi();
183  int qualBX = tsPreviousBX_2->code();
184  if ( abs( (phi >> sh_phi) - (phiBX >> sh_phi) ) <= tolerance &&
185  qualBX > qual ) skipTS = true;
186  }
187 
188  L1MuDTChambPhDigi const * tsNextBX_1 = dttrig->chPhiSegm1(wheel,station,sector,bx+1);
189  if ( tsNextBX_1 ) {
190  int phiBX = tsNextBX_1->phi();
191  int qualBX = tsNextBX_1->code();
192  if ( abs( (phi >> sh_phi) - (phiBX >> sh_phi) ) <= tolerance &&
193  qualBX > qual ) skipTS = true;
194  }
195 
196  L1MuDTChambPhDigi const * tsNextBX_2 = dttrig->chPhiSegm2(wheel,station,sector,bx+1);
197  if ( tsNextBX_2 ) {
198  int phiBX = tsNextBX_2->phi();
199  int qualBX = tsNextBX_2->code();
200  if ( abs( (phi >> sh_phi) - (phiBX >> sh_phi) ) <= tolerance &&
201  qualBX > qual ) skipTS = true;
202  }
203 
204  }
205 
206  if ( !skipTS ) {
207  L1MuDTTrackSegPhi tmpts(wheel,sector,station,phi,phib,
208  static_cast<L1MuDTTrackSegPhi::TSQuality>(qual),
209  tag,bx-bx_offset);
210  m_sp.data()->addTSphi(address-1,tmpts);
211  }
212  }
213 
214  }
215  }
216 
217 }
218 
219 
220 //
221 // receive track segment data from CSC chamber trigger
222 //
224 
225  if ( (m_sp.tf().config()->getCSCTrSInputTag()).label() == "none" ) return;
226 
227  if ( bx < -6 || bx > 6 ) return;
228 
230  e.getByToken(m_CSCTrSToken,csctrig);
231 
232  const int bxCSC = 6;
233 
234  vector<csctf::TrackStub> csc_list;
235  vector<csctf::TrackStub>::const_iterator csc_iter;
236 
237  int station = 1; // only ME13
238  int wheel = m_sp.id().wheel();
239  int side = ( wheel == 3 ) ? 1 : 2;
240  int sector = m_sp.id().sector();
241  int csc_sector = ( sector == 0 ) ? 6 : (sector+1)/2;
242  int subsector = ( sector%2 == 0 ) ? 2 : 1;
243 
244  csc_list = csctrig->get(side,station,csc_sector,subsector,bxCSC+bx);
245  int ncsc = 0;
246  for ( csc_iter = csc_list.begin(); csc_iter != csc_list.end(); csc_iter++ ) {
247  bool etaFlag = ( csc_iter->etaPacked() > 17 );
248  int qualCSC = csc_iter->getQuality();
249 
250  // convert CSC quality code to DTBX quality code
251  unsigned int qual = 7;
252  if ( qualCSC == 2 ) qual = 0;
253  if ( qualCSC == 6 ) qual = 1;
254  if ( qualCSC == 7 ) qual = 2;
255  if ( qualCSC == 8 ) qual = 2;
256  if ( qualCSC == 9 ) qual = 3;
257  if ( qualCSC == 10 ) qual = 3;
258  if ( qualCSC == 11 ) qual = 4;
259  if ( qualCSC == 12 ) qual = 5;
260  if ( qualCSC == 13 ) qual = 5;
261  if ( qualCSC == 14 ) qual = 6;
262  if ( qualCSC == 15 ) qual = 6;
263  if ( qual == 7 ) continue;
264 
265  // convert CSC phi to DTBX phi
266  int phi = csc_iter->phiPacked();
267  if ( phi > 2047 ) phi -= 4096;
268  if ( phi < -2048 || phi > 2047 ) continue;
269 
270  if ( msks->get_inrec_chdis_csc(m_sp.id().wheel(), m_sp.id().sector()) ) continue;
271  if ( qual < pars->get_soc_qual_csc(m_sp.id().wheel(), m_sp.id().sector()) ) continue;
272  if ( pars->get_soc_csc_etacanc(m_sp.id().wheel(), m_sp.id().sector()) && etaFlag ) continue;
273  if ( m_sp.tf().config()->getEtaCanc() && etaFlag ) continue;
274 
275  if ( ncsc < 2 ) {
276  int address = 16 + ncsc;
277  bool tag = (ncsc == 1 ) ? true : false;
278  L1MuDTTrackSegPhi tmpts(wheel,sector,station+2,phi,0,
279  static_cast<L1MuDTTrackSegPhi::TSQuality>(qual),
280  tag,bx,etaFlag);
281  m_sp.data()->addTSphi(address,tmpts);
282  ncsc++;
283  }
284  // else cout << "too many CSC track segments!" << endl;
285  }
286 
287 }
288 
289 
290 
291 //
292 // find the right sector for a given address
293 //
295 
296  int sector = m_sp.id().sector();
297 
298  if ( adr >= 4 && adr <= 7 ) sector = (sector+13)%12; // +1
299  if ( adr >= 8 && adr <= 11 ) sector = (sector+11)%12; // -1
300 
301  return sector;
302 
303 }
304 
305 
306 //
307 // find the right wheel for a given address
308 //
310 
311  int wheel = m_sp.id().locwheel();
312 
313  // for 2, 3, 6, 7, 10, 11
314  if ( (adr/2)%2 == 1 ) wheel = m_sp.id().wheel();
315 
316  return wheel;
317 
318 }
void receiveCSCData(int bx, const edm::Event &e, const edm::EventSetup &c)
receive track segment data from CSC chamber trigger
void receiveDTBXData(int bx, const edm::Event &e, const edm::EventSetup &c)
receive track segment data from DTBX chamber trigger
L1MuDTSectorReceiver(L1MuDTSectorProcessor &, edm::ConsumesCollector &&iC)
constructor
int address2sector(int adr) const
find the right sector for a given address
L1MuDTChambPhDigi const * chPhiSegm1(int wheel, int stat, int sect, int bx) const
void reset()
clear Sector Receiver
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
bool get_soc_nbx_del(int wh, int sc) const
const double tolerance
L1MuDTChambPhDigi const * chPhiSegm2(int wheel, int stat, int sect, int bx) const
bool get_inrec_chdis_st1(int wh, int sc) const
L1MuDTSectorProcessor & m_sp
bool get_inrec_chdis_st3(int wh, int sc) const
edm::InputTag getCSCTrSInputTag() const
Definition: config.py:1
bool getTSOutOfTimeFilter() const
int sector() const
return sector number
const L1MuDTDataBuffer * data() const
return pointer to Data Buffer
bool get_soc_csc_etacanc(int wh, int sc) const
char const * label
bool getEtaCanc() const
int getTSOutOfTimeWindow() const
edm::ESHandle< L1MuDTTFParameters > pars
edm::ESHandle< L1MuDTTFMasks > msks
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool get_inrec_chdis_st2(int wh, int sc) const
const L1MuDTSecProcId & id() const
return Sector Processor identifier
bool overlap() const
int locwheel() const
return physical wheel number (-2,-1,0,+1,+2)
bool get_inrec_chdis_st4(int wh, int sc) const
int address2wheel(int adr) const
find the right wheel for a given address
bool get_inrec_chdis_csc(int wh, int sc) const
static L1MuDTTFConfig * config()
return configuration
void addTSphi(int adr, const L1MuDTTrackSegPhi &)
add new phi track segment to the Data Buffer
T get() const
Definition: EventSetup.h:68
const L1MuDTTrackFinder & tf() const
return reference to barrel MTTF
edm::EDGetTokenT< L1MuDTChambPhContainer > m_DTDigiToken
edm::EDGetTokenT< CSCTriggerContainer< csctf::TrackStub > > m_CSCTrSToken
void run(int bx, const edm::Event &e, const edm::EventSetup &c)
receive track segment data from the DTBX and CSC chamber triggers
int wheel() const
return wheel number
int getNbitsExtPhi() const
bool ovl() const
is it an overlap region Sector Processor?
virtual ~L1MuDTSectorReceiver()
destructor