CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
L1MuBMSectorReceiver Class Reference

#include <L1MuBMSectorReceiver.h>

Public Member Functions

 L1MuBMSectorReceiver (L1MuBMSectorProcessor &, edm::ConsumesCollector &&iC)
 constructor More...
 
void reset ()
 clear Sector Receiver More...
 
void run (int bx, const edm::Event &e, const edm::EventSetup &c)
 receive track segment data from the BBMX and CSC chamber triggers More...
 
virtual ~L1MuBMSectorReceiver ()
 destructor More...
 

Private Member Functions

int address2sector (int adr) const
 find the right sector for a given address More...
 
int address2wheel (int adr) const
 find the right wheel for a given address More...
 
void receiveBBMXData (int bx, const edm::Event &e, const edm::EventSetup &c)
 receive track segment data from BBMX chamber trigger More...
 
void receiveCSCData (int bx, const edm::Event &e, const edm::EventSetup &c)
 receive track segment data from CSC chamber trigger More...
 

Private Attributes

edm::EDGetTokenT
< L1MuDTChambPhContainer
m_DTDigiToken
 
L1MuBMSectorProcessorm_sp
 
edm::ESHandle< L1MuDTTFMasksmsks
 
edm::ESHandle< L1MuDTTFParameterspars
 

Detailed Description

Sector Receiver:

The Sector Receiver receives track segment data from the BBMX and CSC chamber triggers and stores it into the Data Buffer

N. Neumeister CERN EP J. Troconiz UAM Madrid

Definition at line 51 of file L1MuBMSectorReceiver.h.

Constructor & Destructor Documentation

L1MuBMSectorReceiver::L1MuBMSectorReceiver ( L1MuBMSectorProcessor sp,
edm::ConsumesCollector &&  iC 
)

constructor

Definition at line 52 of file L1MuBMSectorReceiver.cc.

52  :
53  m_sp(sp),
55 
56 }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
edm::EDGetTokenT< L1MuDTChambPhContainer > m_DTDigiToken
L1MuBMSectorProcessor & m_sp
static edm::InputTag getBMDigiInputTag()
L1MuBMSectorReceiver::~L1MuBMSectorReceiver ( )
virtual

destructor

Definition at line 62 of file L1MuBMSectorReceiver.cc.

62  {
63 
64 // reset();
65 
66 }

Member Function Documentation

int L1MuBMSectorReceiver::address2sector ( int  adr) const
private

find the right sector for a given address

Definition at line 289 of file L1MuBMSectorReceiver.cc.

References L1MuBMSectorProcessor::id(), m_sp, and L1MuBMSecProcId::sector().

Referenced by receiveBBMXData().

289  {
290 
291  int sector = m_sp.id().sector();
292 
293  if ( adr >= 4 && adr <= 7 ) sector = (sector+13)%12; // +1
294  if ( adr >= 8 && adr <= 11 ) sector = (sector+11)%12; // -1
295 
296  return sector;
297 
298 }
int sector() const
return sector number
L1MuBMSectorProcessor & m_sp
const L1MuBMSecProcId & id() const
return Sector Processor identifier
int L1MuBMSectorReceiver::address2wheel ( int  adr) const
private

find the right wheel for a given address

Definition at line 304 of file L1MuBMSectorReceiver.cc.

References L1MuBMSectorProcessor::id(), L1MuBMSecProcId::locwheel(), m_sp, and L1MuBMSecProcId::wheel().

Referenced by receiveBBMXData().

304  {
305 
306  int wheel = m_sp.id().locwheel();
307 
308  // for 2, 3, 6, 7, 10, 11
309  if ( (adr/2)%2 == 1 ) wheel = m_sp.id().wheel();
310 
311  return wheel;
312 
313 }
L1MuBMSectorProcessor & m_sp
const L1MuBMSecProcId & id() const
return Sector Processor identifier
int locwheel() const
return physical wheel number (-2,-1,0,+1,+2)
int wheel() const
return wheel number
void L1MuBMSectorReceiver::receiveBBMXData ( int  bx,
const edm::Event e,
const edm::EventSetup c 
)
private

receive track segment data from BBMX chamber trigger

Definition at line 103 of file L1MuBMSectorReceiver.cc.

References funct::abs(), address2sector(), address2wheel(), L1MuBMDataBuffer::addTSphi(), L1MuDTChambPhDigi::code(), L1MuBMSectorProcessor::data(), edm::Event::getByToken(), L1MuBMTFConfig::getNbitsExtPhi(), L1MuBMTFConfig::getTSOutOfTimeFilter(), L1MuBMTFConfig::getTSOutOfTimeWindow(), L1MuBMSectorProcessor::id(), m_DTDigiToken, m_sp, msks, pars, phi, L1MuDTChambPhDigi::phi(), L1MuDTChambPhDigi::phiB(), L1MuBMSecProcId::sector(), relativeConstraints::station, GlobalPosition_Frontier_DevDB_cff::tag, and L1MuBMSecProcId::wheel().

Referenced by run().

103  {
105  //e.getByLabel(L1MuBMTFConfig::getBMDigiInputTag(),dttrig);
106  e.getByToken(m_DTDigiToken,dttrig);
107  L1MuDTChambPhDigi const* ts=0;
108 
109  // const int bx_offset = dttrig->correctBX();
110  int bx_offset=0;
111  bx = bx + bx_offset;
112  // get BBMX phi track segments
113  int address = 0;
114  for ( int station = 1; station <= 4; station++ ) {
115  int max_address = (station == 1) ? 2 : 12;
116  for (int reladr =0; reladr < max_address; reladr++) {
117  address++;
118  //if ( m_sp.ovl() && (reladr/2)%2 != 0 ) continue;
119  int wheel = address2wheel(reladr);
120  int sector = address2sector(reladr);
121  if ( reladr%2 == 0 ) ts = dttrig->chPhiSegm1(wheel,station,sector,bx);
122  if ( reladr%2 == 1 ) ts = dttrig->chPhiSegm2(wheel,station,sector,bx);
123  if ( ts ) {
124  int phi = ts->phi();
125  int phib = ts->phiB();
126  int qual = ts->code();
127  bool tag = (reladr%2 == 1) ? true : false;
128 
129  int lwheel = m_sp.id().wheel();
130  lwheel = abs(lwheel)/lwheel*(abs(wheel)+1);
131 
132  if ( station == 1 ) {
133  if ( msks->get_inrec_chdis_st1(lwheel, sector) ) continue;
134  if ( qual < pars->get_inrec_qual_st1(lwheel, sector) ) continue;
135  }
136  else if ( station == 2 ) {
137  if ( msks->get_inrec_chdis_st2(lwheel, sector) ) continue;
138  if ( qual < pars->get_inrec_qual_st2(lwheel, sector) ) continue;
139  }
140  else if ( station == 3 ) {
141  if ( msks->get_inrec_chdis_st3(lwheel, sector) ) continue;
142  if ( qual < pars->get_inrec_qual_st3(lwheel, sector) ) continue;
143  }
144  else if ( station == 4 ) {
145  if ( msks->get_inrec_chdis_st4(lwheel, sector) ) continue;
146  if ( qual < pars->get_inrec_qual_st4(lwheel, sector) ) continue;
147  }
148 
149  if ( reladr/2 == 1 && qual < pars->get_soc_stdis_n(m_sp.id().wheel(), m_sp.id().sector()) ) continue;
150  if ( reladr/2 == 2 && qual < pars->get_soc_stdis_wl(m_sp.id().wheel(), m_sp.id().sector()) ) continue;
151  if ( reladr/2 == 3 && qual < pars->get_soc_stdis_zl(m_sp.id().wheel(), m_sp.id().sector()) ) continue;
152  if ( reladr/2 == 4 && qual < pars->get_soc_stdis_wr(m_sp.id().wheel(), m_sp.id().sector()) ) continue;
153  if ( reladr/2 == 5 && qual < pars->get_soc_stdis_zr(m_sp.id().wheel(), m_sp.id().sector()) ) continue;
154 
155  //
156  // out-of-time TS filter (compare TS at +-1 bx)
157  //
158  bool skipTS = false;
159 
160  bool nbx_del = pars->get_soc_nbx_del(m_sp.id().wheel(), m_sp.id().sector());
161  if ( L1MuBMTFConfig::getTSOutOfTimeFilter() || nbx_del ) {
162 
163  int sh_phi = 12 - L1MuBMTFConfig::getNbitsExtPhi();
164  int tolerance = L1MuBMTFConfig::getTSOutOfTimeWindow();
165 
166  L1MuDTChambPhDigi const * tsPreviousBX_1 = dttrig->chPhiSegm1(wheel,station,sector,bx-1);
167  if ( tsPreviousBX_1 ) {
168  int phiBX = tsPreviousBX_1->phi();
169  int qualBX = tsPreviousBX_1->code();
170  if ( abs( (phi >> sh_phi) - (phiBX >> sh_phi) ) <= tolerance &&
171  qualBX > qual ) skipTS = true;
172  }
173 
174  L1MuDTChambPhDigi const * tsPreviousBX_2 = dttrig->chPhiSegm2(wheel,station,sector,bx-1);
175  if ( tsPreviousBX_2 ) {
176  int phiBX = tsPreviousBX_2->phi();
177  int qualBX = tsPreviousBX_2->code();
178  if ( abs( (phi >> sh_phi) - (phiBX >> sh_phi) ) <= tolerance &&
179  qualBX > qual ) skipTS = true;
180  }
181 
182  L1MuDTChambPhDigi const * tsNextBX_1 = dttrig->chPhiSegm1(wheel,station,sector,bx+1);
183  if ( tsNextBX_1 ) {
184  int phiBX = tsNextBX_1->phi();
185  int qualBX = tsNextBX_1->code();
186  if ( abs( (phi >> sh_phi) - (phiBX >> sh_phi) ) <= tolerance &&
187  qualBX > qual ) skipTS = true;
188  }
189 
190  L1MuDTChambPhDigi const * tsNextBX_2 = dttrig->chPhiSegm2(wheel,station,sector,bx+1);
191  if ( tsNextBX_2 ) {
192  int phiBX = tsNextBX_2->phi();
193  int qualBX = tsNextBX_2->code();
194  if ( abs( (phi >> sh_phi) - (phiBX >> sh_phi) ) <= tolerance &&
195  qualBX > qual ) skipTS = true;
196  }
197 
198  }
199 
200  if ( !skipTS ) {
201  L1MuBMTrackSegPhi tmpts(wheel,sector,station,phi,phib,
202  static_cast<L1MuBMTrackSegPhi::TSQuality>(qual),
203  tag,bx-bx_offset);
204  m_sp.data()->addTSphi(address-1,tmpts);
205  }
206  }
207 
208  }
209  }
210 
211 }
static int getTSOutOfTimeWindow()
int address2wheel(int adr) const
find the right wheel for a given address
static int getNbitsExtPhi()
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
edm::ESHandle< L1MuDTTFParameters > pars
int address2sector(int adr) const
find the right sector for a given address
edm::EDGetTokenT< L1MuDTChambPhContainer > m_DTDigiToken
int sector() const
return sector number
L1MuBMSectorProcessor & m_sp
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
static bool getTSOutOfTimeFilter()
const L1MuBMSecProcId & id() const
return Sector Processor identifier
edm::ESHandle< L1MuDTTFMasks > msks
const L1MuBMDataBuffer * data() const
return pointer to Data Buffer
int wheel() const
return wheel number
void L1MuBMSectorReceiver::receiveCSCData ( int  bx,
const edm::Event e,
const edm::EventSetup c 
)
private

receive track segment data from CSC chamber trigger

void L1MuBMSectorReceiver::reset ( void  )

clear Sector Receiver

Definition at line 95 of file L1MuBMSectorReceiver.cc.

Referenced by L1MuBMSectorProcessor::reset().

95  {
96 
97 }
void L1MuBMSectorReceiver::run ( int  bx,
const edm::Event e,
const edm::EventSetup c 
)

receive track segment data from the BBMX and CSC chamber triggers

Definition at line 76 of file L1MuBMSectorReceiver.cc.

References edm::EventSetup::get(), msks, pars, and receiveBBMXData().

Referenced by L1MuBMSectorProcessor::run().

76  {
77 
78  c.get< L1MuDTTFParametersRcd >().get( pars );
79  c.get< L1MuDTTFMasksRcd >().get( msks );
80 
81  // get track segments from BBMX chamber trigger
82  receiveBBMXData(bx, e, c);
83 
84  // get track segments from CSC chamber trigger
85  //if ( L1MuBMTFConfig::overlap() && m_sp.ovl() ) {
86  //receiveCSCData(bx, e, c);
87  //}
88 
89 }
edm::ESHandle< L1MuDTTFParameters > pars
edm::ESHandle< L1MuDTTFMasks > msks
const T & get() const
Definition: EventSetup.h:56
void receiveBBMXData(int bx, const edm::Event &e, const edm::EventSetup &c)
receive track segment data from BBMX chamber trigger

Member Data Documentation

edm::EDGetTokenT<L1MuDTChambPhContainer> L1MuBMSectorReceiver::m_DTDigiToken
private

Definition at line 87 of file L1MuBMSectorReceiver.h.

Referenced by receiveBBMXData().

L1MuBMSectorProcessor& L1MuBMSectorReceiver::m_sp
private

Definition at line 83 of file L1MuBMSectorReceiver.h.

Referenced by address2sector(), address2wheel(), and receiveBBMXData().

edm::ESHandle< L1MuDTTFMasks > L1MuBMSectorReceiver::msks
private

Definition at line 86 of file L1MuBMSectorReceiver.h.

Referenced by receiveBBMXData(), and run().

edm::ESHandle< L1MuDTTFParameters > L1MuBMSectorReceiver::pars
private

Definition at line 85 of file L1MuBMSectorReceiver.h.

Referenced by receiveBBMXData(), and run().