]> git.uio.no Git - u/mrichter/AliRoot.git/blame - VZERO/VZEROda.cxx
Reset rawReader in MakeRaws (Yves)
[u/mrichter/AliRoot.git] / VZERO / VZEROda.cxx
CommitLineData
bb80ae5b 1/*********************************************************************************
2- Contact: Brigitte Cheynis b.cheynis@ipnl.in2p3.fr
f30630cb 3- Link: http
03e924d7 4- Raw data test file :
5- Reference run number :
bb80ae5b 6- Run Type: PHYSICS
03e924d7 7- DA Type: MON
8- Number of events needed: >=500
f30630cb 9- Input Files: argument list
03e924d7 10- Output Files: local files VZERO_Histos.root, V0_Ped_Width_Gain.dat
bb80ae5b 11 FXS file V0_Ped_Width_Gain.dat
f30630cb 12- Trigger types used: PHYSICS_EVENT
bb80ae5b 13**********************************************************************************/
14
bb80ae5b 15/**********************************************************************************
16* *
17* VZERO Detector Algorithm used for extracting calibration parameters *
18* *
03e924d7 19* This program connects to the DAQ data source passed as argument. *
bb80ae5b 20* It computes calibration parameters, populates local "./V0_Ped_Width_Gain.dat" *
d154146a 21* file, exports it to the FES, and stores it into DAQ DB *
03e924d7 22* The program exits when being asked to shut down (daqDA_checkshutdown) *
c2d54b4c 23* or on End of Run event. *
bb80ae5b 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 ... *
bb80ae5b 26* *
27***********************************************************************************/
f52220c5 28
29// DATE
30#include "event.h"
31#include "monitor.h"
32#include "daqDA.h"
33
34//AliRoot
35#include <AliVZERORawStream.h>
36#include <AliRawReaderDate.h>
37#include <AliRawReader.h>
38#include <AliDAQ.h>
39
40// standard
41#include <stdio.h>
42#include <stdlib.h>
43
44//ROOT
bb80ae5b 45#include "TROOT.h"
46#include "TPluginManager.h"
f52220c5 47#include <TFile.h>
48#include <TH1F.h>
49#include <TMath.h>
50
03e924d7 51/* Main routine --- Arguments: monitoring data source */
f52220c5 52
53int main(int argc, char **argv) {
54
6ee2e288 55/* magic line from Cvetan */
56 gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
57 "*",
58 "TStreamerInfo",
59 "RIO",
60 "TStreamerInfo()");
f52220c5 61 int status;
03e924d7 62 if (argc!=2) {
63 printf("Wrong number of arguments\n");
64 return -1;
65 }
c2d54b4c 66
df6d42e7 67 Double_t ADCmean[128];
03e924d7 68 Double_t ADCsigma[128];
df6d42e7 69 Double_t PEDmean[128];
70 Double_t PEDsigma[128];
c2d54b4c 71
72//___________________________________________________
73// Get cuts from V00DA.config file
74
5a8e0a88 75 Int_t kLowCut; // = 60; low cut on signal distribution - to be tuned
76 Int_t kHighCut; // = 50; high cut on pedestal distribution - to be tuned
77
c2d54b4c 78 status = daqDA_DB_getFile("V00DA.config","./V00DA.config");
79 if (status) {
d154146a 80 printf("Failed to get config file (V00DA.config) from DAQ DB, status=%d\n", status);
c2d54b4c 81 return -1;
82 }
83 /* open the config file and retrieve cuts */
84 FILE *fpConfig = fopen("V00DA.config","r");
85 fscanf(fpConfig,"%d %d",&kLowCut,&kHighCut);
86 fclose(fpConfig);
f52220c5 87
c2d54b4c 88 printf("LowCut on signal = %d ; HighCut on pedestal = %d\n",kLowCut,kHighCut);
89
f52220c5 90//___________________________________________________
91// Book HISTOGRAMS - dynamics of p-p collisions -
92
93 char ADCname[6];
94 char PEDname[6];
df6d42e7 95 TH1F *hADCname[128];
96 TH1F *hPEDname[128];
f52220c5 97
98 char texte[12];
df6d42e7 99 for (Int_t i=0; i<128; i++) {
f52220c5 100 sprintf(ADCname,"hADC%d",i);
101 sprintf(texte,"ADC cell%d",i);
b68c91e6 102 hADCname[i] = new TH1F(ADCname,texte,1024,-0.5, 1023.5);
f52220c5 103 sprintf(PEDname,"hPED%d",i);
104 sprintf(texte,"PED cell%d",i);
b68c91e6 105 hPEDname[i] = new TH1F(PEDname,texte,1024,-0.5, 1023.5);
03e924d7 106 }
f52220c5 107//___________________________________________________
d154146a 108
f52220c5 109 /* open result file to be exported to FES */
110 FILE *fp=NULL;
b68c91e6 111 fp=fopen("./V0_Ped_Width_Gain.dat","w");
f52220c5 112 if (fp==NULL) {
113 printf("Failed to open result file\n");
114 return -1;}
115
03e924d7 116 /* define data source : this is argument 1 */
117 status=monitorSetDataSource( argv[1] );
118 if (status!=0) {
119 printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
120 return -1;
121 }
122
123 /* declare monitoring program */
124 status=monitorDeclareMp( __FILE__ );
125 if (status!=0) {
126 printf("monitorDeclareMp() failed : %s\n",monitorDecodeError(status));
127 return -1;
128 }
129
130 /* define wait event timeout - 1s max */
131 monitorSetNowait();
132 monitorSetNoWaitNetworkTimeout(1000);
133
f52220c5 134 /* init counters on events */
135 int nevents_physics=0;
136 int nevents_total=0;
137
03e924d7 138 /* loop on events (infinite) */
139 for(;;) {
140 struct eventHeaderStruct *event;
141 eventTypeType eventT;
142
143 /* check shutdown condition */
144 if (daqDA_checkShutdown()) {break;}
145
146 /* get next event (blocking call until timeout) */
147 status=monitorGetEventDynamic((void **)&event);
148 if (status==MON_ERR_EOF) {
149 printf ("End of File detected\n");
150 break; /* end of monitoring file has been reached */
151 }
152
153 if (status!=0) {
154 printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
155 break;
156 }
157
158 /* retry if got no event */
159 if (event==NULL) continue;
160
161 /* decode event */
162 eventT=event->eventType;
f52220c5 163
03e924d7 164 switch (event->eventType){
f52220c5 165
03e924d7 166 case START_OF_RUN:
167 break;
f52220c5 168
03e924d7 169 case END_OF_RUN:
170 printf("End Of Run detected\n");
171 break;
f52220c5 172
03e924d7 173 case PHYSICS_EVENT:
174 nevents_physics++;
c2d54b4c 175
03e924d7 176 AliRawReader *rawReader = new AliRawReaderDate((void*)event);
df6d42e7 177
03e924d7 178 AliVZERORawStream* rawStream = new AliVZERORawStream(rawReader);
179 rawStream->Next();
180 for(Int_t i=0; i<64; i++) {
5a8e0a88 181 Int_t nFlag = 0;
182 for(Int_t j=0; j<21; j++) {
183 if((rawStream->GetBBFlag(i,j)) || (rawStream->GetBGFlag(i,j))) nFlag++;
184 }
185 if(nFlag == 0){ // Pedestal
186 for(Int_t j=5;j<16;j++){
187 Int_t Integrator = rawStream->GetIntegratorFlag(i,j);
188 Float_t pedestal = (float)(rawStream->GetPedestal(i,j));
189 hPEDname[i + 64 * Integrator]->Fill(pedestal);
190 }
191 }
192 if((rawStream->GetBBFlag(i,10)) || (rawStream->GetBGFlag(i,10))){ // Charge
193 Int_t Integrator = rawStream->GetIntegratorFlag(i,10);
194 Float_t charge = (float)(rawStream->GetADC(i));
195 hADCname[i + 64 * Integrator]->Fill(charge);
196 }
03e924d7 197 }
198 delete rawStream;
199 rawStream = 0x0;
200 delete rawReader;
201 rawReader = 0x0;
202 } // end of switch on event type
203
204 nevents_total++;
205 /* free resources */
206 free(event);
f52220c5 207
03e924d7 208 /* exit when last event received, no need to wait for TERM signal */
209 if (eventT==END_OF_RUN) {
210 printf("End Of Run event detected\n");
211 break;
212 }
f52220c5 213
03e924d7 214 } // loop over events
c2d54b4c 215
216 printf("%d physics events processed\n",nevents_physics);
f52220c5 217
f52220c5 218//________________________________________________________________________
219// Computes mean values, dumps them into the output text file
220
df6d42e7 221 for(Int_t i=0; i<128; i++) {
f52220c5 222 hPEDname[i]->GetXaxis()->SetRange(0,kHighCut);
223 PEDmean[i] = hPEDname[i]->GetMean();
224 PEDsigma[i] = hPEDname[i]->GetRMS();
5a8e0a88 225 hADCname[i]->GetXaxis()->SetRange(kLowCut,1024);
03e924d7 226 ADCmean[i] = hADCname[i]->GetMean();
227 ADCsigma[i] = hADCname[i]->GetRMS();
228 fprintf(fp," %.3f %.3f %.3f %.3f\n",PEDmean[i],PEDsigma[i],ADCmean[i],ADCsigma[i]);
c2d54b4c 229 }
230
f52220c5 231//________________________________________________________________________
232// Write root file with histos for users further check - just in case -
233
234 TFile *histoFile = new TFile("VZERO_histos.root","RECREATE");
235
df6d42e7 236 for (Int_t i=0; i<128; i++) {
5a8e0a88 237 hADCname[i]->GetXaxis()->SetRange(0,1024);
f52220c5 238 hADCname[i]->Write();
239 hPEDname[i]->Write(); }
240
241 histoFile->Close();
df6d42e7 242 delete histoFile;
f52220c5 243
244//________________________________________________________________________
245
03e924d7 246 /* close result file */
247 fclose(fp);
248
f52220c5 249 /* export result file to FES */
250 status=daqDA_FES_storeFile("./V0_Ped_Width_Gain.dat","V00da_results");
251 if (status) {
252 printf("Failed to export file : %d\n",status);
253 return -1; }
254
d154146a 255 /* store result file into Online DB */
256 status=daqDA_DB_storeFile("./V0_Ped_Width_Gain.dat","V00da_results");
257 if (status) {
258 printf("Failed to store file into Online DB: %d\n",status);
259 return -1; }
260
f52220c5 261 return status;
262}