CMS 3D CMS Logo

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