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