CMS 3D CMS Logo

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