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