]> git.uio.no Git - u/mrichter/AliRoot.git/blob - T0/T0Physda.cxx
Protection in case of wrong time stemp improved
[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, npmtA, npmtC;
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_CF
86     case 'd': {fscanf(inp, "%d", &npmtC ); break;} //number of reference PMTC
87     case 'e': {fscanf(inp, "%d", &npmtA ); break;} //number of reference PMTA
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",ccbx,cclx,ccmx);
128   }
129   TH1F *hVertex = new TH1F("hVertex","Z vertex",ccbx,cclx,ccmx);
130   
131  
132   // Allocation of histograms - end
133
134   Int_t iev=0;
135   /* main loop (infinite) */
136   for(;;) {
137     struct eventHeaderStruct *event;
138     eventTypeType eventT;
139     
140     /* check shutdown condition */
141     if (daqDA_checkShutdown()) {break;}
142     
143     /* get next event (blocking call until timeout) */
144     status=monitorGetEventDynamic((void **)&event);
145     if (status==(int)MON_ERR_EOF) {
146       printf ("End of File detected\n");
147       break; /* end of monitoring file has been reached */
148     }
149     
150     if (status!=0) {
151       printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
152       break;
153     }
154     
155     /* retry if got no event */
156     if (event==NULL) {
157       continue;
158     }
159     
160     /* use event - here, just write event id to result file */
161     eventT=event->eventType;
162     
163     switch (event->eventType){
164       
165       case START_OF_RUN:
166         break;
167         
168     case END_OF_RUN:
169       break;
170       
171       //  case PHYSICS_EVENT:
172           case CALIBRATION_EVENT:
173       iev++;
174       
175       if(iev==1){
176         printf("First event - %i\n",iev);
177       }
178       
179       // Initalize raw-data reading and decoding
180       AliRawReader *reader = new AliRawReaderDate((void*)event);
181       
182       // Enable the following two lines in case of real-data
183       reader->RequireHeader(kTRUE);
184       AliT0RawReader *start = new AliT0RawReader(reader, kTRUE);
185       
186       // Read raw data
187       Int_t allData[105][5];
188       for(Int_t i0=0;i0<105;i0++)
189         for(Int_t j0=0;j0<5;j0++)
190           allData[i0][j0] = 0;
191       
192       if(start->Next()){
193         for (Int_t i=0; i<105; i++) {
194           for(Int_t iHit=0;iHit<5;iHit++){
195             allData[i][iHit]= start->GetData(i,iHit);
196           }
197         }
198       }
199       
200       // Fill the histograms
201       Float_t besttimeA=9999999;
202       Float_t besttimeC=9999999;
203       Float_t time[24]; 
204       Int_t acshift=56;
205       for (Int_t ik = 0; ik<24; ik++)
206         if(allData[ik+1][0]>0 ){
207           if(ik<12 
208              && (allData[ik+13][0]-allData[ik+1][0]) < 530 ){
209             hCFD1minCFD[ik]->Fill(allData[ik+1][0]-allData[npmtC][0]);
210           }
211           if(ik>11
212              && (allData[ik+57][0]-allData[ik+45][0]) <530 ){
213             hCFD1minCFD[ik]->Fill(allData[ik+45][0]-allData[56+npmtA][0]);
214           }
215
216         }
217       delete start;
218       start = 0x0;
219       reader->Reset();
220       // End of fill histograms
221       
222     }
223     
224     /* free resources */
225     free(event);
226     
227     /* exit when last event received, no need to wait for TERM signal */
228     if (eventT==END_OF_RUN) {
229       printf("EOR event detected\n");
230       printf("Number of events processed - %i\n ",iev);         
231       break;
232     }
233   }
234   printf("After loop, before writing histos\n");
235   // write a file with the histograms
236
237   TFile *hist = new TFile(FILE_OUT,"RECREATE");
238
239   for(Int_t j=0;j<24;j++){
240      hCFD1minCFD[j]->Write();
241     }
242   hVertex->Write();
243   hist->Close();
244   delete hist;
245
246   status=0;
247
248   /* export file to FXS */
249   if (daqDA_FES_storeFile(FILE_OUT, "PHYSICS")) {
250     status=-2;
251   }
252
253   return status;
254 }
255
256