]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PMD/PMDGAINda.cxx
feature added to store and fetch from DB for online gain calibration
[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 //#include "daqDA.h"
22
23 #include <Riostream.h>
24 #include <stdio.h>
25 #include <stdlib.h>
26
27 //AliRoot
28 #include "AliRawReaderDate.h"
29 #include "AliPMDCalibPedestal.h"
30 #include "AliPMDCalibGain.h"
31
32 //ROOT
33 #include "TFile.h"
34 #include "TH1F.h"
35 #include "TBenchmark.h"
36 #include "TTree.h"
37 #include "TROOT.h"
38 #include "TPluginManager.h"
39
40
41
42 /* Main routine
43       Arguments: 
44       1- monitoring data source
45 */
46 int main(int argc, char **argv) {
47   
48     /* magic line from Rene */
49     gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
50                                           "*",
51                                           "TStreamerInfo",
52                                           "RIO",
53                                           "TStreamerInfo()");
54
55
56     int status = 0;
57
58
59     Int_t filestatus = -1, totevt = -1;
60     Int_t maxevt = -1;
61
62     // Reads the pedestal file and keep the values in memory for subtraction
63
64     AliPMDCalibGain calibgain;
65
66     // Fetch the pedestal file - PMD_PED.root 
67     status = daqDA_DB_getFile("PMD_PED.root","PMD_PED.root");
68
69     if(!status)
70       {
71         printf("*** Pedestal file retrieved from DB *** \n");
72       }
73     else
74       {
75         printf("*** Pedestal file NOT retrieved from DB *** \n");
76         return -1;
77       }
78
79     Int_t pstatus = calibgain.ExtractPedestal("PMD_PED.root");
80
81     if(pstatus == -3) return -3;
82
83     TTree *ic = NULL;
84
85     // Retrieve the PMD_GAIN_CONFIGFILE
86     status = daqDA_DB_getFile("PMD_GAIN_CONFIGFILE","PMD_GAIN_CONFIGFILE");
87
88     FILE *fp1 = NULL;
89
90     fp1 = fopen("PMD_GAIN_CONFIGFILE","r");
91
92     if (fp1 == NULL)
93       {
94         printf("*** PMD GAIN Configfile doesn't exist,Provide one ***\n");
95         return -1;
96       }
97     else
98       {
99         fscanf(fp1,"%d %d %d\n",&filestatus, &totevt,&maxevt);
100         //printf("%d %d %d\n",filestatus, totevt, maxevt);
101       }
102     fclose(fp1);
103     
104     if (filestatus == 1)
105       {
106         // Retrieve the Temporray ascii file from DB
107         status = daqDA_DB_getFile("pmd_gain_tempfile.dat","pmd_gain_tempfile.dat");
108         if(!status)
109           {
110             calibgain.ReadTempFile("pmd_gain_tempfile.dat");
111           }
112         else
113           {
114             printf("--- pmd_gain_tempfile.dat: not retrieved from DB ---\n");
115           }
116       }
117
118     // decoding the events
119     
120
121     if (argc!=2) {
122         printf("Wrong number of arguments\n");
123         return -1;
124     }
125     
126     
127     /* define data source : this is argument 1 */  
128     status=monitorSetDataSource( argv[1] );
129     if (status!=0) {
130         printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
131         return -1;
132     }
133     
134     /* declare monitoring program */
135     status=monitorDeclareMp( __FILE__ );
136     if (status!=0) {
137         printf("monitorDeclareMp() failed : %s\n",monitorDecodeError(status));
138         return -1;
139     }
140     
141     /* define wait event timeout - 1s max */
142     monitorSetNowait();
143     monitorSetNoWaitNetworkTimeout(1000);
144     
145     /* log start of process */
146     printf("PMD GAIN DA - strted generating the gain of a cell\n");  
147     
148     /* init some counters */
149     int nevents_physics=0;
150     int nevents_total=0;
151     
152     struct eventHeaderStruct *event;
153     eventTypeType eventT = 0;
154
155     Int_t iev=0;
156     
157     /* main loop (infinite) */
158     for(;;) {
159         
160         /* check shutdown condition */
161         if (daqDA_checkShutdown()) {break;}
162         
163         /* get next event (blocking call until timeout) */
164         status=monitorGetEventDynamic((void **)&event);
165         if (status==MON_ERR_EOF) {
166             printf ("End of File detected\n");
167             break; /* end of monitoring file has been reached */
168         }
169         
170         if (status!=0) {
171             printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
172             break;
173         }
174         
175         /* retry if got no event */
176         if (event==NULL) {
177             continue;
178         }
179         
180         iev++; 
181         
182         /* use event - here, just write event id to result file */
183         nevents_total++;
184         eventT=event->eventType;
185         switch (event->eventType){
186       
187             /* START OF RUN */
188             case START_OF_RUN:
189                 break;
190                 /* END START OF RUN */
191                 
192                 /* END OF RUN */
193             case END_OF_RUN:
194                 break;
195                 
196             case PHYSICS_EVENT:
197                 nevents_physics++;
198                 //if(nevents_physics%100 == 0)printf("Physis Events = %d\n",nevents_physics);
199                 AliRawReader *rawReader = new AliRawReaderDate((void*)event);
200                 TObjArray *pmdddlcont = new TObjArray();
201                 calibgain.ProcessEvent(rawReader, pmdddlcont);
202
203                 delete pmdddlcont;
204                 pmdddlcont = 0x0;
205                 delete rawReader;
206                 rawReader = 0x0;
207                 
208         }
209        
210         /* free resources */
211         free(event);
212         
213     }
214
215     /* exit when last event received, no need to wait for TERM signal */
216
217     ic = new TTree("ic","PMD Gain tree");
218
219     totevt += nevents_physics++;
220
221     fp1 = fopen("PMD_GAIN_CONFIGFILE","w+");
222
223     if (totevt < maxevt)
224       {
225         printf("-----------------------------------------------\n");
226         printf("***  Required Number of Events not reached  ***\n");
227         printf("***  Number of Events processed = %d        ***\n",totevt);
228         printf("***  Writing the intermediate ASCII file    ***\n");
229         printf("-----------------------------------------------\n");
230
231         calibgain.WriteTempFile("pmd_gain_tempfile.dat");
232
233         // Store the Intermediate ascii file in the DB
234         status = daqDA_DB_storeFile("pmd_gain_tempfile.dat","pmd_gain_tempfile.dat");
235
236         filestatus = 1;
237         fprintf(fp1,"%d %d %d\n",filestatus,totevt,maxevt);
238
239         // Store the configfile in the DB
240         status = daqDA_DB_storeFile("PMD_GAIN_CONFIGFILE","PMD_GAIN_CONFIGFILE");
241
242       }
243     else if (totevt >= maxevt)
244       {
245         printf("-----------------------------------------------\n");
246         printf("***  Required Number of Events reached = %d ***\n",totevt);
247         printf("***  Writing the PMDGAINS.root file           ***\n");
248         printf("-----------------------------------------------\n");
249
250         calibgain.Analyse(ic);
251
252         TFile * gainRun = new TFile ("PMDGAINS.root","RECREATE"); 
253         ic->Write();
254         gainRun->Close();
255
256         filestatus = 0;
257         totevt     = 0;
258         fprintf(fp1,"%d %d %d\n",filestatus,totevt,maxevt);
259
260         // Store the configfile in the DB
261         status = daqDA_DB_storeFile("PMD_GAIN_CONFIGFILE","PMD_GAIN_CONFIGFILE");
262
263       }
264     fclose(fp1);
265     
266     delete ic;
267     ic = 0;
268     
269
270     /* store the result file on FES */
271  
272     if (filestatus == 0)
273       {
274         printf("root file is created and getting exported\n");
275         status = daqDA_FES_storeFile("PMDGAINS.root","gaincalib");
276       }
277
278     if (status) {
279       status = -2;
280     }
281
282
283     return status;
284 }