]> git.uio.no Git - u/mrichter/AliRoot.git/blame_incremental - PMD/PMDGAINda.cxx
EVE executable
[u/mrichter/AliRoot.git] / PMD / PMDGAINda.cxx
... / ...
CommitLineData
1/*
2PMD DA for online calibration
3
4contact: basanta@iitb.ac.in, Satyajit.Jena@cern.ch
5Link: https://twiki.cern.ch/twiki/bin/view/ALICE/DA
6Reference run:/afs/cern.ch/user/s/sjena/public/run83496.raw
7Run Type: PHYSICS
8DA Type: MON
9Number of events needed: 1 million for PB+PB, 200 milion for p+p
10Input Files: PMD/Calib/Mapping/Run0_999999999_v0_s0.root, PMD_PED.root, PMD_GAIN_CONFIGFILE, pmd_gain_tempfile.dat
11Output Files: PMDGAINS.root, to be exported to the DAQ FES
12Trigger types used: PHYSICS_EVENT
13
14*/
15extern "C" {
16#include <daqDA.h>
17}
18
19#include "event.h"
20#include "monitor.h"
21
22#include <Riostream.h>
23#include <stdio.h>
24#include <stdlib.h>
25
26//AliRoot
27#include "AliRawReaderDate.h"
28#include "AliPMDCalibPedestal.h"
29#include "AliPMDCalibGain.h"
30#include "AliLog.h"
31#include "AliCDBManager.h"
32
33//ROOT
34#include "TFile.h"
35#include "TH1F.h"
36#include "TBenchmark.h"
37#include "TTree.h"
38#include "TROOT.h"
39#include "TPluginManager.h"
40#include "TSystem.h"
41
42
43/* Main routine
44 Arguments:
45 1- monitoring data source
46*/
47int main(int argc, char **argv) {
48
49 /* magic line from Rene */
50 gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
51 "*",
52 "TStreamerInfo",
53 "RIO",
54 "TStreamerInfo()");
55
56
57 int status = 0;
58
59
60 Int_t filestatus = -1, xvar = 5;
61 Int_t totevt = -1, maxevt = -1;
62 Int_t hotevtsize = -1;
63 Bool_t hotfilestatus = false;
64
65 // Reads the pedestal file and keep the values in memory for subtraction
66
67 AliPMDCalibGain *calibgain = new AliPMDCalibGain();
68
69 // Fetch the pedestal file - PMD_PED.root
70
71 status = daqDA_DB_getFile("PMD_PED.root","PMD_PED.root");
72
73 if(!status)
74 {
75 printf("*** Pedestal file retrieved from DB *** \n");
76 }
77 else
78 {
79 printf("*** Pedestal file NOT retrieved from DB *** \n");
80 return -1;
81 }
82
83 Int_t pstatus = calibgain->ExtractPedestal("PMD_PED.root");
84
85 if(pstatus == -3) return -3;
86
87 TTree *ic = NULL;
88 TTree *meanc = NULL;
89
90 // Retrieve the PMD_GAIN_CONFIGFILE
91 status = daqDA_DB_getFile("PMD_GAIN_CONFIGFILE","PMD_GAIN_CONFIGFILE");
92
93 FILE *fp1 = NULL;
94
95 fp1 = fopen("PMD_GAIN_CONFIGFILE","r");
96
97 if (fp1 == NULL)
98 {
99 printf("*** PMD GAIN Configfile doesn't exist,Provide one ***\n");
100 return -1;
101 }
102 else
103 {
104 fscanf(fp1,"%d %d %d %d %d\n",&filestatus, &xvar, &totevt, &maxevt, &hotevtsize);
105 //printf("%d %d %d %d %d\n",filestatus, xvar, totevt, maxevt, hotevtsize);
106 }
107 fclose(fp1);
108
109
110 if (filestatus == 1)
111 {
112 // Retrieve the Temporray ascii file from DB
113 status = daqDA_DB_getFile("pmd_gain_tempfile.dat","pmd_gain_tempfile.dat");
114 if(!status)
115 {
116 calibgain->ReadTempFile("pmd_gain_tempfile.dat");
117 }
118 else
119 {
120 printf("--- pmd_gain_tempfile.dat: not retrieved from DB ---\n");
121 }
122 // Retrieve the hot cell file from DB - PMD_HOT.root
123 status = daqDA_DB_getFile("PMD_HOT.root","PMD_HOT.root");
124 if(!status)
125 {
126 calibgain->ExtractHotChannel("PMD_HOT.root");
127 }
128 else
129 {
130 printf("--- pmd_gain_tempfile.dat: not retrieved from DB ---\n");
131 }
132 }
133
134
135 // decoding the events
136
137
138 if (argc!=2) {
139 printf("Wrong number of arguments\n");
140 return -1;
141 }
142
143
144 /* define data source : this is argument 1 */
145 status=monitorSetDataSource( argv[1] );
146 if (status!=0) {
147 printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
148 return -1;
149 }
150
151 /* declare monitoring program */
152 status=monitorDeclareMp( __FILE__ );
153 if (status!=0) {
154 printf("monitorDeclareMp() failed : %s\n",monitorDecodeError(status));
155 return -1;
156 }
157
158 /* define wait event timeout - 1s max */
159 monitorSetNowait();
160 monitorSetNoWaitNetworkTimeout(1000);
161
162 /* log start of process */
163 printf("PMD GAIN DA - strted generating the gain of a cell\n");
164
165 /* init some counters */
166 int nevents_physics=0;
167 int nevents_total=0;
168
169 struct eventHeaderStruct *event;
170 eventTypeType eventT = 0;
171
172 Int_t iev=0;
173
174 // Get run number
175
176 int runNr = 0;
177
178 if (getenv("DATE_RUN_NUMBER")==0) {
179 printf("DATE_RUN_NUMBER not properly set.\n");
180 printf("Run Number set to Zero \n");
181 //return -1;
182 }
183 else
184 {
185 runNr = atoi(getenv("DATE_RUN_NUMBER"));
186 }
187
188 if (gSystem->AccessPathName("localOCDB/PMD/Calib/Mapping",kFileExists))
189 {
190 if (gSystem->mkdir("localOCDB/PMD/Calib/Mapping",kTRUE) != 0)
191 {
192 printf("Failed to create directory: localOCDB/PMD/Calib/Mapping");
193 return -1;
194 }
195 }
196 status = daqDA_DB_getFile("PMD/Calib/Mapping/Run0_999999999_v0_s0.root","localOCDB/PMD/Calib/Mapping/Run0_999999999_v0_s0.root");
197 if (status)
198 {
199 printf("Failed to get PMD-Mapping file (PMD/Calib/Mapping/Run0_999999999_v0_s0.root) from DAQdetDB, status=%d\n", status);
200 return -1;
201 }
202
203 // Global initializations
204 AliLog::SetGlobalLogLevel(AliLog::kError);
205 AliCDBManager *man = AliCDBManager::Instance();
206 man->SetDefaultStorage("local://localOCDB");
207 man->SetRun(runNr);
208
209
210 /* main loop (infinite) */
211 for(;;) {
212
213 /* check shutdown condition */
214 if (daqDA_checkShutdown()) {break;}
215
216 /* get next event (blocking call until timeout) */
217 status=monitorGetEventDynamic((void **)&event);
218 if (status==MON_ERR_EOF) {
219 printf ("End of File detected\n");
220 break; /* end of monitoring file has been reached */
221 }
222
223 if (status!=0) {
224 printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
225 break;
226 }
227
228 /* retry if got no event */
229 if (event==NULL) {
230 continue;
231 }
232
233 iev++;
234
235 /* use event - here, just write event id to result file */
236 nevents_total++;
237 eventT=event->eventType;
238 switch (event->eventType){
239
240 /* START OF RUN */
241 case START_OF_RUN:
242 break;
243 /* END START OF RUN */
244
245 /* END OF RUN */
246 case END_OF_RUN:
247 break;
248
249 case PHYSICS_EVENT:
250 nevents_physics++;
251 totevt++;
252 //if(nevents_physics%100 == 0)printf("Physis Events = %d\n",nevents_physics);
253 AliRawReader *rawReader = new AliRawReaderDate((void*)event);
254 TObjArray *pmdddlcont = new TObjArray();
255 calibgain->ProcessEvent(rawReader, pmdddlcont);
256
257 if (totevt%hotevtsize == 0) hotfilestatus = true;
258 delete pmdddlcont;
259 pmdddlcont = 0x0;
260 delete rawReader;
261 rawReader = 0x0;
262
263 }
264
265 /* free resources */
266 free(event);
267
268 }
269
270 /* exit when last event received, no need to wait for TERM signal */
271
272
273 ic = new TTree("ic","PMD Gain tree");
274 meanc = new TTree("meanc","PMD Module mean tree");
275
276 if (filestatus == 0)
277 {
278 TFile *hotRun = new TFile ("PMD_HOT.root","RECREATE");
279
280 TTree *hot = new TTree("hot","PMD Hot cell tree");
281
282 calibgain->FindHotCell(hot,xvar);
283
284 hot->Write();
285 hotRun->Close();
286
287 // store the hot cell root file in the DB
288
289 status = daqDA_DB_storeFile("PMD_HOT.root","PMD_HOT.root");
290
291 // store the hot cell root file in the file exchange server
292
293 printf("root file for hot cell is created and getting exported\n");
294 status = daqDA_FES_storeFile("PMD_HOT.root","PMD_HOT.root");
295 }
296
297 if (hotfilestatus)
298 {
299 TFile *hotRun = new TFile ("PMD_HOT.root","RECREATE");
300
301 TTree *hot = new TTree("hot","PMD Hot cell tree");
302
303 calibgain->FindHotCell(hot,xvar);
304
305 hot->Write();
306 hotRun->Close();
307
308 // store the hot cell root file in the DB
309
310 status = daqDA_DB_storeFile("PMD_HOT.root","PMD_HOT.root");
311
312 // store the hot cell root file in the file exchange server
313
314 printf("root file for hot cell is created and getting exported\n");
315 status = daqDA_FES_storeFile("PMD_HOT.root","PMD_HOT.root");
316 }
317
318
319 fp1 = fopen("PMD_GAIN_CONFIGFILE","w+");
320
321 if (totevt < maxevt)
322 {
323 printf("-----------------------------------------------\n");
324 printf("*** Required Number of Events not reached ***\n");
325 printf("*** Number of Events processed = %d ***\n",totevt);
326 printf("*** Writing the intermediate ASCII file ***\n");
327 printf("-----------------------------------------------\n");
328
329 calibgain->WriteTempFile("pmd_gain_tempfile.dat");
330
331 // Store the Intermediate ascii file in the DB
332 status = daqDA_DB_storeFile("pmd_gain_tempfile.dat","pmd_gain_tempfile.dat");
333
334 filestatus = 1;
335 fprintf(fp1,"%d %d %d %d %d\n",filestatus,xvar,totevt,maxevt,hotevtsize);
336 fclose(fp1);
337
338 // Store the configfile in the DB
339 status = daqDA_DB_storeFile("PMD_GAIN_CONFIGFILE","PMD_GAIN_CONFIGFILE");
340
341 }
342 else if (totevt >= maxevt)
343 {
344 printf("-----------------------------------------------\n");
345 printf("*** Required Number of Events reached = %d ***\n",totevt);
346 printf("*** Writing the PMDGAINS.root file ***\n");
347 printf("-----------------------------------------------\n");
348
349 calibgain->Analyse(ic, meanc);
350
351 TFile * gainRun = new TFile ("PMDGAINS.root","RECREATE");
352 ic->Write();
353 gainRun->Close();
354
355 TFile * meanRun = new TFile ("PMD_MEAN_SM.root","RECREATE");
356 meanc->Write();
357 meanRun->Close();
358
359
360 filestatus = 0;
361 totevt = 0;
362 fprintf(fp1,"%d %d %d %d %d\n",filestatus,xvar,totevt,maxevt,hotevtsize);
363 fclose(fp1);
364
365 // Store the configfile in the DB
366 status = daqDA_DB_storeFile("PMD_GAIN_CONFIGFILE","PMD_GAIN_CONFIGFILE");
367 }
368
369 delete ic;
370 ic = 0;
371
372 delete meanc;
373 meanc = 0;
374
375 delete calibgain;
376
377 /* store the result file on FES */
378
379 if (filestatus == 0)
380 {
381 printf("root file for cell gain is created and getting exported\n");
382 status = daqDA_FES_storeFile("PMDGAINS.root","PMDGAINS.root");
383 printf("root file for normalised means of different modules\n");
384 status = daqDA_FES_storeFile("PMD_MEAN_SM.root","PMD_MEAN_SM.root");
385 }
386
387 if (hotfilestatus)
388 {
389 }
390
391
392 if (status) {
393 status = -2;
394 }
395
396
397
398 return status;
399}