CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/EventFilter/DTTFRawToDigi/src/DTTFFEDSim.cc

Go to the documentation of this file.
00001 //-------------------------------------------------
00002 //
00003 //   Class: DTTFFEDSim
00004 //
00005 //   L1 DT Track Finder Digi-to-Raw
00006 //
00007 //
00008 //   $Date: 2009/11/18 13:27:12 $
00009 //   $Revision: 1.13 $
00010 //
00011 //   Author :
00012 //   J. Troconiz  UAM Madrid
00013 //
00014 //--------------------------------------------------
00015 
00016 #include "EventFilter/DTTFRawToDigi/interface/DTTFFEDSim.h"
00017 
00018 #include "DataFormats/Common/interface/Handle.h"
00019 #include "FWCore/Framework/interface/Event.h"
00020 
00021 #include <DataFormats/FEDRawData/interface/FEDRawData.h>
00022 #include <DataFormats/L1DTTrackFinder/interface/L1MuDTChambPhContainer.h>
00023 #include <DataFormats/L1DTTrackFinder/interface/L1MuDTChambThContainer.h>
00024 #include <DataFormats/L1DTTrackFinder/interface/L1MuDTTrackContainer.h>
00025 
00026 #include <iostream>
00027 
00028 using namespace std;
00029 
00030 DTTFFEDSim::DTTFFEDSim(const edm::ParameterSet& pset) : eventNum(0) {
00031 
00032   produces<FEDRawDataCollection>();
00033 
00034   DTDigiInputTag = pset.getParameter<edm::InputTag>("DTDigi_Source");
00035   DTPHTFInputTag = pset.getParameter<edm::InputTag>("DTTracks_Source");
00036 
00037 }
00038 
00039 DTTFFEDSim::~DTTFFEDSim(){}
00040 
00041 void DTTFFEDSim::produce(edm::Event& e, const edm::EventSetup& c) {
00042 
00043   FEDRawDataCollection data;
00044 
00045   if (!fillRawData(e, data)) return;
00046 
00047   auto_ptr<FEDRawDataCollection> fed_product(new FEDRawDataCollection(data));
00048 
00049   e.put(fed_product);
00050 
00051 }
00052 
00053 bool DTTFFEDSim::fillRawData(edm::Event& e,
00054                              FEDRawDataCollection& data) {
00055 
00056   eventNum = e.id().event();
00057 
00058   int lines = 2;
00059 
00060   edm::Handle<L1MuDTChambPhContainer> phtrig;
00061   e.getByLabel(getDTDigiInputTag(),phtrig);
00062   lines += phtrig->bxSize(-1, 1);
00063 
00064   edm::Handle<L1MuDTChambThContainer> thtrig;
00065   e.getByLabel(getDTDigiInputTag(),thtrig);
00066   lines += thtrig->bxSize(-1, 1);
00067 
00068   edm::Handle<L1MuDTTrackContainer>   trtrig;
00069   e.getByLabel(getDTPHTFInputTag(),trtrig);
00070   lines += trtrig->bxSize(-1, 1)*3;
00071 
00072   FEDRawData& dttfdata = data.FEDData(0x30C);
00073   dttfdata.resize(lines*8); // size in bytes
00074   unsigned char* LineFED=dttfdata.data();
00075 
00076   int* dataWord1 = new int;
00077   int* dataWord2 = new int;
00078 
00079   //--> Header
00080 
00081   *dataWord1 = 0x50000000
00082              + (eventNum&0xFFFFFF);
00083   *dataWord2 = 0x00030C00;
00084 
00085   int newCRC =  0xFFFF;
00086   calcCRC(*dataWord1, *dataWord2, newCRC);
00087   *((int*)LineFED)=*dataWord2; 
00088   LineFED+=4;
00089   *((int*)LineFED)=*dataWord1; 
00090 
00091   //--> DTTF data 
00092 
00093   int TS1Id[4], TS2Id[4];   // word identifier for TS #1,#2 for stations 
00094   TS1Id[0] = 0x0E;      
00095   TS2Id[0] = 0x1E;
00096   TS1Id[1] = 0x2E;
00097   TS2Id[1] = 0x3E;
00098   TS1Id[3] = 0x4E;
00099   TS2Id[3] = 0x5E;
00100   TS1Id[2] = 0x8FFF8;
00101   TS2Id[2] = 0x9FFF8;
00102 
00103   //Input
00104   L1MuDTChambPhContainer::Phi_iterator tsphi;
00105 
00106   for ( tsphi =  phtrig->getContainer()->begin();
00107         tsphi != phtrig->getContainer()->end();
00108         tsphi++ ) {
00109     if ( tsphi->code() != 7 ) {
00110 
00111       int wheelID   = tsphi->whNum()+1;
00112       if ( wheelID <= 0 ) wheelID -= 2;
00113       int stationID = tsphi->stNum()-1;
00114       int is2nd     = tsphi->Ts2Tag();
00115 
00116       int channelNr = channel(wheelID, tsphi->scNum(), tsphi->bxNum()-is2nd);
00117       if ( channelNr == 255 ) continue;
00118       int TSId = ( is2nd == 0 ) ? TS1Id[stationID] : TS2Id[stationID];
00119 
00120       *dataWord1 = ((channelNr&0xFF)<<24)
00121                  + 0x00FFFFFF;
00122 
00123       if ( stationID != 2 ){
00124         *dataWord2 = (             (TSId&0x0FF)<<24)
00125                    + (~(tsphi->code()+1)&0x007)
00126                    + (   (~tsphi->phiB()&0x3FF)<<3) 
00127                    + (    (~tsphi->phi()&0xFFF)<<13);
00128       }
00129       else {
00130         *dataWord2 = (             (TSId&0xFFFFF)<<12) 
00131                    + (~(tsphi->code()+1)&0x00007)
00132                    + (    (~tsphi->phi()&0x00FFF)<<3);
00133       }
00134 
00135       calcCRC(*dataWord1, *dataWord2, newCRC);
00136       LineFED+=4;
00137       *((int*)LineFED)=*dataWord2; 
00138       LineFED+=4;
00139       *((int*)LineFED)=*dataWord1; 
00140     }
00141   }
00142   //Input
00143 
00144   //Input
00145   L1MuDTChambThContainer::The_iterator tsthe;
00146 
00147   for ( tsthe =  thtrig->getContainer()->begin();
00148         tsthe != thtrig->getContainer()->end();
00149         tsthe++ ) {
00150 
00151     int wheelTh  = tsthe->whNum();
00152     int sectorID = tsthe->scNum();
00153 
00154     int channelNr = channel(0, sectorID, tsthe->bxNum());
00155     if ( channelNr == 255 ) continue;
00156     int TSId = wheelTh+2;
00157 
00158     *dataWord1 = ((channelNr&0xFF)<<24)
00159                + 0x00FFFFFF;
00160 
00161     *dataWord2 = ((TSId&0x07)<<28)
00162                + 0x0FFFFFFF;
00163 
00164     int stationID = tsthe->stNum()-1;
00165     for ( int bti = 0; bti < 7; bti++ )
00166       if ( wheelTh == -2 || wheelTh == -1 || 
00167            ( wheelTh == 0 && (sectorID == 0 || sectorID == 3 || sectorID == 4 || sectorID == 7 || sectorID == 8 || sectorID == 11) ) )
00168         *dataWord2 -= (tsthe->position(bti)&0x1)<<(stationID*7+bti);
00169       else
00170         *dataWord2 -= (tsthe->position(6-bti)&0x1)<<(stationID*7+bti);
00171 
00172     calcCRC(*dataWord1, *dataWord2, newCRC);
00173     LineFED+=4;
00174     *((int*)LineFED)=*dataWord2; 
00175     LineFED+=4;
00176     *((int*)LineFED)=*dataWord1; 
00177   }
00178   //Input
00179 
00180   //Output
00181   L1MuDTTrackContainer::Trackiterator tstrk;
00182 
00183   for ( tstrk =  trtrig->getContainer()->begin();
00184         tstrk != trtrig->getContainer()->end();
00185         tstrk++ ) {
00186 
00187     int channelNr = channel(tstrk->whNum(), tstrk->scNum(), tstrk->bx());
00188     if ( channelNr == 255 ) continue;
00189     int TSId = ( tstrk->TrkTag() == 0 ) ? 0xAFFF : 0xBFFF;
00190 
00191     *dataWord1 = ((channelNr&0xFF)<<24)
00192                + 0x00FFFFFF;
00193 
00194     *dataWord2 = (           (TSId&0xFFFF)<<16)
00195                + ( tstrk->stNum(4)&0x0000F)
00196                + ((tstrk->stNum(3)&0x0000F)<<4)
00197                + ((tstrk->stNum(2)&0x0000F)<<8)
00198                + ((tstrk->stNum(1)&0x00003)<<12);
00199 
00200     calcCRC(*dataWord1, *dataWord2, newCRC);
00201     LineFED+=4;
00202     *((int*)LineFED)=*dataWord2; 
00203     LineFED+=4;
00204     *((int*)LineFED)=*dataWord1; 
00205 
00206     TSId = ( tstrk->TrkTag() == 0 ) ? 0xCFFE : 0xDFFE;
00207 
00208     *dataWord1 = ((channelNr&0xFF)<<24)
00209                + 0x00FFFFFF;
00210 
00211     *dataWord2 = (                    (TSId&0xFFFE)<<16)
00212                + ( ~tstrk->quality_packed()&0x0007)
00213                + (     (tstrk->phi_packed()&0x00FF)<<3)
00214                + ( (~tstrk->charge_packed()&0x0001)<<11)
00215                + (     (~tstrk->pt_packed()&0x001F)<<12);
00216 
00217     calcCRC(*dataWord1, *dataWord2, newCRC);
00218     LineFED+=4;
00219     *((int*)LineFED)=*dataWord2; 
00220     LineFED+=4;
00221     *((int*)LineFED)=*dataWord1; 
00222 
00223     channelNr = channel(0, tstrk->scNum(), tstrk->bx());
00224     if ( channelNr == 255 ) continue;
00225     TSId = (tstrk->whNum()+3)<<16;
00226     TSId += ( tstrk->whNum() < 0 ) ? 0x8FFFC : 0x7FFFC;
00227 
00228     *dataWord1 = ((channelNr&0xFF)<<24)
00229                + 0x00FFFFFF;
00230 
00231     *dataWord2 = (TSId&0xFFFFC)<<12;
00232 
00233     if ( tstrk->TrkTag() == 0 ) {
00234       *dataWord2 +=                            0x3F80
00235                  + (       tstrk->eta_packed()&0x003F)
00236                  + ((~tstrk->finehalo_packed()&0x0001)<<6);
00237     }
00238     else {
00239       *dataWord2 +=                            0x007F
00240                  + (     ( tstrk->eta_packed()&0x003F)<<7)
00241                  + ((~tstrk->finehalo_packed()&0x0001)<<13);
00242     }
00243 
00244     calcCRC(*dataWord1, *dataWord2, newCRC);
00245     LineFED+=4;
00246     *((int*)LineFED)=*dataWord2; 
00247     LineFED+=4;
00248     *((int*)LineFED)=*dataWord1; 
00249   }
00250   //Output
00251 
00252   //--> Trailer
00253 
00254   *dataWord1 = 0xA0000000
00255              + (lines&0xFFFFFF);
00256   *dataWord2 = 0;
00257 
00258   calcCRC(*dataWord1, *dataWord2&0xFFFF, newCRC);
00259 
00260   *dataWord2 += (newCRC&0xFFFF)<<16;
00261 
00262   LineFED+=4;
00263   *((int*)LineFED)=*dataWord2; 
00264   LineFED+=4;
00265   *((int*)LineFED)=*dataWord1; 
00266 
00267   delete dataWord1;
00268   delete dataWord2;
00269   return true;
00270 }
00271 
00272 int DTTFFEDSim::channel( int wheel, int sector,  int bx ){
00273 
00274   // wheel  :  -3 -2 -1 +1 +2 +3 <=> PHTF's : N2, N1, N0, P0, P1, P2
00275   //                           0 <=> ETTF  
00276   // sector :  0 -> 11
00277   // bx     : -1 -> +1
00278 
00279   int myChannel = 255;
00280 
00281   if ( abs(bx) > 1)               { return myChannel; }
00282   if ( sector < 0 || sector > 11) { return myChannel; }
00283   if ( abs(wheel) > 3)            { return myChannel; }
00284 
00285   myChannel = sector*21 + wheel*3 - bx + 10 ; 
00286 
00287   if (myChannel > 125) myChannel += 2;
00288 
00289   return myChannel;
00290 }
00291 
00292 int DTTFFEDSim::bxNr( int channel ){
00293 
00294   int myChannel = channel;
00295 
00296   if (myChannel > 127) myChannel -= 2;
00297 
00298   if (myChannel < 0 || myChannel > 251 ){ return -999; }
00299 
00300   int myBx = 1-(myChannel%3);
00301 
00302   return myBx;
00303 }
00304 
00305 int DTTFFEDSim::sector( int channel ){
00306 
00307   int myChannel = channel;
00308 
00309   if (myChannel > 127) myChannel -= 2;
00310 
00311   if (myChannel < 0 || myChannel > 251 ){ return -999; }
00312 
00313   return myChannel/21;
00314 }
00315 
00316 int DTTFFEDSim::wheel( int channel ){
00317 
00318   int myChannel = channel;
00319 
00320   if (myChannel > 127) myChannel -= 2;
00321 
00322   if (myChannel < 0 || myChannel > 251 ){ return -999; }
00323 
00324   int myWheel = ((myChannel%21)/3)-3;
00325 
00326   return myWheel;
00327 }
00328 
00329 void DTTFFEDSim::calcCRC(int myD1, int myD2, int &myC){
00330 
00331   int myCRC[16],D[64],C[16];
00332 
00333   for( int i=0; i < 32; i++ ){ D[i]=(myD2>>i)&0x1; }
00334   for( int i=0; i < 32; i++ ){ D[i+32]=(myD1>>i)&0x1; }
00335   for( int i=0; i < 16; i++ ){ C[i]=(myC>>i)&0x1; }
00336 
00337   myCRC[0] = ( D[63] + D[62] + D[61] + D[60] + D[55] + D[54] +
00338                D[53] + D[52] + D[51] + D[50] + D[49] + D[48] +
00339                D[47] + D[46] + D[45] + D[43] + D[41] + D[40] +
00340                D[39] + D[38] + D[37] + D[36] + D[35] + D[34] +
00341                D[33] + D[32] + D[31] + D[30] + D[27] + D[26] +
00342                D[25] + D[24] + D[23] + D[22] + D[21] + D[20] +
00343                D[19] + D[18] + D[17] + D[16] + D[15] + D[13] +
00344                D[12] + D[11] + D[10] + D[9]  + D[8]  + D[7]  +
00345                D[6]  + D[5]  + D[4]  + D[3]  + D[2]  + D[1]  +
00346                D[0]  + C[0]  + C[1]  + C[2]  + C[3]  + C[4]  +
00347                C[5]  + C[6]  + C[7]  + C[12] + C[13] + C[14] +
00348                C[15] )%2;
00349 
00350   myCRC[1] = ( D[63] + D[62] + D[61] + D[56] + D[55] + D[54] +
00351                D[53] + D[52] + D[51] + D[50] + D[49] + D[48] +
00352                D[47] + D[46] + D[44] + D[42] + D[41] + D[40] +
00353                D[39] + D[38] + D[37] + D[36] + D[35] + D[34] +
00354                D[33] + D[32] + D[31] + D[28] + D[27] + D[26] +
00355                D[25] + D[24] + D[23] + D[22] + D[21] + D[20] +
00356                D[19] + D[18] + D[17] + D[16] + D[14] + D[13] +
00357                D[12] + D[11] + D[10] + D[9]  + D[8]  + D[7]  +
00358                D[6]  + D[5]  + D[4]  + D[3]  + D[2]  + D[1]  +
00359                C[0]  + C[1]  + C[2]  + C[3]  + C[4]  + C[5]  +
00360                C[6]  + C[7]  + C[8]  + C[13] + C[14] + C[15] )%2;
00361 
00362   myCRC[2] = ( D[61] + D[60] + D[57] + D[56] + D[46] + D[42] +
00363                D[31] + D[30] + D[29] + D[28] + D[16] + D[14] +
00364                D[1]  + D[0]  + C[8]  + C[9]  + C[12] + C[13] )%2;
00365 
00366   myCRC[3] = ( D[62] + D[61] + D[58] + D[57] + D[47] + D[43] +
00367                D[32] + D[31] + D[30] + D[29] + D[17] + D[15] +
00368                D[2]  + D[1]  + C[9]  + C[10] + C[13] + C[14] )%2;
00369 
00370   myCRC[4] = ( D[63] + D[62] + D[59] + D[58] + D[48] + D[44] +
00371                D[33] + D[32] + D[31] + D[30] + D[18] + D[16] + 
00372                D[3]  + D[2]  + C[0]  + C[10] + C[11] + C[14] +
00373                C[15] )%2;
00374 
00375   myCRC[5] = ( D[63] + D[60] + D[59] + D[49] + D[45] + D[34] +
00376                D[33] + D[32] + D[31] + D[19] + D[17] + D[4]  +
00377                D[3]  + C[1]  + C[11] + C[12] + C[15] )%2;
00378 
00379   myCRC[6] = ( D[61] + D[60] + D[50] + D[46] + D[35] + D[34] +
00380                D[33] + D[32] + D[20] + D[18] + D[5]  + D[4]  +
00381                C[2]  + C[12] + C[13] )%2;
00382 
00383   myCRC[7] = ( D[62] + D[61] + D[51] + D[47] + D[36] + D[35] +
00384                D[34] + D[33] + D[21] + D[19] + D[6]  + D[5]  +
00385                C[3]  + C[13] + C[14] )%2;
00386 
00387   myCRC[8] = ( D[63] + D[62] + D[52] + D[48] + D[37] + D[36] +
00388                D[35] + D[34] + D[22] + D[20] + D[7]  + D[6]  +
00389                C[0]  + C[4]  + C[14] + C[15] )%2;
00390 
00391   myCRC[9] = ( D[63] + D[53] + D[49] + D[38] + D[37] + D[36] +
00392                D[35] + D[23] + D[21] + D[8]  + D[7]  + C[1]  +
00393                C[5]  + C[15] )%2;
00394 
00395   myCRC[10] = ( D[54] + D[50] + D[39] + D[38] + D[37] + D[36] + 
00396                 D[24] + D[22] + D[9]  + D[8]  + C[2]  + C[6] )%2;
00397 
00398   myCRC[11] = ( D[55] + D[51] + D[40] + D[39] + D[38] + D[37] +
00399                 D[25] + D[23] + D[10] + D[9]  + C[3]  + C[7] )%2;
00400 
00401   myCRC[12] = ( D[56] + D[52] + D[41] + D[40] + D[39] + D[38] +
00402                 D[26] + D[24] + D[11] + D[10] + C[4]  + C[8] )%2;
00403 
00404   myCRC[13] = ( D[57] + D[53] + D[42] + D[41] + D[40] + D[39] +
00405                 D[27] + D[25] + D[12] + D[11] + C[5]  + C[9] )%2;
00406 
00407   myCRC[14] = ( D[58] + D[54] + D[43] + D[42] + D[41] + D[40] +
00408                 D[28] + D[26] + D[13] + D[12] + C[6]  + C[10] )%2;
00409 
00410   myCRC[15] = ( D[63] + D[62] + D[61] + D[60] + D[59] + D[54] +
00411                 D[53] + D[52] + D[51] + D[50] + D[49] + D[48] + 
00412                 D[47] + D[46] + D[45] + D[44] + D[42] + D[40] +
00413                 D[39] + D[38] + D[37] + D[36] + D[35] + D[34] + 
00414                 D[33] + D[32] + D[31] + D[30] + D[29] + D[26] +
00415                 D[25] + D[24] + D[23] + D[22] + D[21] + D[20] + 
00416                 D[19] + D[18] + D[17] + D[16] + D[15] + D[14] +
00417                 D[12] + D[11] + D[10] + D[9]  + D[8]  + D[7]  + 
00418                 D[6]  + D[5]  + D[4]  + D[3]  + D[2]  + D[1]  +
00419                 D[0]  + C[0]  + C[1]  + C[2]  + C[3]  + C[4]  + 
00420                 C[5]  + C[6]  + C[11] + C[12] + C[13] + C[14] +
00421                 C[15] )%2;
00422 
00423   int tempC = 0x0;  
00424   for(int i=0; i<16 ; i++){ tempC = tempC + (myCRC[i]<<i); }
00425   myC = tempC ;
00426   return;
00427 }