]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ZDC/ZDCLASERda.cxx
DA for CALIBRATION_MB runs
[u/mrichter/AliRoot.git] / ZDC / ZDCLASERda.cxx
1 /*
2
3 This program reads the DAQ data files passed as argument using the monitoring library.
4
5 The program reports about its processing progress.
6
7 Messages on stdout are exported to DAQ log system.
8
9 DA for ZDC standalone pedestal runs
10
11 Contact: Chiara.Oppedisano@to.infn.it
12 Link: 
13 Run Type: STANDALONE_LASER_RUN
14 DA Type: LDC
15 Number of events needed: no constraint (tipically ~10^3)
16 Input Files: ZDCPedestal.dat
17 Output Files: ZDCLaser.dat
18 Trigger Types Used: Standalone Trigger
19
20 */
21 #define PEDDATA_FILE  "ZDCPedestal.dat"
22 #define MAPDATA_FILE  "ZDCChMapping.dat"
23 #define LASHISTO_FILE "ZDCLaserHisto.root"
24 #define LASDATA_FILE  "ZDCLaserCalib.dat"
25
26 #include <stdio.h>
27 #include <stdlib.h>
28 #include <Riostream.h>
29
30 // DATE
31 #include <event.h>
32 #include <monitor.h>
33 #include <daqDA.h>
34
35 //ROOT
36 #include <TROOT.h>
37 #include <TPluginManager.h>
38 #include <TH1F.h>
39 #include <TF1.h>
40 #include <TFile.h>
41 #include <TFitter.h>
42 #include "TMinuitMinimizer.h"
43
44 //AliRoot
45 #include <AliRawReaderDate.h>
46 #include <AliRawEventHeaderBase.h>
47 #include <AliZDCRawStream.h>
48
49
50 /* Main routine
51       Arguments: list of DATE raw data files
52 */
53 int main(int argc, char **argv) {
54   
55   gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
56                                         "*",
57                                         "TStreamerInfo",
58                                         "RIO",
59                                         "TStreamerInfo()"); 
60
61   TMinuitMinimizer m; 
62   gROOT->GetPluginManager()->AddHandler("ROOT::Math::Minimizer", "Minuit","TMinuitMinimizer",
63       "Minuit", "TMinuitMinimizer(const char *)");
64   TVirtualFitter::SetDefaultFitter("Minuit");
65
66   int status = 0;
67   int const kNChannels = 24;
68   int const kNScChannels = 32;
69   Int_t kFirstADCGeo=0, kLastADCGeo=3;
70   
71   Int_t ich=0;
72   Int_t adcMod[2*kNChannels], adcCh[2*kNChannels], sigCode[2*kNChannels];
73   Int_t det[2*kNChannels], sec[2*kNChannels];
74   for(Int_t y=0; y<2*kNChannels; y++){
75     adcMod[y]=adcCh[y]=sigCode[y]=det[y]=sec[y]=0;
76   }
77
78   /* log start of process */
79   printf("\n ZDC LASER 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   // --- Histograms for LASER runs 
88   //     20 signal channels + 2 reference PTMs
89   //
90   TH1F::AddDirectory(0);
91   // --- Histos for reference PMTs (high gain chains)
92   TH1F *hPMRefChg = new TH1F("hPMRefChg","hPMRefChg", 100,0.,1000.);
93   TH1F *hPMRefAhg = new TH1F("hPMRefAhg","hPMRefAhg", 100,0.,1000.);
94   TH1F *hPMRefClg = new TH1F("hPMRefClg","hPMRefClg", 100,0.,4000.);
95   TH1F *hPMRefAlg = new TH1F("hPMRefAlg","hPMRefAlg", 100,0.,4000.);
96   //
97   // --- Histos for detector PMTs 
98   TH1F *hZNChg[5], *hZPChg[5], *hZNAhg[5], *hZPAhg[5], *hZEMhg[2];
99   TH1F *hZNClg[5], *hZPClg[5], *hZNAlg[5], *hZPAlg[5], *hZEMlg[2];
100   char hnamZNChg[20], hnamZPChg[20], hnamZNAhg[20], hnamZPAhg[20];
101   char hnamZNClg[20], hnamZPClg[20], hnamZNAlg[20], hnamZPAlg[20];
102   char hnamZEMhg[20], hnamZEMlg[20];
103   for(Int_t j=0; j<5; j++){
104     sprintf(hnamZNChg,"ZNChg-tow%d",j);
105     sprintf(hnamZPChg,"ZPChg-tow%d",j);
106     sprintf(hnamZNAhg,"ZNAhg-tow%d",j);
107     sprintf(hnamZPAhg,"ZPAhg-tow%d",j);
108     //
109     hZNChg[j] = new TH1F(hnamZNChg, hnamZNChg, 100, 0., 1000.);
110     hZPChg[j] = new TH1F(hnamZPChg, hnamZPChg, 100, 0., 1000.);
111     hZNAhg[j] = new TH1F(hnamZNAhg, hnamZNAhg, 100, 0., 1000.);
112     hZPAhg[j] = new TH1F(hnamZPAhg, hnamZPAhg, 100, 0., 1000.);
113     //
114     sprintf(hnamZNClg,"ZNClg-tow%d",j);
115     sprintf(hnamZPClg,"ZPClg-tow%d",j);
116     sprintf(hnamZNAlg,"ZNAlg-tow%d",j);
117     sprintf(hnamZPAlg,"ZPAlg-tow%d",j);
118     //
119     hZNClg[j] = new TH1F(hnamZNClg, hnamZNClg, 100, 0., 4000.);
120     hZPClg[j] = new TH1F(hnamZPClg, hnamZPClg, 100, 0., 4000.);
121     hZNAlg[j] = new TH1F(hnamZNAlg, hnamZNAlg, 100, 0., 4000.);
122     hZPAlg[j] = new TH1F(hnamZPAlg, hnamZPAlg, 100, 0., 4000.);
123     //
124     if(j<2){
125       sprintf(hnamZEMhg,"ZEM%dhg",j);
126       sprintf(hnamZEMlg,"ZEM%dlg",j);
127       //
128       hZEMhg[j] = new TH1F(hnamZEMhg, hnamZEMhg, 100, 0., 1000.);      
129       hZEMlg[j] = new TH1F(hnamZEMlg, hnamZEMlg, 100, 0., 4000.);      
130     }
131   }
132
133   /* open result file */
134   FILE *fp=NULL;
135   fp=fopen("./result.txt","a");
136   if (fp==NULL) {
137     printf("Failed to open file\n");
138     return -1;
139   }
140   /* report progress */
141   daqDA_progressReport(10);
142         
143   // *** To analyze LASER events you MUST have a pedestal data file!!!
144   // *** -> check if a pedestal run has been analyzed
145   int read = 0;
146   read = daqDA_DB_getFile(PEDDATA_FILE, PEDDATA_FILE);
147   if(read){
148     printf("\t ERROR!!! ZDCPedestal.dat file NOT FOUND in DAQ db!!!\n");
149     return -1;
150   }
151   else printf("\t ZDCPedestal.dat file retrieved from DAQ db\n");
152   
153   FILE *filePed = fopen(PEDDATA_FILE,"r");
154   if (filePed==NULL) {
155     printf("\t ERROR!!! Can't open ZDCPedestal.dat file!!!\n");
156     return -1;
157   }
158
159   // 144 = 48 in-time + 48 out-of-time + 48 correlations
160   Float_t readValues[2][6*kNChannels];
161   Float_t MeanPedhg[kNChannels], MeanPedlg[kNChannels];
162   Float_t CorrCoeff0[2*kNChannels], CorrCoeff1[2*kNChannels];
163   // ***************************************************
164   //   Unless we have a narrow correlation to fit we
165   //    don't fit and store in-time vs. out-of-time
166   //    histograms -> mean pedstal subtracted!!!!!!
167   // ***************************************************
168   //
169   for(int jj=0; jj<6*kNChannels; jj++){
170     for(int ii=0; ii<2; ii++){
171        fscanf(filePed,"%f",&readValues[ii][jj]);
172     }
173     if(jj<kNChannels){
174       MeanPedhg[jj] = readValues[0][jj];
175       //printf("\t MeanPedhg[%d] = %1.1f\n",jj, MeanPedhg[jj]);
176     }
177     else if(jj>=kNChannels && jj<2*kNChannels){
178       MeanPedlg[jj-kNChannels] = readValues[0][jj];
179       //printf("\t MeanPedlg[%d] = %1.1f\n",jj-kNChannels, MeanPedlg[jj-kNChannels]);
180     }
181     else if(jj>4*kNChannels){
182       CorrCoeff0[jj-4*kNChannels] = readValues[0][jj]; 
183       CorrCoeff1[jj-4*kNChannels] = readValues[1][jj];;
184     }
185   }
186   
187   FILE *mapFile4Shuttle;
188
189   /* report progress */
190   daqDA_progressReport(20);
191
192
193   /* init some counters */
194   int nevents_physics=0;
195   int nevents_total=0;
196
197   struct eventHeaderStruct *event;
198   eventTypeType eventT;
199
200   /* read the data files */
201   int n;
202   for(n=1;n<argc;n++) {
203    
204     status=monitorSetDataSource( argv[n] );
205     if (status!=0) {
206       printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
207       return -1;
208     }
209
210     /* report progress */
211     /* in this example, indexed on the number of files */
212     daqDA_progressReport(20+70*n/argc);
213
214     /* read the file */
215     for(;;) {
216
217       /* get next event */
218       status=monitorGetEventDynamic((void **)&event);
219       if (status==MON_ERR_EOF) break; /* end of monitoring file has been reached */
220       if (status!=0) {
221         printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
222         return -1;
223       }
224
225       /* retry if got no event */
226       if (event==NULL) {
227         break;
228       }
229
230       // Initalize raw-data reading and decoding
231       AliRawReader *reader = new AliRawReaderDate((void*)event);
232       reader->Select("ZDC");
233       // --- Reading event header
234       //UInt_t evtype = reader->GetType();
235       //printf("\n\t ZDCLASERda -> ev. type %d\n",evtype);
236       //printf("\t ZDCLASERda -> run # %d\n",reader->GetRunNumber());
237       //
238       AliZDCRawStream *rawStreamZDC = new AliZDCRawStream(reader);
239         
240
241       /* use event - here, just write event id to result file */
242       eventT=event->eventType;
243             
244       Int_t iScCh=0;
245       Int_t scMod[kNScChannels], scCh[kNScChannels], scSigCode[kNScChannels];
246       Int_t scDet[kNScChannels], scSec[kNScChannels];
247       for(Int_t y=0; y<kNScChannels; y++){
248         scMod[y]=scCh[y]=scSigCode[y]=scDet[y]=scSec[y]=0;
249       }
250       //
251       Int_t modNum=-1, modType=-1;
252       
253       if(eventT==START_OF_DATA){
254                 
255         rawStreamZDC->SetSODReading(kTRUE);
256         
257         // --------------------------------------------------------
258         // --- Writing ascii data file for the Shuttle preprocessor
259         mapFile4Shuttle = fopen(MAPDATA_FILE,"w");
260         if(!rawStreamZDC->Next()) printf(" \t No raw data found!! \n");
261         else{
262           while((rawStreamZDC->Next())){
263             if(rawStreamZDC->IsHeaderMapping()){ // mapping header
264                modNum = rawStreamZDC->GetADCModule();
265                modType = rawStreamZDC->GetModType();
266             }
267             if(rawStreamZDC->IsChMapping()){ 
268               if(modType==1){ // ADC mapping ----------------------
269                 adcMod[ich]  = rawStreamZDC->GetADCModFromMap(ich);
270                 adcCh[ich]   = rawStreamZDC->GetADCChFromMap(ich);
271                 sigCode[ich] = rawStreamZDC->GetADCSignFromMap(ich);
272                 det[ich]     = rawStreamZDC->GetDetectorFromMap(ich);
273                 sec[ich]     = rawStreamZDC->GetTowerFromMap(ich);
274                 //
275                 fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\t%d\t%d\t%d\n",
276                   ich,adcMod[ich],adcCh[ich],sigCode[ich],det[ich],sec[ich]);
277                 //
278                 //printf("  Mapping in DA -> %d ADC: mod %d ch %d, code %d det %d, sec %d\n",
279                 //  ich,adcMod[ich],adcCh[ich],sigCode[ich],det[ich],sec[ich]);
280                 //
281                 ich++;
282               }
283               else if(modType==2){ //VME scaler mapping --------------------
284                 scMod[iScCh]     = rawStreamZDC->GetScalerModFromMap(iScCh);
285                 scCh[iScCh]      = rawStreamZDC->GetScalerChFromMap(iScCh);
286                 scSigCode[iScCh] = rawStreamZDC->GetScalerSignFromMap(iScCh);
287                 scDet[iScCh]     = rawStreamZDC->GetScDetectorFromMap(iScCh);
288                 scSec[iScCh]    = rawStreamZDC->GetScTowerFromMap(iScCh);
289                 //
290                 fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\t%d\t%d\t%d\n",
291                   iScCh,scMod[iScCh],scCh[iScCh],scSigCode[iScCh],scDet[iScCh],scSec[iScCh]);
292                 //
293                 //printf("  Mapping in DA -> %d Scaler: mod %d ch %d, code %d det %d, sec %d\n",
294                 //  iScCh,scMod[iScCh],scCh[iScCh],scSigCode[iScCh],scDet[iScCh],scSec[iScCh]);
295                 //
296                 iScCh++;
297               }
298             }
299           }
300         }
301         fclose(mapFile4Shuttle);
302       }// SOD event
303
304       else if(eventT==PHYSICS_EVENT){
305         // --- Reading data header
306         reader->ReadHeader();
307         const AliRawDataHeader* header = reader->GetDataHeader();
308         if(header) {
309          UChar_t message = header->GetAttributes();
310          if((message & 0x30) == 0x30){ // DEDICATED LASER RUN
311             //printf("\t STANDALONE_LASER_RUN raw data found\n");
312          }
313          else{
314             printf("ZDCLASERda.cxx -> NO STANDALONE_LASER_RUN raw data found\n");
315             return -1;
316          }
317         }
318         else{
319            printf("\t ATTENTION! No Raw Data Header found!!!\n");
320            return -1;
321         }
322
323         rawStreamZDC->SetSODReading(kTRUE);
324
325         if (!rawStreamZDC->Next()) printf(" \t No raw data found!! \n");
326         //
327         // ----- Setting ch. mapping -----
328         for(Int_t jk=0; jk<2*kNChannels; jk++){
329           rawStreamZDC->SetMapADCMod(jk, adcMod[jk]);
330           rawStreamZDC->SetMapADCCh(jk, adcCh[jk]);
331           rawStreamZDC->SetMapADCSig(jk, sigCode[jk]);
332           rawStreamZDC->SetMapDet(jk, det[jk]);
333           rawStreamZDC->SetMapTow(jk, sec[jk]);
334         }
335         //
336         while(rawStreamZDC->Next()){
337           Int_t index=-1;
338           Int_t detector = rawStreamZDC->GetSector(0);
339           Int_t sector = rawStreamZDC->GetSector(1);
340           
341           if(rawStreamZDC->IsADCDataWord() && !(rawStreamZDC->IsUnderflow()) && 
342              !(rawStreamZDC->IsOverflow()) && detector!=-1 &&
343              rawStreamZDC->GetADCModule()>=kFirstADCGeo && rawStreamZDC->GetADCModule()<=kLastADCGeo){
344             
345             if(sector!=5){ // Physics signals
346               if(detector==1)      index = sector;   // *** ZNC
347               else if(detector==2) index = sector+5; // *** ZPC
348               else if(detector==3) index = sector+9; // *** ZEM
349               else if(detector==4) index = sector+12;// *** ZNA
350               else if(detector==5) index = sector+17;// *** ZPA
351             }
352             else{ // Reference PMs
353               index = (detector-1)/3+22;
354             }
355             //
356             if(index==-1) printf("ERROR in ZDCLASERda.cxx -> det %d quad %d res %d index %d ADC %d\n", 
357               detector, sector, rawStreamZDC->GetADCGain(), index, rawStreamZDC->GetADCValue());
358             
359             Float_t Pedestal=0.;
360             if(rawStreamZDC->GetADCGain()==0)      Pedestal = MeanPedhg[index];
361             else if(rawStreamZDC->GetADCGain()==1) Pedestal = MeanPedlg[index];
362             //
363             Float_t CorrADC = rawStreamZDC->GetADCValue() - Pedestal;
364             //
365             //printf("\tdet %d sec %d res %d index %d ped %1.0f ADCcorr %1.0f\n", 
366             //  detector, sector, rawStreamZDC->GetADCGain(), index, Pedestal,CorrADC);
367             
368             // **** Detector PMs
369             if(sector!=5){
370               if(rawStreamZDC->GetADCGain()==0){ // --- High gain chain ---
371                 // ---- side C
372                 if(detector==1)      hZNChg[sector]->Fill(CorrADC);
373                 else if(detector==2) hZPChg[sector]->Fill(CorrADC);
374                 // ---- side A
375                 else if(detector==4) hZNAhg[sector]->Fill(CorrADC);
376                 else if(detector==5) hZPAhg[sector]->Fill(CorrADC);
377                 // ---- ZEM
378                 else if(detector==3){
379                   hZEMhg[sector-1]->Fill(CorrADC);
380                 }
381               }
382               else if(rawStreamZDC->GetADCGain()==1){ // --- Low gain chain ---
383                 // ---- side C
384                 if(detector==1)      hZNClg[sector]->Fill(CorrADC);
385                 else if(detector==2) hZPClg[sector]->Fill(CorrADC);
386                 // ---- side A
387                 else if(detector==4) hZNAlg[sector]->Fill(CorrADC);
388                 else if(detector==5) hZPAlg[sector]->Fill(CorrADC);
389                 // ---- ZEM
390                 else if(detector==3) hZEMlg[sector-1]->Fill(CorrADC);
391               }
392             }
393             // **** Reference PMs
394             else if(sector==5){
395               if(rawStreamZDC->GetADCGain()==0){ // --- High gain chain ---
396                 // ---- PMRef chain side C
397                 if(detector==1) hPMRefChg->Fill(CorrADC);
398                 // ---- PMRef side A
399                 else if(detector==4) hPMRefAhg->Fill(CorrADC);
400               }
401               else if(rawStreamZDC->GetADCGain()==1){ // --- Low gain chain ---
402                 // ---- PMRef chain side C
403                 if(detector==1) hPMRefClg->Fill(CorrADC);
404                 // ---- PMRef side A
405                 else if(detector==4) hPMRefAlg->Fill(CorrADC);
406               }
407             }
408           }//IsADCDataWord()+NOunderflow+NOoverflow
409           //
410          }
411          //
412          nevents_physics++;
413          //
414          delete reader;
415          delete rawStreamZDC;
416
417       }//(if PHYSICS_EVENT) 
418
419       /* exit when last event received, no need to wait for TERM signal */
420       else if(eventT==END_OF_RUN) {
421         printf(" -> EOR event detected\n");
422         break;
423       }
424       
425       
426       nevents_total++;
427     
428     }
429           
430     /* free resources */
431     free(event);
432   }  
433   
434   /* Analysis of the histograms */
435   //
436   Int_t detector[2*kNChannels], quad[2*kNChannels];
437   Int_t maxBin[2*kNChannels], nBin[2*kNChannels];
438   Float_t xMax[2*kNChannels], maxXval[2*kNChannels], xlow[2*kNChannels]; 
439   Float_t mean[2*kNChannels], sigma[2*kNChannels];
440   for(Int_t t=0; t<2*kNChannels; t++){
441     detector[t] = quad[t] = 0;
442     maxBin[t] = nBin[t] = 0;
443     xMax[t] = maxXval[t] = xlow[t] = 0.;
444     mean[t] = sigma[t] = 0.;
445   }
446   TF1 *fun[2*kNChannels];
447   Int_t atLeastOneHisto=0;
448   
449   // ******** High gain chain ********
450   for(Int_t k=0; k<5; k++){
451     // --- ZNC
452     detector[k] = 1;
453     quad[k] = k;
454     maxBin[k] = hZNChg[k]->GetMaximumBin();
455     nBin[k] = (hZNChg[k]->GetXaxis())->GetNbins();
456     xMax[k] = (hZNChg[k]->GetXaxis())->GetXmax();
457     if(nBin[k]!=0) maxXval[k] = maxBin[k]*xMax[k]/nBin[k];
458     if(maxXval[k]-150.<0.) xlow[k]=0.;
459     else xlow[k] = maxXval[k]-150.;
460     // checking if at least one histo is fitted
461     if(hZNChg[k]->GetEntries()!=0 && hZNChg[k]->GetMean()>0){
462       atLeastOneHisto=1;
463       //
464       hZNChg[k]->Fit("gaus","Q","",xlow[k],maxXval[k]+150.);
465       fun[k] = hZNChg[k]->GetFunction("gaus");
466       mean[k]  = (Float_t) (fun[k]->GetParameter(1));
467       sigma[k] = (Float_t) (fun[k]->GetParameter(2));
468     }
469     // --- ZPC
470     detector[k+5] = 2;
471     quad[k+5] = k;
472     maxBin[k+5] = hZPChg[k]->GetMaximumBin();
473     nBin[k+5] = (hZPChg[k]->GetXaxis())->GetNbins();
474     xMax[k+5] = (hZPChg[k]->GetXaxis())->GetXmax();
475     if(nBin[k+5]!=0) maxXval[k+5] = maxBin[k+5]*xMax[k+5]/nBin[k+5];
476     if(maxXval[k+5]-150.<0.) xlow[k+5]=0.;
477     else xlow[k+5] = maxXval[k+5]-150.;
478     if(hZPChg[k]->GetEntries()!=0 && hZPChg[k]->GetMean()>0){
479       atLeastOneHisto=1; 
480       //
481       hZPChg[k]->Fit("gaus","Q","",xlow[k+5],maxXval[k+5]+150.);
482       fun[k+5] = hZPChg[k]->GetFunction("gaus");
483       mean[k+5]  = (Float_t) (fun[k+5]->GetParameter(1));
484       sigma[k+5] = (Float_t) (fun[k+5]->GetParameter(2));
485     }
486     // --- ZEM1
487     if(k<2){
488       detector[k+10] = 3;
489       quad[k+10] = k+1;
490       maxBin[k+10] = hZEMhg[k]->GetMaximumBin();
491       nBin[k+10] = (hZEMhg[k]->GetXaxis())->GetNbins();
492       xMax[k+10] = (hZEMhg[k]->GetXaxis())->GetXmax();
493       if(nBin[k+10]!=0) maxXval[k+10] = maxBin[k+10]*xMax[k+10]/nBin[k+10];
494       if(maxXval[k+10]-150.<0.) xlow[k+10]=0.;
495       else xlow[k+10] = maxXval[k+10]-150.;
496       if(hZEMhg[k]->GetEntries()!=0 && hZEMhg[k]->GetMean()>0){
497         atLeastOneHisto=1; 
498         //
499         hZEMhg[k]->Fit("gaus","Q","",xlow[k+10],maxXval[k+10]+150.);
500         fun[k+10] = hZEMhg[k]->GetFunction("gaus");
501         mean[k+10]  = (Float_t) (fun[k+10]->GetParameter(1));
502         sigma[k+10] = (Float_t) (fun[k+10]->GetParameter(2));
503       }
504     }
505     // --- ZNA
506     detector[k+12] = 4;
507     quad[k+12] = k;
508     maxBin[k+12] = hZNAhg[k]->GetMaximumBin();
509     nBin[k+12] = (hZNAhg[k]->GetXaxis())->GetNbins();
510     xMax[k+12] = (hZNAhg[k]->GetXaxis())->GetXmax();
511     if(nBin[k+12]!=0) maxXval[k+12] = maxBin[k+12]*xMax[k+12]/nBin[k+12];
512     if(maxXval[k+12]-150.<0.) xlow[k+12]=0.;
513     else xlow[k+12] = maxXval[k+12]-150.;
514     if(hZNAhg[k]->GetEntries()!=0 && hZNAhg[k]->GetMean()>0){
515       atLeastOneHisto=1; 
516       //
517       hZNAhg[k]->Fit("gaus","Q","",xlow[k+12],maxXval[k+12]+150.);
518       fun[k+12] = hZNAhg[k]->GetFunction("gaus");
519       mean[k+12]  = (Float_t) (fun[k+12]->GetParameter(1));
520       sigma[k+12] = (Float_t) (fun[k+12]->GetParameter(2));
521     }
522     // --- ZPA
523     detector[k+17] = 4;
524     quad[k+17] = 5;
525     maxBin[k+17] = hZPAhg[k]->GetMaximumBin();
526     nBin[k+17] = (hZPAhg[k]->GetXaxis())->GetNbins();
527     xMax[k+17] = (hZPAhg[k]->GetXaxis())->GetXmax();
528     if(nBin[k+17]!=0) maxXval[k+17] = maxBin[k+17]*xMax[k+17]/nBin[k+17];
529     if(maxXval[k+17]-150.<0.) xlow[k+17]=0.;
530     else xlow[k+17] = maxXval[k+17]-150.;
531     if(hZPAhg[k]->GetEntries()!=0 && hZPAhg[k]->GetMean()>0){
532       atLeastOneHisto=1; 
533       //
534       hZPAhg[k]->Fit("gaus","Q","",xlow[k+17],maxXval[k+17]+150.);
535       fun[k+17] = hZPAhg[k]->GetFunction("gaus");
536       mean[k+17]  = (Float_t) (fun[k+17]->GetParameter(1));
537       sigma[k+17] = (Float_t) (fun[k+17]->GetParameter(2));    
538     }
539   }
540   // ~~~~~~~~ PM Ref side C ~~~~~~~~
541   detector[22] = 1;
542   quad[22] = 5;
543   maxBin[22] = hPMRefChg->GetMaximumBin();
544   nBin[22] = (hPMRefChg->GetXaxis())->GetNbins();
545   xMax[22] = (hPMRefChg->GetXaxis())->GetXmax();
546   if(nBin[22]!=0) maxXval[22] = maxBin[22]*xMax[22]/nBin[22];
547   if(maxXval[22]-150.<0.) xlow[22]=0.;
548   else xlow[22] = maxXval[22]-150.;
549   if(hPMRefChg->GetEntries()!=0 && hPMRefChg->GetMean()>0){
550     atLeastOneHisto=1; 
551     //
552     hPMRefChg->Fit("gaus","Q","",xlow[22],maxXval[22]+150.);
553     fun[22] = hPMRefChg->GetFunction("gaus");
554     mean[22]  = (Float_t) (fun[22]->GetParameter(1));
555     sigma[22] = (Float_t) (fun[22]->GetParameter(2));
556   }
557   // ~~~~~~~~ PM Ref side A ~~~~~~~~
558   detector[23] = 4;
559   quad[23] = 5;
560   maxBin[23] = hPMRefAhg->GetMaximumBin();
561   nBin[23] = (hPMRefAhg->GetXaxis())->GetNbins();
562   xMax[23] = (hPMRefAhg->GetXaxis())->GetXmax();
563   if(nBin[23]!=0) maxXval[23] = maxBin[23]*xMax[23]/nBin[23];
564   if(maxXval[23]-100.<0.) xlow[23]=0.;
565   else xlow[23] = maxXval[23]-150.;
566   if(hPMRefAhg->GetEntries()!=0 && hPMRefAhg->GetMean()>0){
567     atLeastOneHisto=1; 
568     //
569     hPMRefAhg->Fit("gaus","Q","",xlow[23],maxXval[23]+100.);
570     fun[23] = hPMRefAhg->GetFunction("gaus");
571     mean[23]  = (Float_t) (fun[23]->GetParameter(1));
572     sigma[23] = (Float_t) (fun[23]->GetParameter(2));
573   }
574   
575   // ******** Low gain chain ********
576   Int_t kOffset = 24;
577   for(Int_t k=0; k<5; k++){
578     // --- ZNC
579     detector[k+kOffset] = 1;
580     quad[k+kOffset] = k;
581     maxBin[k+kOffset] = hZNClg[k]->GetMaximumBin();
582     nBin[k+kOffset] = (hZNClg[k]->GetXaxis())->GetNbins();
583     xMax[k+kOffset] = (hZNClg[k]->GetXaxis())->GetXmax();
584     if(nBin[k+kOffset]!=0) maxXval[k+kOffset] = maxBin[k+kOffset]*xMax[k+kOffset]/nBin[k+kOffset];
585     if(maxXval[k+kOffset]-150.<0.) xlow[k+kOffset]=0.;
586     else xlow[k+kOffset] = maxXval[k+kOffset]-150.;
587     if(hZNClg[k]->GetEntries()!=0 && hZNClg[k]->GetMean()>0){
588       atLeastOneHisto=1; 
589       //
590       hZNClg[k]->Fit("gaus","Q","",xlow[k+kOffset],maxXval[k+kOffset]+150.);
591       fun[k+kOffset] = hZNClg[k]->GetFunction("gaus");
592       mean[k+kOffset]  = (Float_t) (fun[k+kOffset]->GetParameter(1));
593       sigma[k+kOffset] = (Float_t) (fun[k+kOffset]->GetParameter(2));
594     }
595     // --- ZPC
596     detector[k+kOffset+5] = 2;
597     quad[k+kOffset+5] = k;
598     maxBin[k+kOffset+5] = hZPClg[k]->GetMaximumBin();
599     nBin[k+kOffset+5] = (hZPClg[k]->GetXaxis())->GetNbins();
600     xMax[k+kOffset+5] = (hZPClg[k]->GetXaxis())->GetXmax();
601     if(nBin[k+kOffset+5]!=0) maxXval[k+kOffset+5] = maxBin[k+kOffset+5]*xMax[k+kOffset+5]/nBin[k+kOffset+5];
602     if(maxXval[k+kOffset+5]-150.<0.) xlow[k+kOffset+5]=0.;
603     else xlow[k+kOffset+5] = maxXval[k+kOffset+5]-150.;
604     if(hZPClg[k]->GetEntries()!=0 && hZPClg[k]->GetMean()>0){
605       atLeastOneHisto=1;  
606       //
607       hZPClg[k]->Fit("gaus","Q","",xlow[k+kOffset+5],maxXval[k+kOffset+5]+150.);
608       fun[k+kOffset+5] = hZPClg[k]->GetFunction("gaus");
609       mean[k+kOffset+5]  = (Float_t) (fun[k+kOffset+5]->GetParameter(1));
610       sigma[k+kOffset+5] = (Float_t) (fun[k+kOffset+5]->GetParameter(2));
611     }
612     // --- ZEM1
613     if(k+kOffset<2){
614       detector[k+kOffset+10] = 3;
615       quad[k+kOffset+10] = k+1;
616       maxBin[k+kOffset+10] = hZEMlg[k]->GetMaximumBin();
617       nBin[k+kOffset+10] = (hZEMlg[k]->GetXaxis())->GetNbins();
618       xMax[k+kOffset+10] = (hZEMlg[k]->GetXaxis())->GetXmax();
619       if(nBin[k+kOffset+10]!=0) maxXval[k+kOffset+10] = maxBin[k+kOffset+10]*xMax[k+kOffset+10]/nBin[k+kOffset+10];
620       if(maxXval[k+kOffset+10]-150.<0.) xlow[k+kOffset+10]=0.;
621       else xlow[k+kOffset+10] = maxXval[k+kOffset+10]-150.;
622       if(hZEMlg[k]->GetEntries()!=0 && hZEMlg[k]->GetMean()>0){
623         atLeastOneHisto=1;  
624         //
625         hZEMlg[k]->Fit("gaus","Q","",xlow[k+kOffset+10],maxXval[k+kOffset+10]+150.);
626         fun[k+kOffset+10] = hZEMlg[k]->GetFunction("gaus");
627         mean[k+kOffset+10]  = (Float_t) (fun[k+kOffset+10]->GetParameter(1));
628         sigma[k+kOffset+10] = (Float_t) (fun[k+kOffset+10]->GetParameter(2));
629       }
630     }
631     // --- ZNA
632     detector[k+kOffset+12] = 4;
633     quad[k+kOffset+12] = k;
634     maxBin[k+kOffset+12] = hZNAlg[k]->GetMaximumBin();
635     nBin[k+kOffset+12] = (hZNAlg[k]->GetXaxis())->GetNbins();
636     xMax[k+kOffset+12] = (hZNAlg[k]->GetXaxis())->GetXmax();
637     if(nBin[k+kOffset+12]!=0) maxXval[k+kOffset+12] = maxBin[k+kOffset+12]*xMax[k+kOffset+12]/nBin[k+kOffset+12];
638     if(maxXval[k+kOffset+12]-150.<0.) xlow[k+kOffset+12]=0.;
639     else xlow[k+kOffset+12] = maxXval[k+kOffset+12]-150.;
640     if(hZNAlg[k]->GetEntries()!=0 && hZNAlg[k]->GetMean()>0){
641       atLeastOneHisto=1;
642       //
643       hZNAlg[k]->Fit("gaus","Q","",xlow[k+kOffset+12],maxXval[k+kOffset+12]+150.);
644       fun[k+kOffset+12] = hZNAlg[k]->GetFunction("gaus");
645       mean[k+kOffset+12]  = (Float_t) (fun[k+kOffset+12]->GetParameter(1));
646       sigma[k+kOffset+12] = (Float_t) (fun[k+kOffset+12]->GetParameter(2));
647     }
648     // --- ZPA
649     detector[k+kOffset+17] = 5;
650     quad[k+kOffset+17] = k;
651     maxBin[k+kOffset+17] = hZPAlg[k]->GetMaximumBin();
652     nBin[k+kOffset+17] = (hZPAlg[k]->GetXaxis())->GetNbins();
653     xMax[k+kOffset+17] = (hZPAlg[k]->GetXaxis())->GetXmax();
654     if(nBin[k+kOffset+17]!=0) maxXval[k+kOffset+17] = maxBin[k+kOffset+17]*xMax[k+kOffset+17]/nBin[k+kOffset+17];
655     if(maxXval[k+kOffset+17]-150.<0.) xlow[k+kOffset+17]=0.;
656     else xlow[k+kOffset+17] = maxXval[k+kOffset+17]-150.;
657     if(hZPAlg[k]->GetEntries()!=0 && hZPAlg[k]->GetMean()>0){
658       atLeastOneHisto=1;  
659       //
660       hZPAlg[k]->Fit("gaus","Q","",xlow[k+kOffset+17],maxXval[k+kOffset+17]+150.);
661       fun[k+kOffset+17] = hZPAlg[k]->GetFunction("gaus");
662       mean[k+kOffset+17]  = (Float_t) (fun[k+kOffset+17]->GetParameter(1));
663       sigma[k+kOffset+17] = (Float_t) (fun[k+kOffset+17]->GetParameter(2)); 
664     }   
665   }
666   // ~~~~~~~~ PM Ref side C ~~~~~~~~
667   detector[46] = 1;
668   quad[46] = 5;
669   maxBin[46] = hPMRefClg->GetMaximumBin();
670   nBin[46] = (hPMRefClg->GetXaxis())->GetNbins();
671   xMax[46] = (hPMRefClg->GetXaxis())->GetXmax();
672   if(nBin[46]!=0) maxXval[46] = maxBin[46]*xMax[46]/nBin[46];
673   if(maxXval[46]-150.<0.) xlow[46]=0.;
674   else xlow[46] = maxXval[46]-150.;
675   if(hPMRefClg->GetEntries()!=0 && hPMRefClg->GetMean()>0){
676     atLeastOneHisto=1; 
677     //
678     hPMRefClg->Fit("gaus","Q","",xlow[46],maxXval[46]+150.);
679     fun[46] = hPMRefClg->GetFunction("gaus");
680     mean[46]  = (Float_t) (fun[46]->GetParameter(1));
681     sigma[46] = (Float_t) (fun[46]->GetParameter(2));
682   }
683   // ~~~~~~~~ PM Ref side A ~~~~~~~~
684   detector[47] = 4;
685   quad[47] = 5;
686   maxBin[47] = hPMRefAlg->GetMaximumBin();
687   nBin[47] = (hPMRefAlg->GetXaxis())->GetNbins();
688   xMax[47] = (hPMRefAlg->GetXaxis())->GetXmax();
689   if(nBin[47]!=0) maxXval[47] = maxBin[47]*xMax[47]/nBin[47];
690   if(maxXval[47]-100.<0.) xlow[47]=0.;
691   else xlow[47] = maxXval[47]-150.;
692   if(hPMRefAlg->GetEntries()!=0 && hPMRefAlg->GetMean()>0){
693     atLeastOneHisto=1;  
694     //
695     hPMRefAlg->Fit("gaus","Q","",xlow[47],maxXval[47]+100.);
696     fun[47] = hPMRefAlg->GetFunction("gaus");
697     mean[47]  = (Float_t) (fun[47]->GetParameter(1));
698     sigma[47] = (Float_t) (fun[47]->GetParameter(2));
699   }
700   
701   if(atLeastOneHisto==0){
702     printf("\n WARNING! Empty LASER histos -> ending DA WITHOUT writing output\n\n");
703     return -1;
704   }
705     
706   FILE *fileShuttle;
707   fileShuttle = fopen(LASDATA_FILE,"w");
708   for(Int_t i=0; i<2*kNChannels; i++){
709     fprintf(fileShuttle,"\t%d\t%d\t%f\t%f\n",detector[i],quad[i],mean[i], sigma[i]); 
710   }
711   //                                                   
712   fclose(fileShuttle);
713   /* report progress */
714   daqDA_progressReport(80);
715   //
716   TFile *histofile = new TFile(LASHISTO_FILE,"RECREATE");
717   histofile->cd();
718   for(int j=0; j<5; j++){
719      hZNChg[j]->Write();
720      hZPChg[j]->Write();
721      hZNAhg[j]->Write();
722      hZPAhg[j]->Write();
723      hZNClg[j]->Write();
724      hZPClg[j]->Write();
725      hZNAlg[j]->Write();
726      hZPAlg[j]->Write();  
727      if(j<2){
728        hZEMhg[j]->Write();
729        hZEMlg[j]->Write();
730     }
731   }
732   hPMRefChg->Write();
733   hPMRefAhg->Write();
734   hPMRefClg->Write();
735   hPMRefAlg->Write();  
736   //
737   histofile->Close();
738   //
739   for(Int_t j=0; j<5; j++){
740     delete hZNChg[j];
741     delete hZPChg[j];
742     delete hZNAhg[j];
743     delete hZPAhg[j];
744     delete hZNClg[j];
745     delete hZPClg[j];
746     delete hZNAlg[j];
747     delete hZPAlg[j];
748     if(j<2){
749       delete hZEMhg[j];
750       delete hZEMlg[j];
751     }
752   }
753   delete hPMRefChg;
754   delete hPMRefAhg;
755   delete hPMRefClg;
756   delete hPMRefAlg;
757
758   /* write report */
759   fprintf(fp,"Run #%s, received %d physics events out of %d\n",getenv("DATE_RUN_NUMBER"),nevents_physics,nevents_total);
760
761   /* close result file */
762   fclose(fp);
763   
764   /* report progress */
765   daqDA_progressReport(90);
766   
767   /* store the result file on FES */
768   // [1] File with mapping
769   status = daqDA_FES_storeFile(MAPDATA_FILE, "MAPPING");
770   if(status){
771     printf("Failed to export file : %d\n",status);
772     return -1;
773   }
774   //
775   // [2] File with laser data
776   status = daqDA_FES_storeFile(LASDATA_FILE, "LASERDATA");
777   if(status){
778     printf("Failed to export file : %d\n",status);
779     return -1;
780   }
781   // [3] File with laser histos
782   status = daqDA_FES_storeFile(LASHISTO_FILE, "LASERHISTOS");
783   if(status){
784     printf("Failed to export pedestal histos file to DAQ FES\n");
785     return -1;
786   }
787
788   /* report progress */
789   daqDA_progressReport(100);
790
791   return status;
792 }