]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/TPCQAda.cxx
Fix funny signs in combo boxes
[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 "AliTPCROC.h"
63 #include "AliTPCCalROC.h"
64 #include "AliTPCCalPad.h"
65 #include "TTreeStream.h"
66 #include "AliLog.h"
67 #include "AliTPCConfigDA.h"
68 //
69 //AMORE
70 //
71 #include <AmoreDA.h>
72 //
73 // TPC calibration algorithm includes
74 //
75 #include "AliTPCdataQA.h"
76
77 //functios, implementation below
78 void SendToAmoreDB(TObject *o, unsigned long32 runNb);
79
80 /* Main routine
81       Arguments: list of DATE raw data files
82 */
83 int main(int argc, char **argv) {
84   /* log start of process */
85   printf("TPCQAda: DA started - %s\n",__FILE__);
86   
87   if (argc<2) {
88     printf("TPCQAda: Wrong number of arguments\n");
89     return -1;
90   }
91   
92  gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
93                                          "*",
94                                          "TStreamerInfo",
95                                          "RIO",
96                                          "TStreamerInfo()");
97
98   AliLog::SetClassDebugLevel("AliRawReaderDate",-5);
99   AliLog::SetClassDebugLevel("AliTPCAltroMapping",-5);
100   AliLog::SetModuleDebugLevel("RAW",-5);
101
102   /* declare monitoring program */
103   int status=monitorDeclareMp( __FILE__ );
104   if (status!=0) {
105     printf("TPCQAda: monitorDeclareMp() failed : %s\n",monitorDecodeError(status));
106     return -1;
107   }
108   //Set network timeout
109   monitorSetNowait();
110   monitorSetNoWaitNetworkTimeout(1000);
111   
112   //variables
113   AliTPCmapper *mapping = 0;   // The TPC mapping
114   unsigned long32 runNb=0;      //run number
115   // if  test setup get parameters from $DAQDA_TEST_DIR
116   
117   if (!mapping){
118     /* copy locally the mapping file from daq detector config db */
119     status = daqDA_DB_getFile(MAPPING_FILE,"./tpcMapping.root");
120     if (status) {
121       printf("TPCQAda: Failed to get mapping file (%s) from DAQdetDB, status=%d\n", MAPPING_FILE, status);
122 //       printf("Continue anyway ... maybe it works?\n");              // !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
123       return -1;   // temporarily uncommented for testing on pcald47 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
124     }
125
126     /* open the mapping file and retrieve mapping object */
127     TFile *fileMapping = new TFile(MAPPING_FILE, "read");
128     mapping = (AliTPCmapper*) fileMapping->Get("tpcMapping");
129     delete fileMapping;
130   }
131
132   if (mapping == 0) {
133     printf("TPCQAda: Failed to get mapping object from %s.  ...\n", MAPPING_FILE);
134     //return -1;
135   } else {
136     printf("TPCQAda: Got mapping object from %s\n", MAPPING_FILE);
137   }
138  
139   //
140   // DA configuration from configuration file
141   //
142   //retrieve configuration file
143   char localfile[255];
144   sprintf(localfile,"./%s",CONFIG_FILE);
145   status = daqDA_DB_getFile(CONFIG_FILE,localfile);
146   if (status) {
147     printf("TPCQAda: Failed to get configuration file (%s) from DAQdetDB, status=%d\n", CONFIG_FILE, status);
148     return -1;
149   }
150   AliTPCConfigDA config(CONFIG_FILE);
151
152   // set default configuration options
153   Double_t updateInterval=30; //seconds
154   TString laserTriggerName("C0LSR-ABCE-NOPF-CENT");
155   TString forceLaserTriggerId("-1");
156   TString monitorAttributes=Form("%d",ATTR_ORIGINAL_EVENT);
157   Bool_t  skipAmore=kFALSE;
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("TPCQAda: 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("TPCQAda: Force laser trigger Id: %s.\n",forceLaserTriggerId.Data());
172   }
173   //skip the amore part
174   if ( config.GetConfigurationMap()->GetValue("SkipAmore") ) {
175     skipAmore=((TObjString*)config.GetConfigurationMap()->GetValue("SkipAmore"))->GetString().Atoi();
176     printf("TPCQAda: Skip Amore set in config\n");
177   }
178   //monitoring Attributes
179   if ( config.GetConfigurationMap()->GetValue("MonitorAttributes") ) {
180     monitorAttributes=config.GetConfigurationMap()->GetValue("MonitorAttributes")->GetName();
181     printf("TPCQAda: Monitor attributes set in config: %s\n",monitorAttributes.Data());
182   }
183   
184   //reject laser triggers in a global partition if we have interleaved laser events
185   unsigned char classId=0;
186   int retClassId=daqDA_getClassIdFromName(laserTriggerName.Data(),&classId);
187   //chek if we shall force the laser trigger id. Mainly for test purposes
188   if (forceLaserTriggerId!="-1"){
189     retClassId=0;
190     classId=static_cast<unsigned char>(forceLaserTriggerId.Atoi());
191   }
192   //create trigger mask
193   if (retClassId==0){
194     //interleaved laser in physics runs
195     //reject laser triggered events
196     TString triggerClasses;
197     //TODO
198     //TODO: in the next release of daq put 49 back to 50!!!
199     //TODO
200     for (unsigned char iclassId=0; iclassId<50; ++iclassId){
201       if (iclassId==classId) continue; //exclude laser trigger
202       triggerClasses+=Form("%u|",(unsigned int)iclassId);
203     }
204     triggerClasses.Chop();
205     char *table[5] = {"PHY",
206                       "Y",
207                       const_cast<char*>(monitorAttributes.Data()),
208                       const_cast<char*>(triggerClasses.Data()),NULL};
209     monitorDeclareTableExtended(table);
210     printf("TPCQAda: Using laser trigger class Id: %u\n",(unsigned int)classId);
211     printf("TPCQAda: Accepted trigger class Ids: %s\n",triggerClasses.Data());
212     printf("TPCQAda: Monitor attributes used: %s\n",monitorAttributes.Data());
213   }
214   
215   //
216   // create calibration object
217   //
218   AliTPCdataQA calibQA(config.GetConfigurationMap());   // qa object
219   calibQA.SetAltroMapping(mapping->GetAltroMapping()); // Use altro mapping we got from daqDetDb
220   
221   //timer
222   TStopwatch stopWatch;
223   
224   //===========================//
225   // loop over RAW data files //
226   //==========================//
227   int nevents=0;
228   for(int i=1;i<argc;i++) {
229
230     /* define data source : this is argument i */
231     printf("TPCQAda: Processing file %s\n", argv[i]);
232     status=monitorSetDataSource( argv[i] );
233     if (status!=0) {
234       printf("TPCQAda: monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
235       return -1;
236     }
237
238     /* read until EOF */
239     while (true) {
240       struct eventHeaderStruct *event;
241
242       /* check shutdown condition */
243       if (daqDA_checkShutdown()) {break;}
244
245       /* get next event (blocking call until timeout) */
246       status=monitorGetEventDynamic((void **)&event);
247       if (status==MON_ERR_EOF) {
248         printf ("TPCQAda: End of File %d detected\n",i);
249         break; /* end of monitoring file has been reached */
250       }
251
252       if (status!=0) {
253         printf("TPCQAda: monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
254         break;
255       }
256
257       /* retry if got no event */
258       if (event==NULL) {
259         continue;
260       }
261       nevents++;
262       // get the run number
263       runNb = event->eventRunNb;
264       //  QA
265       AliRawReader *rawReader = new AliRawReaderDate((void*)event);
266       calibQA.ProcessEvent(rawReader);
267       delete rawReader;
268       
269       // sending to AMOREdb
270       if (stopWatch.RealTime()>updateInterval){
271         if (!skipAmore) SendToAmoreDB(&calibQA,runNb);
272         stopWatch.Start();
273       } else {
274         stopWatch.Continue();
275       }
276       
277       /* free resources */
278       free(event);
279     }
280   }
281
282   calibQA.Analyse(); 
283   printf ("TPCQAda: %d events processed\n",nevents);
284
285   TFile * fileTPC = new TFile (RESULT_FILE,"recreate");
286   calibQA.Write("tpcCalibQA");
287   delete fileTPC;
288   printf("TPCQAda: Wrote %s\n",RESULT_FILE);
289
290   /* store the result file on FXS */
291
292   status=daqDA_FES_storeFile(RESULT_FILE,FILE_ID);
293   if (status) {
294     status = -2;
295   }
296   //
297   //Send objects to the AMORE DB
298   //
299   if (!skipAmore){
300     printf ("TPCQAda: AMORE part\n");
301     SendToAmoreDB(&calibQA, runNb);
302   }
303   
304   return status;
305 }
306
307 void SendToAmoreDB(TObject *o, unsigned long32 runNb)
308 {
309   //AMORE
310   const char *amoreDANameorig=gSystem->Getenv("AMORE_DA_NAME");
311   //cheet a little -- temporary solution (hopefully)
312   //
313   //currently amoreDA uses the environment variable AMORE_DA_NAME to create the mysql
314   //table in which the calib objects are stored. This table is dropped each time AmoreDA
315   //is initialised. This of course makes a problem if we would like to store different
316   //calibration entries in the AMORE DB. Therefore in each DA which writes to the AMORE DB
317   //the AMORE_DA_NAME env variable is overwritten.
318   
319   gSystem->Setenv("AMORE_DA_NAME","TPC-dataQA");
320   //
321   // end cheet
322   TDatime time;
323   TObjString info(Form("Run: %u; Date: %s",runNb,time.AsSQLString()));
324   
325   amore::da::AmoreDA amoreDA(amore::da::AmoreDA::kSender);
326   Int_t statusDA=0;
327   statusDA+=amoreDA.Send("DataQA",o);
328   statusDA+=amoreDA.Send("Info",&info);
329   if ( statusDA!=0 )
330     printf("TPCQAda: Waring: Failed to write one of the calib objects to the AMORE database\n");
331   // reset env var
332   if (amoreDANameorig) gSystem->Setenv("AMORE_DA_NAME",amoreDANameorig);
333 }