]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ZDC/ZDCCALIBMBda.cxx
Fixing bug #57328
[u/mrichter/AliRoot.git] / ZDC / ZDCCALIBMBda.cxx
1 /*
2
3 This program reads the DAQ data files passed as argument using the monitoring library.
4
5 It computes the average event size and populates local "./result.txt" file with the 
6 result.
7
8 The program reports about its processing progress.
9
10 Messages on stdout are exported to DAQ log system.
11
12 DA for ZDC standalone CALIBRATION_MB runs
13
14 Contact: Chiara.Oppedisano@to.infn.it
15 Link: 
16 Run Type: CALIBRATION_MB_RUN
17 DA Type: LDC
18 Number of events needed: at least 10^4
19 Input Files: ZDCPedestal.dat
20 Output Files: ZDCMBCalib.root, ZDCChMapping.dat
21 Trigger Types Used: Standalone Trigger
22
23 */
24 #define PEDDATA_FILE  "ZDCPedestal.dat"
25 #define MAPDATA_FILE  "ZDCChMapping.dat"
26 #define MBCALIBDATA_FILE   "ZDCMBCalib.root"
27
28 #include <stdio.h>
29 #include <Riostream.h>
30 #include <Riostream.h>
31
32 // DATE
33 #include <daqDA.h>
34 #include <event.h>
35 #include <monitor.h>
36
37 //ROOT
38 #include <TROOT.h>
39 #include <TPluginManager.h>
40 #include <TH1F.h>
41 #include <TH2F.h>
42 #include <TProfile.h>
43 #include <TF1.h>
44 #include <TFile.h>
45 #include <TFitter.h>
46 #include "TMinuitMinimizer.h"
47
48 //AliRoot
49 #include <AliRawReaderDate.h>
50 #include <AliRawEventHeaderBase.h>
51 #include <AliZDCRawStream.h>
52
53
54 /* Main routine
55       Arguments: 
56       1- monitoring data source
57 */
58 int main(int argc, char **argv) {
59   
60
61   gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
62                                         "*",
63                                         "TStreamerInfo",
64                                         "RIO",
65                                         "TStreamerInfo()"); 
66
67   TMinuitMinimizer m; 
68   gROOT->GetPluginManager()->AddHandler("ROOT::Math::Minimizer", "Minuit","TMinuitMinimizer",
69       "Minuit", "TMinuitMinimizer(const char *)");
70   TVirtualFitter::SetDefaultFitter("Minuit");
71
72   int status = 0;
73   // No. of ZDC cabled ch.
74   int const kNChannels = 24;
75   int const kNScChannels = 32;
76   Int_t kFirstADCGeo=0, kLastADCGeo=3;
77
78   /* log start of process */
79   printf("ZDC EMD program started\n");  
80
81   /* check that we got some arguments = list of files */
82   if (argc<2) {
83     printf("Wrong number of arguments\n");
84     return -1;
85   }
86   
87   // --- Preparing histos for EM dissociation spectra
88   //
89   TH2F * hZDCvsZEM  = new TH2F("hZDCvsZEM","EZDCTot vs. EZEM",100,0.,5.,100,0.,800.);
90   TH2F * hZDCCvsZEM = new TH2F("hZDCCvsZEM","EZDCC vs. EZEM",100,0.,5.,100,0.,400.);
91   TH2F * hZDCAvsZEM = new TH2F("hZDCAvsZEM","EZDCA vs. EZEM",100,0.,5.,100,0.,400.);
92   
93   /* open result file */
94   FILE *fp=NULL;
95   fp=fopen("./result.txt","a");
96   if (fp==NULL) {
97     printf("Failed to open file\n");
98     return -1;
99   }
100   
101   FILE *mapFile4Shuttle;
102
103   // *** To analyze CALIBRATION_MB events a pedestal data file is needed!!!
104   // *** -> check if a pedestal run has been analyzed
105   int read = 0;
106   read = daqDA_DB_getFile(PEDDATA_FILE,PEDDATA_FILE);
107   if(read){
108     printf("\t ERROR!!! ZDCPedestal.dat file NOT FOUND in DAQ db!!!\n");
109     return -1;
110   }
111   else printf("\t ZDCPedestal.dat file retrieved from DAQ db\n");
112   
113   FILE *filePed = fopen(PEDDATA_FILE,"r");
114   if (filePed==NULL) {
115     printf("\t ERROR!!! Can't open ZDCPedestal.dat file!!!\n");
116     return -1;
117   }
118
119   // 144 = 48 in-time + 48 out-of-time + 48 correlations
120   Float_t readValues[2][6*kNChannels];
121   Float_t MeanPedhg[kNChannels], MeanPedlg[kNChannels];
122   Float_t CorrCoeff0[2*kNChannels], CorrCoeff1[2*kNChannels];
123   // ***************************************************
124   //   Unless we have a narrow correlation to fit we
125   //    don't fit and store in-time vs. out-of-time
126   //    histograms -> mean pedstal subtracted!!!!!!
127   // ***************************************************
128   //
129   for(int jj=0; jj<6*kNChannels; jj++){
130     for(int ii=0; ii<2; ii++){
131        fscanf(filePed,"%f",&readValues[ii][jj]);
132     }
133     if(jj<kNChannels){
134       MeanPedhg[jj] = readValues[0][jj];
135       //printf("\t MeanPedhg[%d] = %1.1f\n",jj, MeanPedhg[jj]);
136     }
137     else if(jj>=kNChannels && jj<2*kNChannels){
138       MeanPedlg[jj-kNChannels] = readValues[0][jj];
139       //printf("\t MeanPedlg[%d] = %1.1f\n",jj-kNChannels, MeanPedlg[jj-kNChannels]);
140     }
141     else if(jj>4*kNChannels){
142       CorrCoeff0[jj-4*kNChannels] = readValues[0][jj]; 
143       CorrCoeff1[jj-4*kNChannels] = readValues[1][jj];;
144     }
145   }
146
147   /* report progress */
148   daqDA_progressReport(10);
149
150
151   /* init some counters */
152   int nevents_physics=0;
153   int nevents_total=0;
154
155   struct eventHeaderStruct *event;
156   eventTypeType eventT;
157
158   /* read the data files */
159   int n;
160   for(n=1;n<argc;n++){
161    
162     status=monitorSetDataSource( argv[n] );
163     if (status!=0) {
164       printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
165       return -1;
166     }
167
168     /* report progress */
169     /* in this example, indexed on the number of files */
170     daqDA_progressReport(10+80*n/argc);
171
172     /* read the file */
173     for(;;){
174
175       /* get next event */
176       status=monitorGetEventDynamic((void **)&event);
177       if(status==MON_ERR_EOF) break; /* end of monitoring file has been reached */
178       if(status!=0) {
179         printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
180         return -1;
181       }
182
183       /* retry if got no event */
184       if(event==NULL) {
185         break;
186       }
187       
188       // Initalize raw-data reading and decoding
189       AliRawReader *reader = new AliRawReaderDate((void*)event);
190       reader->Select("ZDC");
191       // --- Reading event header
192       //UInt_t evtype = reader->GetType();
193       //printf("\n\t ZDCEMDda -> ev. type %d\n",evtype);
194       //printf("\t ZDCEMDda -> run # %d\n",reader->GetRunNumber());
195       //
196       AliZDCRawStream *rawStreamZDC = new AliZDCRawStream(reader);
197         
198
199       /* use event - here, just write event id to result file */
200       eventT=event->eventType;
201       
202       Int_t ich=0;
203       Int_t adcMod[2*kNChannels], adcCh[2*kNChannels], sigCode[2*kNChannels];
204       Int_t det[2*kNChannels], sec[2*kNChannels];
205       for(Int_t y=0; y<2*kNChannels; y++){
206         adcMod[y]=adcCh[y]=sigCode[y]=det[y]=sec[y]=0;
207       }
208       
209       Int_t iScCh=0;
210       Int_t scMod[kNScChannels], scCh[kNScChannels], scSigCode[kNScChannels];
211       Int_t scDet[kNScChannels], scSec[kNScChannels];
212       for(Int_t y=0; y<kNScChannels; y++){
213         scMod[y]=scCh[y]=scSigCode[y]=scDet[y]=scSec[y]=0;
214       }
215       //
216       Int_t modNum=-1, modType=-1;
217       
218       if(eventT==START_OF_DATA){
219                 
220         rawStreamZDC->SetSODReading(kTRUE);
221         
222         // --------------------------------------------------------
223         // --- Writing ascii data file for the Shuttle preprocessor
224         mapFile4Shuttle = fopen(MAPDATA_FILE,"w");
225         if(!rawStreamZDC->Next()) printf(" \t No raw data found!! \n");
226         else{
227           while((rawStreamZDC->Next())){
228             if(rawStreamZDC->IsHeaderMapping()){ // mapping header
229                modNum = rawStreamZDC->GetADCModule();
230                modType = rawStreamZDC->GetModType();
231             }
232             if(rawStreamZDC->IsChMapping()){ 
233               if(modType==1){ // ADC mapping ----------------------
234                 adcMod[ich]  = rawStreamZDC->GetADCModFromMap(ich);
235                 adcCh[ich]   = rawStreamZDC->GetADCChFromMap(ich);
236                 sigCode[ich] = rawStreamZDC->GetADCSignFromMap(ich);
237                 det[ich]     = rawStreamZDC->GetDetectorFromMap(ich);
238                 sec[ich]     = rawStreamZDC->GetTowerFromMap(ich);
239                 //
240                 fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\t%d\t%d\t%d\n",
241                   ich,adcMod[ich],adcCh[ich],sigCode[ich],det[ich],sec[ich]);
242                 //
243                 //printf("  Mapping in DA -> %d ADC: mod %d ch %d, code %d det %d, sec %d\n",
244                 //  ich,adcMod[ich],adcCh[ich],sigCode[ich],det[ich],sec[ich]);
245                 //
246                 ich++;
247               }
248               else if(modType==2){ //VME scaler mapping --------------------
249                 scMod[iScCh]     = rawStreamZDC->GetScalerModFromMap(iScCh);
250                 scCh[iScCh]      = rawStreamZDC->GetScalerChFromMap(iScCh);
251                 scSigCode[iScCh] = rawStreamZDC->GetScalerSignFromMap(iScCh);
252                 scDet[iScCh]     = rawStreamZDC->GetScDetectorFromMap(iScCh);
253                 scSec[iScCh]    = rawStreamZDC->GetScTowerFromMap(iScCh);
254                 //
255                 fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\t%d\t%d\t%d\n",
256                   iScCh,scMod[iScCh],scCh[iScCh],scSigCode[iScCh],scDet[iScCh],scSec[iScCh]);
257                 //
258                 //printf("  Mapping in DA -> %d Scaler: mod %d ch %d, code %d det %d, sec %d\n",
259                 //  iScCh,scMod[iScCh],scCh[iScCh],scSigCode[iScCh],scDet[iScCh],scSec[iScCh]);
260                 //
261                 iScCh++;
262               }
263             }
264           }
265         }
266         fclose(mapFile4Shuttle);
267       }// SOD event
268     
269     if(eventT==PHYSICS_EVENT){
270       // --- Reading data header
271       reader->ReadHeader();
272       const AliRawDataHeader* header = reader->GetDataHeader();
273       if(header){
274          UChar_t message = header->GetAttributes();
275          if((message & 0xf0) == 0x60){ // DEDICATED CALIBRATION MB RUN
276             //printf("\t STANDALONE_CALIBRAION_MB_RUN raw data found\n");
277             continue;
278          }
279          else{
280             printf("\t NO STANDALONE_MB_RUN raw data found\n");
281             return -1;
282          }
283       }
284       else{
285          printf("\t ATTENTION! No Raw Data Header found!!!\n");
286          return -1;
287       }
288       
289       rawStreamZDC->SetSODReading(kTRUE);
290
291       if (!rawStreamZDC->Next()) printf(" \t No raw data found!! ");
292       //
293       // ----- Setting ch. mapping -----
294       for(Int_t jk=0; jk<2*kNChannels; jk++){
295         rawStreamZDC->SetMapADCMod(jk, adcMod[jk]);
296         rawStreamZDC->SetMapADCCh(jk, adcCh[jk]);
297         rawStreamZDC->SetMapADCSig(jk, sigCode[jk]);
298         rawStreamZDC->SetMapDet(jk, det[jk]);
299         rawStreamZDC->SetMapTow(jk, sec[jk]);
300       }
301       //
302       Float_t rawZDC=0., rawZEM=0.;
303       Float_t corrZDC=0., corrZEM=0.;
304       Float_t corrZDCTot=0., corrZDCC=0., corrZDCA=0., corrZEMTot=0.;
305       //
306       while(rawStreamZDC->Next()){
307         Int_t detector = rawStreamZDC->GetSector(0);
308         Int_t quad = rawStreamZDC->GetSector(1);
309         Int_t index=-1;
310         
311         if( (rawStreamZDC->IsADCDataWord()) && !(rawStreamZDC->IsUnderflow())
312              && !(rawStreamZDC->IsOverflow()) && (detector!=-1)
313              && ((rawStreamZDC->GetADCGain()) == 0) && // Selecting HIGH RES ch.s
314              (rawStreamZDC->GetADCModule()>=kFirstADCGeo) && (rawStreamZDC->GetADCModule()<=kLastADCGeo)){
315
316             //printf("  IsADCWord %d, IsUnderflow %d, IsOverflow %d\n",
317             //  rawStreamZDC->IsADCDataWord(),rawStreamZDC->IsUnderflow(),rawStreamZDC->IsOverflow());
318           
319             if(quad!=5){ // Physics signals
320               if(detector==1)      index = quad;   // *** ZNC
321               else if(detector==2) index = quad+5; // *** ZPC
322               else if(detector==3) index = quad+9; // *** ZEM
323               else if(detector==4) index = quad+12;// *** ZNA
324               else if(detector==5) index = quad+17;// *** ZPA
325             }
326             else continue;
327             //
328             if(index==-1) printf("ERROR in ZDCCALIBMBda.cxx -> det %d quad %d index %d ADC %d\n", 
329               detector, quad, index, rawStreamZDC->GetADCValue());
330             
331             // Mean pedestal subtraction 
332             Float_t Pedestal = MeanPedhg[index];
333             // Pedestal subtraction from correlation with out-of-time signals
334             //Float_t Pedestal = CorrCoeff0[index]+CorrCoeff1[index]*MeanPedOOT[index];
335             //
336             if(index!=-1 && quad!=5){ 
337               //
338               if(detector==3){
339                rawZEM  = (Float_t) rawStreamZDC->GetADCValue();
340                corrZEM = (rawStreamZDC->GetADCValue()) - Pedestal;
341                corrZEMTot += corrZEM;
342               }
343               else{
344                rawZDC  = (Float_t) rawStreamZDC->GetADCValue();
345                corrZDC = rawZDC - Pedestal;
346                if(detector==1 || detector==2) corrZDCC += corrZDC;
347                else corrZDCA += corrZDC;
348                corrZDCTot += corrZDC;
349               }
350             }
351         }//IsADCDataWord()
352         
353        } // Next()
354         
355        hZDCvsZEM ->Fill(corrZDCTot, corrZEMTot);
356        hZDCCvsZEM->Fill(corrZDCC, corrZEMTot);
357        hZDCAvsZEM->Fill(corrZDCA, corrZEMTot);
358
359        nevents_physics++;
360     }//(if PHYSICS_EVENT)
361     
362     nevents_total++;
363
364    }
365
366    /* free resources */
367    free(event);
368   }
369     
370   //
371   TFile *histofile = new TFile(MBCALIBDATA_FILE,"RECREATE");
372   histofile->cd();
373   hZDCvsZEM ->Write();
374   hZDCCvsZEM->Write();
375   hZDCAvsZEM->Write();
376   histofile->Close();
377   
378   if(hZDCvsZEM)  delete hZDCvsZEM ;
379   if(hZDCCvsZEM) delete hZDCCvsZEM;
380   if(hZDCAvsZEM) delete hZDCAvsZEM;
381   
382   /* write report */
383   fprintf(fp,"Run #%s, received %d physics events out of %d\n",getenv("DATE_RUN_NUMBER"),nevents_physics,nevents_total);
384
385   /* close result file */
386   fclose(fp);
387   
388   /* report progress */
389   daqDA_progressReport(90);
390
391   /* store the result file on FES */
392   status = daqDA_FES_storeFile(MAPDATA_FILE, "MAPPING");
393   if(status){
394     printf("Failed to export file : %d\n",status);
395     return -1;
396   }
397   //
398   status = daqDA_FES_storeFile(MBCALIBDATA_FILE, "MBCALIB");
399   if(status){
400     printf("Failed to export file : %d\n",status);
401     return -1;
402   }
403  
404   /* report progress */
405   daqDA_progressReport(100);
406
407   return status;
408 }