test
CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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),
54  m_DTDigiToken(iC.consumes<L1MuDTChambPhContainer>(L1MuBMTFConfig::getBMDigiInputTag())){
55 
56 }
57 
58 
59 //--------------
60 // Destructor --
61 //--------------
63 
64 // reset();
65 
66 }
67 
68 
69 //--------------
70 // Operations --
71 //--------------
72 
73 //
74 // receive track segment data from the BBMX chamber triggers
75 //
76 void L1MuBMSectorReceiver::run(int bx, const edm::Event& e, const edm::EventSetup& c) {
77 
78  //c.get< L1MuDTTFParametersRcd >().get( pars );
79  //c.get< L1MuDTTFMasksRcd >().get( msks );
80 
81  const L1TMuonBarrelParamsRcd& bmtfParamsRcd = c.get<L1TMuonBarrelParamsRcd>();
82  bmtfParamsRcd.get(bmtfParamsHandle);
83  const L1TMuonBarrelParams& bmtfParams = *bmtfParamsHandle.product();
84  msks = bmtfParams.l1mudttfmasks;
85  pars = bmtfParams.l1mudttfparams;
86  //pars.print();
87  //msks.print();
88 
89  // get track segments from BBMX chamber trigger
90  receiveBBMXData(bx, e, c);
91 
92 }
93 
94 
95 //
96 // clear
97 //
99 
100 }
101 
102 
103 //
104 // receive track segment data from the BBMX chamber trigger
105 //
108  //e.getByLabel(L1MuBMTFConfig::getBMDigiInputTag(),dttrig);
109  e.getByToken(m_DTDigiToken,dttrig);
110  L1MuDTChambPhDigi const* ts=0;
111 
112  // const int bx_offset = dttrig->correctBX();
113  int bx_offset=0;
114  bx = bx + bx_offset;
115  // get BBMX phi track segments
116  int address = 0;
117  for ( int station = 1; station <= 4; station++ ) {
118  int max_address = (station == 1) ? 2 : 12;
119  for (int reladr =0; reladr < max_address; reladr++) {
120  address++;
121  //if ( m_sp.ovl() && (reladr/2)%2 != 0 ) continue;
122  int wheel = address2wheel(reladr);
123  int sector = address2sector(reladr);
124  //if ( (wheel==2 || wheel==-2) && station==1 ) continue;
125 
126  if ( reladr%2 == 0 ) ts = dttrig->chPhiSegm1(wheel,station,sector,bx);
127  if ( reladr%2 == 1 ) ts = dttrig->chPhiSegm2(wheel,station,sector,bx-1);
128  if ( ts ) {
129  int phi = ts->phi();
130 // int phib = ts->phiB();
131  int phib = 0;
132  if(station!=3) phib = ts->phiB();
133 
134  int qual = ts->code();
135  bool tag = (reladr%2 == 1) ? true : false;
136 
137  int lwheel = m_sp.id().wheel();
138  lwheel = abs(lwheel)/lwheel*(abs(wheel)+1);
139 
140  if ( station == 1 ) {
141  if ( msks.get_inrec_chdis_st1(lwheel, sector) ) continue;
142  if ( qual < pars.get_inrec_qual_st1(lwheel, sector) ) continue;
143  }
144  else if ( station == 2 ) {
145  if ( msks.get_inrec_chdis_st2(lwheel, sector) ) continue;
146  if ( qual < pars.get_inrec_qual_st2(lwheel, sector) ) continue;
147  }
148  else if ( station == 3 ) {
149  if ( msks.get_inrec_chdis_st3(lwheel, sector) ) continue;
150  if ( qual < pars.get_inrec_qual_st3(lwheel, sector) ) continue;
151  }
152  else if ( station == 4 ) {
153  if ( msks.get_inrec_chdis_st4(lwheel, sector) ) continue;
154  if ( qual < pars.get_inrec_qual_st4(lwheel, sector) ) continue;
155  }
156 
157  if ( reladr/2 == 1 && qual < pars.get_soc_stdis_n(m_sp.id().wheel(), m_sp.id().sector()) ) continue;
158  if ( reladr/2 == 2 && qual < pars.get_soc_stdis_wl(m_sp.id().wheel(), m_sp.id().sector()) ) continue;
159  if ( reladr/2 == 3 && qual < pars.get_soc_stdis_zl(m_sp.id().wheel(), m_sp.id().sector()) ) continue;
160  if ( reladr/2 == 4 && qual < pars.get_soc_stdis_wr(m_sp.id().wheel(), m_sp.id().sector()) ) continue;
161  if ( reladr/2 == 5 && qual < pars.get_soc_stdis_zr(m_sp.id().wheel(), m_sp.id().sector()) ) 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 
171  int sh_phi = 12 - L1MuBMTFConfig::getNbitsExtPhi();
172  int tolerance = L1MuBMTFConfig::getTSOutOfTimeWindow();
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 &&
179  qualBX > qual ) 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 &&
187  qualBX > qual ) 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 &&
195  qualBX > qual ) 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 &&
203  qualBX > qual ) skipTS = true;
204  }
205 
206  }
207 
208  if ( !skipTS ) {
209 
210  /* if(reladr%2 == 0) {
211  L1MuBMTrackSegPhi tmpts(wheel,sector,station,phi,phib,
212  static_cast<L1MuBMTrackSegPhi::TSQuality>(qual),
213  tag,bx-bx_offset);
214  m_sp.data()->addTSphi(address-1,tmpts);
215  }
216  if(reladr%2 == 1) {
217  L1MuBMTrackSegPhi tmpts(wheel,sector,station,phi,phib,
218  static_cast<L1MuBMTrackSegPhi::TSQuality>(qual),
219  tag,bx+1);
220  m_sp.data()->addTSphi(address-1,tmpts);
221  }*/
222  L1MuBMTrackSegPhi tmpts(wheel,sector,station,phi,phib,
223  static_cast<L1MuBMTrackSegPhi::TSQuality>(qual),
224  tag,bx-bx_offset);
225  m_sp.data()->addTSphi(address-1,tmpts);
226  }
227  }
228 
229  }
230  }
231 
232 }
233 
234 
235 //
236 // find the right sector for a given address
237 //
239 
240  int sector = m_sp.id().sector();
241 
242  if ( adr >= 4 && adr <= 7 ) sector = (sector+13)%12; // +1
243  if ( adr >= 8 && adr <= 11 ) sector = (sector+11)%12; // -1
244 
245  return sector;
246 
247 }
248 
249 
250 //
251 // find the right wheel for a given address
252 //
254 
255  int wheel = m_sp.id().locwheel();
256 
257  // for 2, 3, 6, 7, 10, 11
258  if ( (adr/2)%2 == 1 ) wheel = m_sp.id().wheel();
259 
260  return wheel;
261 
262 }
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
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
bool get_soc_nbx_del(int wh, int sc) 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
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
void get(HolderT &iHolder) 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
const T & get() const
Definition: EventSetup.h:56
T const * product() const
Definition: ESHandle.h:86
void reset()
clear Sector Receiver
const L1MuBMDataBuffer * data() const
return pointer to Data Buffer
unsigned short int get_soc_stdis_wl(int wh, int sc) const
Definition: sp.h:21
edm::ESHandle< L1TMuonBarrelParams > bmtfParamsHandle
unsigned short int get_soc_stdis_n(int wh, int sc) const
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