CMS 3D CMS Logo

L1MuBMSectorReceiver.cc
Go to the documentation of this file.
1 //-------------------------------------------------
2 //
3 // Class: L1MuBMSectorReceiver
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 
37 
42 
43 using namespace std;
44 
45 // --------------------------------
46 // class L1MuBMSectorReceiver
47 //---------------------------------
48 
49 //----------------
50 // Constructors --
51 //----------------
53  : m_sp(sp), m_DTDigiToken(iC.consumes<L1MuDTChambPhContainer>(L1MuBMTFConfig::getBMDigiInputTag())) {}
54 
55 //--------------
56 // Destructor --
57 //--------------
59  // reset();
60 }
61 
62 //--------------
63 // Operations --
64 //--------------
65 
66 //
67 // receive track segment data from the BBMX chamber triggers
68 //
70  //c.get< L1MuDTTFParametersRcd >().get( pars );
71  //c.get< L1MuDTTFMasksRcd >().get( msks );
72 
73  const L1TMuonBarrelParamsRcd& bmtfParamsRcd = c.get<L1TMuonBarrelParamsRcd>();
74  bmtfParamsRcd.get(bmtfParamsHandle);
75  const L1TMuonBarrelParams& bmtfParams = *bmtfParamsHandle.product();
76  msks = bmtfParams.l1mudttfmasks;
77  pars = bmtfParams.l1mudttfparams;
78  //pars.print();
79  //msks.print();
80 
81  // get track segments from BBMX chamber trigger
82  receiveBBMXData(bx, e, c);
83 }
84 
85 //
86 // clear
87 //
89 
90 //
91 // receive track segment data from the BBMX chamber trigger
92 //
95  //e.getByLabel(L1MuBMTFConfig::getBMDigiInputTag(),dttrig);
96  e.getByToken(m_DTDigiToken, dttrig);
97  L1MuDTChambPhDigi const* ts = nullptr;
98 
99  // const int bx_offset = dttrig->correctBX();
100  int bx_offset = 0;
101  bx = bx + bx_offset;
102  // get BBMX phi track segments
103  int address = 0;
104  for (int station = 1; station <= 4; station++) {
105  int max_address = (station == 1) ? 2 : 12;
106  for (int reladr = 0; reladr < max_address; reladr++) {
107  address++;
108  //if ( m_sp.ovl() && (reladr/2)%2 != 0 ) continue;
109  int wheel = address2wheel(reladr);
110  int sector = address2sector(reladr);
111  //if ( (wheel==2 || wheel==-2) && station==1 ) continue;
112 
113  if (reladr % 2 == 0)
114  ts = dttrig->chPhiSegm1(wheel, station, sector, bx);
115  if (reladr % 2 == 1)
116  ts = dttrig->chPhiSegm2(wheel, station, sector, bx - 1);
117  if (ts) {
118  int phi = ts->phi();
119  // int phib = ts->phiB();
120  int phib = 0;
121  if (station != 3)
122  phib = ts->phiB();
123 
124  int qual = ts->code();
125  bool tag = (reladr % 2 == 1) ? true : false;
126 
127  int lwheel = m_sp.id().wheel();
128  lwheel = abs(lwheel) / lwheel * (abs(wheel) + 1);
129 
130  if (station == 1) {
131  if (msks.get_inrec_chdis_st1(lwheel, sector))
132  continue;
133  if (qual < pars.get_inrec_qual_st1(lwheel, sector))
134  continue;
135  } else if (station == 2) {
136  if (msks.get_inrec_chdis_st2(lwheel, sector))
137  continue;
138  if (qual < pars.get_inrec_qual_st2(lwheel, sector))
139  continue;
140  } else if (station == 3) {
141  if (msks.get_inrec_chdis_st3(lwheel, sector))
142  continue;
143  if (qual < pars.get_inrec_qual_st3(lwheel, sector))
144  continue;
145  } else if (station == 4) {
146  if (msks.get_inrec_chdis_st4(lwheel, sector))
147  continue;
148  if (qual < pars.get_inrec_qual_st4(lwheel, sector))
149  continue;
150  }
151 
152  if (reladr / 2 == 1 && qual < pars.get_soc_stdis_n(m_sp.id().wheel(), m_sp.id().sector()))
153  continue;
154  if (reladr / 2 == 2 && qual < pars.get_soc_stdis_wl(m_sp.id().wheel(), m_sp.id().sector()))
155  continue;
156  if (reladr / 2 == 3 && qual < pars.get_soc_stdis_zl(m_sp.id().wheel(), m_sp.id().sector()))
157  continue;
158  if (reladr / 2 == 4 && qual < pars.get_soc_stdis_wr(m_sp.id().wheel(), m_sp.id().sector()))
159  continue;
160  if (reladr / 2 == 5 && qual < pars.get_soc_stdis_zr(m_sp.id().wheel(), m_sp.id().sector()))
161  continue;
162 
163  //
164  // out-of-time TS filter (compare TS at +-1 bx)
165  //
166  bool skipTS = false;
167 
168  bool nbx_del = pars.get_soc_nbx_del(m_sp.id().wheel(), m_sp.id().sector());
169  if (L1MuBMTFConfig::getTSOutOfTimeFilter() || nbx_del) {
170  int sh_phi = 12 - L1MuBMTFConfig::getNbitsExtPhi();
172 
173  L1MuDTChambPhDigi const* tsPreviousBX_1 = dttrig->chPhiSegm1(wheel, station, sector, bx - 1);
174  if (tsPreviousBX_1) {
175  int phiBX = tsPreviousBX_1->phi();
176  int qualBX = tsPreviousBX_1->code();
177  if (abs((phi >> sh_phi) - (phiBX >> sh_phi)) <= tolerance && qualBX > qual)
178  skipTS = true;
179  }
180 
181  L1MuDTChambPhDigi const* tsPreviousBX_2 = dttrig->chPhiSegm2(wheel, station, sector, bx - 1);
182  if (tsPreviousBX_2) {
183  int phiBX = tsPreviousBX_2->phi();
184  int qualBX = tsPreviousBX_2->code();
185  if (abs((phi >> sh_phi) - (phiBX >> sh_phi)) <= tolerance && qualBX > qual)
186  skipTS = true;
187  }
188 
189  L1MuDTChambPhDigi const* tsNextBX_1 = dttrig->chPhiSegm1(wheel, station, sector, bx + 1);
190  if (tsNextBX_1) {
191  int phiBX = tsNextBX_1->phi();
192  int qualBX = tsNextBX_1->code();
193  if (abs((phi >> sh_phi) - (phiBX >> sh_phi)) <= tolerance && qualBX > qual)
194  skipTS = true;
195  }
196 
197  L1MuDTChambPhDigi const* tsNextBX_2 = dttrig->chPhiSegm2(wheel, station, sector, bx + 1);
198  if (tsNextBX_2) {
199  int phiBX = tsNextBX_2->phi();
200  int qualBX = tsNextBX_2->code();
201  if (abs((phi >> sh_phi) - (phiBX >> sh_phi)) <= tolerance && qualBX > qual)
202  skipTS = true;
203  }
204  }
205 
206  if (!skipTS) {
207  /* if(reladr%2 == 0) {
208  L1MuBMTrackSegPhi tmpts(wheel,sector,station,phi,phib,
209  static_cast<L1MuBMTrackSegPhi::TSQuality>(qual),
210  tag,bx-bx_offset);
211  m_sp.data()->addTSphi(address-1,tmpts);
212  }
213  if(reladr%2 == 1) {
214  L1MuBMTrackSegPhi tmpts(wheel,sector,station,phi,phib,
215  static_cast<L1MuBMTrackSegPhi::TSQuality>(qual),
216  tag,bx+1);
217  m_sp.data()->addTSphi(address-1,tmpts);
218  }*/
219  L1MuBMTrackSegPhi tmpts(
220  wheel, sector, station, phi, phib, static_cast<L1MuBMTrackSegPhi::TSQuality>(qual), tag, bx - bx_offset);
221  m_sp.data()->addTSphi(address - 1, tmpts);
222  }
223  }
224  }
225  }
226 }
227 
228 //
229 // find the right sector for a given address
230 //
232  int sector = m_sp.id().sector();
233 
234  if (adr >= 4 && adr <= 7)
235  sector = (sector + 13) % 12; // +1
236  if (adr >= 8 && adr <= 11)
237  sector = (sector + 11) % 12; // -1
238 
239  return sector;
240 }
241 
242 //
243 // find the right wheel for a given address
244 //
246  int wheel = m_sp.id().locwheel();
247 
248  // for 2, 3, 6, 7, 10, 11
249  if ((adr / 2) % 2 == 1)
250  wheel = m_sp.id().wheel();
251 
252  return wheel;
253 }
L1MuDTTFParameters l1mudttfparams
static int getTSOutOfTimeWindow()
L1MuDTTFParameters pars
void run(int bx, const edm::Event &e, const edm::EventSetup &c)
receive track segment data from the BBMX and CSC chamber triggers
int address2wheel(int adr) const
find the right wheel for a given address
static int getNbitsExtPhi()
unsigned short int get_inrec_qual_st4(int wh, int sc) const
void addTSphi(int adr, const L1MuBMTrackSegPhi &)
add new phi track segment to the Data Buffer
L1MuDTChambPhDigi const * chPhiSegm1(int wheel, int stat, int sect, int bx) const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
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
unsigned short int get_soc_stdis_wr(int wh, int sc) const
unsigned short int get_soc_stdis_zr(int wh, int sc) const
bool get_inrec_chdis_st3(int wh, int sc) const
int address2sector(int adr) const
find the right sector for a given address
L1MuBMSectorReceiver(L1MuBMSectorProcessor &, edm::ConsumesCollector &&iC)
constructor
edm::EDGetTokenT< L1MuDTChambPhContainer > m_DTDigiToken
PRODUCT const & get(ESGetToken< PRODUCT, T > const &iToken) const
edm::ESHandle< L1TMuonBarrelParams > bmtfParamsHandle
int sector() const
return sector number
virtual ~L1MuBMSectorReceiver()
destructor
L1MuDTTFMasks l1mudttfmasks
L1MuBMSectorProcessor & m_sp
unsigned short int get_inrec_qual_st3(int wh, int sc) const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool get_inrec_chdis_st2(int wh, int sc) const
static bool getTSOutOfTimeFilter()
unsigned short int get_inrec_qual_st2(int wh, int sc) const
const L1MuBMSecProcId & id() const
return Sector Processor identifier
unsigned short int get_soc_stdis_zl(int wh, int sc) const
bool get_inrec_chdis_st4(int wh, int sc) const
unsigned short int get_inrec_qual_st1(int wh, int sc) const
void reset()
clear Sector Receiver
const L1MuBMDataBuffer * data() const
return pointer to Data Buffer
T get() const
Definition: EventSetup.h:73
unsigned short int get_soc_stdis_wl(int wh, int sc) const
unsigned short int get_soc_stdis_n(int wh, int sc) const
T const * product() const
Definition: ESHandle.h:86
int locwheel() const
return physical wheel number (-2,-1,0,+1,+2)
void receiveBBMXData(int bx, const edm::Event &e, const edm::EventSetup &c)
receive track segment data from BBMX chamber trigger
int wheel() const
return wheel number