Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:20:51

0001 //-------------------------------------------------
0002 //
0003 //   Class: L1MuBMSectorReceiver
0004 //
0005 //   Description: Sector Receiver
0006 //
0007 //
0008 //
0009 //   Author :
0010 //   N. Neumeister            CERN EP
0011 //   J. Troconiz              UAM Madrid
0012 //
0013 //--------------------------------------------------
0014 
0015 //-----------------------
0016 // This Class's Header --
0017 //-----------------------
0018 
0019 #include "L1Trigger/L1TMuonBarrel/src/L1MuBMSectorReceiver.h"
0020 
0021 //---------------
0022 // C++ Headers --
0023 //---------------
0024 
0025 #include <iostream>
0026 #include <cmath>
0027 
0028 //-------------------------------
0029 // Collaborating Class Headers --
0030 //-------------------------------
0031 
0032 #include "L1Trigger/L1TMuonBarrel/src/L1MuBMTFConfig.h"
0033 #include "L1Trigger/L1TMuonBarrel/src/L1MuBMSectorProcessor.h"
0034 #include "L1Trigger/L1TMuonBarrel/src/L1MuBMDataBuffer.h"
0035 #include "DataFormats/L1TMuon/interface/BMTF/L1MuBMTrackSegLoc.h"
0036 #include "DataFormats/L1TMuon/interface/L1MuBMTrackSegPhi.h"
0037 
0038 #include "CondFormats/L1TObjects/interface/L1MuDTTFParameters.h"
0039 #include "CondFormats/DataRecord/interface/L1MuDTTFParametersRcd.h"
0040 #include "CondFormats/L1TObjects/interface/L1MuDTTFMasks.h"
0041 #include "CondFormats/DataRecord/interface/L1MuDTTFMasksRcd.h"
0042 
0043 using namespace std;
0044 
0045 // --------------------------------
0046 //       class L1MuBMSectorReceiver
0047 //---------------------------------
0048 
0049 //----------------
0050 // Constructors --
0051 //----------------
0052 L1MuBMSectorReceiver::L1MuBMSectorReceiver(L1MuBMSectorProcessor& sp, edm::ConsumesCollector&& iC)
0053     : m_sp(sp),
0054       m_bmtfParamsToken(iC.esConsumes()),
0055       m_DTDigiToken(iC.consumes<L1MuDTChambPhContainer>(L1MuBMTFConfig::getBMDigiInputTag())) {}
0056 
0057 //--------------
0058 // Destructor --
0059 //--------------
0060 L1MuBMSectorReceiver::~L1MuBMSectorReceiver() {
0061   //  reset();
0062 }
0063 
0064 //--------------
0065 // Operations --
0066 //--------------
0067 
0068 //
0069 // receive track segment data from the BBMX chamber triggers
0070 //
0071 void L1MuBMSectorReceiver::run(int bx, const edm::Event& e, const edm::EventSetup& c) {
0072   const L1TMuonBarrelParams& bmtfParams = c.getData(m_bmtfParamsToken);
0073   msks = bmtfParams.l1mudttfmasks;
0074   pars = bmtfParams.l1mudttfparams;
0075   //pars.print();
0076   //msks.print();
0077 
0078   // get track segments from BBMX chamber trigger
0079   receiveBBMXData(bx, e);
0080 }
0081 
0082 //
0083 // clear
0084 //
0085 void L1MuBMSectorReceiver::reset() {}
0086 
0087 //
0088 // receive track segment data from the BBMX chamber trigger
0089 //
0090 void L1MuBMSectorReceiver::receiveBBMXData(int bx, const edm::Event& e) {
0091   edm::Handle<L1MuDTChambPhContainer> dttrig;
0092   //e.getByLabel(L1MuBMTFConfig::getBMDigiInputTag(),dttrig);
0093   e.getByToken(m_DTDigiToken, dttrig);
0094   L1MuDTChambPhDigi const* ts = nullptr;
0095 
0096   // const int bx_offset = dttrig->correctBX();
0097   int bx_offset = 0;
0098   bx = bx + bx_offset;
0099   // get BBMX phi track segments
0100   int address = 0;
0101   for (int station = 1; station <= 4; station++) {
0102     int max_address = (station == 1) ? 2 : 12;
0103     for (int reladr = 0; reladr < max_address; reladr++) {
0104       address++;
0105       //if ( m_sp.ovl() && (reladr/2)%2 != 0 ) continue;
0106       int wheel = address2wheel(reladr);
0107       int sector = address2sector(reladr);
0108       //if ( (wheel==2 || wheel==-2) && station==1 ) continue;
0109 
0110       if (reladr % 2 == 0)
0111         ts = dttrig->chPhiSegm1(wheel, station, sector, bx);
0112       if (reladr % 2 == 1)
0113         ts = dttrig->chPhiSegm2(wheel, station, sector, bx - 1);
0114       if (ts) {
0115         int phi = ts->phi();
0116         //        int phib = ts->phiB();
0117         int phib = 0;
0118         if (station != 3)
0119           phib = ts->phiB();
0120 
0121         int qual = ts->code();
0122         bool tag = (reladr % 2 == 1) ? true : false;
0123 
0124         int lwheel = m_sp.id().wheel();
0125         lwheel = abs(lwheel) / lwheel * (abs(wheel) + 1);
0126 
0127         if (station == 1) {
0128           if (msks.get_inrec_chdis_st1(lwheel, sector))
0129             continue;
0130           if (qual < pars.get_inrec_qual_st1(lwheel, sector))
0131             continue;
0132         } else if (station == 2) {
0133           if (msks.get_inrec_chdis_st2(lwheel, sector))
0134             continue;
0135           if (qual < pars.get_inrec_qual_st2(lwheel, sector))
0136             continue;
0137         } else if (station == 3) {
0138           if (msks.get_inrec_chdis_st3(lwheel, sector))
0139             continue;
0140           if (qual < pars.get_inrec_qual_st3(lwheel, sector))
0141             continue;
0142         } else if (station == 4) {
0143           if (msks.get_inrec_chdis_st4(lwheel, sector))
0144             continue;
0145           if (qual < pars.get_inrec_qual_st4(lwheel, sector))
0146             continue;
0147         }
0148 
0149         if (reladr / 2 == 1 && qual < pars.get_soc_stdis_n(m_sp.id().wheel(), m_sp.id().sector()))
0150           continue;
0151         if (reladr / 2 == 2 && qual < pars.get_soc_stdis_wl(m_sp.id().wheel(), m_sp.id().sector()))
0152           continue;
0153         if (reladr / 2 == 3 && qual < pars.get_soc_stdis_zl(m_sp.id().wheel(), m_sp.id().sector()))
0154           continue;
0155         if (reladr / 2 == 4 && qual < pars.get_soc_stdis_wr(m_sp.id().wheel(), m_sp.id().sector()))
0156           continue;
0157         if (reladr / 2 == 5 && qual < pars.get_soc_stdis_zr(m_sp.id().wheel(), m_sp.id().sector()))
0158           continue;
0159 
0160         //
0161         // out-of-time TS filter (compare TS at +-1 bx)
0162         //
0163         bool skipTS = false;
0164 
0165         bool nbx_del = pars.get_soc_nbx_del(m_sp.id().wheel(), m_sp.id().sector());
0166         if (L1MuBMTFConfig::getTSOutOfTimeFilter() || nbx_del) {
0167           int sh_phi = 12 - L1MuBMTFConfig::getNbitsExtPhi();
0168           int tolerance = L1MuBMTFConfig::getTSOutOfTimeWindow();
0169 
0170           L1MuDTChambPhDigi const* tsPreviousBX_1 = dttrig->chPhiSegm1(wheel, station, sector, bx - 1);
0171           if (tsPreviousBX_1) {
0172             int phiBX = tsPreviousBX_1->phi();
0173             int qualBX = tsPreviousBX_1->code();
0174             if (abs((phi >> sh_phi) - (phiBX >> sh_phi)) <= tolerance && qualBX > qual)
0175               skipTS = true;
0176           }
0177 
0178           L1MuDTChambPhDigi const* tsPreviousBX_2 = dttrig->chPhiSegm2(wheel, station, sector, bx - 1);
0179           if (tsPreviousBX_2) {
0180             int phiBX = tsPreviousBX_2->phi();
0181             int qualBX = tsPreviousBX_2->code();
0182             if (abs((phi >> sh_phi) - (phiBX >> sh_phi)) <= tolerance && qualBX > qual)
0183               skipTS = true;
0184           }
0185 
0186           L1MuDTChambPhDigi const* tsNextBX_1 = dttrig->chPhiSegm1(wheel, station, sector, bx + 1);
0187           if (tsNextBX_1) {
0188             int phiBX = tsNextBX_1->phi();
0189             int qualBX = tsNextBX_1->code();
0190             if (abs((phi >> sh_phi) - (phiBX >> sh_phi)) <= tolerance && qualBX > qual)
0191               skipTS = true;
0192           }
0193 
0194           L1MuDTChambPhDigi const* tsNextBX_2 = dttrig->chPhiSegm2(wheel, station, sector, bx + 1);
0195           if (tsNextBX_2) {
0196             int phiBX = tsNextBX_2->phi();
0197             int qualBX = tsNextBX_2->code();
0198             if (abs((phi >> sh_phi) - (phiBX >> sh_phi)) <= tolerance && qualBX > qual)
0199               skipTS = true;
0200           }
0201         }
0202 
0203         if (!skipTS) {
0204           /* if(reladr%2 == 0) {
0205                 L1MuBMTrackSegPhi tmpts(wheel,sector,station,phi,phib,
0206                                   static_cast<L1MuBMTrackSegPhi::TSQuality>(qual),
0207                                   tag,bx-bx_offset);
0208                 m_sp.data()->addTSphi(address-1,tmpts);
0209             }
0210             if(reladr%2 == 1) {
0211                 L1MuBMTrackSegPhi tmpts(wheel,sector,station,phi,phib,
0212                                   static_cast<L1MuBMTrackSegPhi::TSQuality>(qual),
0213                                   tag,bx+1);
0214                 m_sp.data()->addTSphi(address-1,tmpts);
0215             }*/
0216           L1MuBMTrackSegPhi tmpts(
0217               wheel, sector, station, phi, phib, static_cast<L1MuBMTrackSegPhi::TSQuality>(qual), tag, bx - bx_offset);
0218           m_sp.data()->addTSphi(address - 1, tmpts);
0219         }
0220       }
0221     }
0222   }
0223 }
0224 
0225 //
0226 // find the right sector for a given address
0227 //
0228 int L1MuBMSectorReceiver::address2sector(int adr) const {
0229   int sector = m_sp.id().sector();
0230 
0231   if (adr >= 4 && adr <= 7)
0232     sector = (sector + 13) % 12;  // +1
0233   if (adr >= 8 && adr <= 11)
0234     sector = (sector + 11) % 12;  // -1
0235 
0236   return sector;
0237 }
0238 
0239 //
0240 // find the right wheel for a given address
0241 //
0242 int L1MuBMSectorReceiver::address2wheel(int adr) const {
0243   int wheel = m_sp.id().locwheel();
0244 
0245   // for 2, 3, 6, 7, 10, 11
0246   if ((adr / 2) % 2 == 1)
0247     wheel = m_sp.id().wheel();
0248 
0249   return wheel;
0250 }