• Main Page
  • Related Pages
  • Namespaces
  • Data Structures
  • Files
  • File List
  • Globals

beamtime/cern-oct11/go4/TMbsCrateProc.cxx (r4864/r2543)

Go to the documentation of this file.
00001 #include "TMbsCrateProc.h"
00002 
00003 #include <stdlib.h>
00004 
00005 #include "TTimeStamp.h"
00006 #include "TSystem.h"
00007 #include "TROOT.h"
00008 
00009 #include "TGo4Log.h"
00010 
00011 // still need this for subevent types:
00012 #include "roc/Message.h"
00013 #include "roc/Board.h"
00014 
00015 
00016 #define L1182_DEBUG 0
00017 #define L1183_DEBUG 0
00018 #define MADC_DEBUG 0
00019 #define V550_READ_PED 1
00020 
00021 
00022 TMbsCrateProc::TMbsCrateProc(const char* name) : TCBMBeamtimeProc(name),fOutputEvent(0)
00023 {
00024    TGo4Log::Info("TMbsCrateProc: Create instance %s", name);
00025 
00026    fPar = (TMbsCrateParam*) MakeParameter("MbsCratePar", "TMbsCrateParam");
00027 
00028    TString setupmacro="set_MbsPar.C";
00029    if (!gSystem->AccessPathName(setupmacro.Data())) {
00030       cout << "**** Executing MBS parameter setup script " << setupmacro.Data() << endl;
00031       gROOT->ProcessLine(Form(".x %s", setupmacro.Data()));
00032    }
00033 
00034    TString obname, obtitle;   
00035 
00036 
00037    /*********** MADCs *************/
00038    for (int mad = 0; mad < MAX_MADC; mad++) {
00039       obname.Form("MADC%d/MADC%d_sample", mad,mad);
00040       obtitle.Form("Sample of MADC%d readout", mad);
00041       fMadc_trace[mad]=MakeTH1('I', obname.Data(), obtitle.Data(), N_MADC_CHA, 0, N_MADC_CHA,"channel No.");
00042 
00043       for (int ch = 0; ch < N_MADC_CHA; ch++) {
00044          obname.Form("MADC%d/MADC%d_ch%d", mad,mad,ch);
00045          obtitle.Form("MADC %d Channel %d", mad, ch);
00046          fMadc_adc[mad][ch] = MakeTH1('I', obname.Data(), obtitle.Data(), 0x3FF, 0, 0x3FF);
00047 
00048          //          obname.Form("MADC%d/MADC%d_ch%d_pion", mad,mad,ch);
00049          //          obtitle.Form("MADC %d Channel %d pion", mad, ch);
00050          //          fMadc_pion[mad][ch] = MakeTH1('I', obname.Data(), obtitle.Data(), 0x3FF, 0, 0x3FF);
00051          //
00052          //          obname.Form("MADC%d/MADC%d_ch%d_electron", mad,mad,ch);
00053          //          obtitle.Form("MADC %d Channel %d electron", mad, ch);
00054          //          fMadc_electron[mad][ch] = MakeTH1('I', obname.Data(), obtitle.Data(), 0x3FF, 0, 0x3FF);
00055 
00056       }
00057 
00058       obname.Form("MADC%d_overview", mad);
00059       pMadc_overview[mad] = GetPicture(obname.Data());
00060       if (pMadc_overview[mad] == 0) {
00061          obtitle.Form("Channels overview MADC %2d", mad);
00062          pMadc_overview[mad] = new TGo4Picture(obname.Data(), obtitle.Data());
00063          pMadc_overview[mad]->SetDivision(8, 4);
00064          for(int y=0; y<8;++y)
00065             for(int x=0;x<4;++x) {
00066                pMadc_overview[mad]->Pic(y, x)->AddObject(fMadc_adc[mad][4*y+x]);
00067                pMadc_overview[mad]->Pic(y, x)->SetFillAtt(5, 3001);
00068                pMadc_overview[mad]->Pic(y, x)->SetRangeX(0, 200);
00069                pMadc_overview[mad]->Pic(y, x)->SetLogScale(1, 1);
00070             }
00071          AddPicture(pMadc_overview[mad]);
00072       }
00073 
00074       //       obname.Form("MADC%d_pion", mad);
00075       //       if (GetPicture(obname.Data()) == 0) {
00076       //          obtitle.Form("Pion channels overview MADC %2d", mad);
00077       //          TGo4Picture* pic = new TGo4Picture(obname.Data(), obtitle.Data());
00078       //          pic->SetDivision(8, 4);
00079       //          for(int y=0; y<8;++y)
00080       //             for(int x=0;x<4;++x) {
00081       //                pic->Pic(y, x)->AddObject(fMadc_pion[mad][4*y+x]);
00082       //                pic->Pic(y, x)->SetFillAtt(5, 3001);
00083       //                pic->Pic(y, x)->SetRangeX(0, 200);
00084       //                pic->Pic(y, x)->SetLogScale(1, 1);
00085       //             }
00086       //          AddPicture(pic);
00087       //       }
00088       //
00089       //
00090       //       obname.Form("MADC%d_electron", mad);
00091       //       if (GetPicture(obname.Data()) == 0) {
00092       //          obtitle.Form("Electron channels overview MADC %2d", mad);
00093       //          TGo4Picture* pic = new TGo4Picture(obname.Data(), obtitle.Data());
00094       //          pic->SetDivision(8, 4);
00095       //          for(int y=0; y<8;++y)
00096       //             for(int x=0;x<4;++x) {
00097       //                pic->Pic(y, x)->AddObject(fMadc_electron[mad][4*y+x]);
00098       //                pic->Pic(y, x)->SetFillAtt(5, 3001);
00099       //                pic->Pic(y, x)->SetRangeX(0, 200);
00100       //                pic->Pic(y, x)->SetLogScale(1, 1);
00101       //             }
00102       //          AddPicture(pic);
00103       //       }
00104 
00105    }
00106 
00107 
00108    for (int n=0;n<MAX_1290;n++)
00109       fTDC[n].MakeHistos(this, Form("TDC%d",n));
00110 
00111    // 1182 histo
00112 
00113    for (int n=0;n<MAX_1182;n++)
00114       for (int ch=0;ch<NUM_1182_CH;ch++) {
00115          obname.Form("118%d/118%d_channel%d", n+2, n+2, ch);
00116          obtitle.Form("118%d Channel %d", n+2, ch);
00117          f1182h[n][ch] = MakeTH1('I', obname.Data(), obtitle.Data(), 1024, 0, 4096.);
00118       }
00119 
00120    for (int n=0;n<MAX_1182;n++)
00121       for (int ch=0;ch<NUM_1182_CH;ch++) {
00122          obname.Form("118%d/118%d_pulser_channel%d", n+2, n+2, ch);
00123          obtitle.Form("118%d Channel %d under pulser", n+2, ch);
00124          f1182Pulserh[n][ch] = MakeTH1('I', obname.Data(), obtitle.Data(), 1024, 0, 4096.);
00125       }
00126 
00127    /* some histograms for v550 data:*/
00128 
00129    obname.Form("v550/Fifo0_Trace");
00130    obtitle.Form("v550 Fifo 0 channels trace");
00131    f550h_fifo0_trace= MakeTH1('I', obname.Data(), obtitle.Data(),
00132          NUM_V550_CH, 0, NUM_V550_CH,"channel number");
00133    obname.Form("v550/Fifo0_Acc");
00134    obtitle.Form("v550 Fifo 0 channels accumulated");
00135    f550h_fifo0_acc=MakeTH1('I', obname.Data(), obtitle.Data(),
00136          NUM_V550_CH, 0, NUM_V550_CH,"channel number");
00137    obname.Form("v550/Fifo1_Trace");
00138    obtitle.Form("v550 Fifo 1 channels trace");
00139    f550h_fifo1_trace=MakeTH1('I', obname.Data(), obtitle.Data(),
00140          NUM_V550_CH, 0, NUM_V550_CH,"channel number");
00141    obname.Form("v550/Fifo1_Acc");
00142    obtitle.Form("v550 Fifo 1 channels accumulated");
00143 
00144    f550h_fifo1_acc= MakeTH1('I', obname.Data(), obtitle.Data(),
00145          NUM_V550_CH, 0, NUM_V550_CH,"channel number");
00146 
00147    TGo4Log::Info("TMbsCrateProc - Histograms created");
00148 
00149 
00150 //   obname.Form("MADC/SiStrip/Pedestal_Trace");
00151 //   obtitle.Form("SiStrip pedestals trace");
00152 //   fMadc_si_ped_trace= MakeTH1('I', obname.Data(), obtitle.Data(),
00153 //                  N_SISTRIP_CHA, 0, N_SISTRIP_CHA,"channel number");
00154 //   obname.Form("MADC/SiStrip/Pedestal_Acc");
00155 //   obtitle.Form("SiStrip pedestals accumulated");
00156 //   fMadc_si_ped_acc=MakeTH1('I', obname.Data(), obtitle.Data(),
00157 //                  N_SISTRIP_CHA, 0, N_SISTRIP_CHA,"channel number");
00158 
00159 
00160 }
00161 
00162 TMbsCrateProc::~TMbsCrateProc()
00163 {
00164 }
00165 
00166 void TMbsCrateProc::InitEvent(TGo4EventElement* outevnt)
00167 {
00168    TCBMBeamtimeEvent* btevent=dynamic_cast<TCBMBeamtimeEvent*>(outevnt);
00169    if(btevent)
00170    {
00171       // since output event object is never discarded within processor lifetime,
00172       // we just search for subevent by name once to speed up processing
00173       if(fOutputEvent==0)
00174          fOutputEvent=dynamic_cast<TMbsCrateEvent*>(btevent->GetSubEvent("MBSCRATE"));
00175    }
00176    else
00177    {
00178 
00179       fOutputEvent= dynamic_cast<TMbsCrateEvent*>(outevnt);
00180    }
00181    if(fOutputEvent==0) {
00182       GO4_STOP_ANALYSIS_MESSAGE("**** TMbsCrateProc: Fatal error: output event is not a TMbsCrateEvent!!! STOP GO4");
00183    }
00184 }
00185 
00186 void TMbsCrateProc::ProcessSubevent(TGo4MbsSubEvent* subevt)
00187 {
00188    if (fPar->fBucharestSetup) {
00189 
00190       if ((subevt->GetProcid() != 1) || (subevt->GetControl() != 9)) {
00191          cerr << " Unknown subevent tags for Bucharest setup proc = " << subevt->GetProcid() << " control = " << subevt->GetControl() << endl;
00192          return;
00193       }
00194 
00195       if (GetTriggerNumber()>13) {
00196          cout << "Skip trigger number " << GetTriggerNumber() << endl;
00197          return;
00198       }
00199 
00200       //cout << "Process Bucharest event " << GetEventNumber() << "   trigger = " << GetTriggerNumber() << endl;
00201 
00202       const int numtags = 3;
00203       Int_t tags[numtags] = { 0x4d544443,   // TDC
00204             0x4d414443,   // MADC
00205             0x56353630 }; // V560
00206 
00207       Int_t minlen[numtags] = { 1, 1, 16 };
00208 
00209       Int_t *pdata = subevt->GetDataField();
00210       UInt_t lwords = subevt->GetIntLen();
00211 
00212       Int_t ntdc(0), nmadc(0), indx(0);
00213 
00214       while (indx < (Int_t) lwords) {
00215          int tagid = -1;
00216          for (int n=0;n<numtags;n++)
00217             if (tags[n] == pdata[indx]) tagid = n;
00218          if (tagid<0) {
00219             cerr << "Cannot recognize tag id in MBS event from Bucharest" << endl;
00220             return;
00221          }
00222 
00223          // skip tag identifier
00224          indx++;
00225 
00226          Int_t nextindx = indx + minlen[tagid];
00227 
00228          while (nextindx < (Int_t) lwords) {
00229             bool foundnexttag = false;
00230 
00231             for (int n=0;n<numtags;n++)
00232                if (tags[n] == pdata[nextindx]) foundnexttag = true;
00233 
00234             if (foundnexttag) break;
00235 
00236             nextindx++;
00237          }
00238 
00239 //         cout << "   Found TAG = " << tagid << " len = " << nextindx - indx << " x4 bytes" << endl;
00240 
00241          switch (tagid) {
00242             case 0: // TDC
00243                Process1290(-1, pdata + indx, nextindx - indx);
00244                ntdc++;
00245                break;
00246 
00247             case 1: // MADC
00248                ProcessMADC(nmadc++, pdata + indx, nextindx - indx);
00249                break;
00250 
00251             case 2: // V560
00252                Processv560(pdata + indx, nextindx - indx);
00253                break;
00254          }
00255 
00256          indx = nextindx;
00257       }
00258 
00259       fOutputEvent->SetValid(kTRUE);
00260 
00261       return;
00262    }
00263 
00264  // use this to test with old beamtime data JAM:
00265  //  if (subevt->GetProcid() == roc::proc_TRD_MADC) {
00266 
00267    if (subevt->GetProcid() == roc::proc_CERN_Oct11) {
00268       // MBS data from CERN beamtime, October 11
00269       // subevt->GetControl() == 9 from Master MBS, subevt->GetControl() == 10 from Slave MBS
00270 
00271       fOutputEvent->fIsPulser = (GetTriggerNumber()==8);
00272 
00273       Int_t *pdata = subevt->GetDataField();
00274       UInt_t lwords = subevt->GetIntLen();
00275 
00276       // printf("CERN subevent Control %u TOTALLEN = %u\n", subevt->GetControl(), lwords);
00277 
00278       const int numtags = 12;
00279 
00280       const uint64_t tags[numtags] = { 
00281             0x3238313131313832ULL,  // 1182
00282             0x3338313131313833ULL,  // 1183
00283             0x3164616d6d616431ULL,  // mad1
00284             0x3264616d6d616432ULL,  // mad2
00285             0x3364616d6d616433ULL,  // mad3
00286             0x3464616d6d616434ULL,  // mad4
00287             0x3163647474646331ULL,  // tdc1
00288             0x3263647474646332ULL,  // tdc2
00289             0x3363647474646333ULL,  // tdc3
00290             0x3463647474646334ULL,  // tdc4
00291             0x3036357676353630ULL,  // v560
00292             0x3035357676353530ULL,  // v550
00293       };
00294       Int_t* tagpos[numtags];
00295       Int_t taglen[numtags];
00296 
00297       for (int n=0;n<numtags;n++) {
00298          tagpos[n] = 0;
00299          taglen[n] = 0;
00300       }
00301 
00302       int lasttag = -1;
00303 
00304       //   Int_t* p0 = pdata;
00305 
00306       for (UInt_t cur = 0; cur < lwords; ++cur) {
00307          uint64_t* ptag = (uint64_t*) pdata;
00308 
00309          //      printf("DATA 0x%x\n", *pdata);
00310 
00311          bool find = false;
00312 
00313          for (int ntag=0;ntag<numtags;ntag++)
00314             if (*ptag == tags[ntag]) {
00315                if (lasttag>=0) {
00316                   taglen[lasttag] = pdata - tagpos[lasttag];
00317                   lasttag = -1;
00318                }
00319 
00320                if (tagpos[ntag]==0) {
00321                   //               printf("TAG %d  pos %u\n", ntag, pdata - p0);
00322 
00323                   pdata+=2; cur+=1;
00324                   tagpos[ntag] = pdata;
00325                   lasttag = ntag;
00326                } else {
00327                   printf("FORMAT ERROR found tag %lld twice\n", tags[ntag]);
00328                   pdata+=2; cur+=1;
00329                }
00330 
00331                find = true;
00332                break;
00333 
00334             }
00335 
00336          if (!find) pdata++;
00337       }
00338 
00339       if (lasttag>=0)
00340          taglen[lasttag] = pdata - tagpos[lasttag];
00341 
00342       // for (int ntag=0;ntag<numtags;ntag++)
00343       //    if (tagpos[ntag])
00344       //       printf("TAG %d POS %lu LEN %d\n", ntag, tagpos[ntag] - p0, taglen[ntag]);
00345 
00346 
00347 
00348       // SL: in previous version trigger 8 was excluded at all
00349       //if (GetTriggerNumber()==8) {
00350       //   fOutputEvent->SetValid(kFALSE);
00351       //   return;
00352       //}
00353 
00354       // actually process data in subroutines:
00355 
00356       if (tagpos[0] && taglen[0])   Process1182(0, tagpos[0], taglen[0]);  // 1182
00357       if (tagpos[1] && taglen[1])   Process1182(1, tagpos[1], taglen[1]);  // 1183
00358       if (tagpos[2] && taglen[2])   ProcessMADC(0, tagpos[2], taglen[2]);  // mad1
00359       if (tagpos[3] && taglen[3])   ProcessMADC(1, tagpos[3], taglen[3]);  // mad2
00360       if (tagpos[4] && taglen[4])   ProcessMADC(2, tagpos[4], taglen[4]);  // mad3
00361       if (tagpos[5] && taglen[5])   ProcessMADC(3, tagpos[5], taglen[5]);  // mad4
00362       if (tagpos[6] && taglen[6])   Process1290(0, tagpos[6], taglen[6]);  // tdc1
00363       if (tagpos[7] && taglen[7])   Process1290(1, tagpos[7], taglen[7]);  // tdc2
00364       if (tagpos[8] && taglen[8])   Process1290(2, tagpos[8], taglen[8]);  // tdc3
00365       if (tagpos[9] && taglen[9])   Process1290(3, tagpos[9], taglen[9]);  // tdc4
00366       if (tagpos[10] && taglen[10]) Processv560(  tagpos[10], taglen[10]); // v560
00367       if (tagpos[11] && taglen[11]) Processv550(  tagpos[11], taglen[11]); // v550
00368 
00369       fOutputEvent->SetValid(kTRUE);
00370 
00371       return;
00372    }
00373 }
00374 
00375 
00376 
00377 
00378 void TMbsCrateProc::Process1182(int num, int* pdata, int len)
00379 {
00380    if (len < 2 + NUM_1182_CH ) {
00381       cout <<"TBeamMonitorProc:Process1182 data length"<< len << "not sufficient!" <<endl;
00382       return;
00383    }
00384 
00385    if (num>=MAX_1182) {
00386       cout << "Wrong ID=" << num << " of 1182 module" << endl;
00387       return;
00388    }
00389 
00390    // skip event header:
00391    int header = *pdata++;
00392    if (L1182_DEBUG)  printf("1182 header:0x%0x \n",header);
00393 
00394    // skip status
00395    int status = *pdata++;
00396    if (L1182_DEBUG)
00397       printf("1182 status:0x%0x :  \n"
00398             "        conversion complete    (bit0): 0x%x\n"
00399             "        conversion in progress (bit1): 0x%x\n"
00400             "        front(1)/rear(0) panel (bit2): 0x%x\n"
00401             "        Event buffer not full  (bit3): 0x%x\n"
00402             "        Event count       (bits 4..7): 0x%x\n"
00403             "        should be 0            (bit8): 0x%x\n"
00404             "        should be 0            (bit9): 0x%x\n",
00405             status, status & 1, (status >> 1) & 1 , (status >> 2) & 1, (status >> 3) & 1,
00406             (status >> 4) & 0xf,  (status >> 8) & 1, (status >> 9) & 1);
00407 
00408    for(int ch=0;ch<NUM_1182_CH;++ch) {
00409       int val = *pdata++;
00410       fOutputEvent->fData1182[num][ch] = val; // output event to root tree
00411    }
00412 
00413    // do not fill histograms if pulser event - used for later analysis
00414    if (fOutputEvent->IsPulser())
00415       for(int ch=0;ch<NUM_1182_CH;++ch)
00416          f1182Pulserh[num][ch]->Fill(fOutputEvent->fData1182[num][ch]);
00417    else
00418       for(int ch=0;ch<NUM_1182_CH;++ch)
00419          f1182h[num][ch]->Fill(fOutputEvent->fData1182[num][ch]);
00420 }
00421 
00422 void TMbsCrateProc::ProcessMADC(int i, int* pdata, int len) 
00423 {
00424    //      printf("MADC %d len %d  \n", i, len);
00425 
00426    if(i>=MAX_MADC) return;
00427    if(len<2) return;
00428 
00429    fMadc_trace[i]->Reset("");
00430 
00431    int* cursor = pdata;
00432    // skip event header:
00433    unsigned int header=*cursor;
00434    if ((header & 0xff000000) >> 24 != 0x40)
00435       printf("madc(%d) header error: 0x%08x, expected 0x40nnnnnn \n",i,header);
00436    cursor++;
00437    for(int cnt=0;cnt<len-2;cnt++) {
00438       unsigned int val=*cursor++;
00439 
00440       if ((val & 0xff000000) != 0x4000000)
00441          printf("madc(%d) seq %u value 0x%08x error, expected 0x04nnnnnn\n", i,cnt, val);
00442 
00443       unsigned chanid= (val >> 16) & 0x1f;
00444       unsigned adc = val & 0x1fff;
00445 
00446       fMadc_trace[i]->Fill(chanid, adc);
00447       if (!fOutputEvent->IsPulser())
00448          fMadc_adc[i][chanid]->Fill(adc);
00449       fOutputEvent->fMadc[i][chanid] = adc; // to output event, store to root tree
00450    }
00451    int trailer=*cursor;
00452    if ((trailer & 0xf0000000) >> 24 != 0xc0)
00453       printf("madc(%d) trailer error: 0x%08x, expected 0xcnnnnnnnn \n",i, trailer);
00454 
00455    //   printf("madc(%d) trailer: 0x%08x \n",i,trailer);
00456    //#endif
00457 }
00458 
00459 void TMbsCrateProc::Process1290(int num, int* pdata, unsigned int len)
00460 {
00461    if (num >= MAX_1290) {
00462       cout << "Wrong 1290 index " << num << endl;
00463       return;
00464    }
00465 
00466    int expected_geo = 0;
00467 
00468    if (num<0) {
00469       // if number is not specified, try to identify it base on GEO value
00470       expected_geo = T1290Data::FindGeo(pdata, len);
00471       if (expected_geo<=0) {
00472          printf("ERROR: Did not found GEO in the data\n");
00473          return;
00474       }
00475 
00476       switch (expected_geo) {
00477          case 14: num = 0; break;
00478          case 15: num = 1; break;
00479          case 13: num = 2; break;
00480          case 12: num = 3; break;
00481          default:
00482             printf("Unsupported GEO code %d\n", expected_geo);
00483             return;
00484       }
00485 
00486       // printf("Found GEO code %d defined num %d\n", expected_geo, num);
00487 
00488    } else {
00489       switch(num) {
00490          case 0: expected_geo = 14; break;
00491          case 1: expected_geo = 15; break;
00492          case 2: expected_geo = 13; break;
00493          case 3: expected_geo = 12; break;
00494       }
00495    }
00496 
00497    // printf("Process1290 num = %d expected GEO = %d len = %u\n", num, expected_geo, len);
00498 
00499    fOutputEvent->fMtdc[num].Unpack(pdata, len, expected_geo);
00500 
00501    if (!fOutputEvent->IsPulser())
00502       fTDC[num].FillHistos(fOutputEvent->fMtdc[num]);
00503 }
00504 
00505 void TMbsCrateProc::Processv560(int* pdata, unsigned int len)
00506 {
00507    if (len!=NUM_V560_CH) {
00508       cout << "wrong number = " << len << " of longwords in v560 data" << endl;
00509       return;
00510    }
00511 
00512    for (int ch=0;ch<NUM_V560_CH;ch++)
00513       fOutputEvent->f560Scaler[ch] = pdata[ch];
00514 }
00515 
00516 void TMbsCrateProc::Processv550(int* pdata, unsigned int len)
00517 {
00518    UInt_t status=*pdata++;
00519    printf("v550 status:0x%0x \n",status);
00520    UInt_t bothcount = *pdata++;
00521    UShort_t wcount0 = (bothcount >> 16);
00522    UShort_t wcount1 = (bothcount & 0xFFFF);
00523    printf("v550 wc0:0x%0x wc1:0x%0x\n",wcount0, wcount1);
00524    if((int) len < wcount0 + wcount1 + 4)
00525    {
00526       cout << "TMbsCrateProc::Processv550 data length" << len
00527             << "not sufficient!" << endl;
00528       return;
00529    }
00530 
00531    f550h_fifo0_trace->Reset("");
00532    for (UShort_t i=0;i<wcount0;i++)
00533    {
00534       UInt_t val=*pdata++;
00535       UInt_t ch=(val >> 12) & 0x7FF;
00536       UInt_t adc=val & 0xFFF;
00537       printf("si strip FIFO 0: %d has id 0x%0x, adc:0x%0x \n",i,ch,adc);
00538       if(ch >= NUM_V550_CH)
00539       {
00540          cout << "TMbsCrateProc::Processv550 fifo0 - channel " << ch
00541                << "out of range - max=" <<NUM_V550_CH << endl;
00542          continue;
00543       }
00544       fOutputEvent->f550Data0[ch]=adc;
00545       f550h_fifo0_trace->AddBinContent(ch + 1, (Double_t) adc);
00546       f550h_fifo0_acc->AddBinContent(ch + 1, (Double_t) adc);
00547 
00548    }
00549    f550h_fifo1_trace->Reset("");
00550    for (UShort_t i=0;i<wcount1;i++)
00551    {
00552       UInt_t val=*pdata++;
00553       UInt_t ch=(val >> 12) & 0x7FF;
00554       UInt_t adc=val & 0xFFF;
00555       printf("si strip FIFO 1: %d has id 0x%0x, adc:0x%0x \n",i,ch,adc);
00556       if(ch >= NUM_V550_CH)
00557       {
00558          cout << "TMbsCrateProc::Processv550 fifo1 - channel " << ch
00559                << "out of range - max=" <<NUM_V550_CH << endl;
00560          continue;
00561       }
00562       fOutputEvent->f550Data1[ch]=adc;
00563 
00564       if (!fOutputEvent->IsPulser()) {
00565          f550h_fifo1_trace->AddBinContent(ch + 1, (Double_t) adc);
00566          f550h_fifo1_acc->AddBinContent(ch + 1, (Double_t) adc);
00567       }
00568    }
00569 }
00570 

Generated on Tue Dec 10 2013 04:52:18 for ROCsoft by  doxygen 1.7.1