]> git.uio.no Git - u/mrichter/AliRoot.git/blame_incremental - T0/T0Laserda.cxx
Reconstruction with measused calibration
[u/mrichter/AliRoot.git] / T0 / T0Laserda.cxx
... / ...
CommitLineData
1/*
2T0 DA for online calibration
3
4Contact: Michal.Oledzki@cern.ch
5Link: http://users.jyu.fi/~mioledzk/
6Run Type: T0_STANDALONE_LASER
7DA Type: MON
8Number of events needed: 400000
9Input Files: inLaser.dat, external parameters
10Output Files: daLaser.root, to be exported to the DAQ FXS
11Trigger types used: PHYSICS_EVENT
12
13*/
14
15#define FILE_OUT "daLaser.root"
16#define FILE_IN "inLaser.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"
47int cqbx,cqby,clbx,clby,cbx;
48float cqlx,cqmx,cqly,cqmy,cllx,clmx,clly,clmy,clx,cmx;
49/* Main routine
50 Arguments:
51 1- monitoring data source
52*/
53int 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 >>inLaser.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 >>inLaser.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 }
95 }
96 fclose(inp);
97
98 if (argc!=2) {
99 printf("Wrong number of arguments\n");
100 return -1;
101 }
102
103
104 /* define data source : this is argument 1 */
105 status=monitorSetDataSource( argv[1] );
106 if (status!=0) {
107 printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
108 return -1;
109 }
110
111
112 /* declare monitoring program */
113 status=monitorDeclareMp( __FILE__ );
114 if (status!=0) {
115 printf("monitorDeclareMp() failed : %s\n",monitorDecodeError(status));
116 return -1;
117 }
118
119
120 /* define wait event timeout - 1s max */
121 monitorSetNowait();
122 monitorSetNoWaitNetworkTimeout(1000);
123
124
125 /* log start of process */
126 printf("T0 monitoring program started\n");
127
128 // Allocation of histograms - start
129 TH1F *hCFD[24];
130 TH2F *hCFDvsQTC[24];
131 TH2F *hCFDvsLED[24];
132
133 for(Int_t ic=0; ic<24; ic++) {
134 hCFDvsQTC[ic] = new TH2F(Form("CFD_QTC%d",ic+1),"CFD_QTC",cqbx,cqlx,cqmx,cqby,cqly,cqmy);
135 hCFDvsLED[ic] = new TH2F(Form("CFD_LED%d",ic+1),"CFD_LED",clbx,cllx,clmx,clby,clly,clmy);
136 if(ic<12){
137 hCFD[ic] = new TH1F(Form("T0_C_%d_CFD",ic+1),"CFD", cbx,clx,cmx);
138 }
139 else{
140 hCFD[ic] = new TH1F(Form("T0_A_%d_CFD",ic-11),"CFD", cbx,clx,cmx);
141 }
142 }
143
144 // Allocation of histograms - end
145
146 Int_t iev=0;
147 /* main loop (infinite) */
148 for(;;) {
149 struct eventHeaderStruct *event;
150 eventTypeType eventT;
151
152 /* check shutdown condition */
153 if (daqDA_checkShutdown()) {break;}
154
155 /* get next event (blocking call until timeout) */
156 status=monitorGetEventDynamic((void **)&event);
157 if (status==(int)MON_ERR_EOF) {
158 printf ("End of File detected\n");
159 break; /* end of monitoring file has been reached */
160 }
161
162 if (status!=0) {
163 printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
164 break;
165 }
166
167 /* retry if got no event */
168 if (event==NULL) {
169 continue;
170 }
171
172// iev++;
173
174 /* use event - here, just write event id to result file */
175 eventT=event->eventType;
176
177 switch (event->eventType){
178
179 case START_OF_RUN:
180 break;
181
182 case END_OF_RUN:
183 break;
184
185 case CALIBRATION_EVENT:
186// case PHYSICS_EVENT:
187 iev++;
188
189 if(iev==1){
190 printf("First event - %i\n",iev);
191 }
192// if((iev<10) || (iev>11000)){
193// printf(" event number = %i \n",iev);
194// }
195
196
197 // Initalize raw-data reading and decoding
198 AliRawReader *reader = new AliRawReaderDate((void*)event);
199
200 // Enable the following two lines in case of real-data
201 // reader->LoadEquipmentIdsMap("T0map.txt");
202 // reader->RequireHeader(kFALSE);
203 // reader->RequireHeader(kTRUE);
204 AliT0RawReader *start = new AliT0RawReader(reader, kTRUE);
205
206 // Read raw data
207 Int_t allData[105][5];
208 for(Int_t i0=0;i0<105;i0++)
209 for(Int_t j0=0;j0<5;j0++)
210 allData[i0][j0] = 0;
211
212 if(start->Next())
213 for (Int_t i=0; i<105; i++) {
214 for(Int_t iHit=0;iHit<5;iHit++){
215 allData[i][iHit]= start->GetData(i,iHit);
216 }
217 }
218 else
219 printf("No T0 data found!!\n");
220
221 // Fill the histograms
222
223 for (Int_t ik = 0; ik<24; ik+=2)
224 for (Int_t iHt=0; iHt<5; iHt++){
225 Int_t cc = ik/2;
226 if((allData[cc+1][iHt]-allData[0][0]+5000)!=0 && allData[cc+1][iHt]>0){
227 hCFD[cc]->Fill(allData[cc+1][iHt]-allData[0][0]+5000);
228 }
229 if(allData[ik+25][iHt]!=0 && allData[ik+26][iHt]!=0 && allData[cc+1][iHt]!=0){
230 hCFDvsQTC[cc]->Fill((allData[ik+25][iHt]-allData[ik+26][iHt]) , (allData[cc+1][iHt]-allData[0][0]+5000));
231 }
232 if(allData[cc+13][iHt]!=0 && allData[cc+1][iHt]!=0){
233 hCFDvsLED[cc]->Fill(allData[cc+13][iHt]-allData[cc+1][iHt],allData[cc+1][iHt]-allData[0][0]+5000);
234 }
235 }
236
237 for (Int_t ik = 24; ik<48; ik+=2)
238 for (Int_t iHt=0; iHt<5; iHt++){
239 Int_t cc = ik/2;
240 if((allData[cc+45][iHt]-allData[0][0]+5000)!=0 && allData[cc+45][iHt]>0){
241 hCFD[cc]->Fill(allData[cc+45][iHt]-allData[0][0]+5000);
242 }
243 if(allData[ik+57][iHt]!=0 && allData[ik+58][iHt]!=0 && allData[cc+45][iHt]!=0){
244 hCFDvsQTC[cc]->Fill(allData[ik+57][iHt]-allData[ik+58][iHt],allData[cc+45][iHt]-allData[0][0]+5000);
245 }
246 if(allData[cc+57][iHt]!=0 && allData[cc+45][iHt]!=0){
247 hCFDvsLED[cc]->Fill(allData[cc+57][iHt]-allData[cc+45][iHt],allData[cc+45][iHt]-allData[0][0]+5000);
248 }
249 }
250
251
252 delete start;
253 start = 0x0;
254 reader->Reset();
255 // End of fill histograms
256
257 }
258
259
260
261 /* free resources */
262 free(event);
263
264 /* exit when last event received, no need to wait for TERM signal */
265 if (eventT==END_OF_RUN) {
266 printf("EOR event detected\n");
267 printf("Number of events processed - %i\n ",iev);
268 break;
269 }
270 }
271 printf("After loop, before writing histos\n");
272 // write a file with the histograms
273 TFile *hist = new TFile(FILE_OUT,"RECREATE");
274
275 for(Int_t j=0;j<24;j++){
276 hCFDvsQTC[j]->Write();
277 hCFDvsLED[j]->Write();
278 hCFD[j]->Write();
279 }
280 hist->Close();
281 delete hist;
282
283 status=0;
284
285 /* export file to FXS */
286 if (daqDA_FES_storeFile(FILE_OUT, FILE_OUT)) {
287 status=-2;
288 }
289
290 return status;
291}
292