]> git.uio.no Git - u/mrichter/AliRoot.git/blame - VZERO/VZEROda.cxx
TPCNoiseMapComponent included into build (Kelly)
[u/mrichter/AliRoot.git] / VZERO / VZEROda.cxx
CommitLineData
bb80ae5b 1/*********************************************************************************
2- Contact: Brigitte Cheynis b.cheynis@ipnl.in2p3.fr
3- Link: /afs/cern.ch/user/c/cheynis/public/run546.dat
4- Run Type: PHYSICS
5- DA Type: LDC
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
10- Trigger types used:
11**********************************************************************************/
12
13
14/**********************************************************************************
15* *
16* VZERO Detector Algorithm used for extracting calibration parameters *
17* *
18* This program reads the DDL data file passed as argument using the monitoring *
19* library. *
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. *
25* *
26***********************************************************************************/
f52220c5 27
28// DATE
29#include "event.h"
30#include "monitor.h"
31#include "daqDA.h"
32
33//AliRoot
34#include <AliVZERORawStream.h>
35#include <AliRawReaderDate.h>
36#include <AliRawReader.h>
37#include <AliDAQ.h>
38
39// standard
40#include <stdio.h>
41#include <stdlib.h>
42
43//ROOT
bb80ae5b 44#include "TROOT.h"
45#include "TPluginManager.h"
f52220c5 46#include <TFile.h>
47#include <TH1F.h>
48#include <TMath.h>
49
50
51/* Main routine --- Arguments: list of DATE raw data files */
52
53int main(int argc, char **argv) {
54
55 int status;
56
57 printf(" argc = %d, argv = %s \n",argc, &(**argv));
58
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
df6d42e7 61 Double_t ADCmean[128];
62 Double_t PEDmean[128];
63 Double_t PEDsigma[128];
f52220c5 64
65//___________________________________________________
66// Book HISTOGRAMS - dynamics of p-p collisions -
67
68 char ADCname[6];
69 char PEDname[6];
df6d42e7 70 TH1F *hADCname[128];
71 TH1F *hPEDname[128];
f52220c5 72
73 char texte[12];
df6d42e7 74 for (Int_t i=0; i<128; i++) {
f52220c5 75 sprintf(ADCname,"hADC%d",i);
76 sprintf(texte,"ADC cell%d",i);
df6d42e7 77 hADCname[i] = new TH1F(ADCname,texte,1024,0,1023);
f52220c5 78 sprintf(PEDname,"hPED%d",i);
79 sprintf(texte,"PED cell%d",i);
df6d42e7 80 hPEDname[i] = new TH1F(PEDname,texte,1024,0,1023);}
f52220c5 81//___________________________________________________
82
83
84 /* log start of process */
85 printf("VZERO DA program started\n");
86
87 /* check that we got some arguments = list of files */
88 if (argc<2) {
89 printf("Wrong number of arguments\n");
90 return -1;}
91
92 /* open result file to be exported to FES */
93 FILE *fp=NULL;
94 fp=fopen("./V0_Ped_Width_Gain.dat","a");
95 if (fp==NULL) {
96 printf("Failed to open result file\n");
97 return -1;}
98
99 /* open log file to inform user */
100 FILE *flog=NULL;
101 flog=fopen("./V00log.txt","a");
102 if (flog==NULL) {
103 printf("Failed to open log file\n");
104 return -1; }
105
106 /* report progress */
107 daqDA_progressReport(10);
108
109
110 /* init counters on events */
111 int nevents_physics=0;
112 int nevents_total=0;
113
114 /* read the data files, considering n files */
115
116 int n;
117
118 for (n=1;n<argc;n++) {
119
120 status=monitorSetDataSource( argv[n] );
121 if (status!=0) {
122 printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
123 return -1; }
124
125 /* report progress - here indexed on the number of files */
126 daqDA_progressReport(10+80*n/argc);
127
128 /* read the data file */
129 for(;;) {
130 struct eventHeaderStruct *event;
131 eventTypeType eventT;
132
133 /* get next event */
134 status=monitorGetEventDynamic((void **)&event);
135 if (status==MON_ERR_EOF) break; /* end of monitoring file has been reached */
136 if (status!=0) {
137 printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
138 return -1; }
139
140 /* retry if got no event */
141 if (event==NULL) break;
142
143 /* decode event */
144 eventT=event->eventType;
145
146 switch (event->eventType){
147
148 case START_OF_RUN:
149 break;
150
151 case END_OF_RUN:
152 printf("End Of Run detected\n");
153 break;
154
155 case PHYSICS_EVENT:
156 nevents_physics++;
df6d42e7 157
f52220c5 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) );
164
165 AliRawReader *rawReader = new AliRawReaderDate((void*)event);
df6d42e7 166
167 AliVZERORawStream* rawStream = new AliVZERORawStream(rawReader);
f52220c5 168 rawStream->Next();
169 for(Int_t i=0; i<64; i++) {
df6d42e7 170 if(!rawStream->GetIntegratorFlag(i,10))
171 hADCname[i]->Fill(float(rawStream->GetADC(i))); // even integrator - fills 0 to 63
172 else
173 hADCname[i+64]->Fill(float(rawStream->GetADC(i))); // odd integrator - fills 64 to 123
f52220c5 174 for(Int_t j=0; j<21; j++) {
175 if(j==10) continue;
df6d42e7 176 if(!rawStream->GetIntegratorFlag(i,j))
177 { hPEDname[i]->Fill(float(rawStream->GetPedestal(i,j))); } // even integrator
178 else
179 { hPEDname[i+64]->Fill(float(rawStream->GetPedestal(i,j))); } // odd integrator
180 }
f52220c5 181 }
182 delete rawStream;
183 rawStream = 0x0;
184 delete rawReader;
185 rawReader = 0x0;
186 } // end of switch on event type
187
188 nevents_total++;
189 /* free resources */
190 free(event);
191
192 } // loop over events
193
194 } // loop over data files
195//________________________________________________________________________
196// Computes mean values, dumps them into the output text file
197
df6d42e7 198 for(Int_t i=0; i<128; i++) {
f52220c5 199 hPEDname[i]->GetXaxis()->SetRange(0,kHighCut);
200 PEDmean[i] = hPEDname[i]->GetMean();
201 PEDsigma[i] = hPEDname[i]->GetRMS();
df6d42e7 202 hADCname[i]->GetXaxis()->SetRange(kLowCut,1023);
f52220c5 203 ADCmean[i] = hADCname[i]->GetMean() ;
204 fprintf(fp," %.3f %.3f %.3f\n",PEDmean[i],PEDsigma[i],ADCmean[i]);
205 }
206
207//________________________________________________________________________
208// Write root file with histos for users further check - just in case -
209
210 TFile *histoFile = new TFile("VZERO_histos.root","RECREATE");
211
df6d42e7 212 for (Int_t i=0; i<128; i++) {
213 hADCname[i]->GetXaxis()->SetRange(0,1023);
f52220c5 214 hADCname[i]->Write();
215 hPEDname[i]->Write(); }
216
217 histoFile->Close();
df6d42e7 218 delete histoFile;
f52220c5 219
220//________________________________________________________________________
221
222 /* write report */
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);
225
226 /* close result and log files */
227 fclose(fp);
228 fclose(flog);
229
230 /* report progress */
231 daqDA_progressReport(90);
232
233
234 /* export result file to FES */
235 status=daqDA_FES_storeFile("./V0_Ped_Width_Gain.dat","V00da_results");
236 if (status) {
237 printf("Failed to export file : %d\n",status);
238 return -1; }
239
240 /* report progress */
241 daqDA_progressReport(100);
242
243 return status;
244}