]> git.uio.no Git - u/mrichter/AliRoot.git/blob - T0/T0Physda.cxx
coverty warning fixed
[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 #include "TSpectrum.h"
44 #include "TMath.h"
45
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   int  kcbx, kt0bx, knpmtA, knpmtC;
76   float kclx,kcmx, kt0lx, kt0hx;
77   
78   while((c=getc(inp))!=EOF) {
79     switch(c) {
80     case 'a': {fscanf(inp, "%d", &kcbx ); break;} //N of X bins hCFD1_CFD
81     case 'b': {fscanf(inp, "%f", &kclx ); break;} //Low x hCFD1_CFD
82     case 'c': {fscanf(inp, "%f", &kcmx ); break;} //High x hCFD1_CF
83     case 'd': {fscanf(inp, "%d", &knpmtC ); break;} //number of reference PMTC
84     case 'e': {fscanf(inp, "%d", &knpmtA ); break;} //number of reference PMTA
85     case 'f': {fscanf(inp, "%d", &kt0bx ); break;} //N of X bins hT0
86     case 'g': {fscanf(inp, "%f", &kt0lx ); break;} //Low x hT0
87     case 'k': {fscanf(inp, "%f", &kt0hx ); break;} //High x hT0
88     }
89   }
90   fclose(inp);
91
92   if (argc!=2) {
93     printf("Wrong number of arguments\n");
94     return -1;
95   }
96   
97
98   /* define data source : this is argument 1 */  
99   status=monitorSetDataSource( argv[1] );
100   if (status!=0) {
101     printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
102     return -1;
103   }
104   
105   
106   /* declare monitoring program */
107   status=monitorDeclareMp( __FILE__ );
108   if (status!=0) {
109     printf("monitorDeclareMp() failed : %s\n",monitorDecodeError(status));
110     return -1;
111   }
112   
113   
114   /* define wait event timeout - 1s max */
115   monitorSetNowait();
116   monitorSetNoWaitNetworkTimeout(1000);
117   
118   
119   /* log start of process */
120   printf("T0 monitoring program started\n");  
121   
122   // Allocation of histograms - start
123
124   TH1F *hCFD1minCFD[24];  
125   TH1F *hCFD[24];  
126    
127   for(Int_t ic=0; ic<24; ic++) {
128     hCFD1minCFD[ic] = new TH1F(Form("CFD1minCFD%d",ic+1),"CFD-CFD",kcbx,kclx,kcmx);
129     hCFD[ic] = new TH1F(Form("CFD%d",ic+1),"CFD",kt0bx,kt0lx,kt0hx);
130   }
131   TH1F *hVertex = new TH1F("hVertex","T0 time",kt0bx,kt0lx,kt0hx);
132   
133    // Allocation of histograms - end
134
135   Int_t iev=0;
136   /* main loop (infinite) */
137   for(;;) {
138     struct eventHeaderStruct *event;
139     eventTypeType eventT;
140     
141     /* check shutdown condition */
142     if (daqDA_checkShutdown()) {break;}
143     
144     /* get next event (blocking call until timeout) */
145     status=monitorGetEventDynamic((void **)&event);
146     if (status==(int)MON_ERR_EOF) {
147       printf ("End of File detected\n");
148       break; /* end of monitoring file has been reached */
149     }
150     
151     if (status!=0) {
152       printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
153       break;
154     }
155     
156     /* retry if got no event */
157     if (event==NULL) {
158       continue;
159     }
160     
161     /* use event - here, just write event id to result file */
162     eventT=event->eventType;
163     
164     switch (event->eventType){
165       
166     case START_OF_RUN:
167       break;
168         
169     case END_OF_RUN:
170       break;
171       
172     case PHYSICS_EVENT:
173           //    case CALIBRATION_EVENT:
174       iev++;
175       
176       if(iev==1){
177         printf("First event - %i\n",iev);
178       }
179       
180       // Initalize raw-data reading and decoding
181       AliRawReader *reader = new AliRawReaderDate((void*)event);
182       
183       // Enable the following two lines in case of real-data
184       reader->RequireHeader(kTRUE);
185       AliT0RawReader *start = new AliT0RawReader(reader, kTRUE);
186       
187       // Read raw data
188       Int_t allData[110][5];
189       for(Int_t i0=0;i0<107;i0++)
190         for(Int_t j0=0;j0<5;j0++)
191           allData[i0][j0] = 0;
192       
193       if(start->Next()){
194         for (Int_t i=0; i<107; i++) {
195           for(Int_t iHit=0;iHit<5;iHit++){
196             allData[i][iHit]= start->GetData(i,iHit);
197           }
198         }
199       }
200       
201       // Fill the histograms
202       Float_t besttimeA=9999999;
203       Float_t besttimeC=9999999;
204       Float_t time[24]; 
205        Float_t meanShift[24];
206        for (Int_t ik = 0; ik<24; ik++)
207          { 
208            if(ik<12 && allData[ik+1][0]>0 && allData[knpmtC][0]>0 ){
209              hCFD1minCFD[ik]->Fill(allData[ik+1][0]-allData[knpmtC][0]);
210              if(allData[ik+1][0]>0 ) hCFD[ik]->Fill(allData[ik+1][0]);
211            }
212            
213            if(ik>11 && allData[ik+45][0]>0 && allData[56+knpmtA][0]>0 )
214              {
215              hCFD1minCFD[ik]->Fill(allData[ik+45][0]-allData[56+knpmtA][0]);
216              if(allData[ik+1][0]>0 ) hCFD[ik]->Fill(allData[ik+45][0]);
217              }
218            if(iev == 10000) {   
219              meanShift[ik] =  hCFD1minCFD[ik]->GetMean(); 
220              if(ik==knpmtC || ik==(56+knpmtA)) meanShift[ik]=0;
221            }
222          }
223       //fill  mean time _ fast reconstruction
224       if (iev > 10000 )
225         {
226           for (Int_t in=0; in<12; in++)  
227             {
228               time[in] = allData[in+1][0] - meanShift[in]  ;
229               time[in+12] = allData[in+56+1][0] ;
230             }
231           for (Int_t ipmt=0; ipmt<12; ipmt++){
232             if(time[ipmt] > 1 ) {
233               if(time[ipmt]<besttimeC)
234                 besttimeC=time[ipmt]; //timeC
235             }
236           }
237            for ( Int_t ipmt=12; ipmt<24; ipmt++){
238              if(time[ipmt] > 1) {
239                if(time[ipmt]<besttimeA) 
240                  besttimeA=time[ipmt]; //timeA
241              }
242            }
243            if(besttimeA<9999999 &&besttimeC< 9999999) {
244              Float_t t0 =0.001* 24.4 * Float_t( besttimeA+besttimeC)/2.;
245              hVertex->Fill(t0);
246            }
247         }
248            
249      delete start;
250      start = 0x0;
251      delete reader;
252      reader= 0x0;
253       // End of fill histograms
254       
255     }
256     
257     /* free resources */
258     free(event);
259     
260     /* exit when last event received, no need to wait for TERM signal */
261     if (eventT==END_OF_RUN) {
262       printf("EOR event detected\n");
263       printf("Number of events processed - %i\n ",iev);         
264       break;
265     }
266   }
267   printf("After loop, before writing histos\n");
268   // write a file with the histograms
269
270   TFile hist(FILE_OUT,"RECREATE");
271
272   for(Int_t j=0;j<24;j++){
273     hCFD1minCFD[j]->SetDirectory(&hist);
274     hCFD1minCFD[j]->Write();
275     hCFD[j]->Write();
276
277   }
278   hVertex->Write();
279   hist.Close();
280   //delete hist;
281
282   status=0;
283
284   /* export file to FXS */
285   if (daqDA_FES_storeFile(FILE_OUT, "PHYSICS")) {
286     status=-2;
287   }
288
289   return status;
290 }
291
292