]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - ZDC/ZDCEMDda.cxx
Updated SOD to take into account ADD detector + error removal in DA
[u/mrichter/AliRoot.git] / ZDC / ZDCEMDda.cxx
index a3addb4d07c63425da5387ea6e9cc5fd52b922ce..fc93ae10d87bdfb6a4033c251656a40e4bd67fd5 100644 (file)
@@ -23,7 +23,7 @@ Trigger Types Used: Standalone Trigger
 */
 #define PEDDATA_FILE  "ZDCPedestal.dat"
 #define MAPDATA_FILE  "ZDCChMapping.dat"
-#define ENCALIBDATA_FILE   "ZDCEnCalib.dat"
+#define ENCALIBDATA_FILE   "ZDCEnergyCalib.dat"
 #define TOWCALIBDATA_FILE  "ZDCTowerCalib.dat"
 
 #include <stdio.h>
@@ -73,6 +73,15 @@ int main(int argc, char **argv) {
   int status = 0;
   // No. of ZDC cabled ch.
   int const kNChannels = 24;
+  int const kNScChannels = 32;
+  Int_t kFirstADCGeo=0, kLastADCGeo=3;
+      
+  Int_t ich=0;
+  Int_t adcMod[2*kNChannels], adcCh[2*kNChannels], sigCode[2*kNChannels];
+  Int_t det[2*kNChannels], sec[2*kNChannels];
+  for(Int_t y=0; y<2*kNChannels; y++){
+    adcMod[y]=adcCh[y]=sigCode[y]=det[y]=sec[y]=0;
+  }
 
   /* log start of process */
   printf("ZDC EMD program started\n");  
@@ -83,12 +92,11 @@ int main(int argc, char **argv) {
     return -1;
   }
   
-  //
   // --- Preparing histos for EM dissociation spectra
   //
   TH1F* histoEMDRaw[4];
   TH1F* histoEMDCorr[4];
-
+  //
   char namhistr[50], namhistc[50];
   for(Int_t i=0; i<4; i++) {
      if(i==0){
@@ -111,6 +119,25 @@ int main(int argc, char **argv) {
      histoEMDCorr[i] = new TH1F(namhistc,namhistc,100,0.,4000.);
   }
 
+   // --- Preparing histos for tower inter-calibration
+  //
+  TH1F* histZNCtow[4]; TH1F* histZPCtow[4];
+  TH1F* histZNAtow[4]; TH1F* histZPAtow[4];
+  //
+  char namhistznc[50], namhistzpc[50];
+  char namhistzna[50], namhistzpa[50];
+  for(Int_t i=0; i<4; i++) {
+     sprintf(namhistznc,"ZNC-tow%d",i+1);
+     sprintf(namhistzpc,"ZPC-tow%d",i+1);
+     sprintf(namhistzna,"ZNA-tow%d",i+1);
+     sprintf(namhistzpa,"ZPA-tow%d",i+1);
+     //
+     histZNCtow[i] = new TH1F(namhistznc,namhistznc,100,0.,4000.);
+     histZPCtow[i] = new TH1F(namhistzpc,namhistzpc,100,0.,4000.);
+     histZNAtow[i] = new TH1F(namhistzna,namhistzna,100,0.,4000.);
+     histZPAtow[i] = new TH1F(namhistzpa,namhistzpa,100,0.,4000.);
+  }
+  
   /* open result file */
   FILE *fp=NULL;
   fp=fopen("./result.txt","a");
@@ -151,11 +178,11 @@ int main(int argc, char **argv) {
     for(int ii=0; ii<2; ii++){
        fscanf(filePed,"%f",&readValues[ii][jj]);
     }
-    if(jj<kNChannels && jj<2*kNChannels){
+    if(jj<2*kNChannels){
       MeanPed[jj] = readValues[0][jj];
-      //printf("\t MeanPedhg[%d] = %1.1f\n",jj, MeanPedhg[jj]);
+      printf("\t MeanPed[%d] = %1.1f\n",jj, MeanPed[jj]);
     }
-    else if(jj>2*kNChannels && jj>4*kNChannels){
+    else if(jj>2*kNChannels){
       CorrCoeff0[jj-4*kNChannels] = readValues[0][jj]; 
       CorrCoeff1[jj-4*kNChannels] = readValues[1][jj];;
     }
@@ -169,6 +196,9 @@ int main(int argc, char **argv) {
   int nevents_physics=0;
   int nevents_total=0;
 
+  struct eventHeaderStruct *event;
+  eventTypeType eventT;
+
   /* read the data files */
   int n;
   for(n=1;n<argc;n++){
@@ -185,8 +215,6 @@ int main(int argc, char **argv) {
 
     /* read the file */
     for(;;){
-      struct eventHeaderStruct *event;
-      eventTypeType eventT;
 
       /* get next event */
       status=monitorGetEventDynamic((void **)&event);
@@ -215,39 +243,65 @@ int main(int argc, char **argv) {
       /* use event - here, just write event id to result file */
       eventT=event->eventType;
       
-      Int_t ich=0;
-      Int_t adcMod[2*kNChannels], adcCh[2*kNChannels];
-      Int_t sigCode[2*kNChannels], det[2*kNChannels], sec[2*kNChannels];
+      Int_t iScCh=0;
+      Int_t scMod[kNScChannels], scCh[kNScChannels], scSigCode[kNScChannels];
+      Int_t scDet[kNScChannels], scSec[kNScChannels];
+      for(Int_t y=0; y<kNScChannels; y++){
+        scMod[y]=scCh[y]=scSigCode[y]=scDet[y]=scSec[y]=0;
+      }
+      //
+      Int_t modNum=-1, modType=-1;
       
       if(eventT==START_OF_DATA){
-
-       rawStreamZDC->SetSODReading(kTRUE);
                
+       rawStreamZDC->SetSODReading(kTRUE);
+       
        // --------------------------------------------------------
        // --- Writing ascii data file for the Shuttle preprocessor
         mapFile4Shuttle = fopen(MAPDATA_FILE,"w");
        if(!rawStreamZDC->Next()) printf(" \t No raw data found!! \n");
         else{
-         while(rawStreamZDC->Next()){
-            if(rawStreamZDC->IsChMapping()){
-             adcMod[ich] = rawStreamZDC->GetADCModFromMap(ich);
-             adcCh[ich] = rawStreamZDC->GetADCChFromMap(ich);
-             sigCode[ich] = rawStreamZDC->GetADCSignFromMap(ich);
-             det[ich] = rawStreamZDC->GetDetectorFromMap(ich);
-             sec[ich] = rawStreamZDC->GetTowerFromMap(ich);
-             //
-             fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\t%d\t%d\t%d\n",
-               ich,adcMod[ich],adcCh[ich],sigCode[ich],det[ich],sec[ich]);
-             //
-             //printf("ZDCEMDDA.cxx ->  ch.%d mod %d, ch %d, code %d det %d, sec %d\n",
-             //           ich,adcMod[ich],adcCh[ich],sigCode[ich],det[ich],sec[ich]);
-             ich++;
+         while((rawStreamZDC->Next())){
+            if(rawStreamZDC->IsHeaderMapping()){ // mapping header
+              modNum = rawStreamZDC->GetADCModule();
+              modType = rawStreamZDC->GetModType();
+           }
+            if(rawStreamZDC->IsChMapping()){ 
+             if(modType==1){ // ADC mapping ----------------------
+               adcMod[ich]  = rawStreamZDC->GetADCModFromMap(ich);
+               adcCh[ich]   = rawStreamZDC->GetADCChFromMap(ich);
+               sigCode[ich] = rawStreamZDC->GetADCSignFromMap(ich);
+               det[ich]     = rawStreamZDC->GetDetectorFromMap(ich);
+               sec[ich]     = rawStreamZDC->GetTowerFromMap(ich);
+               //
+               fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\t%d\t%d\t%d\n",
+                 ich,adcMod[ich],adcCh[ich],sigCode[ich],det[ich],sec[ich]);
+               //
+               //printf("  Mapping in DA -> %d ADC: mod %d ch %d, code %d det %d, sec %d\n",
+               //  ich,adcMod[ich],adcCh[ich],sigCode[ich],det[ich],sec[ich]);
+               //
+               ich++;
+             }
+             else if(modType==2){ //VME scaler mapping --------------------
+               scMod[iScCh]     = rawStreamZDC->GetScalerModFromMap(iScCh);
+               scCh[iScCh]      = rawStreamZDC->GetScalerChFromMap(iScCh);
+               scSigCode[iScCh] = rawStreamZDC->GetScalerSignFromMap(iScCh);
+               scDet[iScCh]     = rawStreamZDC->GetScDetectorFromMap(iScCh);
+               scSec[iScCh]    = rawStreamZDC->GetScTowerFromMap(iScCh);
+               //
+               fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\t%d\t%d\t%d\n",
+                 iScCh,scMod[iScCh],scCh[iScCh],scSigCode[iScCh],scDet[iScCh],scSec[iScCh]);
+               //
+               //printf("  Mapping in DA -> %d Scaler: mod %d ch %d, code %d det %d, sec %d\n",
+               //  iScCh,scMod[iScCh],scCh[iScCh],scSigCode[iScCh],scDet[iScCh],scSec[iScCh]);
+               //
+               iScCh++;
+             }
            }
          }
        }
         fclose(mapFile4Shuttle);
-      }//SOD event
-      
+      }// SOD event
     
     if(eventT==PHYSICS_EVENT){
       // --- Reading data header
@@ -255,7 +309,7 @@ int main(int argc, char **argv) {
       const AliRawDataHeader* header = reader->GetDataHeader();
       if(header){
          UChar_t message = header->GetAttributes();
-        if(message & 0x70){ // DEDICATED EMD RUN
+        if((message & 0x70) == 0x70){ // DEDICATED EMD RUN
            //printf("\t STANDALONE_EMD_RUN raw data found\n");
            continue;
         }
@@ -293,18 +347,17 @@ int main(int argc, char **argv) {
        Int_t quad = rawStreamZDC->GetSector(1);
         
        if(rawStreamZDC->IsADCDataWord() && !(rawStreamZDC->IsUnderflow())
-            && !(rawStreamZDC->IsOverflow()) && det!=-1
-            && (rawStreamZDC->GetADCGain() == 1)){ // Selecting LOW RES ch.s
-
-         //printf("  IsADCWord %d, IsUnderflow %d, IsOverflow %d\n",
-         //  rawStreamZDC->IsADCDataWord(),rawStreamZDC->IsUnderflow(),rawStreamZDC->IsOverflow());
+            && !(rawStreamZDC->IsOverflow()) && det!=-1 && det!=3 
+            && (rawStreamZDC->GetADCGain() == 1 && // Selecting LOW RES ch.s
+             rawStreamZDC->GetADCModule()>=kFirstADCGeo && rawStreamZDC->GetADCModule()<=kLastADCGeo)){
          
-         // Taking LOW RES channels -> channel+kNChannels !!!!
+         // Taking LOW RES channels -> ch.+kNChannels !!!!
          Int_t DetIndex=999, PedIndex=999;
-         if(det != 3 && quad != 5){ // Not ZEM nor PMRef
+         // Not PMRef
+         if(quad!=5){ 
            if(det == 1){
              DetIndex = det-1;
-             PedIndex = quad+kNChannels;
+             PedIndex = quad+kNChannels; 
            }
            else if(det==2){
              DetIndex = det-1;
@@ -318,73 +371,94 @@ int main(int argc, char **argv) {
              DetIndex = det-2;
              PedIndex = quad+17+kNChannels;
            }
-            //EMD -> LR ADCs
-           if(rawStreamZDC->GetADCGain() == 1 && (DetIndex!=999 || PedIndex!=999)){ 
+           // Mean pedestal subtraction 
+           Float_t Pedestal = MeanPed[PedIndex];
+           // Pedestal subtraction from correlation with out-of-time signals
+           //Float_t Pedestal = CorrCoeff0[PedIndex]+CorrCoeff1[PedIndex]*MeanPedOOT[PedIndex];
+            //
+           if(DetIndex!=999 || PedIndex!=999){ 
              //
              ZDCRawADC[DetIndex] += (Float_t) rawStreamZDC->GetADCValue();
              //
-             // Mean pedestal subtraction 
-             Float_t Pedestal = MeanPed[PedIndex];
-             // Pedestal subtraction from correlation with out-of-time signals
-             //Float_t Pedestal = CorrCoeff0[PedIndex]+CorrCoeff1[PedIndex]*MeanPedOOT[PedIndex];
              //
              ZDCCorrADC[DetIndex] = (rawStreamZDC->GetADCValue()) - Pedestal;
              ZDCCorrADCSum[DetIndex] += ZDCCorrADC[DetIndex];
              //
-             /*printf("\t det %d quad %d res %d pedInd %d detInd %d:"
-                "ADCCorr = %d, ZDCCorrADCSum = %d\n", 
-                det,quad,rawStreamZDC->GetADCGain(),PedIndex,DetIndex, 
-                (Int_t) ZDCCorrADC[DetIndex],(Int_t) ZDCCorrADCSum[DetIndex]);
-             */
+             /*printf("\t det %d quad %d res %d pedInd %d "
+                "Pedestal %1.0f -> ADCCorr = %d ZDCCorrADCSum = %d\n", 
+                det,quad,rawStreamZDC->GetADCGain(),PedIndex,Pedestal, 
+                (Int_t) ZDCCorrADC[DetIndex],(Int_t) ZDCCorrADCSum[DetIndex]);*/
+                     
+           }
+           // Not common PM 
+           if(quad!=0){
+             Float_t corrADCval = (rawStreamZDC->GetADCValue()) - Pedestal;
+             if(det==1)      histZNCtow[quad-1]->Fill(corrADCval);
+             else if(det==2) histZPCtow[quad-1]->Fill(corrADCval);
+             else if(det==4) histZNAtow[quad-1]->Fill(corrADCval);
+             else if(det==5) histZPAtow[quad-1]->Fill(corrADCval);
+             //
+             //printf("\t det %d tow %d fill histo w. value %1.0f\n", 
+             //  det,quad,corrADCval);
            }
+
            if(DetIndex==999 || PedIndex==999) 
-               printf(" WARNING! Detector a/o pedestal index are WRONG!!!\n");
+              printf(" WARNING! Detector a/o pedestal index are WRONG!!!\n");
  
-         }  
+         }//quad!=5
        }//IsADCDataWord()
-        //
+       
        }
        //
        nevents_physics++;
        //
+       delete reader;
+       delete rawStreamZDC;
+       //
        for(Int_t j=0; j<4; j++){
           histoEMDRaw[j]->Fill(ZDCRawADC[j]);
           histoEMDCorr[j]->Fill(ZDCCorrADCSum[j]);
        }
-    }
-    
-    nevents_total++;
-
-    /* free resources */
-    free(event);
-    
+    }//(if PHYSICS_EVENT)
+      
     /* exit when last event received, no need to wait for TERM signal */
-    if (eventT==END_OF_RUN) {
-      printf("EOR event detected\n");
+    else if(eventT==END_OF_RUN) {
+      printf(" -> EOR event detected\n");
       break;
     }
+    
+    nevents_total++;
 
    }
+
+   /* free resources */
+   free(event);
   }
     
   /* Analysis of the histograms */
   //
   FILE *fileShuttle1 = fopen(ENCALIBDATA_FILE,"w");
   //
-  Int_t BinMax[4];
-  Float_t YMax[4];
-  Int_t NBinsx[4];
-  Float_t MeanFitVal[4];
+  Int_t BinMax[4]={0,0,0,0};
+  Float_t YMax[4]={0.,0.,0.,0.};
+  Int_t NBinsx[4]={0,0,0,0};
+  Float_t MeanFitVal[4]={0.,0.,0.,0.};
   TF1 *fitfun[4];
   for(Int_t k=0; k<4; k++){
      if(histoEMDCorr[k]->GetEntries() == 0){
        printf("\n WARNING! Empty histos -> ending DA WITHOUT writing output\n\n");
        return -1;
      } 
+     //
      BinMax[k] = histoEMDCorr[k]->GetMaximumBin();
+     if(BinMax[k]<=6){
+       printf("\n WARNING! Something wrong with det %d histo -> ending DA WITHOUT writing output\n\n", k);
+       return -1;
+     }
+     // 
      YMax[k] = (histoEMDCorr[k]->GetXaxis())->GetXmax();
      NBinsx[k] = (histoEMDCorr[k]->GetXaxis())->GetNbins();
-//     printf("\n\t Det%d -> BinMax = %d, ChXMax = %f\n", k+1, BinMax[k], BinMax[k]*YMax[k]/NBinsx[k]);
+     //printf("\n\t Det%d -> BinMax = %d, ChXMax = %f\n", k+1, BinMax[k], BinMax[k]*YMax[k]/NBinsx[k]);
      histoEMDCorr[k]->Fit("gaus","Q","",BinMax[k]*YMax[k]/NBinsx[k]*0.7,BinMax[k]*YMax[k]/NBinsx[k]*1.25);
      fitfun[k] = histoEMDCorr[k]->GetFunction("gaus");
      MeanFitVal[k] = (Float_t) (fitfun[k]->GetParameter(1));
@@ -392,7 +466,6 @@ int main(int argc, char **argv) {
   }
   //
   Float_t CalibCoeff[6];     
-  Float_t icoeff[5];
   //
   for(Int_t j=0; j<6; j++){
      if(j<4){
@@ -406,16 +479,44 @@ int main(int argc, char **argv) {
      }
   }
   fclose(fileShuttle1);
-  //
   FILE *fileShuttle2 = fopen(TOWCALIBDATA_FILE,"w");
+  //Float_t meanvalznc[4], meanvalzpc[4], meanvalzna[4], meanvalzpa[4];
   for(Int_t j=0; j<4; j++){
+     /*if(histZNCtow[j]->GetEntries() == 0){
+       printf("\n WARNING! Empty histos -> ending DA WITHOUT writing output\n\n");
+       return -1;
+     } 
+     meanvalznc[j] = histZNCtow[j]->GetMean();
+     meanvalzpc[j] = histZPCtow[j]->GetMean();
+     meanvalzna[j] = histZNAtow[j]->GetMean();
+     meanvalzpa[j] = histZPAtow[j]->GetMean();*/
+     
      // Note -> For the moment the inter-calibration coeff. are set to 1 
-     for(Int_t k=0; k<5; k++){  
-       icoeff[k] = 1.;
-       fprintf(fileShuttle2,"\t%f",icoeff[k]);
+     for(Int_t k=0; k<4; k++){  
+       Float_t icoeff = 1.;
+       fprintf(fileShuttle2,"\t%f",icoeff);
        if(k==4) fprintf(fileShuttle2,"\n");
      }
   }
+  //
+  /*if(meanvalznc[1]!=0 && meanvalznc[2]!=0 && meanvalznc[3]!=0 && 
+     meanvalzpc[1]!=0 && meanvalzpc[2]!=0 && meanvalzpc[3]!=0 &&
+     meanvalzna[1]!=0 && meanvalzna[2]!=0 && meanvalzna[3]!=0 &&
+     meanvalzpa[1]!=0 && meanvalzpa[2]!=0 && meanvalzpa[3]!=0){
+    fprintf(fileShuttle2,"\t%f\t%f\t%f\t%f\n",
+       1.0,meanvalznc[0]/meanvalznc[1],meanvalznc[0]/meanvalznc[2],meanvalznc[0]/meanvalznc[3]);
+    fprintf(fileShuttle2,"\t%f\t%f\t%f\t%f\n",
+       1.0,meanvalzpc[0]/meanvalzpc[1],meanvalzpc[0]/meanvalzpc[2],meanvalzpc[0]/meanvalzpc[3]);
+    fprintf(fileShuttle2,"\t%f\t%f\t%f\t%f\n",
+       1.0,meanvalzna[0]/meanvalzna[1],meanvalzpc[0]/meanvalzna[2],meanvalzpc[0]/meanvalzna[3]);
+    fprintf(fileShuttle2,"\t%f\t%f\t%f\t%f\n",
+       1.0,meanvalzpa[0]/meanvalzpa[1],meanvalzpc[0]/meanvalzpa[2],meanvalzpc[0]/meanvalzpa[3]);
+  }
+  else{
+    printf("\n Tower intercalib. coeff. CAN'T be calculated (some mean values are ZERO)!!!\n\n");
+    return -1;
+  }*/
   fclose(fileShuttle2);
   
   for(Int_t ij=0; ij<4; ij++){
@@ -436,19 +537,19 @@ int main(int argc, char **argv) {
   daqDA_progressReport(90);
 
   /* store the result file on FES */
-  status = daqDA_FES_storeFile(MAPDATA_FILE, MAPDATA_FILE);
+  status = daqDA_FES_storeFile(MAPDATA_FILE, "MAPPING");
   if(status){
     printf("Failed to export file : %d\n",status);
     return -1;
   }
   //
-  status = daqDA_FES_storeFile(ENCALIBDATA_FILE, ENCALIBDATA_FILE);
+  status = daqDA_FES_storeFile(ENCALIBDATA_FILE, "EMDENERGYCALIB");
   if(status){
     printf("Failed to export file : %d\n",status);
     return -1;
   }
   //
-  status = daqDA_FES_storeFile(TOWCALIBDATA_FILE, TOWCALIBDATA_FILE);
+  status = daqDA_FES_storeFile(TOWCALIBDATA_FILE, "EMDTOWERCALIB");
   if(status){
     printf("Failed to export file : %d\n",status);
     return -1;