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