ae60351bafdd2e31f1d514a9057a4b529c877d8e
[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 "TObject.h"
40 #include "TBenchmark.h"
41 #include "TString.h"
42 #include "TH1.h"
43
44 int cbx, ccbx, t0bx, npmtA, npmtC;
45 float clx,cmx,cclx,ccmx, t0lx, t0hx;
46
47 /* Main routine
48       Arguments: 
49       1- monitoring data source
50 */
51 int main(int argc, char **argv) {
52 //int main(){
53   int status;
54
55   /* magic line */
56   gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
57                                         "*",
58                                         "TStreamerInfo",
59                                         "RIO",
60                                         "TStreamerInfo()");
61   
62   if(daqDA_DB_getFile(FILE_IN, FILE_IN)){
63     printf("Couldn't get input file >>inPhys.dat<< from DAQ_DB !!!\n");
64     return -1;
65   }
66   
67   
68   FILE *inp;
69   char c;
70   inp = fopen(FILE_IN, "r");
71   if(!inp){
72     printf("Input file >>inPhys.dat<< not found !!!\n");
73     return -1;
74   }
75   
76   while((c=getc(inp))!=EOF) {
77     switch(c) {
78     case 'a': {fscanf(inp, "%d", &ccbx ); break;} //N of X bins hCFD1_CFD
79     case 'b': {fscanf(inp, "%f", &cclx ); break;} //Low x hCFD1_CFD
80     case 'c': {fscanf(inp, "%f", &ccmx ); break;} //High x hCFD1_CF
81     case 'd': {fscanf(inp, "%d", &npmtC ); break;} //number of reference PMTC
82     case 'e': {fscanf(inp, "%d", &npmtA ); break;} //number of reference PMTA
83     case 'f': {fscanf(inp, "%d", &t0bx ); break;} //N of X bins hT0
84     case 'g': {fscanf(inp, "%f", &t0lx ); break;} //Low x hT0
85     case 'k': {fscanf(inp, "%f", &t0hx ); break;} //High x hT0
86     }
87   }
88   fclose(inp);
89
90   if (argc!=2) {
91     printf("Wrong number of arguments\n");
92     return -1;
93   }
94   
95
96   /* define data source : this is argument 1 */  
97   status=monitorSetDataSource( argv[1] );
98   if (status!=0) {
99     printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
100     return -1;
101   }
102   
103   
104   /* declare monitoring program */
105   status=monitorDeclareMp( __FILE__ );
106   if (status!=0) {
107     printf("monitorDeclareMp() failed : %s\n",monitorDecodeError(status));
108     return -1;
109   }
110   
111   
112   /* define wait event timeout - 1s max */
113   monitorSetNowait();
114   monitorSetNoWaitNetworkTimeout(1000);
115   
116   
117   /* log start of process */
118   printf("T0 monitoring program started\n");  
119   
120   // Allocation of histograms - start
121
122   TH1F *hCFD1minCFD[24];  
123    
124   for(Int_t ic=0; ic<24; ic++) {
125     hCFD1minCFD[ic] = new TH1F(Form("CFD1minCFD%d",ic+1),"CFD-CFD",ccbx,cclx,ccmx);
126   }
127   TH1F *hVertex = new TH1F("hVertex","T0 time",t0bx,t0lx,t0hx);
128   
129  
130   // Allocation of histograms - end
131
132   Int_t iev=0;
133   /* main loop (infinite) */
134   for(;;) {
135     struct eventHeaderStruct *event;
136     eventTypeType eventT;
137     
138     /* check shutdown condition */
139     if (daqDA_checkShutdown()) {break;}
140     
141     /* get next event (blocking call until timeout) */
142     status=monitorGetEventDynamic((void **)&event);
143     if (status==(int)MON_ERR_EOF) {
144       printf ("End of File detected\n");
145       break; /* end of monitoring file has been reached */
146     }
147     
148     if (status!=0) {
149       printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
150       break;
151     }
152     
153     /* retry if got no event */
154     if (event==NULL) {
155       continue;
156     }
157     
158     /* use event - here, just write event id to result file */
159     eventT=event->eventType;
160     
161     switch (event->eventType){
162       
163       case START_OF_RUN:
164         break;
165         
166     case END_OF_RUN:
167       break;
168       
169         case PHYSICS_EVENT:
170           //    case CALIBRATION_EVENT:
171       iev++;
172       
173       if(iev==1){
174         printf("First event - %i\n",iev);
175       }
176       
177       // Initalize raw-data reading and decoding
178       AliRawReader *reader = new AliRawReaderDate((void*)event);
179       
180       // Enable the following two lines in case of real-data
181       reader->RequireHeader(kTRUE);
182       AliT0RawReader *start = new AliT0RawReader(reader, kTRUE);
183       
184       // Read raw data
185       Int_t allData[105][5];
186       for(Int_t i0=0;i0<105;i0++)
187         for(Int_t j0=0;j0<5;j0++)
188           allData[i0][j0] = 0;
189       
190       if(start->Next()){
191         for (Int_t i=0; i<105; i++) {
192           for(Int_t iHit=0;iHit<5;iHit++){
193             allData[i][iHit]= start->GetData(i,iHit);
194           }
195         }
196       }
197       
198       // Fill the histograms
199       Float_t besttimeA=9999999;
200       Float_t besttimeC=9999999;
201       Float_t time[24]; 
202        Float_t meanShift[24];
203        for (Int_t ik = 0; ik<24; ik++)
204          { 
205            if(ik<12 && allData[ik+1][0]>0 
206               && (allData[ik+13][0]-allData[ik+1][0]) < 530 ){
207              hCFD1minCFD[ik]->Fill(allData[ik+1][0]-allData[npmtC][0]);
208            }
209            
210            if(ik>11 && allData[ik+45][0]>0 
211               && (allData[ik+57][0]-allData[ik+45][0]) <530 ){
212              hCFD1minCFD[ik]->Fill(allData[ik+45][0]-allData[56+npmtA][0]);
213            }
214            if(iev == 10000) {   
215              meanShift[ik] =  hCFD1minCFD[ik]->GetMean();  
216            }
217          }
218       //fill  mean time _ fast reconstruction
219       if (iev > 10000 )
220         {
221           for (Int_t in=0; in<12; in++)  
222             {
223               time[in] = allData[in+1][0] - meanShift[in]  ;
224               time[in+12] = allData[in+56+1][0] ;
225             }
226           for (Int_t ipmt=0; ipmt<12; ipmt++){
227             if(time[ipmt] > 1 ) {
228               if(time[ipmt]<besttimeC)
229                 besttimeC=time[ipmt]; //timeC
230             }
231           }
232            for ( Int_t ipmt=12; ipmt<24; ipmt++){
233              if(time[ipmt] > 1) {
234                if(time[ipmt]<besttimeA) 
235                  besttimeA=time[ipmt]; //timeA
236              }
237            }
238            if(besttimeA<9999999 &&besttimeC< 9999999) {
239              Float_t t0 =0.001* 24.4 * Float_t( besttimeA+besttimeC)/2.;
240              hVertex->Fill(t0);
241            }
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
264   TFile *hist = new TFile(FILE_OUT,"RECREATE");
265
266   for(Int_t j=0;j<24;j++){
267      hCFD1minCFD[j]->Write();
268     }
269   hVertex->Write();
270   hist->Close();
271   delete hist;
272
273   status=0;
274
275   /* export file to FXS */
276   if (daqDA_FES_storeFile(FILE_OUT, "PHYSICS")) {
277     status=-2;
278   }
279
280   return status;
281 }
282
283