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