00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016 #include "EventFilter/DTTFRawToDigi/interface/DTTFFEDReader.h"
00017
00018 #include "DataFormats/Common/interface/Handle.h"
00019 #include "FWCore/Framework/interface/Event.h"
00020 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00021
00022 #include <DataFormats/FEDRawData/interface/FEDRawData.h>
00023 #include <DataFormats/FEDRawData/interface/FEDRawDataCollection.h>
00024
00025 #include <iostream>
00026
00027 using namespace std;
00028
00029
00030 DTTFFEDReader::DTTFFEDReader(const edm::ParameterSet& pset) {
00031
00032 produces<L1MuDTChambPhContainer>();
00033 produces<L1MuDTChambThContainer>();
00034 produces<L1MuDTTrackContainer>("DATA");
00035
00036 DTTFInputTag = pset.getParameter<edm::InputTag>("DTTF_FED_Source");
00037
00038 verbose_ = pset.getUntrackedParameter<bool>("verbose",false);
00039
00040 }
00041
00042 DTTFFEDReader::~DTTFFEDReader(){}
00043
00044 void DTTFFEDReader::produce(edm::Event& e, const edm::EventSetup& c) {
00045
00046 auto_ptr<L1MuDTChambPhContainer> phi_product(new L1MuDTChambPhContainer);
00047 auto_ptr<L1MuDTChambThContainer> the_product(new L1MuDTChambThContainer);
00048 auto_ptr<L1MuDTTrackContainer> tra_product(new L1MuDTTrackContainer);
00049
00050 L1MuDTChambPhContainer::Phi_Container phi_data;
00051 L1MuDTChambThContainer::The_Container the_data;
00052 L1MuDTTrackContainer::TrackContainer tra_data;
00053
00054 if (!fillRawData(e, phi_data, the_data, tra_data)) return;
00055
00056 phi_product->setContainer(phi_data);
00057 the_product->setContainer(the_data);
00058 tra_product->setContainer(tra_data);
00059
00060 e.put(phi_product);
00061 e.put(the_product);
00062 e.put(tra_product,"DATA");
00063
00064 }
00065
00066 bool DTTFFEDReader::fillRawData(edm::Event& e,
00067 L1MuDTChambPhContainer::Phi_Container& phi_data,
00068 L1MuDTChambThContainer::The_Container& the_data,
00069 L1MuDTTrackContainer::TrackContainer& tra_data) {
00070
00071 analyse(e);
00072
00073 phi_data = p_data();
00074 the_data = t_data();
00075 tra_data = k_data();
00076
00077 return true;
00078 }
00079
00080
00081
00082
00083 void DTTFFEDReader::analyse(edm::Event& e) {
00084 clear();
00085 process(e);
00086 match();
00087 return;
00088 }
00089
00090
00091 void DTTFFEDReader::process(edm::Event& e) {
00092
00093
00094 vector<int> DTTFWordContainer;
00095 vector<int>::iterator DTTFiterator;
00096
00097
00098 int BOEevTy, DTTFId;
00099
00100
00101 int DTTFWord;
00102 int DTTFChan, bitsID;
00103 int addr1[2] = { 3, 3};
00104 int addr2[2] = {15, 15};
00105 int addr3[2] = {15, 15};
00106 int addr4[2] = {15, 15};
00107
00108
00109 int evtLgth , CRC;
00110
00111
00112
00113 edm::Handle<FEDRawDataCollection> data;
00114 e.getByLabel(getDTTFInputTag(),data);
00115 FEDRawData dttfdata = data->FEDData(0x030C);
00116 if ( dttfdata.size() == 0 ) return;
00117
00118 int* dataWord1 = new int;
00119 int* dataWord2 = new int;
00120 unsigned char* LineFED=dttfdata.data();
00121 *dataWord2=*((int*)LineFED);
00122 LineFED+=4;
00123 *dataWord1=*((int*)LineFED);
00124 int lines = 1;
00125
00126 BOEevTy = ((*dataWord1)&0xFF000000)>>24;
00127 DTTFId = ((*dataWord2)&0x000FFF00)>>8;
00128
00129 if( (BOEevTy != 0x50) || ( DTTFId != 0x030C) ){
00130 if ( verbose_ ) edm::LogWarning("dttf_unpacker")
00131 << "Not a DTTF header " << hex << *dataWord1;
00132 delete dataWord1;
00133 delete dataWord2;
00134 return;
00135 }
00136
00137 int newCRC = 0xFFFF;
00138 calcCRC(*dataWord1, *dataWord2, newCRC);
00139
00140
00141
00142
00143 LineFED+=4;
00144 *dataWord2=*((int*)LineFED);
00145 LineFED+=4;
00146 *dataWord1=*((int*)LineFED);
00147 int chkEOE = ((*dataWord1)&0xFFF00000)>>20;
00148 lines++;
00149
00150 while(chkEOE != 0xA00){
00151
00152 calcCRC(*dataWord1, *dataWord2, newCRC);
00153
00154 DTTFWord = *dataWord1;
00155 DTTFWordContainer.push_back(DTTFWord);
00156 DTTFWord = *dataWord2;
00157 DTTFWordContainer.push_back(DTTFWord);
00158
00159 LineFED+=4;
00160 *dataWord2=*((int*)LineFED);
00161 LineFED+=4;
00162 *dataWord1=*((int*)LineFED);
00163 chkEOE = ((*dataWord1)&0xFFF00000)>>20;
00164 lines++;
00165
00166 if(lines > 3026){
00167 if ( verbose_ ) edm::LogWarning("dttf_unpacker")
00168 << "Warning : number of DTTF lines > 3026 ";
00169 delete dataWord1;
00170 delete dataWord2;
00171 return;
00172 }
00173
00174 }
00175
00176
00177
00178
00179 evtLgth = ((*dataWord1)&0x00FFFFFF);
00180 CRC = ((*dataWord2)&0xFFFF0000)>>16;
00181
00182 calcCRC(*dataWord1, (*dataWord2)&0xFFFF, newCRC);
00183 if( newCRC != CRC){
00184 if ( verbose_ ) edm::LogWarning("dttf_unpacker")
00185 << "Calculated CRC " << hex << newCRC << " differs from CRC in trailer " << hex << CRC;
00186 delete dataWord1;
00187 delete dataWord2;
00188 return;
00189 }
00190
00191 if( lines != evtLgth){
00192 if ( verbose_ ) edm::LogWarning("dttf_unpacker")
00193 << "Number of words read != event length " << dec << lines << " " << evtLgth;
00194 delete dataWord1;
00195 delete dataWord2;
00196 return;
00197 }
00198
00199
00200
00201
00202 for( DTTFiterator = DTTFWordContainer.begin();
00203 DTTFiterator != DTTFWordContainer.end();
00204 DTTFiterator++ ){
00205
00206 DTTFChan = ((*DTTFiterator)&0xFF000000)>>24;
00207 DTTFiterator++;
00208 bitsID = ((*DTTFiterator)&0xF0000000)>>28;
00209
00210 int bxID = bxNr(DTTFChan);
00211 if(bxID == -999) continue;
00212 int wheelID = wheel(DTTFChan);
00213 if(wheelID == -999) continue;
00214 int sectorID = sector(DTTFChan);
00215 if(sectorID == -999) continue;
00216
00217
00218 if(wheelID!=0 && bitsID<=0x9){
00219
00220 int wheelPh = (abs(wheelID)-1)*wheelID/abs(wheelID);
00221 int stationID = 0;
00222 int ra = 0;
00223 int ba = 0;
00224 int tsqual = 0;
00225 int ts2tag = 0;
00226
00227 if ( ( bitsID >> 1 ) == 0 ){ stationID = 1;}
00228 if ( ( bitsID >> 1 ) == 1 ){ stationID = 2;}
00229 if ( ( bitsID >> 1 ) == 4 ){ stationID = 3;}
00230 if ( ( bitsID >> 1 ) == 2 ){ stationID = 4;}
00231
00232 if(stationID != 3){
00233
00234 ts2tag = (bitsID)&0x1;
00235 tsqual = (~(*DTTFiterator)&0x07)-1;
00236 ba = (~(*DTTFiterator)&0x1FF8)>>3;
00237 if( ba>0x1FF) ba-=0x400;
00238 ra = (~(*DTTFiterator)&0x1FFE000)>>13;
00239 if( ra>0x7FF) ra-=0x1000;
00240 }
00241 else{
00242
00243 ts2tag = (bitsID)&0x1;
00244 tsqual = (~(*DTTFiterator)&0x07)-1;
00245 ra = (~(*DTTFiterator)&0x7FF8)>>3;
00246 if( ra>0x7FF) ra-=0x1000;
00247 }
00248
00249 if(tsqual!=7 && wheelID!=-1){
00250 phiSegments.push_back(
00251 L1MuDTChambPhDigi( bxID+ts2tag, wheelPh, sectorID, stationID,
00252 ra, ba, tsqual, ts2tag, 0) );
00253 }
00254 }
00255
00256
00257
00258 if(wheelID==0 && bitsID<=0x4){
00259
00260 int wheelTh = bitsID-2;
00261
00262 int posALL, posBTI[7];
00263
00264 if ( wheelTh == -2 || wheelTh == -1 ||
00265 ( wheelTh == 0 && (sectorID == 0 || sectorID == 3 || sectorID == 4 || sectorID == 7 || sectorID == 8 || sectorID == 11) ) ) {
00266
00267 posALL = ~(*DTTFiterator)&0x7F;
00268 posBTI[0] = ~(*DTTFiterator)&0x01;
00269 posBTI[1] = (~(*DTTFiterator)&0x02)>>1;
00270 posBTI[2] = (~(*DTTFiterator)&0x04)>>2;
00271 posBTI[3] = (~(*DTTFiterator)&0x08)>>3;
00272 posBTI[4] = (~(*DTTFiterator)&0x10)>>4;
00273 posBTI[5] = (~(*DTTFiterator)&0x20)>>5;
00274 posBTI[6] = (~(*DTTFiterator)&0x40)>>6;
00275
00276 if(posALL){
00277 theSegments.push_back(
00278 L1MuDTChambThDigi( bxID, wheelTh, sectorID, 1, posBTI) );
00279 }
00280
00281 posALL = ~(*DTTFiterator)&0x3F80;
00282 posBTI[0] = (~(*DTTFiterator)&0x0080)>>7;
00283 posBTI[1] = (~(*DTTFiterator)&0x0100)>>8;
00284 posBTI[2] = (~(*DTTFiterator)&0x0200)>>9;
00285 posBTI[3] = (~(*DTTFiterator)&0x0400)>>10;
00286 posBTI[4] = (~(*DTTFiterator)&0x0800)>>11;
00287 posBTI[5] = (~(*DTTFiterator)&0x1000)>>12;
00288 posBTI[6] = (~(*DTTFiterator)&0x2000)>>13;
00289
00290 if(posALL){
00291 theSegments.push_back(
00292 L1MuDTChambThDigi( bxID, wheelTh, sectorID, 2, posBTI) );
00293 }
00294
00295 posALL = ~(*DTTFiterator)&0x1FC000;
00296 posBTI[0] = (~(*DTTFiterator)&0x004000)>>14;
00297 posBTI[1] = (~(*DTTFiterator)&0x008000)>>15;
00298 posBTI[2] = (~(*DTTFiterator)&0x010000)>>16;
00299 posBTI[3] = (~(*DTTFiterator)&0x020000)>>17;
00300 posBTI[4] = (~(*DTTFiterator)&0x040000)>>18;
00301 posBTI[5] = (~(*DTTFiterator)&0x080000)>>19;
00302 posBTI[6] = (~(*DTTFiterator)&0x100000)>>20;
00303
00304 if(posALL){
00305 theSegments.push_back(
00306 L1MuDTChambThDigi( bxID, wheelTh, sectorID, 3, posBTI) );
00307 }
00308 }
00309
00310 else {
00311
00312 posALL = ~(*DTTFiterator)&0x7F;
00313 posBTI[6] = ~(*DTTFiterator)&0x01;
00314 posBTI[5] = (~(*DTTFiterator)&0x02)>>1;
00315 posBTI[4] = (~(*DTTFiterator)&0x04)>>2;
00316 posBTI[3] = (~(*DTTFiterator)&0x08)>>3;
00317 posBTI[2] = (~(*DTTFiterator)&0x10)>>4;
00318 posBTI[1] = (~(*DTTFiterator)&0x20)>>5;
00319 posBTI[0] = (~(*DTTFiterator)&0x40)>>6;
00320
00321 if(posALL){
00322 theSegments.push_back(
00323 L1MuDTChambThDigi( bxID, wheelTh, sectorID, 1, posBTI) );
00324 }
00325
00326 posALL = ~(*DTTFiterator)&0x3F80;
00327 posBTI[6] = (~(*DTTFiterator)&0x0080)>>7;
00328 posBTI[5] = (~(*DTTFiterator)&0x0100)>>8;
00329 posBTI[4] = (~(*DTTFiterator)&0x0200)>>9;
00330 posBTI[3] = (~(*DTTFiterator)&0x0400)>>10;
00331 posBTI[2] = (~(*DTTFiterator)&0x0800)>>11;
00332 posBTI[1] = (~(*DTTFiterator)&0x1000)>>12;
00333 posBTI[0] = (~(*DTTFiterator)&0x2000)>>13;
00334
00335 if(posALL){
00336 theSegments.push_back(
00337 L1MuDTChambThDigi( bxID, wheelTh, sectorID, 2, posBTI) );
00338 }
00339
00340 posALL = ~(*DTTFiterator)&0x1FC000;
00341 posBTI[6] = (~(*DTTFiterator)&0x004000)>>14;
00342 posBTI[5] = (~(*DTTFiterator)&0x008000)>>15;
00343 posBTI[4] = (~(*DTTFiterator)&0x010000)>>16;
00344 posBTI[3] = (~(*DTTFiterator)&0x020000)>>17;
00345 posBTI[2] = (~(*DTTFiterator)&0x040000)>>18;
00346 posBTI[1] = (~(*DTTFiterator)&0x080000)>>19;
00347 posBTI[0] = (~(*DTTFiterator)&0x100000)>>20;
00348
00349 if(posALL){
00350 theSegments.push_back(
00351 L1MuDTChambThDigi( bxID, wheelTh, sectorID, 3, posBTI) );
00352 }
00353 }
00354 }
00355
00356
00357
00358 if(wheelID!=0 && bitsID>=0xA && bitsID<=0xB){
00359
00360 int candID = bitsID - 0xA;
00361
00362 addr4[candID] = ((*DTTFiterator)&0x0F);
00363 addr3[candID] = ((*DTTFiterator)&0xF0)>>4;
00364 addr2[candID] = ((*DTTFiterator)&0xF00)>>8;
00365 addr1[candID] = ((*DTTFiterator)&0x3000)>>12;
00366 }
00367
00368
00369
00370 if(wheelID!=0 && bitsID>=0xC){
00371
00372 int muonID = 0;
00373 int pt = 0;
00374 int ch = 0;
00375 int phi = 0;
00376 int qual = 0;
00377
00378 muonID = (bitsID&0x1);
00379 qual = (~(*DTTFiterator)&0x07);
00380 phi = ((*DTTFiterator)&0x7F8)>>3;
00381 ch = (~(*DTTFiterator)&0x800)>>11;
00382 pt = (~(*DTTFiterator)&0x1F000)>>12;
00383
00384 if(qual!=0){
00385 dtTracks.push_back(
00386 L1MuDTTrackCand( 0, phi, 0, pt, ch, 1, 0, qual, bxID, wheelID, sectorID,
00387 muonID, addr1[muonID], addr2[muonID], addr3[muonID], addr4[muonID] ) );
00388 }
00389 }
00390
00391
00392
00393 if(wheelID==0 && bitsID>=0x8){
00394
00395 int wheelTh = bitsID&0x7;
00396
00397 int etaALL;
00398
00399 etaALL = ~(*DTTFiterator)&0x007F;
00400 if (etaALL) {
00401 etTrack[bxID+1][sectorID][wheelTh][0] = (*DTTFiterator)&0x003F;
00402 efTrack[bxID+1][sectorID][wheelTh][0] = (~(*DTTFiterator)&0x0040)>>6;
00403 }
00404
00405 etaALL = (~(*DTTFiterator)&0x3F80)>>7;
00406 if (etaALL) {
00407 etTrack[bxID+1][sectorID][wheelTh][1]= ( (*DTTFiterator)&0x1F80)>>7;
00408 efTrack[bxID+1][sectorID][wheelTh][1]= (~(*DTTFiterator)&0x2000)>>13;
00409 }
00410 }
00411
00412
00413 }
00414
00415 delete dataWord1;
00416 delete dataWord2;
00417 return;
00418 }
00419
00420 void DTTFFEDReader::match() {
00421
00422 for ( L1MuDTTrackContainer::TrackIterator i = dtTracks.begin();
00423 i != dtTracks.end();
00424 i++ ) {
00425 int bxTh = i->bx()+1;
00426 int sectorTh = i->scNum();
00427 int wheelTh = i->whNum()+3;
00428 if(wheelTh > 3) wheelTh-=1;
00429 int muonTh = i->TrkTag();
00430
00431 i->setEtaPacked(etTrack[bxTh][sectorTh][wheelTh][muonTh]);
00432 i->setFineHaloPacked(efTrack[bxTh][sectorTh][wheelTh][muonTh]);
00433 }
00434
00435 return;
00436 }
00437
00438
00439 const L1MuDTChambPhContainer::Phi_Container& DTTFFEDReader::p_data() {
00440 return phiSegments;
00441 }
00442
00443 const L1MuDTChambThContainer::The_Container& DTTFFEDReader::t_data() {
00444 return theSegments;
00445 }
00446
00447 const L1MuDTTrackContainer::TrackContainer& DTTFFEDReader::k_data() {
00448 return dtTracks;
00449 }
00450
00451 void DTTFFEDReader::clear() {
00452 phiSegments.clear();
00453 theSegments.clear();
00454 dtTracks.clear();
00455
00456 for(int i=0; i<3; i++){
00457 for(int j=0; j<12; j++){
00458 for(int k=0; k<6; k++){
00459 for(int l=0; l<2; l++){
00460 etTrack[i][j][k][l] = 0;
00461 efTrack[i][j][k][l] = 0;
00462 }
00463 }
00464 }
00465 }
00466
00467 return;
00468 }
00469
00470 int DTTFFEDReader::channel( int wheel, int sector, int bx ){
00471
00472
00473
00474
00475
00476
00477 int myChannel = 255;
00478
00479 if ( abs(bx) > 1) { return myChannel; }
00480 if ( sector < 0 || sector > 11) { return myChannel; }
00481 if ( abs(wheel) > 3) { return myChannel; }
00482
00483 myChannel = sector*21 + wheel*3 - bx + 10 ;
00484
00485 if (myChannel > 125) myChannel += 2;
00486
00487 return myChannel;
00488 }
00489
00490 int DTTFFEDReader::bxNr( int channel ){
00491
00492 int myChannel = channel;
00493
00494 if (myChannel > 127) myChannel -= 2;
00495
00496 if (myChannel < 0 || myChannel > 251 ){ return -999; }
00497
00498 int myBx = 1-(myChannel%3);
00499
00500 return myBx;
00501 }
00502
00503 int DTTFFEDReader::sector( int channel ){
00504
00505 int myChannel = channel;
00506
00507 if (myChannel > 127) myChannel -= 2;
00508
00509 if (myChannel < 0 || myChannel > 251 ){ return -999; }
00510
00511 return myChannel/21;
00512 }
00513
00514 int DTTFFEDReader::wheel( int channel ){
00515
00516 int myChannel = channel;
00517
00518 if (myChannel > 127) myChannel -= 2;
00519
00520 if (myChannel < 0 || myChannel > 251 ){ return -999; }
00521
00522 int myWheel = ((myChannel%21)/3)-3;
00523
00524 return myWheel;
00525 }
00526
00527 void DTTFFEDReader::calcCRC(int myD1, int myD2, int &myC){
00528
00529 int myCRC[16],D[64],C[16];
00530
00531 for( int i=0; i < 32; i++ ){ D[i]=(myD2>>i)&0x1; }
00532 for( int i=0; i < 32; i++ ){ D[i+32]=(myD1>>i)&0x1; }
00533 for( int i=0; i < 16; i++ ){ C[i]=(myC>>i)&0x1; }
00534
00535 myCRC[0] = ( D[63] + D[62] + D[61] + D[60] + D[55] + D[54] +
00536 D[53] + D[52] + D[51] + D[50] + D[49] + D[48] +
00537 D[47] + D[46] + D[45] + D[43] + D[41] + D[40] +
00538 D[39] + D[38] + D[37] + D[36] + D[35] + D[34] +
00539 D[33] + D[32] + D[31] + D[30] + D[27] + D[26] +
00540 D[25] + D[24] + D[23] + D[22] + D[21] + D[20] +
00541 D[19] + D[18] + D[17] + D[16] + D[15] + D[13] +
00542 D[12] + D[11] + D[10] + D[9] + D[8] + D[7] +
00543 D[6] + D[5] + D[4] + D[3] + D[2] + D[1] +
00544 D[0] + C[0] + C[1] + C[2] + C[3] + C[4] +
00545 C[5] + C[6] + C[7] + C[12] + C[13] + C[14] +
00546 C[15] )%2;
00547
00548 myCRC[1] = ( D[63] + D[62] + D[61] + D[56] + D[55] + D[54] +
00549 D[53] + D[52] + D[51] + D[50] + D[49] + D[48] +
00550 D[47] + D[46] + D[44] + D[42] + D[41] + D[40] +
00551 D[39] + D[38] + D[37] + D[36] + D[35] + D[34] +
00552 D[33] + D[32] + D[31] + D[28] + D[27] + D[26] +
00553 D[25] + D[24] + D[23] + D[22] + D[21] + D[20] +
00554 D[19] + D[18] + D[17] + D[16] + D[14] + D[13] +
00555 D[12] + D[11] + D[10] + D[9] + D[8] + D[7] +
00556 D[6] + D[5] + D[4] + D[3] + D[2] + D[1] +
00557 C[0] + C[1] + C[2] + C[3] + C[4] + C[5] +
00558 C[6] + C[7] + C[8] + C[13] + C[14] + C[15] )%2;
00559
00560 myCRC[2] = ( D[61] + D[60] + D[57] + D[56] + D[46] + D[42] +
00561 D[31] + D[30] + D[29] + D[28] + D[16] + D[14] +
00562 D[1] + D[0] + C[8] + C[9] + C[12] + C[13] )%2;
00563
00564 myCRC[3] = ( D[62] + D[61] + D[58] + D[57] + D[47] + D[43] +
00565 D[32] + D[31] + D[30] + D[29] + D[17] + D[15] +
00566 D[2] + D[1] + C[9] + C[10] + C[13] + C[14] )%2;
00567
00568 myCRC[4] = ( D[63] + D[62] + D[59] + D[58] + D[48] + D[44] +
00569 D[33] + D[32] + D[31] + D[30] + D[18] + D[16] +
00570 D[3] + D[2] + C[0] + C[10] + C[11] + C[14] +
00571 C[15] )%2;
00572
00573 myCRC[5] = ( D[63] + D[60] + D[59] + D[49] + D[45] + D[34] +
00574 D[33] + D[32] + D[31] + D[19] + D[17] + D[4] +
00575 D[3] + C[1] + C[11] + C[12] + C[15] )%2;
00576
00577 myCRC[6] = ( D[61] + D[60] + D[50] + D[46] + D[35] + D[34] +
00578 D[33] + D[32] + D[20] + D[18] + D[5] + D[4] +
00579 C[2] + C[12] + C[13] )%2;
00580
00581 myCRC[7] = ( D[62] + D[61] + D[51] + D[47] + D[36] + D[35] +
00582 D[34] + D[33] + D[21] + D[19] + D[6] + D[5] +
00583 C[3] + C[13] + C[14] )%2;
00584
00585 myCRC[8] = ( D[63] + D[62] + D[52] + D[48] + D[37] + D[36] +
00586 D[35] + D[34] + D[22] + D[20] + D[7] + D[6] +
00587 C[0] + C[4] + C[14] + C[15] )%2;
00588
00589 myCRC[9] = ( D[63] + D[53] + D[49] + D[38] + D[37] + D[36] +
00590 D[35] + D[23] + D[21] + D[8] + D[7] + C[1] +
00591 C[5] + C[15] )%2;
00592
00593 myCRC[10] = ( D[54] + D[50] + D[39] + D[38] + D[37] + D[36] +
00594 D[24] + D[22] + D[9] + D[8] + C[2] + C[6] )%2;
00595
00596 myCRC[11] = ( D[55] + D[51] + D[40] + D[39] + D[38] + D[37] +
00597 D[25] + D[23] + D[10] + D[9] + C[3] + C[7] )%2;
00598
00599 myCRC[12] = ( D[56] + D[52] + D[41] + D[40] + D[39] + D[38] +
00600 D[26] + D[24] + D[11] + D[10] + C[4] + C[8] )%2;
00601
00602 myCRC[13] = ( D[57] + D[53] + D[42] + D[41] + D[40] + D[39] +
00603 D[27] + D[25] + D[12] + D[11] + C[5] + C[9] )%2;
00604
00605 myCRC[14] = ( D[58] + D[54] + D[43] + D[42] + D[41] + D[40] +
00606 D[28] + D[26] + D[13] + D[12] + C[6] + C[10] )%2;
00607
00608 myCRC[15] = ( D[63] + D[62] + D[61] + D[60] + D[59] + D[54] +
00609 D[53] + D[52] + D[51] + D[50] + D[49] + D[48] +
00610 D[47] + D[46] + D[45] + D[44] + D[42] + D[40] +
00611 D[39] + D[38] + D[37] + D[36] + D[35] + D[34] +
00612 D[33] + D[32] + D[31] + D[30] + D[29] + D[26] +
00613 D[25] + D[24] + D[23] + D[22] + D[21] + D[20] +
00614 D[19] + D[18] + D[17] + D[16] + D[15] + D[14] +
00615 D[12] + D[11] + D[10] + D[9] + D[8] + D[7] +
00616 D[6] + D[5] + D[4] + D[3] + D[2] + D[1] +
00617 D[0] + C[0] + C[1] + C[2] + C[3] + C[4] +
00618 C[5] + C[6] + C[11] + C[12] + C[13] + C[14] +
00619 C[15] )%2;
00620
00621 int tempC = 0x0;
00622 for(int i=0; i<16 ; i++){ tempC = tempC + (myCRC[i]<<i); }
00623 myC = tempC ;
00624 return;
00625 }