]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/TPCCEnewda.cxx
bug corrected
[u/mrichter/AliRoot.git] / TPC / TPCCEnewda.cxx
1 /*
2 TPC DA for online calibration
3
4 Contact: Haavard.Helstrup@cern.ch
5 Link: 
6 Run Type: PHYSICS STANDALONE DAQ
7 DA Type: MON
8 Number of events needed: 500
9 Input Files: 
10 Output Files: tpcCE.root, to be exported to the DAQ FXS
11 fileId:   CEnew
12 Trigger types used: PHYSICS_EVENT
13
14 */
15
16 /*
17
18 TPCCEda.cxx - calibration algorithm for TPC Central Electrode events
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 19/09/2008  J.Wiechula@gsi.de:           Added support for configuration files.
23
24 contact: marian.ivanov@cern.ch
25
26
27 This process reads RAW data from the files provided as command line arguments
28 and save results in a file (named from RESULT_FILE define - see below).
29
30 */
31
32 #define RESULT_FILE "tpcCEnew.root"
33 #define FILE_ID "CEnew"
34 #define MAPPING_FILE "tpcMapping.root"
35 #define CONFIG_FILE "TPCCEda.conf"
36 #define LASER_TRACK_FILE "LaserTracks.root"
37
38 #include <daqDA.h>
39 #include "event.h"
40 #include "monitor.h"
41 #include <stdio.h>
42 #include <stdlib.h>
43 #include <vector>
44 //
45 //Root includes
46 //
47 #include <TFile.h>
48 #include "TROOT.h"
49 #include "TPluginManager.h"
50 #include "TSystem.h"
51 #include "TString.h"
52 #include "TObjString.h"
53 #include "TDatime.h"
54 #include "TStopwatch.h"
55 #include "TMap.h"
56 #include "TGraph.h"
57 #include "TMath.h"
58 //
59 //AliRoot includes
60 //
61 #include "AliRawReader.h"
62 #include "AliRawReaderDate.h"
63 #include "AliTPCmapper.h"
64 #include "AliTPCRawStream.h"
65 #include "AliTPCROC.h"
66 #include "AliTPCCalROC.h"
67 #include "AliTPCCalPad.h"
68 #include "AliMathBase.h"
69 #include "TTreeStream.h"
70 #include "AliLog.h"
71 #include "AliTPCConfigDA.h"
72 //
73 //AMORE
74 //
75 #include <AmoreDA.h>
76 //
77 // TPC calibration algorithm includes
78 //
79 #include "AliTPCCalibCE.h"
80
81
82 //functios, implementation below
83 void SendToAmoreDB(AliTPCCalibCE &calibCE, unsigned long32 runNb);
84 //for threaded processing
85
86
87 /* Main routine
88       Arguments: list of DATE raw data files
89 */
90 int main(int argc, char **argv) {
91   /* log start of process */
92   printf("TPCCEda: DA started - %s\n",__FILE__);
93   
94   if (argc<2) {
95     printf("TPCCEda: Wrong number of arguments\n");
96     return -1;
97   }
98   
99   AliLog::SetClassDebugLevel("AliTPCRawStream",-5);
100   AliLog::SetClassDebugLevel("AliRawReaderDate",-5);
101   AliLog::SetClassDebugLevel("AliTPCAltroMapping",-5);
102   AliLog::SetModuleDebugLevel("RAW",-5);
103   
104   gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
105                                         "*",
106                                         "TStreamerInfo",
107                                         "RIO",
108                                         "TStreamerInfo()");
109   
110  /* declare monitoring program */
111   int i,status;
112   status=monitorDeclareMp( __FILE__ );
113   if (status!=0) {
114     printf("TPCCEnewda: monitorDeclareMp() failed : %s\n",monitorDecodeError(status));
115     return -1;
116   }
117   monitorSetNowait();
118   monitorSetNoWaitNetworkTimeout(1000);
119
120   //variables
121   AliTPCmapper *mapping = 0;   // The TPC mapping
122   char localfile[255];
123   unsigned long32 runNb=0;      //run number
124   
125  
126   //
127   // DA configuration from configuration file
128   //
129   //retrieve configuration file
130   sprintf(localfile,"./%s",CONFIG_FILE);
131   status = daqDA_DB_getFile(CONFIG_FILE,localfile);
132   if (status) {
133     printf("TPCCEda: Failed to get configuration file (%s) from DAQdetDB, status=%d\n", CONFIG_FILE, status);
134     return -1;
135   }
136   AliTPCConfigDA config(CONFIG_FILE);
137   // check configuration options
138   TString laserTriggerName("C0LSR-ABCE-NOPF-CENT");
139   TString monitoringType("YES");
140   TString cdb("local:///local/cdb");
141   Int_t   forceTriggerId=-1;
142   Bool_t  skipAmore=kFALSE;
143   
144   if ( config.GetConfigurationMap()->GetValue("LaserTriggerName") ) {
145     laserTriggerName=config.GetConfigurationMap()->GetValue("LaserTriggerName")->GetName();
146     printf("TPCCEnewda: Laser trigger class name set to: %s.\n",laserTriggerName.Data());
147   }
148
149   if ( config.GetConfigurationMap()->GetValue("MonitoringType") ) {
150     monitoringType=config.GetConfigurationMap()->GetValue("MonitoringType")->GetName();
151     printf("TPCCEnewda: Monitoring type set to: %s.\n",monitoringType.Data());
152   }
153
154   if ( config.GetConfigurationMap()->GetValue("ForceLaserTriggerId") ) {
155     forceTriggerId=TMath::Nint(config.GetValue("ForceLaserTriggerId"));
156     printf("TPCCEnewda: Only processing triggers with Id: %d.\n",forceTriggerId);
157   }
158
159   if ( config.GetConfigurationMap()->GetValue("SkipAmore") ) {
160     skipAmore=((TObjString*)config.GetConfigurationMap()->GetValue("SkipAmore"))->GetString().Atoi();
161     printf("TPCCENEWda: Skip Amore set in config\n");
162   }
163
164   //subsribe to laser triggers only in physics partition
165   //if the trigger class is not available the return value is -1
166   //in this case we are most probably running as a standalone
167   //  laser run and should request all events
168   unsigned char classIdptr=0;
169   int retClassId=daqDA_getClassIdFromName(laserTriggerName.Data(),&classIdptr);
170   if (retClassId==0){
171     //interleaved laser in physics runs
172     //select proper trigger class id
173     char c[5];
174     snprintf(c,sizeof(c),"%u",(unsigned int)classIdptr);
175     char *table[5] = {"PHY",const_cast<char*>(monitoringType.Data()),"*",c,NULL};
176     monitorDeclareTableExtended(table);
177     printf("TPCCEnewda: Using monitoring table: (PHY, %s, *, %s)\n",monitoringType.Data(),c);
178   } else if (retClassId==-1){
179     //global partition without laser triggered events
180     //the DA should exit properly without processing
181     printf("TPCCEnewda: Laser trigger class '%s' was not found among trigger class names. Will stop processing.\n",laserTriggerName.Data());
182     return 0;
183   } else if (retClassId==-2){
184     //standalone case, accept all physics events
185     char *table[5] = {"PHY","Y","*","*",NULL};
186     monitorDeclareTableExtended(table);
187     printf("TPCCEnewda: Using all trigger class Ids\n");
188   } else {
189     printf("TPCCEnewda: Unknown return value of 'daqDA_getClassIdFromName': %d\n",retClassId);
190     return -2;
191   }
192
193   //see if we should force the trigger id
194   if (forceTriggerId>-1){
195     char c[5];
196     sprintf(c,"%d",forceTriggerId);
197     char *table[5] = {"PHY","Y","*",c,NULL};
198     monitorDeclareTableExtended(table);
199   }
200   
201   
202   // if  test setup get parameters from $DAQDA_TEST_DIR
203   if (!mapping){
204     /* copy locally the mapping file from daq detector config db */
205     sprintf(localfile,"./%s",MAPPING_FILE);
206     status = daqDA_DB_getFile(MAPPING_FILE,localfile);
207     if (status) {
208       printf("TPCCEnewda: Failed to get mapping file (%s) from DAQdetDB, status=%d\n", MAPPING_FILE, status);
209       return -1;
210     }
211     
212     /* open the mapping file and retrieve mapping object */
213     TFile *fileMapping = new TFile(MAPPING_FILE, "read");
214     mapping = (AliTPCmapper*) fileMapping->Get("tpcMapping");
215     delete fileMapping;
216   }
217   
218   if (mapping == 0) {
219     printf("TPCCEnewda: Failed to get mapping object from %s.  ...\n", MAPPING_FILE);
220     return -1;
221   } else {
222     printf("TPCCEnewda: Got mapping object from %s\n", MAPPING_FILE);
223   }
224   
225   //========== Laser Track data base
226   sprintf(localfile,"./%s",LASER_TRACK_FILE);
227   status = daqDA_DB_getFile(LASER_TRACK_FILE,localfile);
228   if (status) {
229     printf("TPCCEnewda: Failed to get the laser track file (%s) from DAQdetDB, status=%d\n", MAPPING_FILE, status);
230     return -1;
231   }
232
233   
234   //create calibration object with the new algorithm
235   AliTPCCalibCE calibCE(config.GetConfigurationMap());   // central electrode calibration
236   calibCE.SetAltroMapping(mapping->GetAltroMapping()); // Use altro mapping we got from daqDetDb
237   calibCE.SetProcessOld(kFALSE);
238   calibCE.SetProcessNew();
239   
240   //amore update interval
241   Double_t updateInterval=300; //seconds
242   Double_t valConf=config.GetValue("AmoreUpdateInterval");
243   if ( valConf>0 ) updateInterval=valConf;
244   //timer
245   TStopwatch stopWatch;
246   
247   //===========================//
248   // loop over RAW data files //
249   //==========================//
250   int nevents=0;
251   int neventsOld=0;
252   size_t counter=0;
253   for ( i=1; i<argc; i++) {
254     
255     /* define data source : this is argument i */
256     printf("TPCCEnewda: Processing file %s\n", argv[i]);
257     status=monitorSetDataSource( argv[i] );
258     if (status!=0) {
259       printf("TPCCEnewda: monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
260       return -1;
261     }
262
263     
264     /* read until EOF */
265     while (true) {
266       struct eventHeaderStruct *event;
267       
268       /* check shutdown condition */
269       if (daqDA_checkShutdown()) {break;}
270         
271       /* get next event (blocking call until timeout) */
272       status=monitorGetEventDynamic((void **)&event);
273       if (status==MON_ERR_EOF) {
274         printf ("TPCCEnewda: End of File %d detected\n",i);
275         break; /* end of monitoring file has been reached */
276       }
277       
278       if (status!=0) {
279         printf("TPCCEnewda: monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
280         break;
281       }
282       
283         /* retry if got no event */
284       if (event==NULL){
285         //use time in between bursts to
286         // send the data to AMOREdb
287         if (stopWatch.RealTime()>updateInterval){
288           calibCE.Analyse();
289           if (!skipAmore) SendToAmoreDB(calibCE,runNb);
290           stopWatch.Start();
291         } else {
292           stopWatch.Continue();
293         }
294         //debug output
295         if (nevents>neventsOld){
296           printf ("TPCCEnewda: %d events processed, %d used\n",nevents,calibCE.GetNeventsProcessed());
297           neventsOld=nevents;
298         }
299         
300         continue;
301       }
302       
303       /* skip start/end of run events */
304       if ( (event->eventType != physicsEvent) && (event->eventType != calibrationEvent) ){
305         free(event);
306         continue;
307       }
308       
309       
310       // get the run number
311       runNb = event->eventRunNb;
312       
313       // CE calibration
314       calibCE.ProcessEvent(event);
315       
316       /* free resources */
317       free(event);
318       ++nevents;
319     }
320   }
321   
322   //
323   // Analyse CE data and write them to rootfile
324   //
325   calibCE.Analyse();
326   printf ("TPCCEnewda: %d events processed, %d used\n",nevents,calibCE.GetNeventsProcessed());
327   
328 //   TFile * fileTPC = new TFile (RESULT_FILE,"recreate");
329 //   calibCE.Write("tpcCalibCE");
330 //   delete fileTPC;
331   
332   calibCE.DumpToFile(RESULT_FILE,"name=tpcCalibCE,type=3");
333   printf("TPCCEnewda: Wrote %s\n",RESULT_FILE);
334   
335   /* store the result file on FES */
336
337   status=daqDA_FES_storeFile(RESULT_FILE,FILE_ID);
338   if (status) {
339     status = -2;
340   }
341
342   if (!skipAmore) {
343     printf("TPCCEnewda: Amore part %s\n",RESULT_FILE);
344     SendToAmoreDB(calibCE,runNb);
345   }
346   
347   return status;
348 }
349
350
351 void SendToAmoreDB(AliTPCCalibCE &calibCE, unsigned long32 runNb)
352 {
353   //AMORE
354 //   printf ("AMORE part\n");
355   const char *amoreDANameorig=gSystem->Getenv("AMORE_DA_NAME");
356   //cheet a little -- temporary solution (hopefully)
357   //
358   //currently amoreDA uses the environment variable AMORE_DA_NAME to create the mysql
359   //table in which the calib objects are stored. This table is dropped each time AmoreDA
360   //is initialised. This of course makes a problem if we would like to store different
361   //calibration entries in the AMORE DB. Therefore in each DA which writes to the AMORE DB
362   //the AMORE_DA_NAME env variable is overwritten.
363   gSystem->Setenv("AMORE_DA_NAME",Form("TPC-%s",FILE_ID));
364   //
365   // end cheet
366   TGraph *grA=calibCE.MakeGraphTimeCE(-1,0,2);
367   TGraph *grC=calibCE.MakeGraphTimeCE(-2,0,2);
368   TDatime time;
369   TObjString info(Form("Run: %u; Date: %s",runNb,time.AsSQLString()));
370   amore::da::AmoreDA amoreDA(amore::da::AmoreDA::kSender);
371   Int_t statusDA=0;
372 //   statusDA+=amoreDA.Send("CET0",calibCE.GetCalPadT0());
373 //   statusDA+=amoreDA.Send("CEQ",calibCE.GetCalPadQ());
374 //   statusDA+=amoreDA.Send("CERMS",calibCE.GetCalPadRMS());
375 //   statusDA+=amoreDA.Send("DriftA",grA);
376 //   statusDA+=amoreDA.Send("DriftC",grC);
377   statusDA+=amoreDA.Send("FitGraphs",calibCE.GetArrFitGraphs());
378   statusDA+=amoreDA.Send("Info",&info);
379   if ( statusDA!=0 )
380     printf("TPCCEnewda: Waring: Failed to write one of the calib objects to the AMORE database\n");
381   // reset env var
382   if (amoreDANameorig) gSystem->Setenv("AMORE_DA_NAME",amoreDANameorig);
383   if (grA) delete grA;
384   if (grC) delete grC;
385 }
386
387