]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/TOFda.old.cxx
FindROOT + dictionary generating macro
[u/mrichter/AliRoot.git] / TOF / TOFda.old.cxx
1 /*
2
3 TOF DA for online calibration
4
5 Contact: Chiara.Zampolli@bo.infn.it
6 Link: www.bo.infn.it/~zampolli
7 Run Type: PHYSICS
8 DA Type: MON
9 Number of events needed: depending on the run, being run-level
10 Input Files: TOFdaTotal.root, to be updated if existing
11 Output Files: TOFdaRun.root, TOFdaTotal.root, both to be exported to the DAQ FXS
12 Trigger types used: PHYSICS_EVENT
13
14 */
15
16 #define FILE_TOTAL "TOFdaTotal.root"
17 #define FILE_RUN "TOFdaRun.root"
18
19 // DATE
20 #include "event.h"
21 #include "monitor.h"
22 #include "daqDA.h"
23
24 #include <stdio.h>
25 #include <stdlib.h>
26
27 //AliRoot
28 #include <AliTOFRawStream.h>
29 #include <AliRawReaderDate.h>
30 #include <AliRawReader.h>
31 #include <AliTOFGeometry.h>
32 #include <AliT0RawReader.h>
33 #include <AliDAQ.h>
34 #include <AliTOFHitData.h>
35 #include <AliTOFHitDataBuffer.h>
36 #include <AliTOFDaConfigHandler.h>
37
38 //ROOT
39 #include <TFile.h>
40 #include <TKey.h>
41 #include <TH2S.h>
42 #include <TObject.h>
43 #include <TMath.h>
44 #include <TSystem.h>
45 #include "TROOT.h"
46 #include "TPluginManager.h"
47 #include "TSAXParser.h"
48
49 /* Main routine
50       Arguments: 
51       1- monitoring data source
52 */
53 int main(int argc, char **argv) {
54
55   /* magic line from Rene */
56   gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
57                     "*",
58                     "TStreamerInfo",
59                     "RIO",
60                     "TStreamerInfo()");
61
62   AliTOFGeometry * geom = new AliTOFGeometry();
63
64   static const Int_t size = AliTOFGeometry::NPadXSector()*AliTOFGeometry::NSectors();
65   static const Int_t nbins = 500;
66   static const Int_t binmin = -20;
67   const Float_t c = 2.99792458E10; //speed of light [cm/s]
68   TH1F::AddDirectory(0);
69   TH2S * htofPartial = new TH2S("htof","histo with delays",
70                                 size,-0.5,size*1.-0.5,
71                                 nbins,binmin-0.5,nbins*1.+binmin-0.5);
72   
73   int status;
74
75   /* log start of process */
76   printf("TOF DA started\n");  
77
78   /* check that we got some arguments = list of files */
79   if (argc!=2) {
80     printf("Wrong number of arguments\n");
81     return -1;
82   }
83
84   /* retrieve config file */
85   int getConfigFile = daqDA_DB_getFile("TOFPhysicsConfig.xml","TOFPhysicsConfig.xml");
86   if (getConfigFile != 0){
87     printf("Failed to retrieve config file from DB! returning...\n");
88     return -1;
89   }
90
91   AliTOFDaConfigHandler* tofHandler = new AliTOFDaConfigHandler();
92   TSAXParser *parser = new TSAXParser();
93   parser->ConnectToHandler("AliTOFDaConfigHandler", tofHandler);  
94   if (parser->ParseFile("./TOFPhysicsConfig.xml") != 0) {
95     printf("Failed parsing config file! retunring... \n");
96     return -1;
97   }
98
99   Int_t debugFlag = tofHandler->GetDebugFlag();
100   printf("The debug flag is %i\n",debugFlag);
101   Int_t t0Flag = tofHandler->GetT0Flag();
102   printf("The T0 flag is %i\n\n",t0Flag);
103   if (t0Flag) {
104           printf("The T0 time will be subtracted from the measured TOF time. So, in TDC bins: \n");
105           printf("tof = tofRaw - (rawReaderT0->GetData(51,0)+rawReaderT0->GetData(52,0))/2) \n\n");
106   }
107   else {
108           printf("The T0 time will not be used.\n");
109           printf("tof = tofRaw \n\n");
110   }
111
112   /* init some counters */
113   int nevents_physics=0;
114   int nevents_total=0;
115
116   Int_t iev=0;
117
118   Int_t nPDBEntriesToT = 0;
119   Int_t nDBEntriesToT = 0;
120   AliTOFHitData *HitData;
121   Int_t dummy = -1;
122   Int_t Volume[5];
123   for (Int_t i=0;i<5;i++) Volume[i]=-1;
124   AliTOFRawStream *rawStreamTOF = new AliTOFRawStream();
125   AliTOFHitDataBuffer *DataBuffer;
126   AliTOFHitDataBuffer *PackedDataBuffer;
127   Int_t nDBEntries = 0;
128   Int_t nPDBEntries = 0;
129   
130   struct eventHeaderStruct *event;
131   eventTypeType eventT;
132
133   /* define data source : this is argument 1 */  
134   status=monitorSetDataSource( argv[1] );
135   if (status!=0) {
136     printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
137     return -1;
138   }
139
140   /* declare monitoring program */
141   status=monitorDeclareMp( __FILE__ );
142   if (status!=0) {
143     printf("monitorDeclareMp() failed : %s\n",monitorDecodeError(status));
144     return -1;
145   }
146
147   /* define wait event timeout - 1s max */
148   monitorSetNowait();
149   monitorSetNoWaitNetworkTimeout(1000);
150   
151   /* main loop (infinite) */
152   for(;;) {
153     
154     /* check shutdown condition */
155     if (daqDA_checkShutdown()) break;
156     
157     /* get next event (blocking call until timeout) */
158     status=monitorGetEventDynamic((void **)&event);
159     if (status==MON_ERR_EOF) {
160       printf ("End of File detected\n");
161       break; /* end of monitoring file has been reached */
162     }
163     
164     if (status!=0) {
165       printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
166       break;
167     }
168     
169     /* retry if got no event */
170     if (event==NULL) continue;
171
172     iev++; 
173
174    /* use event - here, just write event id to result file */
175     nevents_total++;
176     eventT=event->eventType;
177     switch (event->eventType) {
178       
179       /* START OF RUN */
180     case START_OF_RUN:
181       break;
182       /* END START OF RUN */
183       
184       /* END OF RUN */
185     case END_OF_RUN:
186       break;
187       /* END END OF RUN */
188       
189     case PHYSICS_EVENT:
190       nevents_physics++;
191       AliRawReader *rawReader = new AliRawReaderDate((void*)event);
192       //rawReader->RequireHeader(kFALSE);
193
194       //T0 event
195       Int_t meantime = 0;     
196       AliT0RawReader *rawReaderT0 = new AliT0RawReader(rawReader,kTRUE);
197       if (!rawReaderT0->Next()) {
198         printf("T0: no raw data found!\n");
199       } 
200       else {
201         /*
202         Int_t allData[105][5];
203         for (Int_t i=0; i<105; i++) {
204           allData[i][0]=rawReaderT0->GetData(i,0);
205         }
206         meantime = allData[49][0];
207         */
208         //meantime = rawReaderT0->GetData(49,0); //OLD
209         meantime = (Int_t)((rawReaderT0->GetData(51,0)+rawReaderT0->GetData(52,0))/2.); //Alla
210         if (debugFlag > 0) {
211                 printf("\nT0 for the current event:\n");   // debugging purpose
212                 printf("time zero = %i (TDC bin)\n", meantime);   // debugging purpose
213                 printf("time zero = %f (ns)\n\n", (Float_t)(meantime)*24.4*1E-3);   // debugging purpose
214         }
215       }
216       
217       delete rawReaderT0;
218       rawReaderT0 = 0x0;
219       rawReader->Reset();
220       
221       //TOF event
222       dummy = -1;
223       for (Int_t ii=0; ii<5; ii++) Volume[ii]=-1;
224       rawStreamTOF->SetRawReader(rawReader);
225       //rawReader->ReadHeader();
226       rawStreamTOF->ResetBuffers();
227       rawStreamTOF->DecodeDDL(0, AliDAQ::NumberOfDdls("TOF") - 1,0);
228       nPDBEntriesToT = 0;
229       nDBEntriesToT = 0;
230       for (Int_t iDDL = 0; iDDL < AliDAQ::NumberOfDdls("TOF"); iDDL++) {
231
232         /* read decoded data */
233         DataBuffer = rawStreamTOF->GetDataBuffer(iDDL);
234         PackedDataBuffer = rawStreamTOF->GetPackedDataBuffer(iDDL);
235         
236         /* get buffer entries */
237         nDBEntries = DataBuffer->GetEntries();
238         nPDBEntries = PackedDataBuffer->GetEntries();
239         nPDBEntriesToT+=nPDBEntries;
240         nDBEntriesToT+=nDBEntries;
241
242         //for (Int_t iHit = 0; iHit < nDBEntries; iHit++) {
243         // HitData = DataBuffer->GetHit(iHit);
244           /* store volume information */
245         // rawStreamTOF->EquipmentId2VolumeId(HitData, Volume);
246         //}
247         /* reset buffer */
248         DataBuffer->Reset();
249
250         /* read data buffer hits */
251         for (Int_t iHit = 0; iHit < nPDBEntries; iHit++) {
252           HitData = PackedDataBuffer->GetHit(iHit);
253           /* add volume information */
254           HitData->SetDDLID(iDDL);
255           rawStreamTOF->EquipmentId2VolumeId(HitData, Volume);
256           if (Volume[0]==-1 ||
257               Volume[1]==-1 ||
258               Volume[2]==-1 ||
259               Volume[3]==-1 ||
260               Volume[4]==-1) continue;
261           else {
262             dummy = Volume[3];
263             Volume[3] = Volume[4];
264             Volume[4] = dummy;
265             Int_t tofRaw = (Int_t)((Double_t)HitData->GetTime()*1E3/AliTOFGeometry::TdcBinWidth());
266             Int_t tof;
267             if (!t0Flag) tof = tofRaw;
268             else tof = tofRaw - meantime;
269             Int_t index = geom->GetIndex(Volume);
270             Float_t pos[3];
271             geom->GetPosPar(Volume,pos);
272             Float_t texp = TMath::Sqrt(pos[0]*pos[0]+pos[1]*pos[1]+pos[2]*pos[2])/c*1E9; //expected time in ns
273             Float_t texpBin = texp*1E3/AliTOFGeometry::TdcBinWidth(); //expected time in number of TDC bin
274             Int_t deltabin = tof-TMath::Nint(texpBin);   //to be used with real data; rounding expected time to Int_t
275             htofPartial->Fill(index,deltabin); //channel index start from 0, bin index from 1
276             //debugging printings
277             if (debugFlag > 1) {
278                     printf("tofRaw = %i, tof = %i \n",tofRaw,tof);
279             }
280             if (debugFlag > 2) {
281                     printf("sector %2d, plate %1d, strip %2d, padz %1d, padx %2d \n",Volume[0],Volume[1],Volume[2],Volume[3],Volume[4]); // too verbose
282                     printf("pos x = %f, pos y = %f, pos z = %f \n",pos[0],pos[1],pos[2]); // too verbose
283                     printf("expected time = %f (ns)\n",texp); // too verbose
284                     printf("expected time bin = %f (TDC bin)\n",texpBin); // too verbose
285                     printf("measured time bin = %i (TDC bin) with %f (ns) and ACQ bit = %i \n",tof, HitData->GetTime(), HitData->GetACQ()); // too verbose
286                     printf("index = %6d, deltabin = %d , filling index = %6d, and bin = %d\n",index, deltabin, index, deltabin); // too verbose
287             }
288
289           }
290           /* reset buffer */
291           PackedDataBuffer->Reset();
292         }
293       }
294       //if (debugFlag) {
295       //  printf(" Packed Hit Buffer Entries = %i \n",nPDBEntriesToT); // too verbose
296       //  printf(" Hit Buffer Entries = %i \n",nDBEntriesToT); // too verbose
297       //}
298
299       delete rawReader;
300       rawReader = 0x0;
301     }
302     
303     /* free resources */
304     free(event);
305     
306     /* exit when last event received, no need to wait for TERM signal */
307     if (eventT==END_OF_RUN) {
308       printf("EOR event detected\n");
309       break;
310     }
311
312   }
313
314   delete rawStreamTOF;
315   rawStreamTOF = 0x0;
316
317   delete geom;
318   geom = 0x0;
319
320   //write the Run level file   
321   TFile * fileRun = new TFile (FILE_RUN,"RECREATE"); 
322   htofPartial->Write();
323   fileRun->Close();
324
325   //write the Total file
326   TH2S *htoftot = 0x0;
327   TFile * filetot = 0x0;
328   Bool_t isThere=kFALSE;
329   const char *dirname = "./";
330   TString filename = FILE_TOTAL;
331   if((gSystem->FindFile(dirname,filename))!=NULL){
332     isThere=kTRUE;
333     printf("%s found \n",FILE_TOTAL);
334   }
335   if (isThere) {
336
337     TFile * filetot1 = new TFile (FILE_TOTAL,"READ"); 
338     //look for the file
339     if (!filetot1->IsZombie()){
340       printf("updating file %s \n",FILE_TOTAL);
341       TIter next(filetot1->GetListOfKeys());
342       TKey *key;
343       //look for the histogram
344       while ((key=(TKey*)next())){
345         const char * namekey = key->GetName();
346         if (strcmp(namekey,"htoftot")==0) {
347           printf(" histo found \n");
348           htoftot = (TH2S*) filetot1->Get("htoftot");
349           htoftot->AddDirectory(0);
350           htoftot->Add(htofPartial);
351           break;
352         }
353       }
354     }
355     filetot1->Close();
356     delete filetot1;
357     filetot1=0x0;
358   }
359   else {
360     printf(" no %s file found \n",FILE_TOTAL);
361     htoftot = new TH2S(*htofPartial);
362     htoftot->SetName("htoftot");
363     htoftot->AddDirectory(0);
364   }
365
366   filetot = new TFile (FILE_TOTAL,"RECREATE");
367   filetot->cd();
368   htoftot->Write();
369   filetot->Close();
370   
371   delete fileRun;
372   delete filetot;
373   delete htofPartial;
374   delete htoftot;
375
376   fileRun = 0x0;
377   filetot = 0x0;
378   htofPartial = 0x0;
379   htoftot = 0x0;
380
381   /* write report */
382   printf("Run #%s, received %d physics events out of %d\n",
383          getenv("DATE_RUN_NUMBER"),nevents_physics,nevents_total);
384
385   status = 0;
386
387   /* export file to FXS */
388   if (daqDA_FES_storeFile(FILE_RUN, "RUNLevel"))
389     status=-2;
390   if (daqDA_FES_storeFile(FILE_TOTAL, "DELAYS"))
391     status=-2;
392
393   return status;
394 }