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