]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/TOFda.cxx
AliHLTDataTypes.h included
[u/mrichter/AliRoot.git] / TOF / TOFda.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 <AliTOFNoiseConfigHandler.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("TOFNoiseConfig.xml","TOFNoiseConfig.xml");
86   if (getConfigFile != 0){
87     printf("Failed to retrieve config file from DB! returning...\n");
88     return -1;
89   }
90
91   AliTOFNoiseConfigHandler* tofHandler = new AliTOFNoiseConfigHandler();
92   TSAXParser *parser = new TSAXParser();
93   parser->ConnectToHandler("AliTOFNoiseConfigHandler", tofHandler);  
94   if (parser->ParseFile("./TOFNoiseConfig.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
102   /* init some counters */
103   int nevents_physics=0;
104   int nevents_total=0;
105
106   Int_t iev=0;
107
108   Int_t nPDBEntriesToT = 0;
109   Int_t nDBEntriesToT = 0;
110   AliTOFHitData *HitData;
111   Int_t dummy = -1;
112   Int_t Volume[5];
113   for (Int_t i=0;i<5;i++) Volume[i]=-1;
114   AliTOFRawStream *rawStreamTOF = new AliTOFRawStream();
115   AliTOFHitDataBuffer DataBuffer;
116   AliTOFHitDataBuffer PackedDataBuffer;
117   Int_t nDBEntries = 0;
118   Int_t nPDBEntries = 0;
119   
120   struct eventHeaderStruct *event;
121   eventTypeType eventT;
122
123   /* define data source : this is argument 1 */  
124   status=monitorSetDataSource( argv[1] );
125   if (status!=0) {
126     printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
127     return -1;
128   }
129
130   /* declare monitoring program */
131   status=monitorDeclareMp( __FILE__ );
132   if (status!=0) {
133     printf("monitorDeclareMp() failed : %s\n",monitorDecodeError(status));
134     return -1;
135   }
136
137   /* define wait event timeout - 1s max */
138   monitorSetNowait();
139   monitorSetNoWaitNetworkTimeout(1000);
140   
141   /* main loop (infinite) */
142   for(;;) {
143     
144     /* check shutdown condition */
145     if (daqDA_checkShutdown()) break;
146     
147     /* get next event (blocking call until timeout) */
148     status=monitorGetEventDynamic((void **)&event);
149     if (status==MON_ERR_EOF) {
150       printf ("End of File detected\n");
151       break; /* end of monitoring file has been reached */
152     }
153     
154     if (status!=0) {
155       printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
156       break;
157     }
158     
159     /* retry if got no event */
160     if (event==NULL) continue;
161
162     iev++; 
163
164    /* use event - here, just write event id to result file */
165     nevents_total++;
166     eventT=event->eventType;
167     switch (event->eventType) {
168       
169       /* START OF RUN */
170     case START_OF_RUN:
171       break;
172       /* END START OF RUN */
173       
174       /* END OF RUN */
175     case END_OF_RUN:
176       break;
177       /* END END OF RUN */
178       
179     case PHYSICS_EVENT:
180       nevents_physics++;
181       AliRawReader *rawReader = new AliRawReaderDate((void*)event);
182       //rawReader->RequireHeader(kFALSE);
183
184       //T0 event
185       Int_t meantime = 0;     
186       AliT0RawReader *rawReaderT0 = new AliT0RawReader(rawReader,kTRUE);
187       if (!rawReaderT0->Next()) {
188         printf("T0: no raw data found!\n");
189       } 
190       else {
191         /*
192         Int_t allData[105][5];
193         for (Int_t i=0; i<105; i++) {
194           allData[i][0]=rawReaderT0->GetData(i,0);
195         }
196         meantime = allData[49][0];
197         */
198         //meantime = rawReaderT0->GetData(49,0); //OLD
199         meantime = (Int_t)((rawReaderT0->GetData(51,0)+rawReaderT0->GetData(52,0))/2.); //Alla
200         //        printf("time zero (ns) = %i (%f) \n", meantime, (meantime*24.4-200)*1E-3);   // debugging purpose
201       }
202       
203       delete rawReaderT0;
204       rawReaderT0 = 0x0;
205       rawReader->Reset();
206       
207       //TOF event
208       dummy = -1;
209       for (Int_t ii=0; ii<5; ii++) Volume[ii]=-1;
210       rawStreamTOF->SetRawReader(rawReader);
211       //rawReader->ReadHeader();
212       rawStreamTOF->ResetBuffers();
213       rawStreamTOF->DecodeDDL(0, AliDAQ::NumberOfDdls("TOF") - 1,0);
214       nPDBEntriesToT = 0;
215       nDBEntriesToT = 0;
216       for (Int_t iDDL = 0; iDDL < AliDAQ::NumberOfDdls("TOF"); iDDL++) {
217
218         /* read decoded data */
219         DataBuffer = rawStreamTOF->GetDataBuffer(iDDL);
220         PackedDataBuffer = rawStreamTOF->GetPackedDataBuffer(iDDL);
221         
222         /* get buffer entries */
223         nDBEntries = DataBuffer.GetEntries();
224         nPDBEntries = PackedDataBuffer.GetEntries();
225         nPDBEntriesToT+=nPDBEntries;
226         nDBEntriesToT+=nDBEntries;
227
228         //for (Int_t iHit = 0; iHit < nDBEntries; iHit++) {
229         // HitData = DataBuffer->GetHit(iHit);
230           /* store volume information */
231         // rawStreamTOF->EquipmentId2VolumeId(HitData, Volume);
232         //}
233         /* reset buffer */
234         DataBuffer.Reset();
235
236         /* read data buffer hits */
237         for (Int_t iHit = 0; iHit < nPDBEntries; iHit++) {
238           HitData = PackedDataBuffer.GetHit(iHit);
239           /* add volume information */
240           HitData->SetDDLID(iDDL);
241           rawStreamTOF->EquipmentId2VolumeId(HitData, Volume);
242           if (Volume[0]==-1 ||
243               Volume[1]==-1 ||
244               Volume[2]==-1 ||
245               Volume[3]==-1 ||
246               Volume[4]==-1) continue;
247           else {
248             dummy = Volume[3];
249             Volume[3] = Volume[4];
250             Volume[4] = dummy;
251             Int_t tof = (Int_t)((Double_t)HitData->GetTime()*1E3/AliTOFGeometry::TdcBinWidth());
252             Int_t index = geom->GetIndex(Volume);
253             Float_t pos[3];
254             geom->GetPosPar(Volume,pos);
255             Float_t texp = TMath::Sqrt(pos[0]*pos[0]+pos[1]*pos[1]+pos[2]*pos[2])/c*1E9; //expected time in ns
256             Float_t texpBin = texp*1E3/AliTOFGeometry::TdcBinWidth(); //expected time in number of TDC bin
257             Int_t deltabin = tof-TMath::Nint(texpBin);   //to be used with real data; rounding expected time to Int_t
258             htofPartial->Fill(index,deltabin); //channel index start from 0, bin index from 1
259             //debugging printings
260             //if (debugFlag) {
261             //  printf("sector %2d, plate %1d, strip %2d, padz %1d, padx %2d \n",Volume[0],Volume[1],Volume[2],Volume[3],Volume[4]); // too verbose
262             //  printf("pos x = %f, pos y = %f, pos z = %f \n",pos[0],pos[1],pos[2]); // too verbose
263             //  printf("expected time = %f (ns)\n",texp); // too verbose
264             //  printf("expected time bin = %f (TDC bin)\n",texpBin); // too verbose
265             //  printf("measured time bin = %i (TDC bin) with %f (ns) and ACQ bit = %i \n",tof, HitData->GetTime(), HitData->GetACQ()); // too verbose
266             //  printf("index = %6d, deltabin = %d , filling index = %6d, and bin = %d\n",index, deltabin, index, deltabin); // too verbose
267             //}
268
269           }
270           /* reset buffer */
271           PackedDataBuffer.Reset();
272         }
273       }
274       //if (debugFlag) {
275       //  printf(" Packed Hit Buffer Entries = %i \n",nPDBEntriesToT); // too verbose
276       //  printf(" Hit Buffer Entries = %i \n",nDBEntriesToT); // too verbose
277       //}
278
279       delete rawReader;
280       rawReader = 0x0;
281     }
282     
283     /* free resources */
284     free(event);
285     
286     /* exit when last event received, no need to wait for TERM signal */
287     if (eventT==END_OF_RUN) {
288       printf("EOR event detected\n");
289       break;
290     }
291
292   }
293
294   delete rawStreamTOF;
295   rawStreamTOF = 0x0;
296
297   delete geom;
298   geom = 0x0;
299
300   //write the Run level file   
301   TFile * fileRun = new TFile (FILE_RUN,"RECREATE"); 
302   htofPartial->Write();
303   fileRun->Close();
304
305   //write the Total file
306   TH2S *htoftot = 0x0;
307   TFile * filetot = 0x0;
308   Bool_t isThere=kFALSE;
309   const char *dirname = "./";
310   TString filename = FILE_TOTAL;
311   if((gSystem->FindFile(dirname,filename))!=NULL){
312     isThere=kTRUE;
313     printf("%s found \n",FILE_TOTAL);
314   }
315   if (isThere) {
316
317     TFile * filetot1 = new TFile (FILE_TOTAL,"READ"); 
318     //look for the file
319     if (!filetot1->IsZombie()){
320       printf("updating file %s \n",FILE_TOTAL);
321       TIter next(filetot1->GetListOfKeys());
322       TKey *key;
323       //look for the histogram
324       while ((key=(TKey*)next())){
325         const char * namekey = key->GetName();
326         if (strcmp(namekey,"htoftot")==0) {
327           printf(" histo found \n");
328           htoftot = (TH2S*) filetot1->Get("htoftot");
329           htoftot->AddDirectory(0);
330           htoftot->Add(htofPartial);
331           break;
332         }
333       }
334     }
335     filetot1->Close();
336     delete filetot1;
337     filetot1=0x0;
338   }
339   else {
340     printf(" no %s file found \n",FILE_TOTAL);
341     htoftot = new TH2S(*htofPartial);
342     htoftot->SetName("htoftot");
343     htoftot->AddDirectory(0);
344   }
345
346   filetot = new TFile (FILE_TOTAL,"RECREATE");
347   filetot->cd();
348   htoftot->Write();
349   filetot->Close();
350   
351   delete fileRun;
352   delete filetot;
353   delete htofPartial;
354   delete htoftot;
355
356   fileRun = 0x0;
357   filetot = 0x0;
358   htofPartial = 0x0;
359   htoftot = 0x0;
360
361   /* write report */
362   printf("Run #%s, received %d physics events out of %d\n",
363          getenv("DATE_RUN_NUMBER"),nevents_physics,nevents_total);
364
365   status = 0;
366
367   /* export file to FXS */
368   if (daqDA_FES_storeFile(FILE_RUN, "RUNLevel"))
369     status=-2;
370   if (daqDA_FES_storeFile(FILE_TOTAL, "DELAYS"))
371     status=-2;
372
373   return status;
374 }