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