]> git.uio.no Git - u/mrichter/AliRoot.git/blame - T0/T0Laserda.cxx
lookup table in txt
[u/mrichter/AliRoot.git] / T0 / T0Laserda.cxx
CommitLineData
99194e07 1//extern "C"
2
3#include <daqDA.h>
4
5#include <event.h>
6#include <monitor.h>
7
8#include <Riostream.h>
9#include <stdio.h>
10#include <stdlib.h>
11
12//AliRoot
13#include <AliRawReaderDate.h>
14#include <AliRawReader.h>
15#include <AliT0RawReader.h>
16#include <AliCDBManager.h>
17
18//ROOT
19#include "TFile.h"
20#include "TKey.h"
21#include "TH2S.h"
22#include "TObject.h"
23#include "TBenchmark.h"
24#include "TRandom.h"
25#include "TMath.h"
26#include "TCanvas.h"
27#include "TString.h"
28#include "TH1.h"
29#include "TF1.h"
30#include "TSpectrum.h"
31#include "TVirtualFitter.h"
32#include "TProfile.h"
33
34/* Main routine
35 Arguments:
36 1- monitoring data source
37*/
38int main(int argc, char **argv) {
39 int status;
40
41
42 if (argc!=2) {
43 printf("Wrong number of arguments\n");
44 return -1;
45 }
46
47
48 /* define data source : this is argument 1 */
49 status=monitorSetDataSource( argv[1] );
50 if (status!=0) {
51 printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
52 return -1;
53 }
54
55
56 /* declare monitoring program */
57 status=monitorDeclareMp( __FILE__ );
58 if (status!=0) {
59 printf("monitorDeclareMp() failed : %s\n",monitorDecodeError(status));
60 return -1;
61 }
62
63
64 /* define wait event timeout - 1s max */
65 monitorSetNowait();
66 monitorSetNoWaitNetworkTimeout(1000);
67
68
69 /* log start of process */
70 printf("T0 monitoring program started\n");
71
72 // Allocation of histograms - start
73 TH1F *hCFD[24];
74 TH2F *hCFD_QTC[24]; TH2F *hCFD_LED[24];
75 Char_t buf1[10], buf2[10],buf3[10],buf4[10];
76
77 for(Int_t ic=0; ic<24; ic++) {
78 sprintf(buf1,"CFD_QTC%i",ic+1);
79 sprintf(buf2,"CFD_LED%i",ic+1);
80
81 hCFD_QTC[ic] = new TH2F(buf1,"CFD_QTC",8000,0.5,8000.5,2000,10000.5,20000.5);
82 hCFD_LED[ic] = new TH2F(buf2,"CFD_LED",1000,-500.0,500.0,100,14600.0,14700.0);
83 if(ic<12){
84 sprintf(buf3,"T0_C_%i_CFD",ic+1);
85 hCFD[ic] = new TH1F(buf3,"CFD", 30000,0., 30000.);
86 }
87 else{
88 sprintf(buf4,"T0_A_%i_CFD",ic-11);
89 hCFD[ic] = new TH1F(buf4,"CFD", 30000,0., 30000.);
90 }
91 }
92
93 // Allocation of histograms - end
94
95 Int_t iev=0;
96 /* main loop (infinite) */
97 for(;;) {
98 struct eventHeaderStruct *event;
99 eventTypeType eventT;
100
101 /* check shutdown condition */
102 if (daqDA_checkShutdown()) {break;}
103
104 /* get next event (blocking call until timeout) */
105 status=monitorGetEventDynamic((void **)&event);
106 if (status==(int)MON_ERR_EOF) {
107 printf ("End of File detected\n");
108 break; /* end of monitoring file has been reached */
109 }
110
111 if (status!=0) {
112 printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
113 break;
114 }
115
116 /* retry if got no event */
117 if (event==NULL) {
118 continue;
119 }
120
121 iev++;
122
123 /* use event - here, just write event id to result file */
124 eventT=event->eventType;
125
126 switch (event->eventType){
127
128 case START_OF_RUN:
129 break;
130
131 case END_OF_RUN:
132 break;
133
134// case CALIBRATION_EVENT:
135 case PHYSICS_EVENT:
136 printf(" event number = %i \n",iev);
137
138
139 // Initalize raw-data reading and decoding
140 AliRawReader *reader = new AliRawReaderDate((void*)event);
141
142 // Enable the following two lines in case of real-data
143 // reader->LoadEquipmentIdsMap("T0map.txt");
144 // reader->RequireHeader(kFALSE);
145 // reader->RequireHeader(kTRUE);
146 AliT0RawReader *start = new AliT0RawReader(reader, kTRUE);
147
148 // Read raw data
149 Int_t allData[105][5];
150 for(Int_t i0=0;i0<105;i0++)
151 for(Int_t j0=0;j0<5;j0++)
152 allData[i0][j0] = 0;
153
154 if(start->Next())
155 for (Int_t i=0; i<105; i++) {
156 for(Int_t iHit=0;iHit<5;iHit++){
157 allData[i][iHit]= start->GetData(i,iHit);
158 }
159 }
160 else
161 printf("No T0 data found!!\n");
162
163 // Fill the histograms
164
165 for (Int_t ik = 0; ik<24; ik+=2)
166 for (Int_t iHt=0; iHt<5; iHt++){
167 Int_t cc = ik/2;
168 if((allData[cc+1][iHt]-allData[0][0])>0){
169 hCFD[cc]->Fill(allData[cc+1][iHt]-allData[0][0]);
170 }
171 if(allData[ik+25][iHt]!=0 && allData[ik+26][iHt]!=0 && allData[cc+1][iHt]!=0){
172 hCFD_QTC[cc]->Fill((allData[ik+25][iHt]-allData[ik+26][iHt]) , (allData[cc+1][iHt]-allData[0][0]));
173 }
174 if(allData[cc+13][iHt]!=0 && allData[cc+1][iHt]!=0){
175 hCFD_LED[cc]->Fill(allData[cc+13][iHt]-allData[cc+1][iHt],allData[cc+1][iHt]-allData[0][0]);
176 }
177 }
178
179 for (Int_t ik = 24; ik<48; ik+=2)
180 for (Int_t iHt=0; iHt<5; iHt++){
181 Int_t cc = ik/2;
182 if((allData[cc+45][iHt]-allData[0][0])>0){
183 hCFD[cc]->Fill(allData[cc+45][iHt]-allData[0][0]);
184 }
185
186 if(allData[ik+57][iHt]!=0 && allData[ik+58][iHt]!=0){
187 if(allData[cc+45][iHt]!=0)
188 hCFD_QTC[cc]->Fill(allData[ik+57][iHt]-allData[ik+58][iHt],allData[cc+45][iHt]-allData[0][0]);
189 }
190 if(allData[cc+57][iHt]!=0 && allData[cc+45][iHt]!=0){
191 hCFD_LED[cc]->Fill(allData[cc+57][iHt]-allData[cc+45][iHt],allData[cc+45][iHt]-allData[0][0]);
192 }
193 }
194
195
196 delete start;
197 start = 0x0;
198 reader->Reset();
199 // End of fill histograms
200
201 }
202
203
204
205 /* free resources */
206 free(event);
207
208 /* exit when last event received, no need to wait for TERM signal */
209 if (eventT==END_OF_RUN) {
210 printf("EOR event detected\n");
211 break;
212 }
213 }
214printf("After loop, before writing histos\n");
215 // write a file with the histograms
216 Char_t filehist[21];
217 sprintf(filehist,"daLaser.root");
218 TFile *hist = new TFile(filehist,"RECREATE");
219
220 for(Int_t j=0;j<24;j++){
221 hCFD_QTC[j]->Write();
222 hCFD_LED[j]->Write();
223 hCFD[j]->Write();
224 }
225 hist->Close();
226 delete hist;
227
228 return status;
229}
230