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