]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/TOFda.cxx
Very preliminary version of TOF-T0 online calibration algorithm. TO be used as an...
[u/mrichter/AliRoot.git] / TOF / TOFda.cxx
1 /*
2
3 DAcase2.c
4
5 This program connects to the DAQ data source passed as argument
6 and populates local "./result.txt" file with the ids of events received
7 during the run.
8
9 The program exits when being asked to shut down (daqDA_checkshutdown)
10 or End of Run event.
11
12 Messages on stdout are exported to DAQ log system.
13
14 contact: alice-datesupport@cern.ch
15
16 */
17 extern "C" {
18 #include <daqDA.h>
19 }
20
21 #include "event.h"
22 #include "monitor.h"
23 //#include "daqDA.h"
24
25 #include <Riostream.h>
26 #include <stdio.h>
27 #include <stdlib.h>
28
29 //AliRoot
30 #include "AliTOFRawStream.h"
31 #include "AliRawReaderDate.h"
32 #include "AliTOFGeometry.h"
33 #include "AliTOFGeometryV5.h"
34 #include "AliT0digit.h"
35 #include "AliT0RawReader.h"
36
37 //ROOT
38 #include "TFile.h"
39 #include "TKey.h"
40 #include "TH2S.h"
41 #include "TObject.h"
42 #include "TBenchmark.h"
43 #include "TRandom.h"
44
45
46 /* Main routine
47       Arguments: 
48       1- monitoring data source
49 */
50 int main(int argc, char **argv) {
51   
52   AliTOFGeometry * geomV5 = new AliTOFGeometryV5();
53   AliTOFGeometry * geom = new AliTOFGeometry();
54
55   static const Int_t size = geomV5->NPadXSector()*geomV5->NSectors();
56   static const Int_t nbins = 500;
57   static const Int_t binmin = -20;
58   const Float_t c = 2.99792458E10; //speed of light
59   TH1F::AddDirectory(0);
60   TH2S * htofPartial = new TH2S("htof","histo with delays", size,-0.5,size*1.-0.5,nbins,binmin-0.5,nbins*1.+binmin-0.5);
61   
62   TRandom *rand = new TRandom();  //to be used for testing with cosmic data
63   rand->SetSeed(0);               //to be used for testing with cosmic data
64
65   //TTree *tree = 0x0;  // tree for T0 decoder
66       
67   // decoding the events
68   
69   int status;
70   bool decodeStatus;
71   if (argc!=2) {
72     printf("Wrong number of arguments\n");
73     return -1;
74   }
75
76   /* open result file */
77   FILE *fp=NULL;
78   fp=fopen("./result.txt","a");
79   if (fp==NULL) {
80     printf("Failed to open file\n");
81     return -1;
82   }
83
84   /* define data source : this is argument 1 */  
85   status=monitorSetDataSource( argv[1] );
86   if (status!=0) {
87     printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
88     return -1;
89   }
90
91   /* declare monitoring program */
92   status=monitorDeclareMp( __FILE__ );
93   if (status!=0) {
94     printf("monitorDeclareMp() failed : %s\n",monitorDecodeError(status));
95     return -1;
96   }
97
98   /* define wait event timeout - 1s max */
99   monitorSetNowait();
100   monitorSetNoWaitNetworkTimeout(1000);
101   
102   /* log start of process */
103   printf("DA example case2 monitoring program started\n");  
104
105   /* init some counters */
106   int nevents_physics=0;
107   int nevents_total=0;
108
109   struct equipmentStruct *equipment;
110   int *eventEnd;
111   int *eventData;
112   int *equipmentEnd;
113   int *equipmentData;
114   int *equipmentID;
115   struct eventHeaderStruct *event;
116   eventTypeType eventT;
117   Int_t iev=0;
118
119   /* main loop (infinite) */
120   for(;;) {
121     
122     /* check shutdown condition */
123     if (daqDA_checkShutdown()) {break;}
124     
125     /* get next event (blocking call until timeout) */
126     status=monitorGetEventDynamic((void **)&event);
127     if (status==MON_ERR_EOF) {
128       printf ("End of File detected\n");
129       break; /* end of monitoring file has been reached */
130     }
131     
132     if (status!=0) {
133       printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
134       break;
135     }
136     
137     /* retry if got no event */
138     if (event==NULL) {
139       continue;
140     }
141
142     iev++; 
143
144    /* use event - here, just write event id to result file */
145     eventT=event->eventType;
146     switch (event->eventType){
147       
148       /* START OF RUN */
149     case START_OF_RUN:
150       break;
151       /* END START OF RUN */
152       
153     /* END OF RUN */
154     case END_OF_RUN:
155       break;
156       
157     case PHYSICS_EVENT:
158       printf(" event number = %i \n",iev);
159       //if (iev%10000 ==0) printf(" event number = %i \n",iev);
160       //if (iev > 50000) break;
161       AliRawReader *rawReader = new AliRawReaderDate((void*)event);
162       //      rawReader->RequireHeader(kFALSE);
163       //rawReader->LoadEquipmentIdsMap("TOFmap.txt");  //to be used if no exact mapping implemented 
164       Int_t meantime = 0;     
165
166             
167       //      TTree *tree = new TTree();  // tree for T0 decoder
168       AliT0RawReader *rawReaderT0 = new AliT0RawReader(rawReader);
169       //      printf("rawReaderT0 = %p \n", rawReaderT0);
170       //printf("rawReader = %p \n", rawReader);
171       //printf("event = %p \n", event);
172       //      AliT0digit *digit = 0x0;
173       //tree->SetBranchAddress("T0",&digit);
174       
175       if (!rawReaderT0->Next())
176        printf(" no raw data found!! %i", rawReaderT0->Next());
177       Int_t allData[110][5];
178       for (Int_t i=0; i<110; i++) {
179         allData[i][0]=rawReaderT0->GetData(i,0);
180       }
181       meantime = allData[49][0];
182       printf("time zero (ns) = %i (%f) \n", meantime, meantime*25*1E-3-200);
183       
184       delete rawReaderT0;
185       rawReaderT0 = 0x0;
186       rawReader->Reset();
187       
188       AliTOFRawStream *rawStreamTOF = new AliTOFRawStream(rawReader);
189       
190       //loop the event data               
191       Int_t detectorIndex[5] = {-1, -1, -1, -1, -1};
192       
193       printf("before T0 \n");
194       
195       while(rawStreamTOF->Next()){
196         for (Int_t ii=0; ii<5; ii++) detectorIndex[ii] = -1;
197         
198         detectorIndex[0] = (Int_t)rawStreamTOF->GetSector();
199         detectorIndex[1] = (Int_t)rawStreamTOF->GetPlate();
200         detectorIndex[2] = (Int_t)rawStreamTOF->GetStrip();
201         detectorIndex[3] = (Int_t)rawStreamTOF->GetPadZ();
202         detectorIndex[4] = (Int_t)rawStreamTOF->GetPadX();
203         
204         if (detectorIndex[0]==-1 ||
205             detectorIndex[1]==-1 ||
206             detectorIndex[2]==-1 ||
207             detectorIndex[3]==-1 ||
208             detectorIndex[4]==-1) continue;
209         else {
210           Float_t tdc = (Float_t)rawStreamTOF->GetTDC();
211           Int_t tof = (Int_t)rawStreamTOF->GetTofBin();
212           Int_t detID[5]={0,0,0,0,0};
213           detectorIndex[0] = rawStreamTOF->GetSector(); 
214           detectorIndex[1] = rawStreamTOF->GetPlate();
215           detectorIndex[2] = rawStreamTOF->GetStrip(); 
216           detectorIndex[3] = rawStreamTOF->GetPadZ(); 
217           detectorIndex[4] = rawStreamTOF->GetPadX(); 
218           Int_t index = rawStreamTOF->GetIndex(detectorIndex);
219           Float_t pos[3];
220           geomV5->GetPosPar(detectorIndex,pos);
221           Float_t texp=TMath::Sqrt(pos[0]*pos[0]+pos[1]*pos[1]+pos[2]*pos[2])/c*1E9; //expected time in ns
222           Float_t texpBin=(texp*1E3-32)/AliTOFGeometry::TdcBinWidth(); //expected time in number of TDC bin
223           //to be uded with cosmics
224           //Float_t tsim = (Float_t)rand->Gaus(texp,1E-1); //TDC measured time in ns, simulated to be used for testing with cosmic data  
225           //Float_t delta = tsim-texp;  
226           //Int_t deltabin = (Int_t)((delta*1E3-32)/AliTOFGeometry::TdcBinWidth());    //to be used for testing with cosmics data
227           Int_t deltabin = tof-TMath::Nint(texpBin);   //to be used with real data; rounding expected time to Int_t
228           htofPartial->Fill(index,deltabin); //channel index start from 0, bin index from 1
229           //debugging printings
230           printf("sector %i, plate %i, strip %i, padz %i, padx %i \n",detectorIndex[0],detectorIndex[1],detectorIndex[2],detectorIndex[3],detectorIndex[4]);
231           //printf("pos x = %f, pos y = %f, pos z = %f \n",pos[0],pos[1],pos[2]);
232           //printf ("expected time = %f (ns)\n",texp);
233           //printf ("expected time bin = %f (TDC bin)\n",texpBin);
234           printf ("measured time bin = %i (TDC bin)\n",tof);
235           //      printf("index = %i, deltabin = %i , filling index = %i, and bin = % i\n",index, deltabin, index, deltabin);
236         }
237       }
238       
239       delete rawStreamTOF;
240       rawStreamTOF = 0x0;
241  
242       delete rawReader;
243       rawReader = 0x0;
244
245     }
246        
247     /* free resources */
248     free(event);
249     
250     /* exit when last event received, no need to wait for TERM signal */
251     if (eventT==END_OF_RUN) {
252       printf("EOR event detected\n");
253       break;
254     }
255   }
256   
257   delete geomV5;
258   geomV5 = 0x0;
259   delete geom;
260   geom = 0x0;
261
262   //write the Run level file   
263   TFile * fileRun = new TFile ("outciccioTOFT0daRun.root","RECREATE"); 
264   TBenchmark *bench = new TBenchmark();
265   bench->Start("t0");
266   htofPartial->Write();
267   bench->Stop("t0");
268   bench->Print("t0");
269   fileRun->Close();
270
271   TFile * filetot = new TFile ("outciccioTOFT0daTotal.root","READ"); 
272   printf("dopo aver aperto il file in modalita' lettura \n");
273   TH2S *htoftot = 0x0;
274   //look for the file
275   if (!filetot->IsZombie()){
276     printf("il file non e' zombie \n");
277     TIter next(filetot->GetListOfKeys());
278     TKey *key;
279     //look for the histogram
280     while ((key=(TKey*)next())){
281       const char * namekey = key->GetName();
282       if (strcmp(namekey,"htoftot")==0){
283         printf(" histo found \n");
284         htoftot = (TH2S*) filetot->Get("htoftot");
285         htoftot->AddDirectory(0);
286         htoftot->Add(htofPartial);
287         break;
288       }
289     }
290   filetot->Close();
291   }
292   else {
293     printf(" no file found \n");
294     htoftot = new TH2S(*htofPartial);
295     htoftot->SetName("htoftot");
296     htoftot->AddDirectory(0);
297   }
298   filetot=0x0;
299   filetot = new TFile ("outciccioTOFT0daTotal.root","RECREATE");
300   filetot->cd();
301   htoftot->Write();
302   filetot->Close();
303   
304   delete fileRun;
305   delete filetot;
306   delete bench;
307   delete htofPartial;
308   delete htoftot;
309
310   fileRun = 0x0;
311   filetot = 0x0;
312   bench = 0x0;
313   htofPartial = 0x0;
314   htoftot = 0x0;
315
316   /* write report */
317   fprintf(fp,"Run #%s, received %d physics events out of %d\n",getenv("DATE_RUN_NUMBER"),nevents_physics,nevents_total);
318
319   /* close result file */
320   fclose(fp);
321
322
323   return status;
324 }