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