CMS 3D CMS Logo

L1MuDTEtaProcessor.cc

Go to the documentation of this file.
00001 //-------------------------------------------------
00002 //
00003 //   Class: L1MuDTEtaProcessor
00004 //
00005 //   Description: Eta Processor
00006 //
00007 //                An Eta Processor consists of :
00008 //                a receiver unit,
00009 //                one Eta Track Finder (ETF) and 
00010 //                one Eta Matching Unit (EMU) 
00011 //
00012 //   $Date: 2008/05/12 15:00:09 $
00013 //   $Revision: 1.11 $
00014 //
00015 //   Author :
00016 //   N. Neumeister            CERN EP
00017 //   J. Troconiz              UAM Madrid
00018 //
00019 //--------------------------------------------------
00020 
00021 //-----------------------
00022 // This Class's Header --
00023 //-----------------------
00024 
00025 #include "L1Trigger/DTTrackFinder/src/L1MuDTEtaProcessor.h"
00026 
00027 //---------------
00028 // C++ Headers --
00029 //---------------
00030 
00031 #include <iostream>
00032 #include <iomanip>
00033 #include <bitset>
00034 
00035 //-------------------------------
00036 // Collaborating Class Headers --
00037 //-------------------------------
00038 
00039 #include "L1Trigger/DTTrackFinder/src/L1MuDTTFConfig.h"
00040 #include "L1Trigger/DTTrackFinder/src/L1MuDTTrackSegEta.h"
00041 #include "L1Trigger/DTTrackFinder/src/L1MuDTSecProcId.h"
00042 #include "L1Trigger/DTTrackFinder/src/L1MuDTSectorProcessor.h"
00043 #include "L1Trigger/DTTrackFinder/interface/L1MuDTTrackFinder.h"
00044 #include "L1Trigger/DTTrackFinder/interface/L1MuDTTrack.h"
00045 #include "CondFormats/L1TObjects/interface/L1MuDTEtaPattern.h"
00046 #include "CondFormats/L1TObjects/interface/L1MuDTEtaPatternLut.h"
00047 #include "CondFormats/DataRecord/interface/L1MuDTEtaPatternLutRcd.h"
00048 #include "CondFormats/L1TObjects/interface/L1MuDTQualPatternLut.h"
00049 #include "CondFormats/DataRecord/interface/L1MuDTQualPatternLutRcd.h"
00050 #include "DataFormats/L1DTTrackFinder/interface/L1MuDTChambThDigi.h"
00051 #include "DataFormats/L1DTTrackFinder/interface/L1MuDTChambThContainer.h"
00052 
00053 using namespace std;
00054 
00055 // --------------------------------
00056 //       class L1MuDTEtaProcessor
00057 //---------------------------------
00058 
00059 //----------------
00060 // Constructors --
00061 //----------------
00062 
00063 L1MuDTEtaProcessor::L1MuDTEtaProcessor(const L1MuDTTrackFinder& tf, int id) :
00064       m_tf(tf), m_epid(id), m_foundPattern(0), m_tseta(15) {
00065 
00066   m_tseta.reserve(15);
00067   
00068 }
00069 
00070 
00071 //--------------
00072 // Destructor --
00073 //--------------
00074 
00075 L1MuDTEtaProcessor::~L1MuDTEtaProcessor() {}
00076 
00077 
00078 //--------------
00079 // Operations --
00080 //--------------
00081 
00082 //
00083 // run Eta Processor
00084 //
00085 void L1MuDTEtaProcessor::run(int bx, const edm::Event& e, const edm::EventSetup& c) {
00086 
00087   if ( L1MuDTTFConfig::getEtaTF() ) {
00088     receiveData(bx,e);
00089     runEtaTrackFinder(c);
00090   }
00091 
00092   receiveAddresses();
00093   runEtaMatchingUnit(c);
00094 
00095   assign();
00096 
00097 }
00098 
00099 
00100 //
00101 // reset Eta Processor
00102 //
00103 void L1MuDTEtaProcessor::reset() {
00104 
00105   vector<const L1MuDTTrackSegEta*>::iterator iter = m_tseta.begin();
00106   while ( iter != m_tseta.end() ) {
00107     if ( *iter ) {
00108       delete *iter;
00109       *iter = 0;
00110     }
00111     iter++;
00112   }
00113 
00114   m_tseta.clear();
00115   
00116   for ( int i = 0; i < 12; i++ ) {
00117     m_eta[i] = 99;
00118     m_fine[i] = false;
00119     m_pattern[i] = 0;
00120     m_address[i] = 0;
00121     m_TrackCand[i] = 0;
00122     m_TracKCand[i] = 0;
00123   }
00124 
00125   m_foundPattern.clear();
00126 
00127 } 
00128 
00129 
00130 //
00131 // print track candidates found in Eta Processor
00132 //
00133 void L1MuDTEtaProcessor::print() const {
00134 
00135   bool empty1 = true;
00136   for ( int i = 0; i < 15; i++ ) {
00137     empty1 &= ( m_tseta[i] == 0 || m_tseta[i]->empty() );
00138   }
00139 
00140   bool empty2 = true;
00141   for ( int i = 0; i < 12; i++ ) {  
00142     empty2 &= ( m_address[i] == 0 );
00143   }
00144 
00145   if ( !empty1 || !empty2 ) {
00146     cout << "Eta processor " << m_epid << " : " << endl;
00147   
00148    // print local pattern
00149    if ( !empty1 ) {
00150      cout << "Local pattern : " << endl;
00151      for ( int i = 0; i < 15; i++ ) {
00152        if ( (i+5)%5 == 0 ) cout << "station " << m_tseta[i]->station() << " : ";
00153        bitset<7> pos(m_tseta[i]->position());
00154        bitset<7> qua(m_tseta[i]->quality());
00155        for ( int j = 6; j >= 0; j-- ) {
00156          cout << pos[j]+qua[j];
00157        }
00158        cout << " ";
00159        if ( (i+1)%5 == 0 ) cout << endl;
00160      }
00161      cout << "Found patterns :" << endl;
00162      vector<int>::const_iterator iter;
00163      for ( iter = m_foundPattern.begin(); iter != m_foundPattern.end(); iter++ ) {
00164         const L1MuDTEtaPattern p = theEtaPatternLUT->getPattern(*iter);
00165         int qualitycode = p.quality();
00166         cout << "ID = " << setw(4) << p.id() << "  "
00167              << "eta = " << setw(3) << p.eta() << "  "
00168              << "quality = " << setw(2) << qualitycode << " ("
00169              << quality(qualitycode,1) << " "
00170              << quality(qualitycode,2) << " " 
00171              << quality(qualitycode,3) << ")";
00172         for ( int i = 0; i < 12; i++ ) { 
00173           if ( m_pattern[i] ==  p.id() ) cout << " <--";
00174         }
00175         cout << endl;     
00176       }
00177     }
00178     
00179     cout << "Received addresses : " << endl;
00180     for ( int i = 0; i < 12; i++ ) cout << setw(3) << m_address[i] << " ";
00181     cout << endl;
00182     
00183     if ( !empty1 ) {
00184       cout << "Matched patterns : " << endl;
00185       for ( int i = 0; i < 12; i++ ) {      
00186         if ( m_fine[i] ) {
00187           const L1MuDTEtaPattern p = theEtaPatternLUT->getPattern(m_pattern[i]);
00188           int fineeta = p.eta();
00189           int coarseeta = theQualPatternLUT->getCoarseEta(i/2+1,m_address[i]);
00190           cout << "Index = " << setw(2) << i << ", "
00191                << "address = " << setw(2) << m_address[i] << " --> " 
00192                << "pattern = " << setw(4) << m_pattern[i] << " "
00193                << "eta (coarse) = " << setw(3) << coarseeta << " "
00194                << "eta (fine) = " << setw(3) << fineeta << " "
00195                << "quality = " << setw(2) << p.quality() << endl;
00196         }
00197       }
00198     }
00199 
00200     cout << "Eta values and fine bits : " << endl;
00201     for ( int i = 0; i < 12; i++ ) cout << setw(3) << m_eta[i] << " ";
00202     cout << endl;
00203     for ( int i = 0; i < 12; i++ ) cout << setw(3) << m_fine[i] << " ";
00204     cout << endl;
00205   }
00206 
00207 }
00208 
00209 
00210 //
00211 // receive data ( 15*3 DTBX eta trigger primitives )
00212 //
00213 void L1MuDTEtaProcessor::receiveData(int bx, const edm::Event& e) {
00214 
00215   edm::Handle<L1MuDTChambThContainer> dttrig;
00216   e.getByLabel(L1MuDTTFConfig::getDTDigiInputTag(),dttrig);
00217 
00218   // const int bx_offset = dttrig->correctBX();
00219   int bx_offset=0;
00220   bx = bx + bx_offset;
00221   
00222   //
00223   // get 5*3 eta track segments
00224   //
00225   int sector = m_epid;
00226   for ( int stat = 1; stat <= 3; stat++ ) {
00227     for ( int wheel = -2; wheel <= 2; wheel++ ) {
00228       L1MuDTChambThDigi* tseta = dttrig->chThetaSegm(wheel,stat,sector,bx);
00229       bitset<7> pos;
00230       bitset<7> qual;
00231 
00232       if ( tseta ) {
00233         if ( wheel == -2 || wheel == -1 || 
00234              ( wheel == 0 && (sector == 0 || sector == 3 || sector == 4 || sector == 7 || sector == 8 || sector == 11) ) ) {
00235           for ( int i = 0; i < 7; i++ ) {
00236             if ( tseta->position(i) ) pos.set(6-i);
00237             if ( tseta->quality(i) ) qual.set(6-i);
00238           }
00239         } else {
00240           for ( int i = 0; i < 7; i++ ) {
00241             if ( tseta->position(i) ) pos.set(i);
00242             if ( tseta->quality(i) ) qual.set(i);
00243           }
00244         }
00245       }
00246 
00247       const L1MuDTTrackSegEta* tmpts = new L1MuDTTrackSegEta(wheel,sector,stat,pos.to_ulong(),qual.to_ulong(),bx-bx_offset);
00248       m_tseta.push_back(tmpts);
00249     } 
00250   }
00251   
00252 }
00253 
00254 
00255 //
00256 // receive track addresses from 6 Sector Processors
00257 //
00258 void L1MuDTEtaProcessor::receiveAddresses() {
00259 
00260   // get track address code of all track segments
00261   // 6*2 times 5 bits; valid range [1-22]
00262 
00263   int sector = m_epid;
00264   
00265   int i = 0;
00266   for ( int wheel = -3; wheel <= 3; wheel++ ) {
00267     if ( wheel == 0 ) continue;
00268     L1MuDTSecProcId tmpspid(wheel,sector);
00269     for ( int number = 0; number < 2; number++ ) { 
00270       const L1MuDTTrack* cand = m_tf.sp(tmpspid)->track(number);
00271       const L1MuDTTrack* canD = m_tf.sp(tmpspid)->tracK(number);
00272       if ( cand ) {
00273         m_address[i] = cand->address().trackAddressCode();
00274         if ( !cand->empty() ) {
00275           m_TrackCand[i] = const_cast<L1MuDTTrack*>(cand);
00276           m_TracKCand[i] = const_cast<L1MuDTTrack*>(canD);
00277         }
00278         i++;
00279       }
00280     }
00281   }
00282   
00283 }  
00284 
00285 
00286 //
00287 // run Eta Track Finder (ETF)
00288 //
00289 void L1MuDTEtaProcessor::runEtaTrackFinder(const edm::EventSetup& c) {
00290 
00291   c.get< L1MuDTEtaPatternLutRcd >().get( theEtaPatternLUT );
00292 
00293   // check if there are any data
00294   bool empty = true;
00295   for ( int i = 0; i < 15; i++ ) {
00296     empty &= m_tseta[i]->empty();
00297   }
00298   if ( empty ) return;
00299 
00300   // Pattern comparator: 
00301   // loop over all patterns and compare with local chamber pattern
00302   // result : list of valid pattern IDs ( m_foundPattern )
00303   L1MuDTEtaPatternLut::ETFLut_iter it = theEtaPatternLUT->begin();
00304   while ( it != theEtaPatternLUT->end() ) {
00305   
00306     const L1MuDTEtaPattern pattern = (*it).second;
00307     int qualitycode = pattern.quality();
00308 
00309     bool good = true;
00310 
00311     for ( int station = 0; station < 3; station++) {
00312       int q = quality(qualitycode,station+1);
00313       int wheel = pattern.wheel(station+1);
00314       int bin = pattern.position(station+1);
00315       if ( bin == 0 ) continue;
00316       bitset<7> pos  = m_tseta[wheel+2 + station*5]->position();
00317       bitset<7> qual = m_tseta[wheel+2 + station*5]->quality();
00318       // compare position
00319       good &= pos.test(bin-1);
00320       // compare quality
00321       if ( q == 2 ) good &= qual.test(bin-1);  
00322     }
00323 
00324     if ( good ) m_foundPattern.push_back(pattern.id());
00325 
00326     it++;
00327     
00328   }
00329 
00330 }
00331 
00332 
00333 //
00334 // run Eta Matching Unit (EMU)
00335 //
00336 void L1MuDTEtaProcessor::runEtaMatchingUnit(const edm::EventSetup& c) {
00337   
00338   c.get< L1MuDTQualPatternLutRcd >().get( theQualPatternLUT );
00339 
00340   // loop over all addresses
00341   for ( int i = 0; i < 12; i++ ) {
00342   
00343     int adr = m_address[i];
00344     if ( adr == 0 ) continue;
00345     int sp = i/2 + 1;       //sector processor [1,6]
00346     
00347     // assign coarse eta value
00348     m_eta[i] = theQualPatternLUT->getCoarseEta(sp,adr);
00349     
00350     if ( m_foundPattern.empty() ) continue;
00351     
00352     // get list of qualified patterns ordered by quality 
00353     // and compare with found patterns
00354     const vector<short>& qualifiedPatterns = theQualPatternLUT->getQualifiedPatterns(sp,adr);
00355     vector<short>::const_iterator iter;
00356     vector<int>::const_iterator f_iter;
00357     for ( iter = qualifiedPatterns.begin(); iter != qualifiedPatterns.end(); iter++ ) {
00358       f_iter = find(m_foundPattern.begin(),m_foundPattern.end(),(*iter));
00359       // found
00360       if ( f_iter != m_foundPattern.end() ) {
00361         const L1MuDTEtaPattern p = theEtaPatternLUT->getPattern(*f_iter);
00362         // assign fine eta value     
00363         m_fine[i] = true;
00364         m_eta[i]  = p.eta();  // improved eta
00365         m_pattern[i] = (*f_iter);
00366         break;
00367       }
00368     }
00369     
00370   }
00371   
00372   // if both tracks from one sector processor deliver the same track address
00373   // the second track gets only a coarse eta value!  
00374   
00375   // loop over sector processors
00376   for ( int i = 0; i < 6; i++ ) {
00377     int idx1 = 2*i;
00378     int idx2 = 2*i+1;
00379     int adr1 = m_address[idx1];
00380     int adr2 = m_address[idx2];
00381     if ( adr1 == 0 || adr2 == 0 ) continue;
00382     if ( adr1 == adr2 ) {
00383       // second track gets coarse (default) eta value
00384       m_eta[idx2]  = theQualPatternLUT->getCoarseEta(i+1,adr2);
00385       m_pattern[idx2] = 0;
00386       m_fine[idx2] = false; 
00387     }  
00388   }
00389 
00390 }  
00391 
00392 
00393 //
00394 // assign values to track candidates
00395 //
00396 void L1MuDTEtaProcessor::assign() {
00397 
00398   for ( int i = 0; i < 12; i++ ) {
00399     if ( m_TrackCand[i] ) {
00400       if ( m_eta[i] != 99 ) { 
00401         m_TrackCand[i]->setEta(m_eta[i]);
00402         m_TracKCand[i]->setEta(m_eta[i]);
00403       }
00404       else {  
00405         if ( i/2 != 2 ) cerr << "L1MuDTEtaProcessor: assign invalid eta" << " " << m_address[i] << endl;
00406       }
00407       if ( m_fine[i] ) {
00408         m_TrackCand[i]->setFineEtaBit();
00409         m_TracKCand[i]->setFineEtaBit();
00410         // find all contributing track segments 
00411         const L1MuDTEtaPattern p = theEtaPatternLUT->getPattern(m_pattern[i]);
00412         vector<const L1MuDTTrackSegEta*> TSeta;
00413         const L1MuDTTrackSegEta* ts = 0;
00414         for ( int stat = 0; stat < 3; stat++ ) {
00415           int wh = p.wheel(stat+1);
00416           int pos = p.position(stat+1);
00417           if ( pos == 0 ) continue;
00418           ts = m_tseta[wh+2 + stat*5];
00419           TSeta.push_back(ts);
00420         }
00421         m_TrackCand[i]->setTSeta(TSeta);
00422         m_TracKCand[i]->setTSeta(TSeta);
00423       }  
00424     }
00425   }
00426 
00427 }
00428   
00429 
00430 //
00431 // get quality:  id [0,26], stat [1,3]
00432 //    
00433 int L1MuDTEtaProcessor::quality(int id, int stat) {
00434 
00435   // quality codes as defined in CMS Note 2001/027
00436   // This QualityPatterns are used to have a defined Quality-Identifier for
00437   // all possible found tracks.
00438   // Therefore three integer numbers ( Station 1, 2, 3 from left to right )
00439   // can have numbers like 0, 1 or 2
00440   //    0 ... no hit is given
00441   //    1 ... a LTRG is given
00442   //    2 ... a HTRG is given
00443 
00444   const int qualcode[27][3] = { {0,0,0},{1,0,0},{0,1,0},{0,0,1},{2,0,0},
00445                                 {0,2,0},{0,0,2},{1,1,0},{1,0,1},{0,1,1},
00446                                 {2,1,0},{1,2,0},{2,0,1},{1,0,2},{0,2,1},
00447                                 {0,1,2},{2,2,0},{2,0,2},{0,2,2},{1,1,1},
00448                                 {2,1,1},{1,2,1},{1,1,2},{2,2,1},{2,1,2},
00449                                 {1,2,2},{2,2,2} };
00450 
00451   return qualcode[id][stat-1];
00452 
00453 }

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