1 /********************************************************************************
3 * VZERO Detector Algorithm used for extracting calibration parameters *
5 * This program reads the DDL data file passed as argument using the monitoring *
7 * It computes calibration parameters, populates local "./V0_Ped_Width_Gain.dat" *
8 * file and exports it to the FES. *
9 * We have 128 channels instead of 64 as expected for V0 due to the two sets of *
10 * charge integrators which are used by the FEE ... *
11 * The program reports about its processing progress. *
13 *********************************************************************************/
21 #include <AliVZERORawStream.h>
22 #include <AliRawReaderDate.h>
23 #include <AliRawReader.h>
37 /* Main routine --- Arguments: list of DATE raw data files */
39 int main(int argc, char **argv) {
43 printf(" argc = %d, argv = %s \n",argc, &(**argv));
45 Int_t kHighCut = 50; // high cut on pedestal distribution - to be tuned
46 Int_t kLowCut = 30; // low cut on signal distribution - to be tuned
47 Double_t ADCmean[128];
48 Double_t PEDmean[128];
49 Double_t PEDsigma[128];
51 //___________________________________________________
52 // Book HISTOGRAMS - dynamics of p-p collisions -
60 for (Int_t i=0; i<128; i++) {
61 sprintf(ADCname,"hADC%d",i);
62 sprintf(texte,"ADC cell%d",i);
63 hADCname[i] = new TH1F(ADCname,texte,1024,0,1023);
64 sprintf(PEDname,"hPED%d",i);
65 sprintf(texte,"PED cell%d",i);
66 hPEDname[i] = new TH1F(PEDname,texte,1024,0,1023);}
67 //___________________________________________________
70 /* log start of process */
71 printf("VZERO DA program started\n");
73 /* check that we got some arguments = list of files */
75 printf("Wrong number of arguments\n");
78 /* open result file to be exported to FES */
80 fp=fopen("./V0_Ped_Width_Gain.dat","a");
82 printf("Failed to open result file\n");
85 /* open log file to inform user */
87 flog=fopen("./V00log.txt","a");
89 printf("Failed to open log file\n");
93 daqDA_progressReport(10);
96 /* init counters on events */
97 int nevents_physics=0;
100 /* read the data files, considering n files */
104 for (n=1;n<argc;n++) {
106 status=monitorSetDataSource( argv[n] );
108 printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
111 /* report progress - here indexed on the number of files */
112 daqDA_progressReport(10+80*n/argc);
114 /* read the data file */
116 struct eventHeaderStruct *event;
117 eventTypeType eventT;
120 status=monitorGetEventDynamic((void **)&event);
121 if (status==MON_ERR_EOF) break; /* end of monitoring file has been reached */
123 printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
126 /* retry if got no event */
127 if (event==NULL) break;
130 eventT=event->eventType;
132 switch (event->eventType){
138 printf("End Of Run detected\n");
144 fprintf(flog,"Run #%lu, event size: %lu, BC:%u, Orbit:%u, Period:%u\n",
145 (unsigned long)event->eventRunNb,
146 (unsigned long)event->eventSize,
147 EVENT_ID_GET_BUNCH_CROSSING(event->eventId),
148 EVENT_ID_GET_ORBIT(event->eventId),
149 EVENT_ID_GET_PERIOD(event->eventId) );
151 AliRawReader *rawReader = new AliRawReaderDate((void*)event);
153 AliVZERORawStream* rawStream = new AliVZERORawStream(rawReader);
155 for(Int_t i=0; i<64; i++) {
156 if(!rawStream->GetIntegratorFlag(i,10))
157 hADCname[i]->Fill(float(rawStream->GetADC(i))); // even integrator - fills 0 to 63
159 hADCname[i+64]->Fill(float(rawStream->GetADC(i))); // odd integrator - fills 64 to 123
160 for(Int_t j=0; j<21; j++) {
162 if(!rawStream->GetIntegratorFlag(i,j))
163 { hPEDname[i]->Fill(float(rawStream->GetPedestal(i,j))); } // even integrator
165 { hPEDname[i+64]->Fill(float(rawStream->GetPedestal(i,j))); } // odd integrator
172 } // end of switch on event type
178 } // loop over events
180 } // loop over data files
181 //________________________________________________________________________
182 // Computes mean values, dumps them into the output text file
184 for(Int_t i=0; i<128; i++) {
185 hPEDname[i]->GetXaxis()->SetRange(0,kHighCut);
186 PEDmean[i] = hPEDname[i]->GetMean();
187 PEDsigma[i] = hPEDname[i]->GetRMS();
188 hADCname[i]->GetXaxis()->SetRange(kLowCut,1023);
189 ADCmean[i] = hADCname[i]->GetMean() ;
190 fprintf(fp," %.3f %.3f %.3f\n",PEDmean[i],PEDsigma[i],ADCmean[i]);
193 //________________________________________________________________________
194 // Write root file with histos for users further check - just in case -
196 TFile *histoFile = new TFile("VZERO_histos.root","RECREATE");
198 for (Int_t i=0; i<128; i++) {
199 hADCname[i]->GetXaxis()->SetRange(0,1023);
200 hADCname[i]->Write();
201 hPEDname[i]->Write(); }
206 //________________________________________________________________________
209 fprintf(flog,"Run #%s, received %d physics events out of %d\n",getenv("DATE_RUN_NUMBER"),nevents_physics,nevents_total);
210 printf("Run #%s, received %d physics events out of %d\n",getenv("DATE_RUN_NUMBER"),nevents_physics,nevents_total);
212 /* close result and log files */
216 /* report progress */
217 daqDA_progressReport(90);
220 /* export result file to FES */
221 status=daqDA_FES_storeFile("./V0_Ped_Width_Gain.dat","V00da_results");
223 printf("Failed to export file : %d\n",status);
226 /* report progress */
227 daqDA_progressReport(100);