]>
Commit | Line | Data |
---|---|---|
a84a9142 | 1 | /* |
2 | T0 DA for online calibration | |
271b5bd3 | 3 | |
a84a9142 | 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) { | |
271b5bd3 | 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 | |
a84a9142 | 82 | } |
83 | } | |
84 | fclose(inp); | |
85 | ||
86 | if (argc!=2) { | |
87 | printf("Wrong number of arguments\n"); | |
88 | return -1; | |
89 | } | |
90 | ||
91 | ||
92 | /* define data source : this is argument 1 */ | |
93 | status=monitorSetDataSource( argv[1] ); | |
94 | if (status!=0) { | |
95 | printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status)); | |
96 | return -1; | |
97 | } | |
98 | ||
99 | ||
100 | /* declare monitoring program */ | |
101 | status=monitorDeclareMp( __FILE__ ); | |
102 | if (status!=0) { | |
103 | printf("monitorDeclareMp() failed : %s\n",monitorDecodeError(status)); | |
104 | return -1; | |
105 | } | |
106 | ||
107 | ||
108 | /* define wait event timeout - 1s max */ | |
109 | monitorSetNowait(); | |
110 | monitorSetNoWaitNetworkTimeout(1000); | |
111 | ||
112 | ||
113 | /* log start of process */ | |
114 | printf("T0 monitoring program started\n"); | |
115 | ||
116 | // Allocation of histograms - start | |
a84a9142 | 117 | TH1F *hCFD1minCFD[24]; |
a84a9142 | 118 | |
119 | for(Int_t ic=0; ic<24; ic++) { | |
a84a9142 | 120 | hCFD1minCFD[ic] = new TH1F(Form("CFD1-CFD%d",ic+1),"CFD-CFD",ccbx,cclx,ccmx); |
a84a9142 | 121 | } |
122 | ||
123 | // Allocation of histograms - end | |
124 | ||
125 | Int_t iev=0; | |
126 | /* main loop (infinite) */ | |
127 | for(;;) { | |
128 | struct eventHeaderStruct *event; | |
129 | eventTypeType eventT; | |
130 | ||
131 | /* check shutdown condition */ | |
132 | if (daqDA_checkShutdown()) {break;} | |
133 | ||
134 | /* get next event (blocking call until timeout) */ | |
135 | status=monitorGetEventDynamic((void **)&event); | |
136 | if (status==(int)MON_ERR_EOF) { | |
137 | printf ("End of File detected\n"); | |
138 | break; /* end of monitoring file has been reached */ | |
139 | } | |
140 | ||
141 | if (status!=0) { | |
142 | printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status)); | |
143 | break; | |
144 | } | |
145 | ||
146 | /* retry if got no event */ | |
147 | if (event==NULL) { | |
148 | continue; | |
149 | } | |
150 | ||
151 | // iev++; | |
152 | ||
153 | /* use event - here, just write event id to result file */ | |
154 | eventT=event->eventType; | |
155 | ||
156 | switch (event->eventType){ | |
157 | ||
158 | case START_OF_RUN: | |
159 | break; | |
160 | ||
161 | case END_OF_RUN: | |
162 | break; | |
163 | ||
164 | case CALIBRATION_EVENT: | |
165 | // case PHYSICS_EVENT: | |
166 | iev++; | |
167 | ||
168 | if(iev==1){ | |
169 | printf("First event - %i\n",iev); | |
170 | } | |
171 | ||
172 | // Initalize raw-data reading and decoding | |
173 | AliRawReader *reader = new AliRawReaderDate((void*)event); | |
174 | ||
175 | // Enable the following two lines in case of real-data | |
176 | reader->RequireHeader(kTRUE); | |
177 | AliT0RawReader *start = new AliT0RawReader(reader, kTRUE); | |
178 | ||
179 | // Read raw data | |
180 | Int_t allData[105][5]; | |
181 | for(Int_t i0=0;i0<105;i0++) | |
182 | for(Int_t j0=0;j0<5;j0++) | |
183 | allData[i0][j0] = 0; | |
271b5bd3 | 184 | |
185 | if(start->Next()){ | |
186 | for (Int_t i=0; i<105; i++) { | |
a84a9142 | 187 | for(Int_t iHit=0;iHit<5;iHit++){ |
188 | allData[i][iHit]= start->GetData(i,iHit); | |
189 | } | |
271b5bd3 | 190 | } |
a84a9142 | 191 | } |
271b5bd3 | 192 | else |
193 | printf("No T0 data found!!!\n"); | |
a84a9142 | 194 | |
195 | // Fill the histograms | |
196 | ||
271b5bd3 | 197 | for (Int_t ik = 0; ik<24; ik++) |
a84a9142 | 198 | for (Int_t iHt=0; iHt<5; iHt++){ |
271b5bd3 | 199 | if(allData[ik+1][iHt]!=0 ){ |
200 | if(ik<12){ | |
201 | hCFD1minCFD[ik]->Fill(allData[ik+1][iHt]-allData[1][iHt]); | |
202 | } | |
203 | if(ik>11){ | |
204 | hCFD1minCFD[ik]->Fill(allData[ik+45][iHt]-allData[57][iHt]); | |
205 | } | |
a84a9142 | 206 | } |
a84a9142 | 207 | } |
208 | ||
209 | delete start; | |
210 | start = 0x0; | |
211 | reader->Reset(); | |
212 | // End of fill histograms | |
213 | ||
214 | } | |
215 | ||
216 | /* free resources */ | |
217 | free(event); | |
218 | ||
219 | /* exit when last event received, no need to wait for TERM signal */ | |
220 | if (eventT==END_OF_RUN) { | |
221 | printf("EOR event detected\n"); | |
222 | printf("Number of events processed - %i\n ",iev); | |
223 | break; | |
224 | } | |
225 | } | |
226 | printf("After loop, before writing histos\n"); | |
227 | // write a file with the histograms | |
228 | TFile *hist = new TFile(FILE_OUT,"RECREATE"); | |
229 | ||
230 | for(Int_t j=0;j<24;j++){ | |
a84a9142 | 231 | hCFD1minCFD[j]->Write(); |
232 | } | |
233 | hist->Close(); | |
234 | delete hist; | |
235 | ||
236 | status=0; | |
237 | ||
238 | /* export file to FXS */ | |
239 | if (daqDA_FES_storeFile(FILE_OUT, "COSMIC")) { | |
240 | status=-2; | |
241 | } | |
242 | ||
243 | return status; | |
244 | } | |
245 |