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