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