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