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