CMS 3D CMS Logo

L1MuDTSectorReceiver.cc

Go to the documentation of this file.
00001 //-------------------------------------------------
00002 //
00003 //   Class: L1MuDTSectorReceiver
00004 //
00005 //   Description: Sector Receiver 
00006 //
00007 //
00008 //   $Date: 2008/10/01 14:41:05 $
00009 //   $Revision: 1.13 $
00010 //
00011 //   Author :
00012 //   N. Neumeister            CERN EP
00013 //   J. Troconiz              UAM Madrid
00014 //
00015 //--------------------------------------------------
00016 
00017 //-----------------------
00018 // This Class's Header --
00019 //-----------------------
00020 
00021 #include "L1Trigger/DTTrackFinder/src/L1MuDTSectorReceiver.h"
00022 
00023 //---------------
00024 // C++ Headers --
00025 //---------------
00026 
00027 #include <iostream>
00028 #include <cmath>
00029 
00030 //-------------------------------
00031 // Collaborating Class Headers --
00032 //-------------------------------
00033 
00034 #include "L1Trigger/DTTrackFinder/src/L1MuDTTFConfig.h"
00035 #include "L1Trigger/DTTrackFinder/src/L1MuDTSectorProcessor.h"
00036 #include "L1Trigger/DTTrackFinder/src/L1MuDTDataBuffer.h"
00037 #include "L1Trigger/DTTrackFinder/src/L1MuDTTrackSegLoc.h"
00038 #include "L1Trigger/DTTrackFinder/src/L1MuDTTrackSegPhi.h"
00039 #include "L1Trigger/CSCCommonTrigger/interface/CSCConstants.h"
00040 #include "DataFormats/L1DTTrackFinder/interface/L1MuDTChambPhDigi.h"
00041 #include "DataFormats/L1DTTrackFinder/interface/L1MuDTChambPhContainer.h"
00042 #include "DataFormats/L1CSCTrackFinder/interface/TrackStub.h"
00043 #include "DataFormats/L1CSCTrackFinder/interface/CSCTriggerContainer.h"
00044 
00045 using namespace std;
00046 
00047 // --------------------------------
00048 //       class L1MuDTSectorReceiver
00049 //---------------------------------
00050 
00051 //----------------
00052 // Constructors --
00053 //----------------
00054 L1MuDTSectorReceiver::L1MuDTSectorReceiver(L1MuDTSectorProcessor& sp) : 
00055         m_sp(sp) {
00056 
00057 }
00058 
00059 
00060 //--------------
00061 // Destructor --
00062 //--------------
00063 L1MuDTSectorReceiver::~L1MuDTSectorReceiver() { 
00064 
00065 //  reset();
00066   
00067 }
00068 
00069 
00070 //--------------
00071 // Operations --
00072 //--------------
00073 
00074 //
00075 // receive track segment data from the DTBX and CSC chamber triggers
00076 //
00077 void L1MuDTSectorReceiver::run(int bx, const edm::Event& e) {
00078 
00079   // get track segments from DTBX chamber trigger
00080   receiveDTBXData(bx,e);
00081   
00082   // get track segments from CSC chamber trigger
00083   if ( L1MuDTTFConfig::overlap() && m_sp.ovl() ) { 
00084     receiveCSCData(bx,e);
00085   }
00086 
00087 }
00088 
00089 
00090 //
00091 // clear
00092 //
00093 void L1MuDTSectorReceiver::reset() {
00094 
00095 }
00096 
00097 
00098 //
00099 // receive track segment data from the DTBX chamber trigger
00100 //
00101 void L1MuDTSectorReceiver::receiveDTBXData(int bx, const edm::Event& e) {
00102 
00103   edm::Handle<L1MuDTChambPhContainer> dttrig;
00104   e.getByLabel(L1MuDTTFConfig::getDTDigiInputTag(),dttrig);
00105 
00106   L1MuDTChambPhDigi* ts=0;
00107 
00108   // const int bx_offset = dttrig->correctBX();
00109   int bx_offset=0;
00110   bx = bx + bx_offset;
00111 
00112   // get DTBX phi track segments  
00113   int address = 0;
00114   for ( int station = 1; station <= 4; station++ ) {
00115     int max_address = (station == 1) ? 2 : 12;
00116     for (int reladr =0; reladr < max_address; reladr++) {
00117       address++;
00118       if ( m_sp.ovl() && (reladr/2)%2 != 0 ) continue;
00119       int wheel  = address2wheel(reladr);
00120       int sector = address2sector(reladr);     
00121       if ( reladr%2 == 0 ) ts = dttrig->chPhiSegm1(wheel,station,sector,bx);
00122       if ( reladr%2 == 1 ) ts = dttrig->chPhiSegm2(wheel,station,sector,bx);
00123       if ( ts ) {
00124         int phi  = ts->phi();
00125         int phib = ts->phiB();
00126         int qual = ts->code();
00127         bool tag = (reladr%2 == 1) ? true : false;
00128         //
00129         // out-of-time TS filter (compare TS at +-1 bx)
00130         // 
00131         bool skipTS = false;
00132 
00133         if ( L1MuDTTFConfig::getTSOutOfTimeFilter() ) {
00134  
00135           int sh_phi = 12 - L1MuDTTFConfig::getNbitsExtPhi();
00136           int tolerance = L1MuDTTFConfig::getTSOutOfTimeWindow();
00137 
00138           L1MuDTChambPhDigi* tsPreviousBX_1 = dttrig->chPhiSegm1(wheel,station,sector,bx-1);
00139           if ( tsPreviousBX_1 ) {
00140             int phiBX  = tsPreviousBX_1->phi();
00141             int qualBX = tsPreviousBX_1->code();
00142             if ( abs( (phi >> sh_phi) - (phiBX >> sh_phi) ) <= tolerance &&
00143                  qualBX > qual ) skipTS = true;
00144           }
00145           
00146           L1MuDTChambPhDigi* tsPreviousBX_2 = dttrig->chPhiSegm2(wheel,station,sector,bx-1);
00147           if ( tsPreviousBX_2 ) {
00148             int phiBX  = tsPreviousBX_2->phi();
00149             int qualBX = tsPreviousBX_2->code();
00150             if ( abs( (phi >> sh_phi) - (phiBX >> sh_phi) ) <= tolerance &&
00151                  qualBX > qual ) skipTS = true;
00152           }
00153      
00154           L1MuDTChambPhDigi* tsNextBX_1 = dttrig->chPhiSegm1(wheel,station,sector,bx+1);
00155           if ( tsNextBX_1 ) {
00156             int phiBX  = tsNextBX_1->phi();
00157             int qualBX = tsNextBX_1->code();
00158             if ( abs( (phi >> sh_phi) - (phiBX >> sh_phi) ) <= tolerance &&
00159                  qualBX > qual ) skipTS = true;
00160           }
00161 
00162           L1MuDTChambPhDigi* tsNextBX_2 = dttrig->chPhiSegm2(wheel,station,sector,bx+1);
00163           if ( tsNextBX_2 ) {
00164             int phiBX  = tsNextBX_2->phi();
00165             int qualBX = tsNextBX_2->code();
00166             if ( abs( (phi >> sh_phi) - (phiBX >> sh_phi) ) <= tolerance &&
00167                  qualBX > qual ) skipTS = true;
00168           }
00169 
00170         }
00171 
00172         if ( !skipTS ) { 
00173           L1MuDTTrackSegPhi tmpts(wheel,sector,station,phi,phib,
00174                                   static_cast<L1MuDTTrackSegPhi::TSQuality>(qual),
00175                                   tag,bx-bx_offset);
00176           m_sp.data()->addTSphi(address-1,tmpts);
00177         }
00178       }
00179 
00180     }
00181   }
00182 
00183 }
00184  
00185 
00186 //
00187 // receive track segment data from CSC chamber trigger
00188 //
00189 void L1MuDTSectorReceiver::receiveCSCData(int bx, const edm::Event& e) {
00190   
00191   if ( (L1MuDTTFConfig::getCSCTrSInputTag()).label() == "none" ) return;
00192 
00193   if ( bx < CSCConstants::MIN_BUNCH || bx > CSCConstants::MAX_BUNCH ) return;
00194 
00195   edm::Handle<CSCTriggerContainer<csctf::TrackStub> > csctrig;
00196   e.getByLabel(L1MuDTTFConfig::getCSCTrSInputTag(),csctrig);
00197 
00198   const int bxCSC = CSCConstants::TIME_OFFSET;
00199   
00200   vector<csctf::TrackStub> csc_list;
00201   vector<csctf::TrackStub>::const_iterator csc_iter;  
00202   
00203   int station = 1; // only ME13
00204   int wheel = m_sp.id().wheel();
00205   int side = ( wheel == 3 ) ? 1 : 2;
00206   int sector = m_sp.id().sector();
00207   int csc_sector = ( sector == 0 ) ? 6 : (sector+1)/2;
00208   int subsector = ( sector%2 == 0 ) ? 2 : 1;
00209 
00210   csc_list = csctrig->get(side,station,csc_sector,subsector,bxCSC+bx);
00211   int ncsc = 0;
00212   for ( csc_iter = csc_list.begin(); csc_iter != csc_list.end(); csc_iter++ ) {
00213     if ( csc_iter->etaPacked() > 17 ) continue;
00214     bool etaFlag = ( csc_iter->etaPacked() > 17 ); 
00215     int phiCSC = csc_iter->phiPacked();
00216     int qualCSC = csc_iter->getQuality();
00217            
00218     // convert CSC quality code to DTBX quality code
00219     unsigned int qual = 7;
00220     if ( qualCSC ==  2 ) qual = 0;
00221     if ( qualCSC ==  6 ) qual = 1;
00222     if ( qualCSC ==  7 ) qual = 2;
00223     if ( qualCSC ==  8 ) qual = 2;
00224     if ( qualCSC ==  9 ) qual = 3;
00225     if ( qualCSC == 10 ) qual = 3;
00226     if ( qualCSC == 11 ) qual = 4;
00227     if ( qualCSC == 12 ) qual = 5;
00228     if ( qualCSC == 13 ) qual = 5;
00229     if ( qualCSC == 14 ) qual = 6;
00230     if ( qualCSC == 15 ) qual = 6;
00231     if ( qual == 7 ) continue;
00232 
00233     if ( subsector == 1 && phiCSC >= 2048 ) continue;
00234     if ( subsector == 2 && phiCSC <  2048 ) continue;
00235         
00236     // convert CSC phi to DTBX phi
00237     double dphi = ((double) phiCSC) - 2048.*16./31.;
00238     if ( sector%2 == 0 ) dphi = dphi  - 2048.*30./ 31.;
00239     dphi = dphi*62.*M_PI/180.;
00240     int phi = static_cast<int>(floor( dphi ));
00241     if ( phi < -2048 || phi > 2047 ) continue; 
00242 
00243     if ( ncsc < 2 ) {
00244       int address = 16 + ncsc;
00245       bool tag = (ncsc == 1 ) ? true : false;
00246       L1MuDTTrackSegPhi tmpts(wheel,sector,station+2,phi,0,
00247                               static_cast<L1MuDTTrackSegPhi::TSQuality>(qual),
00248                               tag,bx,etaFlag);
00249       m_sp.data()->addTSphi(address,tmpts);
00250       ncsc++;
00251     }
00252     else cout << "too many CSC track segments!" << endl;
00253   }  
00254 
00255 }
00256 
00257 
00258 
00259 //
00260 // find the right sector for a given address
00261 //
00262 int L1MuDTSectorReceiver::address2sector(int adr) const {
00263 
00264   int sector = m_sp.id().sector();
00265 
00266   if ( adr >= 4 && adr <= 7  ) sector = (sector+13)%12; // +1
00267   if ( adr >= 8 && adr <= 11 ) sector = (sector+11)%12; // -1
00268 
00269   return sector;
00270 
00271 }
00272 
00273 
00274 //
00275 // find the right wheel for a given address
00276 //
00277 int L1MuDTSectorReceiver::address2wheel(int adr) const {
00278 
00279   int wheel = m_sp.id().locwheel();
00280 
00281   // for 2, 3, 6, 7, 10, 11
00282   if ( (adr/2)%2 == 1 ) wheel = m_sp.id().wheel();
00283 
00284   return wheel;
00285 
00286 }

Generated on Tue Jun 9 17:40:01 2009 for CMSSW by  doxygen 1.5.4