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