]> git.uio.no Git - u/mrichter/AliRoot.git/blob - T0/T0Physda.cxx
expantion of recoParam (Markus) and fix bug in HLT mode
[u/mrichter/AliRoot.git] / T0 / T0Physda.cxx
1 /*
2 T0 DA for online calibration
3  
4 Contact: Michal.Oledzki@cern.ch
5 Link: http://users.jyu.fi/~mioledzk/
6 Run Type: PHYSICS
7 DA Type: MON
8 Number of events needed: 500000 
9 Input Files: inPhys.dat, external parameters
10 Output Files: daPhys.root, to be exported to the DAQ FXS
11 Trigger types used: PHYSICS_EVENT
12 ------------------------------Alla
13 Now trigger type changed to CALICRATION_EVENT
14 to have data during test.
15 SOULD BE CHANGED BACK BEFORE BEAM
16 ------------------------------- Alla
17 */
18
19 #define FILE_OUT "daPhys.root"
20 #define FILE_IN "inPhys.dat"
21 #include <daqDA.h>
22 #include <event.h>
23 #include <monitor.h>
24  
25 #include <Riostream.h>
26 #include <stdio.h>
27 #include <stdlib.h>
28
29 //AliRoot
30 #include <AliRawReaderDate.h>
31 #include <AliRawReader.h>
32 #include <AliT0RawReader.h>
33
34 //ROOT
35 #include "TROOT.h"
36 #include "TPluginManager.h"
37 #include "TFile.h"
38 #include "TKey.h"
39 #include "TH2S.h"
40 #include "TObject.h"
41 #include "TBenchmark.h"
42 #include "TRandom.h"
43 #include "TCanvas.h"
44 #include "TString.h"
45 #include "TH1.h"
46 #include "TF1.h"
47 #include "TSpectrum.h"
48 #include "TVirtualFitter.h"
49 int cbx, ccbx;
50 float clx,cmx,cclx,ccmx;
51
52 /* Main routine
53       Arguments: 
54       1- monitoring data source
55 */
56 int main(int argc, char **argv) {
57 //int main(){
58   int status;
59
60   /* magic line */
61   gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
62                                         "*",
63                                         "TStreamerInfo",
64                                         "RIO",
65                                         "TStreamerInfo()");
66   
67   if(daqDA_DB_getFile(FILE_IN, FILE_IN)){
68     printf("Couldn't get input file >>inPhys.dat<< from DAQ_DB !!!\n");
69     return -1;
70   }
71   
72   
73   FILE *inp;
74   char c;
75   inp = fopen(FILE_IN, "r");
76   if(!inp){
77     printf("Input file >>inPhys.dat<< not found !!!\n");
78     return -1;
79   }
80   
81   while((c=getc(inp))!=EOF) {
82     switch(c) {
83     case 'a': {fscanf(inp, "%d", &ccbx ); break;} //N of X bins hCFD1_CFD
84     case 'b': {fscanf(inp, "%f", &cclx ); break;} //Low x hCFD1_CFD
85     case 'c': {fscanf(inp, "%f", &ccmx ); break;} //High x hCFD1_CFD
86       //      case 'd': {fscanf(inp, "%d", &cbx ); break;} //N of X bins hCFD
87       //      case 'e': {fscanf(inp, "%f", &clx ); break;} //Low x hCFD
88 //      case 'f': {fscanf(inp, "%f", &cmx ); break;} //High x hCFD
89     }
90   }
91   fclose(inp);
92
93   if (argc!=2) {
94     printf("Wrong number of arguments\n");
95     return -1;
96   }
97   
98
99   /* define data source : this is argument 1 */  
100   status=monitorSetDataSource( argv[1] );
101   if (status!=0) {
102     printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
103     return -1;
104   }
105   
106   
107   /* declare monitoring program */
108   status=monitorDeclareMp( __FILE__ );
109   if (status!=0) {
110     printf("monitorDeclareMp() failed : %s\n",monitorDecodeError(status));
111     return -1;
112   }
113   
114   
115   /* define wait event timeout - 1s max */
116   monitorSetNowait();
117   monitorSetNoWaitNetworkTimeout(1000);
118   
119   
120   /* log start of process */
121   printf("T0 monitoring program started\n");  
122   
123   // Allocation of histograms - start
124
125   TH1F *hCFD1minCFD[24];  
126    
127   for(Int_t ic=0; ic<24; ic++) {
128     hCFD1minCFD[ic] = new TH1F(Form("CFD1-CFD%d",ic+1),"CFD-CFD",ccbx,cclx,ccmx);
129   }
130   TH1F *hVertex = new TH1F("hVertex","Z vertex",ccbx,cclx,ccmx);
131   
132   Float_t meanShift[24];
133   
134   // Allocation of histograms - end
135
136   Int_t iev=0;
137   /* main loop (infinite) */
138   for(;;) {
139     struct eventHeaderStruct *event;
140     eventTypeType eventT;
141     
142     /* check shutdown condition */
143     if (daqDA_checkShutdown()) {break;}
144     
145     /* get next event (blocking call until timeout) */
146     status=monitorGetEventDynamic((void **)&event);
147     if (status==(int)MON_ERR_EOF) {
148       printf ("End of File detected\n");
149       break; /* end of monitoring file has been reached */
150     }
151     
152     if (status!=0) {
153       printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
154       break;
155     }
156     
157     /* retry if got no event */
158     if (event==NULL) {
159       continue;
160     }
161     
162     /* use event - here, just write event id to result file */
163     eventT=event->eventType;
164     
165     switch (event->eventType){
166       
167       case START_OF_RUN:
168         break;
169         
170     case END_OF_RUN:
171       break;
172       
173     case PHYSICS_EVENT:
174       //    case CALIBRATION_EVENT:
175       iev++;
176       
177       if(iev==1){
178         printf("First event - %i\n",iev);
179       }
180       
181       // Initalize raw-data reading and decoding
182       AliRawReader *reader = new AliRawReaderDate((void*)event);
183       
184       // Enable the following two lines in case of real-data
185       reader->RequireHeader(kTRUE);
186       AliT0RawReader *start = new AliT0RawReader(reader, kTRUE);
187       
188       // Read raw data
189       Int_t allData[105][5];
190       for(Int_t i0=0;i0<105;i0++)
191         for(Int_t j0=0;j0<5;j0++)
192           allData[i0][j0] = 0;
193       
194       if(start->Next()){
195         for (Int_t i=0; i<105; i++) {
196           for(Int_t iHit=0;iHit<5;iHit++){
197             allData[i][iHit]= start->GetData(i,iHit);
198           }
199         }
200       }
201       else 
202         printf("No T0 data found!!!\n");
203       
204       // Fill the histograms
205       Float_t besttimeA=9999999;
206       Float_t besttimeC=9999999;
207       Float_t time[24]; 
208       for (Int_t ik = 0; ik<24; ik++)
209         if(allData[ik+1][0]!=0 ){
210           if(ik<12){
211              hCFD1minCFD[ik]->Fill(allData[ik+1][0]-allData[1][0]);
212              if(iev == 20000)   meanShift[ik] =  hCFD1minCFD[ik]->GetMean();  
213           }
214           if(ik>11){
215             hCFD1minCFD[ik]->Fill(allData[ik+45][0]-allData[57][0]);
216             if(iev == 20000)    
217               meanShift[ik] =  hCFD1minCFD[ik]->GetMean();  
218           }
219         }
220       //fill vertex & mean time _ fast reconstruction
221       if (iev > 20000 && iev <50000)
222         {
223           for (Int_t in=0; in<12; in++)  
224             {
225               time[in] = allData[in+1][0] - meanShift[in] + 5000 - allData[0][0] ;
226               time[in+12] = allData[in+56+1][0] - meanShift[in+12] + 5000 - allData[0][0];
227             }
228           for (Int_t ipmt=0; ipmt<12; ipmt++){
229             if(time[ipmt] > 1 ) {
230               if(time[ipmt]<besttimeC)
231                 besttimeC=time[ipmt]; //timeC
232             }
233           }
234            for ( Int_t ipmt=12; ipmt<24; ipmt++){
235              if(time[ipmt] > 1) {
236                if(time[ipmt]<besttimeA) 
237                  besttimeA=time[ipmt]; //timeA
238              }
239            }
240            Float_t vertex = 0.0299792 *(besttimeC - besttimeA)*24.4/2.; 
241            hVertex->Fill(vertex);
242            
243         }
244       delete start;
245       start = 0x0;
246       reader->Reset();
247       // End of fill histograms
248       
249     }
250     
251     /* free resources */
252     free(event);
253     
254     /* exit when last event received, no need to wait for TERM signal */
255     if (eventT==END_OF_RUN) {
256       printf("EOR event detected\n");
257       printf("Number of events processed - %i\n ",iev);         
258       break;
259     }
260   }
261   printf("After loop, before writing histos\n");
262   // write a file with the histograms
263   TFile *hist = new TFile(FILE_OUT,"RECREATE");
264
265   for(Int_t j=0;j<24;j++){
266      hCFD1minCFD[j]->Write();
267     }
268   hVertex->Write();
269   hist->Close();
270   delete hist;
271
272   status=0;
273
274   /* export file to FXS */
275   if (daqDA_FES_storeFile(FILE_OUT, "PHYSICS")) {
276     status=-2;
277   }
278
279   return status;
280 }
281
282