]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HMPID/HMPIDda.cxx
b01727e6a2ed5e8e76a2f33dbae061ded0fd146b
[u/mrichter/AliRoot.git] / HMPID / HMPIDda.cxx
1 /*
2
3 HMPID DA for online calibration
4
5 Contact: Levente.Molnar@cern.ch, Giacomo.Volpe@ba.infn.it
6 Link for validation file from AliEn: /alice/data/2008/LHC08c_HMPID/000047350/raw/08000047350020.10.root 
7 Run Type: PEDESTAL -- but we select on the PHYSICS_EVENTS in th HMPIDda.cxx
8 DA Type: LDC
9 Number of events needed: ~ 1000 events
10 Input Files: Raw pedestal file, EXTERNAL config files: HmpidSigmaCut.txt, HmpDeadChannelMap.txt on all HMPID LDCs
11 Output Files: 2 x 14 txt files including pedestal and error values
12 Trigger types used: PEDESTAL RUN (but selecting on PHYSICS_EVENT)
13
14 */
15
16 extern "C" {
17 #include <daqDA.h>
18 }
19
20 #include "event.h"
21 #include "monitor.h"
22
23 #include <Riostream.h>
24 #include "fstream.h"
25 #include <stdio.h>
26 #include <stdlib.h>
27
28 //AliRoot
29 #include "AliHMPIDRawStream.h"
30 #include "AliHMPIDCalib.h"
31 #include "AliRawReaderDate.h"
32 #include "AliBitPacking.h"
33
34 //ROOT
35 #include "TROOT.h"
36 #include "TFile.h"
37 #include "TSystem.h"
38 #include "TString.h"
39 #include "TObject.h"
40 #include "TPluginManager.h"
41
42 //AMORE
43 #include <AmoreDA.h>
44
45
46 int main(int argc, char **argv){ 
47
48   int status;
49   const Char_t         *hmpConfigFile        = "HmpDaqDaConfig.txt"; 
50   const Char_t         *hmpDeadChannelMapFile = "HmpDeadChannelMap.txt"; 
51   const Int_t               ddlOffset = 1536;
52         TString                 hmpIn,hmpIn2;
53         TString                 feeIn;
54
55   
56   /* log start of process */
57   printf("HMPID DA program started\n");  
58
59   gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo","*","TStreamerInfo","RIO","TStreamerInfo()");
60   
61   /* check that we got some arguments = list of files */
62   if (argc<2) {
63     printf("Wrong number of arguments\n");
64     return -1;
65   }
66
67   /* copy locally a file from daq detector config db */
68   
69   hmpIn=Form("./%s",hmpConfigFile);
70   status+=daqDA_DB_getFile(hmpConfigFile,hmpIn.Data());
71   if (status) { printf("Failed to get HMPID config file status: %d\n",status);return -1; }
72   hmpIn2=Form("./%s",hmpDeadChannelMapFile);
73   status+=daqDA_DB_getFile(hmpDeadChannelMapFile,hmpIn.Data());
74   if (status) { printf("Failed to get HMPID dead channel file status: %d\n",status);return -1; }
75   
76
77   /* report progress */
78   daqDA_progressReport(10);
79  
80   /* define wait event timeout - 1s max */
81   monitorSetNowait();
82   monitorSetNoWaitNetworkTimeout(1000);
83
84   /* set local storage of ped files for Fe2C */
85
86   
87   /* init the pedestal calculation */
88   AliHMPIDCalib *pCal=new AliHMPIDCalib();
89   /* Set the number of sigma cuts inside the file HmpidSigmaCut.txt on all LDCs! */
90   /* If the file is NOT present then the default cut 3 will be used!*/
91   pCal->SetSigCutFromFile(hmpIn);
92   pCal->SetDeadChannelMapFromFile(hmpIn2);  
93   
94   /* init event counter */
95   Int_t firstEvt=0;
96   Int_t iEvtNcal=firstEvt;                      //Start from 1 not 0!                                                 
97   ULong_t runNum=0;
98   ULong_t ldcId=0;
99
100   int n;
101   for (n=1;n<argc;n++) {
102
103     status=monitorSetDataSource( argv[n] );
104     if (status!=0) {
105       printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
106       return -1;
107     }
108
109     /* report progress */
110     /* in this example, indexed on the number of files */
111     daqDA_progressReport(10+80*n/argc);
112
113     for(;;) { // infinite loop 
114       
115        /* check shutdown condition */
116     if (daqDA_checkShutdown()) {break;}
117     
118       struct eventHeaderStruct *event;
119       eventTypeType eventT;
120
121       /* get next event */
122       status=monitorGetEventDynamic((void **)&event);
123       if (status==MON_ERR_EOF)                                              /* end of monitoring file has been reached */
124       {
125         printf("End of monitoring file has been reached! \n");
126         break;
127         }
128       if (status!=0) {
129         printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
130         return -1;
131       }
132
133       /* retry if got no event */
134       if (event==NULL) {
135         //break;
136         continue;
137       }
138
139       /* use event - here, just write event id to result file */
140       eventT=event->eventType;
141
142       if (eventT==PHYSICS_EVENT || eventT==CALIBRATION_EVENT  ) {                 //updated: 18/02/2008 based on http://alice-ecs.web.cern.ch/alice-ecs/runtypes_3.16.html
143         runNum=(unsigned long)event->eventRunNb;                                  //assuming that only one run is processed at a time
144         ldcId=(unsigned long)event->eventLdcId;
145         
146         iEvtNcal++;        
147         AliRawReader *reader = new AliRawReaderDate((void*)event);
148         AliHMPIDRawStream stream(reader);stream.SetTurbo(kTRUE);                  //raw data decoding without error checks SetTurbo(kTRUE)
149         while(stream.Next())
150           {
151              for(Int_t iPad=0;iPad<stream.GetNPads();iPad++) {
152              pCal->FillPedestal(stream.GetPadArray()[iPad],stream.GetChargeArray()[iPad]);
153               } //pads
154           }//while -- loop on det load in one event
155           
156          for(Int_t iddl=0;iddl<stream.GetNDDL();iddl++){   
157                  pCal->FillDDLCnt(iddl,stream.GetnDDLInStream()[iddl],stream.GetnDDLOutStream()[iddl]);                     
158            for(Int_t ierr=0; ierr < stream.GetNErrors(); ierr++) {
159                pCal->FillErrors(iddl,ierr,stream.GetErrors(iddl,ierr));
160                }
161           }//err   
162           
163         pCal->SetRunParams(runNum,stream.GetTimeStamp(),stream.GetLDCNumber());   //Get the last TimeStam read and the LDC ID
164         stream.Delete();            
165       }// if CALIBRATION_EVENT
166
167       /* exit when last event received, no need to wait for TERM signal */
168       if (eventT==END_OF_RUN) {
169         printf("EOR event detected\n");
170         break;    
171     
172       } // events loop   
173
174       free(event);
175     }
176
177   }//arg
178
179   /* write report */
180   printf("HMPID DA processed RUN #%s on LDC#%d, with %d calibration events\n",getenv("DATE_RUN_NUMBER"),ldcId,iEvtNcal);
181
182   if (!iEvtNcal) {
183     printf("No calibration events have been read. Exiting\n");
184     return -1;
185   }
186
187   /* report progress */
188   daqDA_progressReport(90);
189   /* send pedestal and error files to DAQ FXS */
190   for(Int_t nDDL=0; nDDL < AliHMPIDRawStream::kNDDL; nDDL++) {   
191     feeIn=Form("thr%d.dat",ddlOffset+nDDL);         
192     /* Calculate pedestal for the given ddl, if there is no ddl go t next */
193     if(pCal->CalcPedestal(nDDL,Form("HmpidPedDdl%02i.txt",nDDL),Form("%s",feeIn.Data()),iEvtNcal)) {
194       status=daqDA_DB_storeFile(feeIn.Data(),feeIn.Data());                                               //store a single threshold file for a DDL in DAQ DB  
195       if (status) { printf("Failed to store file %s in DAQ DB, status: %d\n",feeIn.Data(),status); }
196       status=daqDA_FES_storeFile(Form("HmpidPedDdl%02i.txt",nDDL),Form("HmpidPedDdl%02i.txt",nDDL));      //store a single pedestal file for a DDL
197       if (status) { printf("Failed to export file on FES: %d\n",status); }
198     }//pedestals and thresholds
199     if(pCal->WriteErrors(nDDL,Form("HmpidErrorsDdl%02i.txt",nDDL),iEvtNcal)) {
200       status=daqDA_FES_storeFile(Form("HmpidErrorsDdl%02i.txt",nDDL),Form("HmpidErrorsDdl%02i.txt",nDDL));
201       if (status) { printf("Failed to export file : %d\n",status); }
202     }//errors
203   }//nDDL
204
205   /* send files to AMORE DB */
206   daqDA_progressReport(95);
207   Int_t statusAmoreDA=0; 
208   amore::da::AmoreDA amoreDA(amore::da::AmoreDA::kSender);
209   for(Int_t iCh=0; iCh < AliHMPIDParam::kMaxCh; iCh++) {  
210   statusAmoreDA+=amoreDA.Send(Form("fPedMeanMap%d",iCh), pCal->GetPedMeanMap(iCh));  
211   statusAmoreDA+=amoreDA.Send(Form("fPedSigMap%d",iCh),  pCal->GetPedSigMap(iCh));  
212   }
213   for(Int_t iCh=0;iCh<=AliHMPIDParam::kMaxCh;iCh++)
214    {
215     for(Int_t iFee=0;iFee<6;iFee++)
216      {
217          statusAmoreDA+=amoreDA.Send(Form("f1DPedMean_Ch%d_FEE_%d",iCh,iFee),pCal->GetPedMean(6*iCh+iFee));
218          statusAmoreDA+=amoreDA.Send(Form("f1DPedSigma_Ch%d_FEE_%d",iCh,iFee),pCal->GetPedSigma(6*iCh+iFee));
219        }
220      }
221          
222   delete pCal;  
223   if(status) return -1;
224   if(statusAmoreDA) return -1;
225   
226   /* report progress */
227   daqDA_progressReport(100);
228   
229   
230   return status;
231 }
232