]> git.uio.no Git - u/mrichter/AliRoot.git/blob - T0/T0Physda.cxx
c++ bug fixed Thanks to Yves
[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
13 */
14
15 #define FILE_OUT "daPhys.root"
16 #define FILE_IN "inPhys.dat"
17 #include <daqDA.h>
18 #include <event.h>
19 #include <monitor.h>
20  
21 #include <Riostream.h>
22 #include <stdio.h>
23 #include <stdlib.h>
24
25 //AliRoot
26 #include <AliRawReaderDate.h>
27 #include <AliRawReader.h>
28 #include <AliT0RawReader.h>
29
30 //ROOT
31 #include "TROOT.h"
32 #include "TPluginManager.h"
33 #include "TFile.h"
34 #include "TKey.h"
35 #include "TH2S.h"
36 #include "TObject.h"
37 #include "TBenchmark.h"
38 #include "TRandom.h"
39 #include "TCanvas.h"
40 #include "TString.h"
41 #include "TH1.h"
42 #include "TF1.h"
43 #include "TSpectrum.h"
44 #include "TVirtualFitter.h"
45 int cbx, ccbx;
46 float clx,cmx,cclx,ccmx;
47
48 /* Main routine
49       Arguments: 
50       1- monitoring data source
51 */
52 int main(int argc, char **argv) {
53 //int main(){
54   int status;
55
56   /* magic line */
57   gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
58                                         "*",
59                                         "TStreamerInfo",
60                                         "RIO",
61                                         "TStreamerInfo()");
62
63   if(daqDA_DB_getFile(FILE_IN, FILE_IN)){
64      printf("Couldn't get input file >>inPhys.dat<< from DAQ_DB !!!\n");
65      return -1;
66   }
67
68   
69   FILE *inp;
70   char c;
71   inp = fopen(FILE_IN, "r");
72   if(!inp){
73         printf("Input file >>inPhys.dat<< not found !!!\n");
74         return -1;
75   }
76
77   while((c=getc(inp))!=EOF) {
78     switch(c) {
79       case 'a': {fscanf(inp, "%d", &ccbx ); break;} //N of X bins hCFD1_CFD
80       case 'b': {fscanf(inp, "%f", &cclx ); break;} //Low x hCFD1_CFD
81       case 'c': {fscanf(inp, "%f", &ccmx ); break;} //High x hCFD1_CFD
82 //      case 'd': {fscanf(inp, "%d", &cbx ); break;} //N of X bins hCFD
83 //      case 'e': {fscanf(inp, "%f", &clx ); break;} //Low x hCFD
84 //      case 'f': {fscanf(inp, "%f", &cmx ); break;} //High x hCFD
85     }
86   }
87   fclose(inp);
88
89   if (argc!=2) {
90     printf("Wrong number of arguments\n");
91     return -1;
92   }
93
94
95   /* define data source : this is argument 1 */  
96   status=monitorSetDataSource( argv[1] );
97   if (status!=0) {
98     printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
99     return -1;
100   }
101
102
103   /* declare monitoring program */
104   status=monitorDeclareMp( __FILE__ );
105   if (status!=0) {
106     printf("monitorDeclareMp() failed : %s\n",monitorDecodeError(status));
107     return -1;
108   }
109
110
111   /* define wait event timeout - 1s max */
112   monitorSetNowait();
113   monitorSetNoWaitNetworkTimeout(1000);
114   
115
116   /* log start of process */
117   printf("T0 monitoring program started\n");  
118
119   // Allocation of histograms - start
120
121   TH1F *hCFD1minCFD[24];  
122    
123    for(Int_t ic=0; ic<24; ic++) {
124         hCFD1minCFD[ic] = new TH1F(Form("CFD1-CFD%d",ic+1),"CFD-CFD",ccbx,cclx,ccmx);
125     }
126
127   // Allocation of histograms - end
128
129   Int_t iev=0;
130   /* main loop (infinite) */
131   for(;;) {
132     struct eventHeaderStruct *event;
133     eventTypeType eventT;
134   
135     /* check shutdown condition */
136     if (daqDA_checkShutdown()) {break;}
137     
138     /* get next event (blocking call until timeout) */
139     status=monitorGetEventDynamic((void **)&event);
140     if (status==(int)MON_ERR_EOF) {
141       printf ("End of File detected\n");
142       break; /* end of monitoring file has been reached */
143     }
144     
145     if (status!=0) {
146       printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
147       break;
148     }
149
150     /* retry if got no event */
151     if (event==NULL) {
152       continue;
153     }
154
155     /* use event - here, just write event id to result file */
156     eventT=event->eventType;
157    
158     switch (event->eventType){
159
160       case START_OF_RUN:
161         break;
162
163       case END_OF_RUN:
164         break;
165
166       case PHYSICS_EVENT:
167 //      case CALIBRATION_EVENT:
168       iev++;
169
170       if(iev==1){
171            printf("First event - %i\n",iev);
172       }
173
174       // Initalize raw-data reading and decoding
175       AliRawReader *reader = new AliRawReaderDate((void*)event);
176           
177       // Enable the following two lines in case of real-data
178           reader->RequireHeader(kTRUE);
179       AliT0RawReader *start = new AliT0RawReader(reader, kTRUE);
180
181       // Read raw data
182       Int_t allData[105][5];
183       for(Int_t i0=0;i0<105;i0++)
184         for(Int_t j0=0;j0<5;j0++)
185                 allData[i0][j0] = 0;
186
187        if(start->Next()){
188        for (Int_t i=0; i<105; i++) {
189         for(Int_t iHit=0;iHit<5;iHit++){
190           allData[i][iHit]= start->GetData(i,iHit);
191         }
192        }
193       }
194         else 
195         printf("No T0 data found!!!\n");
196
197       // Fill the histograms
198         
199       for (Int_t ik = 0; ik<24; ik++)
200          for (Int_t iHt=0; iHt<5; iHt++){
201                 if(allData[ik+1][iHt]!=0 ){
202                   if(ik<12){
203                          hCFD1minCFD[ik]->Fill(allData[ik+1][iHt]-allData[1][iHt]);
204                   }
205                   if(ik>11){
206                          hCFD1minCFD[ik]->Fill(allData[ik+45][iHt]-allData[57][iHt]);
207                   }
208                 }
209         }
210
211      delete start;
212         start = 0x0;
213      reader->Reset();
214       // End of fill histograms
215
216     }
217
218     /* free resources */
219     free(event);
220     
221     /* exit when last event received, no need to wait for TERM signal */
222     if (eventT==END_OF_RUN) {
223       printf("EOR event detected\n");
224       printf("Number of events processed - %i\n ",iev);         
225       break;
226     }
227   }
228   printf("After loop, before writing histos\n");
229   // write a file with the histograms
230   TFile *hist = new TFile(FILE_OUT,"RECREATE");
231
232   for(Int_t j=0;j<24;j++){
233      hCFD1minCFD[j]->Write();
234     }
235   hist->Close();
236   delete hist;
237
238   status=0;
239
240   /* export file to FXS */
241   if (daqDA_FES_storeFile(FILE_OUT, "PHYSICS")) {
242     status=-2;
243   }
244
245   return status;
246 }
247
248