]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/TPCCEda.cxx
correcting baryon mass subtraction for visible energy in MC
[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   
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   
157   //subsribe to laser triggers only in physics partition
158   //if the trigger class is not available the return value is -1
159   //in this case we are most probably running as a standalone
160   //  laser run and should request all events
161   unsigned char classIdptr=0;
162   int retClassId=daqDA_getClassIdFromName(laserTriggerName.Data(),&classIdptr);
163   if (retClassId==0){
164     //interleaved laser in physics runs
165     //select proper trigger class id
166     char c[5];
167     snprintf(c,sizeof(c),"%u",(unsigned int)classIdptr);
168     char *table[5] = {"PHY",const_cast<char*>(monitoringType.Data()),"*",c,NULL};
169     monitorDeclareTableExtended(table);
170     printf("TPCCEda: Using monitoring table: (PHY, %s, *, %s)\n",monitoringType.Data(),c);
171   } else if (retClassId==-1){
172     //global partition without laser triggered events
173     //the DA should exit properly without processing
174     printf("TPCCEda: Laser trigger class '%s' was not found among trigger class names. Will stop processing.\n",laserTriggerName.Data());
175     return 0;
176   } else if (retClassId==-2){
177     //standalone case, accept all physics events
178     char *table[5] = {"PHY","Y","*","*",NULL};
179     monitorDeclareTableExtended(table);
180     printf("TPCCEda: Using all trigger class Ids\n");
181   } else {
182     printf("TPCCEda: Unknown return value of 'daqDA_getClassIdFromName': %d\n",retClassId);
183     return -2;
184   }
185
186   //see if we should force the trigger id
187   if (forceTriggerId>-1){
188     char c[5];
189     sprintf(c,"%d",forceTriggerId);
190     char *table[5] = {"PHY","Y","*",c,NULL};
191     monitorDeclareTableExtended(table);
192   }
193   
194   
195   // if  test setup get parameters from $DAQDA_TEST_DIR
196   if (!mapping){
197     /* copy locally the mapping file from daq detector config db */
198     sprintf(localfile,"./%s",MAPPING_FILE);
199     status = daqDA_DB_getFile(MAPPING_FILE,localfile);
200     if (status) {
201       printf("TPCCEda: Failed to get mapping file (%s) from DAQdetDB, status=%d\n", MAPPING_FILE, status);
202       return -1;
203     }
204     
205     /* open the mapping file and retrieve mapping object */
206     TFile *fileMapping = new TFile(MAPPING_FILE, "read");
207     mapping = (AliTPCmapper*) fileMapping->Get("tpcMapping");
208     delete fileMapping;
209   }
210   
211   if (mapping == 0) {
212     printf("TPCCEda: Failed to get mapping object from %s.  ...\n", MAPPING_FILE);
213     return -1;
214   } else {
215     printf("TPCCEda: Got mapping object from %s\n", MAPPING_FILE);
216   }
217   
218     
219   //create calibration object
220   AliTPCCalibCE calibCE(config.GetConfigurationMap());   // central electrode calibration
221   calibCE.SetAltroMapping(mapping->GetAltroMapping()); // Use altro mapping we got from daqDetDb
222
223   //amore update interval
224   Double_t updateInterval=300; //seconds
225   Double_t valConf=config.GetValue("AmoreUpdateInterval");
226   if ( valConf>0 ) updateInterval=valConf;
227   //timer
228   TStopwatch stopWatch;
229   
230   //===========================//
231   // loop over RAW data files //
232   //==========================//
233   int nevents=0;
234   int neventsOld=0;
235   size_t counter=0;
236   for ( i=1; i<argc; i++) {
237     
238     /* define data source : this is argument i */
239     printf("TPCCEda: Processing file %s\n", argv[i]);
240     status=monitorSetDataSource( argv[i] );
241     if (status!=0) {
242       printf("TPCCEda: monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
243       return -1;
244     }
245
246     
247     /* read until EOF */
248     while (true) {
249       struct eventHeaderStruct *event;
250       
251       /* check shutdown condition */
252       if (daqDA_checkShutdown()) {break;}
253         
254       /* get next event (blocking call until timeout) */
255       status=monitorGetEventDynamic((void **)&event);
256       if (status==MON_ERR_EOF) {
257         printf ("TPCCEda: End of File %d detected\n",i);
258         break; /* end of monitoring file has been reached */
259       }
260       
261       if (status!=0) {
262         printf("TPCCEda: monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
263         break;
264       }
265       
266         /* retry if got no event */
267       if (event==NULL){
268         //use time in between bursts to
269         // send the data to AMOREdb
270         if (stopWatch.RealTime()>updateInterval){
271           calibCE.Analyse();
272           SendToAmoreDB(calibCE,runNb);
273           stopWatch.Start();
274         } else {
275           stopWatch.Continue();
276         }
277         //debug output
278         if (nevents>neventsOld){
279           printf ("TPCCEda: %d events processed, %d used\n",nevents,calibCE.GetNeventsProcessed());
280           neventsOld=nevents;
281         }
282         
283         continue;
284       }
285       
286       /* skip start/end of run events */
287       if ( (event->eventType != physicsEvent) && (event->eventType != calibrationEvent) ){
288         free(event);
289         continue;
290       }
291       
292       
293       // get the run number
294       runNb = event->eventRunNb;
295       
296       // CE calibration
297       calibCE.ProcessEvent(event);
298       
299       /* free resources */
300       free(event);
301       ++nevents;
302     }
303   }
304   
305   //
306   // Analyse CE data and write them to rootfile
307   //
308   calibCE.Analyse();
309   printf ("TPCCEda: %d events processed, %d used\n",nevents,calibCE.GetNeventsProcessed());
310   
311   TFile * fileTPC = new TFile (RESULT_FILE,"recreate");
312   calibCE.Write("tpcCalibCE");
313   delete fileTPC;
314   printf("TPCCEda: Wrote %s\n",RESULT_FILE);
315   
316   /* store the result file on FES */
317   
318   status=daqDA_FES_storeFile(RESULT_FILE,FILE_ID);
319   if (status) {
320     status = -2;
321   }
322   
323   SendToAmoreDB(calibCE,runNb);
324   
325   return status;
326 }
327
328
329 void SendToAmoreDB(AliTPCCalibCE &calibCE, unsigned long32 runNb)
330 {
331   //AMORE
332 //   printf ("AMORE part\n");
333   const char *amoreDANameorig=gSystem->Getenv("AMORE_DA_NAME");
334   //cheet a little -- temporary solution (hopefully)
335   //
336   //currently amoreDA uses the environment variable AMORE_DA_NAME to create the mysql
337   //table in which the calib objects are stored. This table is dropped each time AmoreDA
338   //is initialised. This of course makes a problem if we would like to store different
339   //calibration entries in the AMORE DB. Therefore in each DA which writes to the AMORE DB
340   //the AMORE_DA_NAME env variable is overwritten.
341   gSystem->Setenv("AMORE_DA_NAME",Form("TPC-%s",FILE_ID));
342   //
343   // end cheet
344   TGraph *grA=calibCE.MakeGraphTimeCE(-1,0,2);
345   TGraph *grC=calibCE.MakeGraphTimeCE(-2,0,2);
346   TDatime time;
347   TObjString info(Form("Run: %u; Date: %s",runNb,time.AsSQLString()));
348   amore::da::AmoreDA amoreDA(amore::da::AmoreDA::kSender);
349   Int_t statusDA=0;
350   statusDA+=amoreDA.Send("CET0",calibCE.GetCalPadT0());
351   statusDA+=amoreDA.Send("CEQ",calibCE.GetCalPadQ());
352   statusDA+=amoreDA.Send("CERMS",calibCE.GetCalPadRMS());
353   statusDA+=amoreDA.Send("DriftA",grA);
354   statusDA+=amoreDA.Send("DriftC",grC);
355   statusDA+=amoreDA.Send("Info",&info);
356   if ( statusDA!=0 )
357     printf("TPCCEda: Waring: Failed to write one of the calib objects to the AMORE database\n");
358   // reset env var
359   if (amoreDANameorig) gSystem->Setenv("AMORE_DA_NAME",amoreDANameorig);
360   if (grA) delete grA;
361   if (grC) delete grC;
362 }
363
364