]> git.uio.no Git - u/mrichter/AliRoot.git/blob - T0/T0Physda.cxx
Added fit macro from M. Putis
[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    
126   for(Int_t ic=0; ic<24; ic++) {
127     hCFD1minCFD[ic] = new TH1F(Form("CFD1minCFD%d",ic+1),"CFD-CFD",kcbx,kclx,kcmx);
128   }
129   TH1F *hVertex = new TH1F("hVertex","T0 time",kt0bx,kt0lx,kt0hx);
130   
131    // Allocation of histograms - end
132
133   Int_t iev=0;
134   /* main loop (infinite) */
135   for(;;) {
136     struct eventHeaderStruct *event;
137     eventTypeType eventT;
138     
139     /* check shutdown condition */
140     if (daqDA_checkShutdown()) {break;}
141     
142     /* get next event (blocking call until timeout) */
143     status=monitorGetEventDynamic((void **)&event);
144     if (status==(int)MON_ERR_EOF) {
145       printf ("End of File detected\n");
146       break; /* end of monitoring file has been reached */
147     }
148     
149     if (status!=0) {
150       printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
151       break;
152     }
153     
154     /* retry if got no event */
155     if (event==NULL) {
156       continue;
157     }
158     
159     /* use event - here, just write event id to result file */
160     eventT=event->eventType;
161     
162     switch (event->eventType){
163       
164     case START_OF_RUN:
165       break;
166         
167     case END_OF_RUN:
168       break;
169       
170     case PHYSICS_EVENT:
171           //    case CALIBRATION_EVENT:
172       iev++;
173       
174       if(iev==1){
175         printf("First event - %i\n",iev);
176       }
177       
178       // Initalize raw-data reading and decoding
179       AliRawReader *reader = new AliRawReaderDate((void*)event);
180       
181       // Enable the following two lines in case of real-data
182       reader->RequireHeader(kTRUE);
183       AliT0RawReader *start = new AliT0RawReader(reader, kTRUE);
184       
185       // Read raw data
186       Int_t allData[105][5];
187       for(Int_t i0=0;i0<105;i0++)
188         for(Int_t j0=0;j0<5;j0++)
189           allData[i0][j0] = 0;
190       
191       if(start->Next()){
192         for (Int_t i=0; i<105; i++) {
193           for(Int_t iHit=0;iHit<5;iHit++){
194             allData[i][iHit]= start->GetData(i,iHit);
195           }
196         }
197       }
198       
199       // Fill the histograms
200       Float_t besttimeA=9999999;
201       Float_t besttimeC=9999999;
202       Float_t time[24]; 
203        Float_t meanShift[24];
204        for (Int_t ik = 0; ik<24; ik++)
205          { 
206            if(ik<12 && allData[ik+1][0]>0 && allData[knpmtC][0]>0 ){
207              hCFD1minCFD[ik]->Fill(allData[ik+1][0]-allData[knpmtC][0]);
208            }
209            
210            if(ik>11 && allData[ik+45][0]>0 && allData[56+knpmtA][0]>0 )
211              {
212              hCFD1minCFD[ik]->Fill(allData[ik+45][0]-allData[56+knpmtA][0]);
213              }
214            if(iev == 10000) {   
215              meanShift[ik] =  hCFD1minCFD[ik]->GetMean(); 
216              if(ik==knpmtC || ik==(56+knpmtA)) meanShift[ik]=0;
217            }
218          }
219       //fill  mean time _ fast reconstruction
220       if (iev > 10000 )
221         {
222           for (Int_t in=0; in<12; in++)  
223             {
224               time[in] = allData[in+1][0] - meanShift[in]  ;
225               time[in+12] = allData[in+56+1][0] ;
226             }
227           for (Int_t ipmt=0; ipmt<12; ipmt++){
228             if(time[ipmt] > 1 ) {
229               if(time[ipmt]<besttimeC)
230                 besttimeC=time[ipmt]; //timeC
231             }
232           }
233            for ( Int_t ipmt=12; ipmt<24; ipmt++){
234              if(time[ipmt] > 1) {
235                if(time[ipmt]<besttimeA) 
236                  besttimeA=time[ipmt]; //timeA
237              }
238            }
239            if(besttimeA<9999999 &&besttimeC< 9999999) {
240              Float_t t0 =0.001* 24.4 * Float_t( besttimeA+besttimeC)/2.;
241              hVertex->Fill(t0);
242            }
243         }
244            
245      delete start;
246      start = 0x0;
247      delete reader;
248      reader= 0x0;
249       // End of fill histograms
250       
251     }
252     
253     /* free resources */
254     free(event);
255     
256     /* exit when last event received, no need to wait for TERM signal */
257     if (eventT==END_OF_RUN) {
258       printf("EOR event detected\n");
259       printf("Number of events processed - %i\n ",iev);         
260       break;
261     }
262   }
263   printf("After loop, before writing histos\n");
264   // write a file with the histograms
265
266   TFile hist(FILE_OUT,"RECREATE");
267
268   for(Int_t j=0;j<24;j++){
269     hCFD1minCFD[j]->SetDirectory(&hist);
270     hCFD1minCFD[j]->Write();
271   }
272   hVertex->Write();
273   hist.Close();
274   //delete hist;
275
276   status=0;
277
278   /* export file to FXS */
279   if (daqDA_FES_storeFile(FILE_OUT, "PHYSICS")) {
280     status=-2;
281   }
282
283   return status;
284 }
285
286