]> git.uio.no Git - u/mrichter/AliRoot.git/blob - T0/T0Physda.cxx
DA : 1st hit only collected
[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
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         else 
199         printf("No T0 data found!!!\n");
200
201       // Fill the histograms
202         
203       for (Int_t ik = 0; ik<24; ik++)
204          for (Int_t iHt=0; iHt<1; iHt++){
205                 if(allData[ik+1][iHt]!=0 ){
206                   if(ik<12){
207                          hCFD1minCFD[ik]->Fill(allData[ik+1][iHt]-allData[1][iHt]);
208                   }
209                   if(ik>11){
210                          hCFD1minCFD[ik]->Fill(allData[ik+45][iHt]-allData[57][iHt]);
211                   }
212                 }
213         }
214
215      delete start;
216         start = 0x0;
217      reader->Reset();
218       // End of fill histograms
219
220     }
221
222     /* free resources */
223     free(event);
224     
225     /* exit when last event received, no need to wait for TERM signal */
226     if (eventT==END_OF_RUN) {
227       printf("EOR event detected\n");
228       printf("Number of events processed - %i\n ",iev);         
229       break;
230     }
231   }
232   printf("After loop, before writing histos\n");
233   // write a file with the histograms
234   TFile *hist = new TFile(FILE_OUT,"RECREATE");
235
236   for(Int_t j=0;j<24;j++){
237      hCFD1minCFD[j]->Write();
238     }
239   hist->Close();
240   delete hist;
241
242   status=0;
243
244   /* export file to FXS */
245   if (daqDA_FES_storeFile(FILE_OUT, "PHYSICS")) {
246     status=-2;
247   }
248
249   return status;
250 }
251
252