]> git.uio.no Git - u/mrichter/AliRoot.git/blob - T0/T0Cosmicda.cxx
Add new type oh histograms and analysis
[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 #define FILE_IN "inCosmic.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 "TMath.h"
40 #include "TCanvas.h"
41 #include "TString.h"
42 #include "TH1.h"
43 #include "TF1.h"
44 #include "TSpectrum.h"
45 #include "TVirtualFitter.h"
46 //#include "TProfile.h"
47 int cqbx,cqby,clbx,clby,cbx,ccbx;
48 float cqlx,cqmx,cqly,cqmy,cllx,clmx,clly,clmy,clx,cmx,cclx,ccmx;
49 /* Main routine
50       Arguments: 
51       1- monitoring data source
52 */
53 int main(int argc, char **argv) {
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 >>inCosmic.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 >>inCosmic.dat<< not found !!!\n");
74         return -1;
75   }
76
77   while((c=getc(inp))!=EOF) {
78     switch(c) {
79       case 'a': {fscanf(inp, "%d", &cqbx ); break;} //N of X bins hCFD_QTC
80       case 'b': {fscanf(inp, "%f", &cqlx ); break;} //Low x hCFD_QTC
81       case 'c': {fscanf(inp, "%f", &cqmx ); break;} //High x hCFD_QTC
82       case 'd': {fscanf(inp, "%d", &cqby ); break;} //N of Y bins hCFD_QTC
83       case 'e': {fscanf(inp, "%f", &cqly ); break;} //Low y hCFD_QTC
84       case 'f': {fscanf(inp, "%f", &cqmy ); break;} //High y hCFD_QTC
85       case 'g': {fscanf(inp, "%d", &clbx ); break;} //N of X bins hCFD_LED
86       case 'h': {fscanf(inp, "%f", &cllx ); break;} //Low x hCFD_LED
87       case 'i': {fscanf(inp, "%f", &clmx ); break;} //High x hCFD_LED
88       case 'j': {fscanf(inp, "%d", &clby ); break;} //N of Y bins hCFD_LED
89       case 'k': {fscanf(inp, "%f", &clly ); break;} //Low y hCFD_LED
90       case 'l': {fscanf(inp, "%f", &clmy ); break;} //High y hCFD_LED
91       case 'm': {fscanf(inp, "%d", &cbx ); break;}  //N of Y bins hCFD 
92       case 'n': {fscanf(inp, "%f", &clx ); break;}  //Low x hCFD
93       case 'o': {fscanf(inp, "%f", &cmx ); break;}  //High x hCFD
94       case 'p': {fscanf(inp, "%d", &ccbx ); break;} //N of X bins hCFD1_CFD
95       case 'r': {fscanf(inp, "%f", &cclx ); break;} //Low x hCFD1_CFD
96       case 's': {fscanf(inp, "%f", &ccmx ); break;} //High x hCFD1_CFD
97     }
98   }
99   fclose(inp);
100
101   if (argc!=2) {
102     printf("Wrong number of arguments\n");
103     return -1;
104   }
105
106
107   /* define data source : this is argument 1 */  
108   status=monitorSetDataSource( argv[1] );
109   if (status!=0) {
110     printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
111     return -1;
112   }
113
114
115   /* declare monitoring program */
116   status=monitorDeclareMp( __FILE__ );
117   if (status!=0) {
118     printf("monitorDeclareMp() failed : %s\n",monitorDecodeError(status));
119     return -1;
120   }
121
122
123   /* define wait event timeout - 1s max */
124   monitorSetNowait();
125   monitorSetNoWaitNetworkTimeout(1000);
126   
127
128   /* log start of process */
129   printf("T0 monitoring program started\n");  
130
131   // Allocation of histograms - start
132   TH1F *hCFD[24];
133   TH1F *hCFD1minCFD[24]; 
134   TH2F *hCFDvsQTC[24];
135   TH2F *hCFDvsLED[24]; 
136
137    for(Int_t ic=0; ic<24; ic++) {
138       hCFDvsQTC[ic] = new TH2F(Form("CFD_QTC%d",ic+1),"CFD_QTC",cqbx,cqlx,cqmx,cqby,cqly,cqmy);
139       hCFDvsLED[ic] = new TH2F(Form("CFD_LED%d",ic+1),"CFD_LED",clbx,cllx,clmx,clby,clly,clmy);
140       hCFD1minCFD[ic] = new TH1F(Form("CFD1-CFD%d",ic+1),"CFD-CFD",ccbx,cclx,ccmx);
141       if(ic<12){
142         hCFD[ic] = new TH1F(Form("T0_C_%d_CFD",ic+1),"CFD", cbx,clx,cmx);       
143       }
144       else{
145         hCFD[ic] = new TH1F(Form("T0_A_%d_CFD",ic-11),"CFD", cbx,clx,cmx); 
146       }
147     }
148
149   // Allocation of histograms - end
150
151   Int_t iev=0;
152   /* main loop (infinite) */
153   for(;;) {
154     struct eventHeaderStruct *event;
155     eventTypeType eventT;
156   
157     /* check shutdown condition */
158     if (daqDA_checkShutdown()) {break;}
159     
160     /* get next event (blocking call until timeout) */
161     status=monitorGetEventDynamic((void **)&event);
162     if (status==(int)MON_ERR_EOF) {
163       printf ("End of File detected\n");
164       break; /* end of monitoring file has been reached */
165     }
166     
167     if (status!=0) {
168       printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
169       break;
170     }
171
172     /* retry if got no event */
173     if (event==NULL) {
174       continue;
175     }
176  
177 //    iev++; 
178
179     /* use event - here, just write event id to result file */
180     eventT=event->eventType;
181    
182     switch (event->eventType){
183
184       case START_OF_RUN:
185         break;
186
187       case END_OF_RUN:
188         break;
189
190       case CALIBRATION_EVENT:
191 //      case PHYSICS_EVENT:
192       iev++; 
193
194       if(iev==1){
195            printf("First event - %i\n",iev);
196       }
197
198       // Initalize raw-data reading and decoding
199       AliRawReader *reader = new AliRawReaderDate((void*)event);
200           
201       // Enable the following two lines in case of real-data
202           reader->RequireHeader(kTRUE);
203       AliT0RawReader *start = new AliT0RawReader(reader, kTRUE);
204
205       // Read raw data
206       Int_t allData[105][5];
207       for(Int_t i0=0;i0<105;i0++)
208         for(Int_t j0=0;j0<5;j0++)
209                 allData[i0][j0] = 0;
210      
211       if(start->Next())
212       for (Int_t i=0; i<105; i++) {
213         for(Int_t iHit=0;iHit<5;iHit++){
214           allData[i][iHit]= start->GetData(i,iHit);
215         }
216       }
217         else 
218         printf("No T0 data found!!\n");
219
220       // Fill the histograms
221         
222       for (Int_t ik = 0; ik<24; ik+=2)
223          for (Int_t iHt=0; iHt<5; iHt++){
224                  Int_t cc = ik/2;
225                 if((allData[cc+1][iHt]-allData[0][0]+5000)!=0 && allData[cc+1][iHt]>0){
226                  hCFD[cc]->Fill(allData[cc+1][iHt]-allData[0][0]+5000);
227                 }
228                 if((allData[cc+1][iHt]!=0) && (allData[1][iHt]!=0)){
229                  hCFD1minCFD[cc]->Fill(allData[cc+1][iHt]-allData[1][iHt]);
230                 }
231                 if(allData[ik+25][iHt]!=0 && allData[ik+26][iHt]!=0 && allData[cc+1][iHt]!=0){
232                  hCFDvsQTC[cc]->Fill((allData[ik+25][iHt]-allData[ik+26][iHt]) , (allData[cc+1][iHt]-allData[0][0]+5000));
233                 } 
234                 if(allData[cc+13][iHt]!=0 && allData[cc+1][iHt]!=0){
235                  hCFDvsLED[cc]->Fill(allData[cc+13][iHt]-allData[cc+1][iHt],allData[cc+1][iHt]-allData[0][0]+5000);
236                 }
237         }
238
239       for (Int_t ik = 24; ik<48; ik+=2)
240          for (Int_t iHt=0; iHt<5; iHt++){
241                  Int_t cc = ik/2;
242                 if((allData[cc+45][iHt]-allData[0][0]+5000)!=0 && allData[cc+45][iHt]>0){
243                  hCFD[cc]->Fill(allData[cc+45][iHt]-allData[0][0]+5000);
244                 }
245                 if((allData[cc+45][iHt]!=0) && (allData[57][iHt]!=0)){
246                  hCFD1minCFD[cc]->Fill(allData[cc+45][iHt]-allData[57][iHt]);
247                 }
248                 if(allData[ik+57][iHt]!=0 && allData[ik+58][iHt]!=0 && allData[cc+45][iHt]!=0){
249                  hCFDvsQTC[cc]->Fill(allData[ik+57][iHt]-allData[ik+58][iHt],allData[cc+45][iHt]-allData[0][0]+5000);
250                 }
251                 if(allData[cc+57][iHt]!=0 && allData[cc+45][iHt]!=0){
252                  hCFDvsLED[cc]->Fill(allData[cc+57][iHt]-allData[cc+45][iHt],allData[cc+45][iHt]-allData[0][0]+5000);
253                 }
254         }
255         
256      delete start;
257         start = 0x0;
258      reader->Reset();
259       // End of fill histograms
260
261     }
262
263     /* free resources */
264     free(event);
265     
266     /* exit when last event received, no need to wait for TERM signal */
267     if (eventT==END_OF_RUN) {
268       printf("EOR event detected\n");
269       printf("Number of events processed - %i\n ",iev);
270       break;
271     }
272   }
273   printf("After loop, before writing histos\n");
274   // write a file with the histograms
275   TFile *hist = new TFile(FILE_OUT,"RECREATE");
276
277   for(Int_t j=0;j<24;j++){
278      hCFDvsQTC[j]->Write();
279      hCFDvsLED[j]->Write();
280      hCFD[j]->Write();
281      hCFD1minCFD[j]->Write();
282     }
283   hist->Close();
284   delete hist;
285
286   status=0;
287
288   /* export file to FXS */
289   if (daqDA_FES_storeFile(FILE_OUT, "COSMIC")) {
290     status=-2;
291   }
292
293   return status;
294 }
295