1 /*********************************************************************************
2 - Contact: Brigitte Cheynis b.cheynis@ipnl.in2p3.fr
5 /afs/cern.ch/user/c/cheynis/public/run546.dat
8 - Number of events needed: >=100
9 - Input Files: argument list
10 - Output Files: local files V00Log.txt, VZERO_Histos.root, V0_Ped_Width_Gain.dat
11 FXS file V0_Ped_Width_Gain.dat
12 - Trigger types used: PHYSICS_EVENT
13 **********************************************************************************/
16 /**********************************************************************************
18 * VZERO Detector Algorithm used for extracting calibration parameters *
20 * This program reads the DDL data file passed as argument using the monitoring *
22 * It computes calibration parameters, populates local "./V0_Ped_Width_Gain.dat" *
23 * file and exports it to the FES. *
24 * We have 128 channels instead of 64 as expected for V0 due to the two sets of *
25 * charge integrators which are used by the FEE ... *
26 * The program reports about its processing progress. *
28 ***********************************************************************************/
36 #include <AliVZERORawStream.h>
37 #include <AliRawReaderDate.h>
38 #include <AliRawReader.h>
47 #include "TPluginManager.h"
53 /* Main routine --- Arguments: list of DATE raw data files */
55 int main(int argc, char **argv) {
59 printf(" argc = %d, argv = %s \n",argc, &(**argv));
61 Int_t kHighCut = 50; // high cut on pedestal distribution - to be tuned
62 Int_t kLowCut = 30; // low cut on signal distribution - to be tuned
63 Double_t ADCmean[128];
64 Double_t PEDmean[128];
65 Double_t PEDsigma[128];
67 //___________________________________________________
68 // Book HISTOGRAMS - dynamics of p-p collisions -
76 for (Int_t i=0; i<128; i++) {
77 sprintf(ADCname,"hADC%d",i);
78 sprintf(texte,"ADC cell%d",i);
79 hADCname[i] = new TH1F(ADCname,texte,1024,0,1023);
80 sprintf(PEDname,"hPED%d",i);
81 sprintf(texte,"PED cell%d",i);
82 hPEDname[i] = new TH1F(PEDname,texte,1024,0,1023);}
83 //___________________________________________________
86 /* log start of process */
87 printf("VZERO DA program started\n");
89 /* check that we got some arguments = list of files */
91 printf("Wrong number of arguments\n");
94 /* open result file to be exported to FES */
96 fp=fopen("./V0_Ped_Width_Gain.dat","a");
98 printf("Failed to open result file\n");
101 /* open log file to inform user */
103 flog=fopen("./V00log.txt","a");
105 printf("Failed to open log file\n");
108 /* report progress */
109 daqDA_progressReport(10);
112 /* init counters on events */
113 int nevents_physics=0;
116 /* read the data files, considering n files */
120 for (n=1;n<argc;n++) {
122 status=monitorSetDataSource( argv[n] );
124 printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
127 /* report progress - here indexed on the number of files */
128 daqDA_progressReport(10+80*n/argc);
130 /* read the data file */
132 struct eventHeaderStruct *event;
133 eventTypeType eventT;
136 status=monitorGetEventDynamic((void **)&event);
137 if (status==MON_ERR_EOF) break; /* end of monitoring file has been reached */
139 printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
142 /* retry if got no event */
143 if (event==NULL) break;
146 eventT=event->eventType;
148 switch (event->eventType){
154 printf("End Of Run detected\n");
160 fprintf(flog,"Run #%lu, event size: %lu, BC:%u, Orbit:%u, Period:%u\n",
161 (unsigned long)event->eventRunNb,
162 (unsigned long)event->eventSize,
163 EVENT_ID_GET_BUNCH_CROSSING(event->eventId),
164 EVENT_ID_GET_ORBIT(event->eventId),
165 EVENT_ID_GET_PERIOD(event->eventId) );
167 AliRawReader *rawReader = new AliRawReaderDate((void*)event);
169 AliVZERORawStream* rawStream = new AliVZERORawStream(rawReader);
171 for(Int_t i=0; i<64; i++) {
172 if(!rawStream->GetIntegratorFlag(i,10))
173 hADCname[i]->Fill(float(rawStream->GetADC(i))); // even integrator - fills 0 to 63
175 hADCname[i+64]->Fill(float(rawStream->GetADC(i))); // odd integrator - fills 64 to 123
176 for(Int_t j=0; j<21; j++) {
178 if(!rawStream->GetIntegratorFlag(i,j))
179 { hPEDname[i]->Fill(float(rawStream->GetPedestal(i,j))); } // even integrator
181 { hPEDname[i+64]->Fill(float(rawStream->GetPedestal(i,j))); } // odd integrator
188 } // end of switch on event type
194 } // loop over events
196 } // loop over data files
197 //________________________________________________________________________
198 // Computes mean values, dumps them into the output text file
200 for(Int_t i=0; i<128; i++) {
201 hPEDname[i]->GetXaxis()->SetRange(0,kHighCut);
202 PEDmean[i] = hPEDname[i]->GetMean();
203 PEDsigma[i] = hPEDname[i]->GetRMS();
204 hADCname[i]->GetXaxis()->SetRange(kLowCut,1023);
205 ADCmean[i] = hADCname[i]->GetMean() ;
206 fprintf(fp," %.3f %.3f %.3f\n",PEDmean[i],PEDsigma[i],ADCmean[i]);
209 //________________________________________________________________________
210 // Write root file with histos for users further check - just in case -
212 TFile *histoFile = new TFile("VZERO_histos.root","RECREATE");
214 for (Int_t i=0; i<128; i++) {
215 hADCname[i]->GetXaxis()->SetRange(0,1023);
216 hADCname[i]->Write();
217 hPEDname[i]->Write(); }
222 //________________________________________________________________________
225 fprintf(flog,"Run #%s, received %d physics events out of %d\n",getenv("DATE_RUN_NUMBER"),nevents_physics,nevents_total);
226 printf("Run #%s, received %d physics events out of %d\n",getenv("DATE_RUN_NUMBER"),nevents_physics,nevents_total);
228 /* close result and log files */
232 /* report progress */
233 daqDA_progressReport(90);
236 /* export result file to FES */
237 status=daqDA_FES_storeFile("./V0_Ped_Width_Gain.dat","V00da_results");
239 printf("Failed to export file : %d\n",status);
242 /* report progress */
243 daqDA_progressReport(100);