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