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