1 /*********************************************************************************
2 - Contact: Brigitte Cheynis b.cheynis@ipnl.in2p3.fr
3 - Link: /afs/cern.ch/user/c/cheynis/public/run546.dat
6 - Number of events needed: >=100
7 - Input Files: raw_data_file_on_LDC
8 - Output Files: local files V00Log.txt, VZERO_Histos.root, V0_Ped_Width_Gain.dat
9 FXS file V0_Ped_Width_Gain.dat
11 **********************************************************************************/
14 /**********************************************************************************
16 * VZERO Detector Algorithm used for extracting calibration parameters *
18 * This program reads the DDL data file passed as argument using the monitoring *
20 * It computes calibration parameters, populates local "./V0_Ped_Width_Gain.dat" *
21 * file and exports it to the FES. *
22 * We have 128 channels instead of 64 as expected for V0 due to the two sets of *
23 * charge integrators which are used by the FEE ... *
24 * The program reports about its processing progress. *
26 ***********************************************************************************/
34 #include <AliVZERORawStream.h>
35 #include <AliRawReaderDate.h>
36 #include <AliRawReader.h>
45 #include "TPluginManager.h"
51 /* Main routine --- Arguments: list of DATE raw data files */
53 int main(int argc, char **argv) {
57 printf(" argc = %d, argv = %s \n",argc, &(**argv));
59 Int_t kHighCut = 50; // high cut on pedestal distribution - to be tuned
60 Int_t kLowCut = 30; // low cut on signal distribution - to be tuned
61 Double_t ADCmean[128];
62 Double_t PEDmean[128];
63 Double_t PEDsigma[128];
65 //___________________________________________________
66 // Book HISTOGRAMS - dynamics of p-p collisions -
74 for (Int_t i=0; i<128; i++) {
75 sprintf(ADCname,"hADC%d",i);
76 sprintf(texte,"ADC cell%d",i);
77 hADCname[i] = new TH1F(ADCname,texte,1024,0,1023);
78 sprintf(PEDname,"hPED%d",i);
79 sprintf(texte,"PED cell%d",i);
80 hPEDname[i] = new TH1F(PEDname,texte,1024,0,1023);}
81 //___________________________________________________
84 /* log start of process */
85 printf("VZERO DA program started\n");
87 /* check that we got some arguments = list of files */
89 printf("Wrong number of arguments\n");
92 /* open result file to be exported to FES */
94 fp=fopen("./V0_Ped_Width_Gain.dat","a");
96 printf("Failed to open result file\n");
99 /* open log file to inform user */
101 flog=fopen("./V00log.txt","a");
103 printf("Failed to open log file\n");
106 /* report progress */
107 daqDA_progressReport(10);
110 /* init counters on events */
111 int nevents_physics=0;
114 /* read the data files, considering n files */
118 for (n=1;n<argc;n++) {
120 status=monitorSetDataSource( argv[n] );
122 printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
125 /* report progress - here indexed on the number of files */
126 daqDA_progressReport(10+80*n/argc);
128 /* read the data file */
130 struct eventHeaderStruct *event;
131 eventTypeType eventT;
134 status=monitorGetEventDynamic((void **)&event);
135 if (status==MON_ERR_EOF) break; /* end of monitoring file has been reached */
137 printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
140 /* retry if got no event */
141 if (event==NULL) break;
144 eventT=event->eventType;
146 switch (event->eventType){
152 printf("End Of Run detected\n");
158 fprintf(flog,"Run #%lu, event size: %lu, BC:%u, Orbit:%u, Period:%u\n",
159 (unsigned long)event->eventRunNb,
160 (unsigned long)event->eventSize,
161 EVENT_ID_GET_BUNCH_CROSSING(event->eventId),
162 EVENT_ID_GET_ORBIT(event->eventId),
163 EVENT_ID_GET_PERIOD(event->eventId) );
165 AliRawReader *rawReader = new AliRawReaderDate((void*)event);
167 AliVZERORawStream* rawStream = new AliVZERORawStream(rawReader);
169 for(Int_t i=0; i<64; i++) {
170 if(!rawStream->GetIntegratorFlag(i,10))
171 hADCname[i]->Fill(float(rawStream->GetADC(i))); // even integrator - fills 0 to 63
173 hADCname[i+64]->Fill(float(rawStream->GetADC(i))); // odd integrator - fills 64 to 123
174 for(Int_t j=0; j<21; j++) {
176 if(!rawStream->GetIntegratorFlag(i,j))
177 { hPEDname[i]->Fill(float(rawStream->GetPedestal(i,j))); } // even integrator
179 { hPEDname[i+64]->Fill(float(rawStream->GetPedestal(i,j))); } // odd integrator
186 } // end of switch on event type
192 } // loop over events
194 } // loop over data files
195 //________________________________________________________________________
196 // Computes mean values, dumps them into the output text file
198 for(Int_t i=0; i<128; i++) {
199 hPEDname[i]->GetXaxis()->SetRange(0,kHighCut);
200 PEDmean[i] = hPEDname[i]->GetMean();
201 PEDsigma[i] = hPEDname[i]->GetRMS();
202 hADCname[i]->GetXaxis()->SetRange(kLowCut,1023);
203 ADCmean[i] = hADCname[i]->GetMean() ;
204 fprintf(fp," %.3f %.3f %.3f\n",PEDmean[i],PEDsigma[i],ADCmean[i]);
207 //________________________________________________________________________
208 // Write root file with histos for users further check - just in case -
210 TFile *histoFile = new TFile("VZERO_histos.root","RECREATE");
212 for (Int_t i=0; i<128; i++) {
213 hADCname[i]->GetXaxis()->SetRange(0,1023);
214 hADCname[i]->Write();
215 hPEDname[i]->Write(); }
220 //________________________________________________________________________
223 fprintf(flog,"Run #%s, received %d physics events out of %d\n",getenv("DATE_RUN_NUMBER"),nevents_physics,nevents_total);
224 printf("Run #%s, received %d physics events out of %d\n",getenv("DATE_RUN_NUMBER"),nevents_physics,nevents_total);
226 /* close result and log files */
230 /* report progress */
231 daqDA_progressReport(90);
234 /* export result file to FES */
235 status=daqDA_FES_storeFile("./V0_Ped_Width_Gain.dat","V00da_results");
237 printf("Failed to export file : %d\n",status);
240 /* report progress */
241 daqDA_progressReport(100);