]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/TPCQAda.cxx
added fake OCDB object needed for QA
[u/mrichter/AliRoot.git] / TPC / TPCQAda.cxx
1 /*
2 TPC DA for online calibration
3
4 Contact: Haavard.Helstrup@cern.ch, peter.christiansen@hep.lu.se
5 Link: 
6 Run Type: PHYSICS STANDALONE DAQ
7 DA Type: MON
8 Number of events needed: 500
9 Input Files: /castor/cern.ch/alice/raw/global/2009/08/22/11/09000080958023.30.root
10 Output Files: tpcQA.root, to be exported to the DAQ FXS
11 fileId:   QA
12 Trigger types used: PHYSICS_EVENT
13
14 */
15
16 /*
17
18 TPCQAda.cxx - algorithm for TPC RAW QA
19
20 10/06/2007  sylvain.chapeland@cern.ch :  first version - clean skeleton based on DAQ DA case1
21 06/12/2007  haavard.helstrup@cern.ch  :  created CE DA based on pulser code
22 09/06/2008  peter.christiansen@hep.lu.se and haavard.helstrup@cern.ch  :  created QA DA based on AliTPCdataQA code
23
24 10/09/2009  Jens.Wiechula@cern.ch:     Export object to AMOREdb after a defined update interval for QA
25 26/01/2010  Jens.Wiechula@cern.ch:     Exclude laser triggers when running in a global partition
26
27 contact: marian.ivanov@cern.ch, peter.christiansen@hep.lu.se
28
29 This process reads RAW data from the files provided as command line arguments
30 and save results in a file (named from RESULT_FILE define - see below).
31
32 */
33
34 #define RESULT_FILE "tpcQA.root"
35 #define FILE_ID "QA"
36 #define MAPPING_FILE "tpcMapping.root"
37 #define CONFIG_FILE "TPCQAda.conf"
38
39
40 #include <daqDA.h>
41 #include "event.h"
42 #include "monitor.h"
43 #include <stdio.h>
44 #include <stdlib.h>
45
46 //
47 //Root includes
48 //
49 #include <TFile.h>
50 #include <TROOT.h>
51 #include <TPluginManager.h>
52 #include <TSystem.h>
53 #include <TStopwatch.h>
54 #include <TObject.h>
55 #include <TMap.h>
56 //
57 //AliRoot includes
58 //
59 #include "AliRawReader.h"
60 #include "AliRawReaderDate.h"
61 #include "AliTPCmapper.h"
62 #include "AliTPCRawStream.h"
63 #include "AliTPCROC.h"
64 #include "AliTPCCalROC.h"
65 #include "AliTPCCalPad.h"
66 #include "TTreeStream.h"
67 #include "AliLog.h"
68 #include "AliTPCConfigDA.h"
69 //
70 //AMORE
71 //
72 #include <AmoreDA.h>
73 //
74 // TPC calibration algorithm includes
75 //
76 #include "AliTPCdataQA.h"
77
78 //functios, implementation below
79 void SendToAmoreDB(TObject *o, unsigned long32 runNb);
80
81 /* Main routine
82       Arguments: list of DATE raw data files
83 */
84 int main(int argc, char **argv) {
85   /* log start of process */
86   printf("TPCQAda: DA started - %s\n",__FILE__);
87   
88   if (argc<2) {
89     printf("TPCQAda: Wrong number of arguments\n");
90     return -1;
91   }
92   
93  gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
94                                          "*",
95                                          "TStreamerInfo",
96                                          "RIO",
97                                          "TStreamerInfo()");
98
99   AliLog::SetClassDebugLevel("AliTPCRawStream",-5);
100   AliLog::SetClassDebugLevel("AliRawReaderDate",-5);
101   AliLog::SetClassDebugLevel("AliTPCAltroMapping",-5);
102   AliLog::SetModuleDebugLevel("RAW",-5);
103
104   /* declare monitoring program */
105   int status=monitorDeclareMp( __FILE__ );
106   if (status!=0) {
107     printf("TPCQAda: monitorDeclareMp() failed : %s\n",monitorDecodeError(status));
108     return -1;
109   }
110   //Set network timeout
111   monitorSetNowait();
112   monitorSetNoWaitNetworkTimeout(1000);
113   
114   //variables
115   AliTPCmapper *mapping = 0;   // The TPC mapping
116   unsigned long32 runNb=0;      //run number
117   // if  test setup get parameters from $DAQDA_TEST_DIR
118   
119   if (!mapping){
120     /* copy locally the mapping file from daq detector config db */
121     status = daqDA_DB_getFile(MAPPING_FILE,"./tpcMapping.root");
122     if (status) {
123       printf("TPCQAda: Failed to get mapping file (%s) from DAQdetDB, status=%d\n", MAPPING_FILE, status);
124 //       printf("Continue anyway ... maybe it works?\n");              // !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
125       return -1;   // temporarily uncommented for testing on pcald47 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
126     }
127
128     /* open the mapping file and retrieve mapping object */
129     TFile *fileMapping = new TFile(MAPPING_FILE, "read");
130     mapping = (AliTPCmapper*) fileMapping->Get("tpcMapping");
131     delete fileMapping;
132   }
133
134   if (mapping == 0) {
135     printf("TPCQAda: Failed to get mapping object from %s.  ...\n", MAPPING_FILE);
136     //return -1;
137   } else {
138     printf("TPCQAda: Got mapping object from %s\n", MAPPING_FILE);
139   }
140  
141   //
142   // DA configuration from configuration file
143   //
144   //retrieve configuration file
145   char localfile[255];
146   sprintf(localfile,"./%s",CONFIG_FILE);
147   status = daqDA_DB_getFile(CONFIG_FILE,localfile);
148   if (status) {
149     printf("TPCQAda: Failed to get configuration file (%s) from DAQdetDB, status=%d\n", CONFIG_FILE, status);
150     return -1;
151   }
152   AliTPCConfigDA config(CONFIG_FILE);
153
154   // set default configuration options
155   Double_t updateInterval=30; //seconds
156   TString laserTriggerName("C0LSR-ABCE-NOPF-CENT");
157   TString forceLaserTriggerId("-1");
158   
159   //amore update interval
160   Double_t valConf=config.GetValue("AmoreUpdateInterval");
161   if ( valConf>0 ) updateInterval=valConf;
162   
163   //laser trigger class name
164   if ( config.GetConfigurationMap()->GetValue("LaserTriggerName") ) {
165     laserTriggerName=config.GetConfigurationMap()->GetValue("LaserTriggerName")->GetName();
166     printf("TPCRAWda: Laser trigger class name set to: %s.\n",laserTriggerName.Data());
167   }
168   //force laser trigger id
169   if ( config.GetConfigurationMap()->GetValue("ForceLaserTriggerId") ) {
170     forceLaserTriggerId=config.GetConfigurationMap()->GetValue("ForceLaserTriggerId")->GetName();
171     printf("TPCRAWda: Force laser trigger Id: %s.\n",forceLaserTriggerId.Data());
172   }
173   
174   
175   //reject laser triggers in a global partition if we have interleaved laser events
176   unsigned char classId=0;
177   int retClassId=daqDA_getClassIdFromName(laserTriggerName.Data(),&classId);
178   //chek if we shall force the laser trigger id. Mainly for test purposes
179   if (forceLaserTriggerId!="-1"){
180     retClassId=0;
181     classId=static_cast<unsigned char>(forceLaserTriggerId.Atoi());
182   }
183   //create trigger mask
184   if (retClassId==0){
185     //interleaved laser in physics runs
186     //reject laser triggered events
187     TString triggerClasses;
188     //TODO
189     //TODO: in the next release of daq put 49 back to 50!!!
190     //TODO
191     for (unsigned char iclassId=0; iclassId<49; ++iclassId){
192       if (iclassId==classId) continue; //exclude laser trigger
193       triggerClasses+=Form("%u|",(unsigned int)iclassId);
194     }
195     triggerClasses.Chop();
196     char *table[5] = {"PHY","Y","*",const_cast<char*>(triggerClasses.Data()),NULL};
197     monitorDeclareTableExtended(table);
198     printf("TPCRAWda: Using laser trigger class Id: %u\n",(unsigned int)classId);
199     printf("TPCRAWda: Accepted trigger class Ids: %s\n",triggerClasses.Data());
200   }
201   
202   //
203   // create calibration object
204   //
205   AliTPCdataQA calibQA(config.GetConfigurationMap());   // qa object
206   calibQA.SetAltroMapping(mapping->GetAltroMapping()); // Use altro mapping we got from daqDetDb
207   
208   //timer
209   TStopwatch stopWatch;
210   
211   //===========================//
212   // loop over RAW data files //
213   //==========================//
214   int nevents=0;
215   for(int i=1;i<argc;i++) {
216
217     /* define data source : this is argument i */
218     printf("TPCQAda: Processing file %s\n", argv[i]);
219     status=monitorSetDataSource( argv[i] );
220     if (status!=0) {
221       printf("TPCQAda: monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
222       return -1;
223     }
224
225     /* read until EOF */
226     while (true) {
227       struct eventHeaderStruct *event;
228
229       /* check shutdown condition */
230       if (daqDA_checkShutdown()) {break;}
231
232       /* get next event (blocking call until timeout) */
233       status=monitorGetEventDynamic((void **)&event);
234       if (status==MON_ERR_EOF) {
235         printf ("TPCQAda: End of File %d detected\n",i);
236         break; /* end of monitoring file has been reached */
237       }
238
239       if (status!=0) {
240         printf("TPCQAda: monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
241         break;
242       }
243
244       /* retry if got no event */
245       if (event==NULL) {
246         continue;
247       }
248       nevents++;
249       // get the run number
250       runNb = event->eventRunNb;
251       //  QA
252       AliRawReader *rawReader = new AliRawReaderDate((void*)event);
253       calibQA.ProcessEvent(rawReader);
254       delete rawReader;
255       
256       // sending to AMOREdb
257       if (stopWatch.RealTime()>updateInterval){
258         SendToAmoreDB(&calibQA,runNb);
259         stopWatch.Start();
260       } else {
261         stopWatch.Continue();
262       }
263       
264       /* free resources */
265       free(event);
266     }
267   }
268
269   calibQA.Analyse(); 
270   printf ("TPCQAda: %d events processed\n",nevents);
271
272   TFile * fileTPC = new TFile (RESULT_FILE,"recreate");
273   calibQA.Write("tpcCalibQA");
274   delete fileTPC;
275   printf("TPCQAda: Wrote %s\n",RESULT_FILE);
276
277   /* store the result file on FXS */
278
279   status=daqDA_FES_storeFile(RESULT_FILE,FILE_ID);
280   if (status) {
281     status = -2;
282   }
283   //
284   //Send objects to the AMORE DB
285   //
286   printf ("TPCQAda: AMORE part\n");
287   SendToAmoreDB(&calibQA, runNb);
288   
289   return status;
290 }
291
292 void SendToAmoreDB(TObject *o, unsigned long32 runNb)
293 {
294   //AMORE
295   const char *amoreDANameorig=gSystem->Getenv("AMORE_DA_NAME");
296   //cheet a little -- temporary solution (hopefully)
297   //
298   //currently amoreDA uses the environment variable AMORE_DA_NAME to create the mysql
299   //table in which the calib objects are stored. This table is dropped each time AmoreDA
300   //is initialised. This of course makes a problem if we would like to store different
301   //calibration entries in the AMORE DB. Therefore in each DA which writes to the AMORE DB
302   //the AMORE_DA_NAME env variable is overwritten.
303   
304   gSystem->Setenv("AMORE_DA_NAME","TPC-dataQA");
305   //
306   // end cheet
307   TDatime time;
308   TObjString info(Form("Run: %u; Date: %s",runNb,time.AsSQLString()));
309   
310   amore::da::AmoreDA amoreDA(amore::da::AmoreDA::kSender);
311   Int_t statusDA=0;
312   statusDA+=amoreDA.Send("DataQA",o);
313   statusDA+=amoreDA.Send("Info",&info);
314   if ( statusDA!=0 )
315     printf("TPCQAda: Waring: Failed to write one of the calib objects to the AMORE database\n");
316   // reset env var
317   if (amoreDANameorig) gSystem->Setenv("AMORE_DA_NAME",amoreDANameorig);
318 }