]> git.uio.no Git - u/mrichter/AliRoot.git/blame - VZERO/VZEROda.cxx
Channel mapping applied to output pedestal file sent to FXS
[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
26206101 51Int_t GetOfflineChannel(Int_t channel);
52
03e924d7 53/* Main routine --- Arguments: monitoring data source */
f52220c5 54
55int main(int argc, char **argv) {
56
6ee2e288 57/* magic line from Cvetan */
58 gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
59 "*",
60 "TStreamerInfo",
61 "RIO",
62 "TStreamerInfo()");
f52220c5 63 int status;
03e924d7 64 if (argc!=2) {
65 printf("Wrong number of arguments\n");
66 return -1;
67 }
c2d54b4c 68
26206101 69 // Online values (using FEE channel numbering):
df6d42e7 70 Double_t ADCmean[128];
03e924d7 71 Double_t ADCsigma[128];
df6d42e7 72 Double_t PEDmean[128];
73 Double_t PEDsigma[128];
26206101 74 // Offline values(same but ordered as in aliroot for offliners):
75 Double_t ADCmean_Off[128];
76 Double_t ADCsigma_Off[128];
77 Double_t PEDmean_Off[128];
78 Double_t PEDsigma_Off[128];
79
c2d54b4c 80//___________________________________________________
81// Get cuts from V00DA.config file
82
5a8e0a88 83 Int_t kLowCut; // = 60; low cut on signal distribution - to be tuned
84 Int_t kHighCut; // = 50; high cut on pedestal distribution - to be tuned
85
c2d54b4c 86 status = daqDA_DB_getFile("V00DA.config","./V00DA.config");
87 if (status) {
d154146a 88 printf("Failed to get config file (V00DA.config) from DAQ DB, status=%d\n", status);
c2d54b4c 89 return -1;
90 }
91 /* open the config file and retrieve cuts */
92 FILE *fpConfig = fopen("V00DA.config","r");
93 fscanf(fpConfig,"%d %d",&kLowCut,&kHighCut);
94 fclose(fpConfig);
f52220c5 95
c2d54b4c 96 printf("LowCut on signal = %d ; HighCut on pedestal = %d\n",kLowCut,kHighCut);
97
f52220c5 98//___________________________________________________
99// Book HISTOGRAMS - dynamics of p-p collisions -
100
101 char ADCname[6];
102 char PEDname[6];
df6d42e7 103 TH1F *hADCname[128];
104 TH1F *hPEDname[128];
f52220c5 105
106 char texte[12];
df6d42e7 107 for (Int_t i=0; i<128; i++) {
f52220c5 108 sprintf(ADCname,"hADC%d",i);
109 sprintf(texte,"ADC cell%d",i);
b68c91e6 110 hADCname[i] = new TH1F(ADCname,texte,1024,-0.5, 1023.5);
f52220c5 111 sprintf(PEDname,"hPED%d",i);
112 sprintf(texte,"PED cell%d",i);
b68c91e6 113 hPEDname[i] = new TH1F(PEDname,texte,1024,-0.5, 1023.5);
03e924d7 114 }
f52220c5 115//___________________________________________________
d154146a 116
f52220c5 117 /* open result file to be exported to FES */
118 FILE *fp=NULL;
b68c91e6 119 fp=fopen("./V0_Ped_Width_Gain.dat","w");
f52220c5 120 if (fp==NULL) {
121 printf("Failed to open result file\n");
122 return -1;}
123
03e924d7 124 /* define data source : this is argument 1 */
125 status=monitorSetDataSource( argv[1] );
126 if (status!=0) {
127 printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
128 return -1;
129 }
130
131 /* declare monitoring program */
132 status=monitorDeclareMp( __FILE__ );
133 if (status!=0) {
134 printf("monitorDeclareMp() failed : %s\n",monitorDecodeError(status));
135 return -1;
136 }
137
138 /* define wait event timeout - 1s max */
139 monitorSetNowait();
140 monitorSetNoWaitNetworkTimeout(1000);
141
f52220c5 142 /* init counters on events */
143 int nevents_physics=0;
144 int nevents_total=0;
145
03e924d7 146 /* loop on events (infinite) */
147 for(;;) {
148 struct eventHeaderStruct *event;
149 eventTypeType eventT;
150
151 /* check shutdown condition */
152 if (daqDA_checkShutdown()) {break;}
153
154 /* get next event (blocking call until timeout) */
155 status=monitorGetEventDynamic((void **)&event);
156 if (status==MON_ERR_EOF) {
157 printf ("End of File detected\n");
158 break; /* end of monitoring file has been reached */
159 }
160
161 if (status!=0) {
162 printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
163 break;
164 }
165
166 /* retry if got no event */
167 if (event==NULL) continue;
168
169 /* decode event */
170 eventT=event->eventType;
f52220c5 171
03e924d7 172 switch (event->eventType){
f52220c5 173
03e924d7 174 case START_OF_RUN:
175 break;
f52220c5 176
03e924d7 177 case END_OF_RUN:
178 printf("End Of Run detected\n");
179 break;
f52220c5 180
03e924d7 181 case PHYSICS_EVENT:
182 nevents_physics++;
c2d54b4c 183
03e924d7 184 AliRawReader *rawReader = new AliRawReaderDate((void*)event);
df6d42e7 185
03e924d7 186 AliVZERORawStream* rawStream = new AliVZERORawStream(rawReader);
187 rawStream->Next();
188 for(Int_t i=0; i<64; i++) {
5a8e0a88 189 Int_t nFlag = 0;
190 for(Int_t j=0; j<21; j++) {
191 if((rawStream->GetBBFlag(i,j)) || (rawStream->GetBGFlag(i,j))) nFlag++;
192 }
193 if(nFlag == 0){ // Pedestal
194 for(Int_t j=5;j<16;j++){
195 Int_t Integrator = rawStream->GetIntegratorFlag(i,j);
196 Float_t pedestal = (float)(rawStream->GetPedestal(i,j));
197 hPEDname[i + 64 * Integrator]->Fill(pedestal);
198 }
199 }
200 if((rawStream->GetBBFlag(i,10)) || (rawStream->GetBGFlag(i,10))){ // Charge
201 Int_t Integrator = rawStream->GetIntegratorFlag(i,10);
202 Float_t charge = (float)(rawStream->GetADC(i));
203 hADCname[i + 64 * Integrator]->Fill(charge);
204 }
03e924d7 205 }
206 delete rawStream;
207 rawStream = 0x0;
208 delete rawReader;
209 rawReader = 0x0;
210 } // end of switch on event type
211
212 nevents_total++;
213 /* free resources */
214 free(event);
f52220c5 215
03e924d7 216 /* exit when last event received, no need to wait for TERM signal */
217 if (eventT==END_OF_RUN) {
218 printf("End Of Run event detected\n");
219 break;
220 }
f52220c5 221
03e924d7 222 } // loop over events
c2d54b4c 223
224 printf("%d physics events processed\n",nevents_physics);
f52220c5 225
26206101 226//___________________________________________________________________________
227// Computes mean values, converts FEE channels into Offline AliRoot channels
228// and dumps the ordered values into the output text file for SHUTTLE
f52220c5 229
df6d42e7 230 for(Int_t i=0; i<128; i++) {
f52220c5 231 hPEDname[i]->GetXaxis()->SetRange(0,kHighCut);
232 PEDmean[i] = hPEDname[i]->GetMean();
233 PEDsigma[i] = hPEDname[i]->GetRMS();
5a8e0a88 234 hADCname[i]->GetXaxis()->SetRange(kLowCut,1024);
03e924d7 235 ADCmean[i] = hADCname[i]->GetMean();
236 ADCsigma[i] = hADCname[i]->GetRMS();
26206101 237// printf(" i = %d, %.3f %.3f %.3f %.3f\n",i,PEDmean[i],PEDsigma[i],ADCmean[i],ADCsigma[i]);
238 if (i < 64) {
239 Int_t j = GetOfflineChannel(i);
240 PEDmean_Off[j] = PEDmean[i];
241 PEDsigma_Off[j] = PEDsigma[i];
242 ADCmean_Off[j] = ADCmean[i];
243 ADCsigma_Off[j] = ADCsigma[i]; }
244 else{
245 Int_t j = GetOfflineChannel(i-64);
246 PEDmean_Off[j+64] = PEDmean[i];
247 PEDsigma_Off[j+64] = PEDsigma[i];
248 ADCmean_Off[j+64] = ADCmean[i];
249 ADCsigma_Off[j+64] = ADCsigma[i];
250 }
251 }
252
253 for(Int_t j=0; j<128; j++) {
254// printf(" j = %d, %.3f %.3f %.3f %.3f\n",j,PEDmean_Off[j],PEDsigma_Off[j],
255// ADCmean_Off[j],ADCsigma_Off[j]);
256 fprintf(fp," %.3f %.3f %.3f %.3f\n",PEDmean_Off[j],PEDsigma_Off[j],
257 ADCmean_Off[j],ADCsigma_Off[j]);
c2d54b4c 258 }
259
f52220c5 260//________________________________________________________________________
261// Write root file with histos for users further check - just in case -
262
263 TFile *histoFile = new TFile("VZERO_histos.root","RECREATE");
264
df6d42e7 265 for (Int_t i=0; i<128; i++) {
5a8e0a88 266 hADCname[i]->GetXaxis()->SetRange(0,1024);
f52220c5 267 hADCname[i]->Write();
268 hPEDname[i]->Write(); }
269
270 histoFile->Close();
df6d42e7 271 delete histoFile;
f52220c5 272
273//________________________________________________________________________
274
03e924d7 275 /* close result file */
276 fclose(fp);
277
f52220c5 278 /* export result file to FES */
279 status=daqDA_FES_storeFile("./V0_Ped_Width_Gain.dat","V00da_results");
280 if (status) {
281 printf("Failed to export file : %d\n",status);
282 return -1; }
283
d154146a 284 /* store result file into Online DB */
285 status=daqDA_DB_storeFile("./V0_Ped_Width_Gain.dat","V00da_results");
286 if (status) {
287 printf("Failed to store file into Online DB: %d\n",status);
288 return -1; }
289
f52220c5 290 return status;
291}
26206101 292
293 Int_t GetOfflineChannel(Int_t channel) {
294
295// Channel mapping Online - Offline:
296
297 Int_t fOfflineChannel[64] = {39, 38, 37, 36, 35, 34, 33, 32,
298 47, 46, 45, 44, 43, 42, 41, 40,
299 55, 54, 53, 52, 51, 50, 49, 48,
300 63, 62, 61, 60, 59, 58, 57, 56,
301 7, 6, 5, 4, 3, 2, 1, 0,
302 15, 14, 13, 12, 11, 10, 9, 8,
303 23, 22, 21, 20, 19, 18, 17, 16,
304 31, 30, 29, 28, 27, 26, 25, 24};
305 return fOfflineChannel[channel];
306}