]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - VZERO/VZEROda.cxx
Smalle changes required by HLT (Jochen)
[u/mrichter/AliRoot.git] / VZERO / VZEROda.cxx
index 2cc1255aabd9d6fea8eefdb47432ca11874f99c6..cda344a88a83a74c4dbc316ef01dff9aa31f7c5e 100755 (executable)
@@ -7,8 +7,8 @@
 - DA Type:    MON
 - Number of events needed: >=500
 - Input Files:  argument list
-- Output Files: local files  VZERO_Histos.root, V0_Ped_Width_Gain.dat
-                FXS file     V0_Ped_Width_Gain.dat
+- Output Files: local files  VZERO_Histos.root, V0_Pedestals.dat (Online mapping)
+                FXS file     V0_Ped_Width_Gain.dat (Offline mapping)
 - Trigger types used: PHYSICS_EVENT
 **********************************************************************************/
 
@@ -18,7 +18,7 @@
 *                                                                                 *
 * This program connects to the DAQ data source passed as argument.                *
 * It computes calibration parameters, populates local "./V0_Ped_Width_Gain.dat"   *            
-* file and exports it to the FES.                                                 *
+* file, exports it to the FES, and stores it into DAQ DB                          *
 * The program exits when being asked to shut down (daqDA_checkshutdown)           *
 * or on End of Run event.                                                         *
 * We have 128 channels instead of 64 as expected for V0 due to the two sets of    *
@@ -48,6 +48,8 @@
 #include <TH1F.h>
 #include <TMath.h>
 
+Int_t GetOfflineChannel(Int_t channel);
+
 /* Main routine --- Arguments: monitoring data source */
       
 int main(int argc, char **argv) {
@@ -64,28 +66,51 @@ int main(int argc, char **argv) {
      return -1;
   }
   
+  // Online values (using FEE channel numbering), 
+  // stored into local V0_Pedestals.dat:
   Double_t ADCmean[128];
   Double_t ADCsigma[128];
   Double_t PEDmean[128];
   Double_t PEDsigma[128];
-     
+  // Offline values(same but ordered as in aliroot for offliners)
+  // stored into V0_Ped_Width_Gain.dat: 
+  Double_t ADCmean_Off[128];
+  Double_t ADCsigma_Off[128];
+  Double_t PEDmean_Off[128];
+  Double_t PEDsigma_Off[128];
+       
 //___________________________________________________
 // Get cuts from V00DA.config file
 
-  Int_t    kHighCut;    // = 50; high cut on pedestal distribution - to be tuned
-  Int_t    kLowCut;     // = 30; low cut on signal distribution - to be tuned
-
+  Int_t    kClockMin;   // = 16;   LHC Clock Min for pedestal calculation
+  Int_t    kClockMax;   // = 19;   LHC Clock Max for pedestal calculation
+  Int_t    kLowCut;     // = 60;   low cut on signal distribution - to be tuned
+  Int_t    kHighCut;    // = 50;   high cut on pedestal distribution - to be tuned
+  
   status = daqDA_DB_getFile("V00DA.config","./V00DA.config");
   if (status) {
-      printf("Failed to get config file (V00DA.config) from DAQdetDB, status=%d\n", status);
-      return -1;   
+      printf("Failed to get Config file (V00DA.config) from DAQ DB, status=%d\n", status);
+      printf("Take default values of parameters for pedestal calculation \n");
+      kClockMin  =  16; 
+      kClockMax  =  19; 
+      kLowCut    =  60;   
+      kHighCut   =  50;  
+  } else {
+      /* open the config file and retrieve cuts */
+      FILE *fpConfig = fopen("V00DA.config","r");
+      int res = fscanf(fpConfig,"%d %d %d %d ",&kClockMin,&kClockMax,&kLowCut,&kHighCut);
+      if(res!=4) {
+           printf("Failed to get values from Config file (V00DA.config): wrong file format - 4 integers are expected - \n");
+           kClockMin  =  16; 
+            kClockMax  =  19; 
+            kLowCut    =  60;   
+            kHighCut   =  50; 
+      }
+      fclose(fpConfig);
   }
-  /* open the config file and retrieve cuts */
-  FILE *fpConfig = fopen("V00DA.config","r");
-  fscanf(fpConfig,"%d %d",&kLowCut,&kHighCut);
-  fclose(fpConfig);
   
-  printf("LowCut on signal = %d ; HighCut on pedestal = %d\n",kLowCut,kHighCut);
+  printf("LHC Clock Min for pedestal calculation = %d; LHC Clock Max for pedestal calculation = %d; LowCut on signal = %d ; HighCut on pedestal = %d\n",
+          kClockMin, kClockMax, kLowCut, kHighCut);
 
 //___________________________________________________
 // Book HISTOGRAMS - dynamics of p-p collisions -
@@ -99,17 +124,23 @@ int main(int argc, char **argv) {
   for (Int_t i=0; i<128; i++) {
        sprintf(ADCname,"hADC%d",i);
        sprintf(texte,"ADC cell%d",i);
-       hADCname[i]  = new TH1F(ADCname,texte,1024,0,1023);
+       hADCname[i]  = new TH1F(ADCname,texte,1024,-0.5, 1023.5);
        sprintf(PEDname,"hPED%d",i);
        sprintf(texte,"PED cell%d",i);
-       hPEDname[i]  = new TH1F(PEDname,texte,1024,0,1023);
+       hPEDname[i]  = new TH1F(PEDname,texte,1024,-0.5, 1023.5);
   }
 //___________________________________________________ 
-  
+
+  /* open result file to be exported to FES */
+  FILE *fpLocal=NULL;
+  fpLocal=fopen("./V0_Pedestals.dat","w");
+  if (fpLocal==NULL) {
+      printf("Failed to open local result file\n");
+      return -1;}
+   
   /* open result file to be exported to FES */
   FILE *fp=NULL;
-  fp=fopen("./V0_Ped_Width_Gain.dat","a");
+  fp=fopen("./V0_Ped_Width_Gain.dat","w");
   if (fp==NULL) {
       printf("Failed to open result file\n");
       return -1;}
@@ -177,20 +208,26 @@ int main(int argc, char **argv) {
           AliRawReader *rawReader = new AliRawReaderDate((void*)event);
   
           AliVZERORawStream* rawStream  = new AliVZERORawStream(rawReader); 
-          rawStream->Next();   
+          if (rawStream->Next()) {     
            for(Int_t i=0; i<64; i++) {
-              if(!rawStream->GetIntegratorFlag(i,10))
-                   hADCname[i]->Fill(float(rawStream->GetADC(i)));       // even integrator - fills 0 to 63
-              else 
-                  hADCname[i+64]->Fill(float(rawStream->GetADC(i)));    // odd integrator  - fills 64 to 123
-              for(Int_t j=0; j<21; j++) {
-                  if(j==10) continue;
-                  if(!rawStream->GetIntegratorFlag(i,j))
-                       { hPEDname[i]->Fill(float(rawStream->GetPedestal(i,j))); }     // even integrator
-                  else 
-                      { hPEDname[i+64]->Fill(float(rawStream->GetPedestal(i,j))); }  // odd integrator 
-              }
-           }    
+               Int_t nFlag = 0;
+               for(Int_t j=kClockMin; j <= kClockMax; j++) {  // Check flags on clock range used for pedestal calculation
+                  if((rawStream->GetBBFlag(i,j)) || (rawStream->GetBGFlag(i,j))) nFlag++; 
+               }
+               if(nFlag == 0){       // Fill 64*2 pedestal histograms  - 2 integrators -
+                  for(Int_t j=kClockMin;j <= kClockMax;j++){
+                      Int_t Integrator = rawStream->GetIntegratorFlag(i,j);
+                      Float_t pedestal = (float)(rawStream->GetPedestal(i,j));
+                      hPEDname[i + 64 * Integrator]->Fill(pedestal);
+                  }    
+               } 
+               if((rawStream->GetBBFlag(i,10)) || (rawStream->GetBGFlag(i,10))){ // Charge
+                   Int_t Integrator = rawStream->GetIntegratorFlag(i,10);
+                   Float_t charge = (float)(rawStream->GetADC(i));   // Fill 64*2 ADCmax histograms 
+                   hADCname[i + 64 * Integrator]->Fill(charge);
+               }                          
+             } 
+          }   
            delete rawStream;
            rawStream = 0x0;      
            delete rawReader;
@@ -211,17 +248,40 @@ int main(int argc, char **argv) {
   
   printf("%d physics events processed\n",nevents_physics);
     
-//________________________________________________________________________
-//  Computes mean values, dumps them into the output text file
+//___________________________________________________________________________
+//  Computes mean values, converts FEE channels into Offline AliRoot channels
+//  and dumps the ordered values into the output text file for SHUTTLE
        
   for(Int_t i=0; i<128; i++) {
       hPEDname[i]->GetXaxis()->SetRange(0,kHighCut);
       PEDmean[i]  = hPEDname[i]->GetMean(); 
       PEDsigma[i] = hPEDname[i]->GetRMS(); 
-      hADCname[i]->GetXaxis()->SetRange(kLowCut,1023);
+      hADCname[i]->GetXaxis()->SetRange(kLowCut,1024);
       ADCmean[i]  = hADCname[i]->GetMean();
       ADCsigma[i] = hADCname[i]->GetRMS(); 
-      fprintf(fp," %.3f %.3f %.3f %.3f\n",PEDmean[i],PEDsigma[i],ADCmean[i],ADCsigma[i]);
+//      printf(" i = %d, %.3f %.3f %.3f %.3f\n",i,PEDmean[i],PEDsigma[i],ADCmean[i],ADCsigma[i]);
+      fprintf(fpLocal," %.3f %.3f %.3f %.3f\n",PEDmean[i],PEDsigma[i],
+                                               ADCmean[i],ADCsigma[i]);
+      if (i < 64) {
+          Int_t j = GetOfflineChannel(i);     
+          PEDmean_Off[j]  = PEDmean[i];
+          PEDsigma_Off[j] = PEDsigma[i];
+          ADCmean_Off[j]  = ADCmean[i];
+          ADCsigma_Off[j] = ADCsigma[i]; }
+      else{
+          Int_t j = GetOfflineChannel(i-64);     
+          PEDmean_Off[j+64]  = PEDmean[i];
+          PEDsigma_Off[j+64] = PEDsigma[i];
+          ADCmean_Off[j+64]  = ADCmean[i];
+          ADCsigma_Off[j+64] = ADCsigma[i];
+      }
+  }
+  
+  for(Int_t j=0; j<128; j++) {
+//      printf(" j = %d, %.3f %.3f %.3f %.3f\n",j,PEDmean_Off[j],PEDsigma_Off[j],
+//                                                ADCmean_Off[j],ADCsigma_Off[j]);
+      fprintf(fp," %.3f %.3f %.3f %.3f\n",PEDmean_Off[j],PEDsigma_Off[j],
+                                          ADCmean_Off[j],ADCsigma_Off[j]);                                    
   }
    
 //________________________________________________________________________
@@ -230,7 +290,7 @@ int main(int argc, char **argv) {
   TFile *histoFile = new TFile("VZERO_histos.root","RECREATE");
     
   for (Int_t i=0; i<128; i++) {
-       hADCname[i]->GetXaxis()->SetRange(0,1023);
+       hADCname[i]->GetXaxis()->SetRange(0,1024);
        hADCname[i]->Write(); 
        hPEDname[i]->Write(); }
 
@@ -239,8 +299,8 @@ int main(int argc, char **argv) {
   
 //________________________________________________________________________
    
-
-  /* close result file */
+  /* close local result file and FXS result file*/
+  fclose(fpLocal);
   fclose(fp);
  
   /* export result file to FES */
@@ -249,5 +309,26 @@ int main(int argc, char **argv) {
       printf("Failed to export file : %d\n",status);
       return -1; }
 
+  /* store result file into Online DB */
+  status=daqDA_DB_storeFile("./V0_Pedestals.dat","V00da_results");
+  if (status)    {
+      printf("Failed to store file into Online DB: %d\n",status);
+      return -1; }
+      
   return status;
 }
+
+ Int_t GetOfflineChannel(Int_t channel) {
+
+// Channel mapping Online - Offline:
+ Int_t fOfflineChannel[64] = {39, 38, 37, 36, 35, 34, 33, 32, 
+                              47, 46, 45, 44, 43, 42, 41, 40, 
+                             55, 54, 53, 52, 51, 50, 49, 48, 
+                             63, 62, 61, 60, 59, 58, 57, 56,
+                              7,  6,  5,  4,  3,  2,  1,  0, 
+                             15, 14, 13, 12, 11, 10,  9,  8,
+                             23, 22, 21, 20, 19, 18, 17, 16, 
+                             31, 30, 29, 28, 27, 26, 25, 24};
+ return        fOfflineChannel[channel];                     
+}