]> git.uio.no Git - u/mrichter/AliRoot.git/blame_incremental - TOF/TOFda.cxx
Protection added against no entries to compute the mean for pulser run.
[u/mrichter/AliRoot.git] / TOF / TOFda.cxx
... / ...
CommitLineData
1/*
2
3TOF DA for online calibration
4
5Contact: Chiara.Zampolli@bo.infn.it
6Link: www.bo.infn.it/~zampolli
7Run Type: PHYSICS
8DA Type: MON
9Number of events needed: depending on the run, being run-level
10Input Files: TOFdaTotal.root, to be updated if existing
11Output Files: TOFdaRun.root, TOFdaTotal.root, both to be exported to the DAQ FXS
12Trigger types used: PHYSICS_EVENT
13
14*/
15
16#define FILE_TOTAL "TOFdaTotal.root"
17#define FILE_RUN "TOFdaRun.root"
18
19// DATE
20#include <daqDA.h>
21#include <event.h>
22#include <monitor.h>
23
24//AliRoot
25#include <AliTOFRawStream.h>
26#include <AliRawReaderDate.h>
27#include <AliRawReader.h>
28#include <AliTOFGeometry.h>
29#include <AliT0RawReader.h>
30#include <AliDAQ.h>
31#include <AliTOFHitData.h>
32#include <AliTOFHitDataBuffer.h>
33
34//ROOT
35#include <TFile.h>
36#include <TKey.h>
37#include <TH2S.h>
38#include <TObject.h>
39#include <TMath.h>
40#include <TSystem.h>
41#include "TROOT.h"
42#include "TPluginManager.h"
43
44
45/* Main routine
46 Arguments:
47 1- monitoring data source
48*/
49
50int main(int argc, char **argv) {
51
52 /* magic line from Rene */
53 gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
54 "*",
55 "TStreamerInfo",
56 "RIO",
57 "TStreamerInfo()");
58
59 AliTOFGeometry * geom = new AliTOFGeometry();
60
61 static const Int_t size = AliTOFGeometry::NPadXSector()*AliTOFGeometry::NSectors();
62 static const Int_t nbins = 500;
63 static const Int_t binmin = -20;
64 const Float_t c = 2.99792458E10; //speed of light
65 TH1F::AddDirectory(0);
66 TH2S * htofPartial = new TH2S("htof","histo with delays", size,-0.5,size*1.-0.5,nbins,binmin-0.5,nbins*1.+binmin-0.5);
67
68 // decoding the events
69
70 int status;
71 if (argc!=2) {
72 printf("Wrong number of arguments\n");
73 return -1;
74 }
75
76 /* define data source : this is argument 1 */
77 status=monitorSetDataSource( argv[1] );
78 if (status!=0) {
79 printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
80 return -1;
81 }
82
83 /* declare monitoring program */
84 status=monitorDeclareMp( __FILE__ );
85 if (status!=0) {
86 printf("monitorDeclareMp() failed : %s\n",monitorDecodeError(status));
87 return -1;
88 }
89
90 /* define wait event timeout - 1s max */
91 monitorSetNowait();
92 monitorSetNoWaitNetworkTimeout(1000);
93
94 /* log start of process */
95 printf("TOF DA started\n");
96
97 /* init some counters */
98 int nevents_physics=0;
99 int nevents_total=0;
100
101 struct eventHeaderStruct *event;
102 eventTypeType eventT;
103 Int_t iev=0;
104
105 /* main loop (infinite) */
106 for(;;) {
107
108 /* check shutdown condition */
109 if (daqDA_checkShutdown()) {break;}
110
111 /* get next event (blocking call until timeout) */
112 status=monitorGetEventDynamic((void **)&event);
113 if (status==MON_ERR_EOF) {
114 printf ("End of File detected\n");
115 break; /* end of monitoring file has been reached */
116 }
117
118 if (status!=0) {
119 printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
120 break;
121 }
122
123 /* retry if got no event */
124 if (event==NULL) {
125 continue;
126 }
127
128 iev++;
129
130 /* use event - here, just write event id to result file */
131 nevents_total++;
132 eventT=event->eventType;
133 switch (event->eventType){
134
135 /* START OF RUN */
136 case START_OF_RUN:
137 break;
138 /* END START OF RUN */
139
140 /* END OF RUN */
141 case END_OF_RUN:
142 break;
143
144 case PHYSICS_EVENT:
145 nevents_physics++;
146 AliRawReader *rawReader = new AliRawReaderDate((void*)event);
147 //rawReader->RequireHeader(kFALSE);
148
149 //T0 event
150 Int_t meantime = 0;
151 AliT0RawReader *rawReaderT0 = new AliT0RawReader(rawReader,kTRUE);
152 if (!rawReaderT0->Next()) {
153 printf("T0: no raw data found!\n");
154 }
155 else {
156 /*
157 Int_t allData[105][5];
158 for (Int_t i=0; i<105; i++) {
159 allData[i][0]=rawReaderT0->GetData(i,0);
160 }
161 meantime = allData[49][0];
162 */
163 meantime = rawReaderT0->GetData(49,0);
164 // printf("time zero (ns) = %i (%f) \n", meantime, (meantime*24.4-200)*1E-3); // debugging purpose
165 }
166
167 delete rawReaderT0;
168 rawReaderT0 = 0x0;
169 rawReader->Reset();
170
171 //TOF event
172 Int_t dummy = -1;
173 Int_t Volume[5];
174 AliTOFHitData *HitData;
175 AliTOFHitDataBuffer *DataBuffer;
176 AliTOFHitDataBuffer *PackedDataBuffer;
177 AliTOFRawStream *rawStreamTOF = new AliTOFRawStream(rawReader);
178 // rawReader->ReadHeader();
179 rawStreamTOF->ResetBuffers();
180 rawStreamTOF->DecodeDDL(0, AliDAQ::NumberOfDdls("TOF") - 1,0);
181 Int_t nPDBEntriesToT = 0;
182 Int_t nDBEntriesToT = 0;
183 for (Int_t iDDL = 0; iDDL < AliDAQ::NumberOfDdls("TOF"); iDDL++){
184
185 /* read decoded data */
186 DataBuffer = rawStreamTOF->GetDataBuffer(iDDL);
187 PackedDataBuffer = rawStreamTOF->GetPackedDataBuffer(iDDL);
188
189 /* get buffer entries */
190 Int_t nDBEntries = DataBuffer->GetEntries();
191 Int_t nPDBEntries = PackedDataBuffer->GetEntries();
192 nPDBEntriesToT+=nPDBEntries;
193 nDBEntriesToT+=nDBEntries;
194
195 //for (Int_t iHit = 0; iHit < nDBEntries; iHit++){
196 // HitData = DataBuffer->GetHit(iHit);
197 /* store volume information */
198 // rawStreamTOF->EquipmentId2VolumeId(HitData, Volume);
199 //}
200 /* reset buffer */
201 DataBuffer->Reset();
202
203 /* read data buffer hits */
204 for (Int_t iHit = 0; iHit < nPDBEntries; iHit++){
205 HitData = PackedDataBuffer->GetHit(iHit);
206 /* add volume information */
207 HitData->SetDDLID(iDDL);
208 rawStreamTOF->EquipmentId2VolumeId(HitData, Volume);
209 if (Volume[0]==-1 ||
210 Volume[1]==-1 ||
211 Volume[2]==-1 ||
212 Volume[3]==-1 ||
213 Volume[4]==-1) continue;
214 else {
215 dummy = Volume[3];
216 Volume[3] = Volume[4];
217 Volume[4] = dummy;
218 Int_t tof = (Int_t)((Double_t)HitData->GetTime()*1E3/AliTOFGeometry::TdcBinWidth());
219 Int_t index = geom->GetIndex(Volume);
220 Float_t pos[3];
221 geom->GetPosPar(Volume,pos);
222 Float_t texp=TMath::Sqrt(pos[0]*pos[0]+pos[1]*pos[1]+pos[2]*pos[2])/c*1E9; //expected time in ns
223 Float_t texpBin=(texp*1E3-32)/AliTOFGeometry::TdcBinWidth(); //expected time in number of TDC bin
224 Int_t deltabin = tof-TMath::Nint(texpBin); //to be used with real data; rounding expected time to Int_t
225 htofPartial->Fill(index,deltabin); //channel index start from 0, bin index from 1
226 //debugging printings
227 //printf("sector %i, plate %i, strip %i, padz %i, padx %i \n",Volume[0],Volume[1],Volume[2],Volume[3],Volume[4]);
228 //printf("pos x = %f, pos y = %f, pos z = %f \n",pos[0],pos[1],pos[2]);
229 //printf ("expected time = %f (ns)\n",texp);
230 //printf ("expected time bin = %f (TDC bin)\n",texpBin);
231 //printf ("measured time bin = %i (TDC bin) with %f (ns) and ACQ bit = %i \n",tof, HitData->GetTime(), HitData->GetACQ());
232 // printf("index = %i, deltabin = %i , filling index = %i, and bin = % i\n",index, deltabin, index, deltabin);
233
234 }
235 /* reset buffer */
236 PackedDataBuffer->Reset();
237 }
238 }
239 //printf(" Packed Hit Buffer Entries = %i \n",nPDBEntriesToT);
240 //printf(" Hit Buffer Entries = %i \n",nDBEntriesToT);
241
242 delete rawStreamTOF;
243 rawStreamTOF = 0x0;
244
245 delete rawReader;
246 rawReader = 0x0;
247 }
248
249 /* free resources */
250 free(event);
251
252 /* exit when last event received, no need to wait for TERM signal */
253 if (eventT==END_OF_RUN) {
254 printf("EOR event detected\n");
255 break;
256 }
257 }
258
259 delete geom;
260 geom = 0x0;
261
262 //write the Run level file
263 TFile * fileRun = new TFile (FILE_RUN,"RECREATE");
264 htofPartial->Write();
265 fileRun->Close();
266
267 //write the Total file
268 TH2S *htoftot = 0x0;
269 TFile * filetot = 0x0;
270 Bool_t isThere=kFALSE;
271 const char *dirname = "./";
272 TString filename = FILE_TOTAL;
273 if((gSystem->FindFile(dirname,filename))!=NULL){
274 isThere=kTRUE;
275 printf("%s found \n",FILE_TOTAL);
276 }
277 if(isThere){
278
279 TFile * filetot1 = new TFile (FILE_TOTAL,"READ");
280 //look for the file
281 if (!filetot1->IsZombie()){
282 printf("updating file %s \n",FILE_TOTAL);
283 TIter next(filetot1->GetListOfKeys());
284 TKey *key;
285 //look for the histogram
286 while ((key=(TKey*)next())){
287 const char * namekey = key->GetName();
288 if (strcmp(namekey,"htoftot")==0){
289 printf(" histo found \n");
290 htoftot = (TH2S*) filetot1->Get("htoftot");
291 htoftot->AddDirectory(0);
292 htoftot->Add(htofPartial);
293 break;
294 }
295 }
296 }
297 filetot1->Close();
298 delete filetot1;
299 filetot1=0x0;
300 }
301 else {
302 printf(" no %s file found \n",FILE_TOTAL);
303 htoftot = new TH2S(*htofPartial);
304 htoftot->SetName("htoftot");
305 htoftot->AddDirectory(0);
306 }
307
308 filetot = new TFile (FILE_TOTAL,"RECREATE");
309 filetot->cd();
310 htoftot->Write();
311 filetot->Close();
312
313 delete fileRun;
314 delete filetot;
315 delete htofPartial;
316 delete htoftot;
317
318 fileRun = 0x0;
319 filetot = 0x0;
320 htofPartial = 0x0;
321 htoftot = 0x0;
322
323 /* write report */
324 printf("Run #%s, received %d physics events out of %d\n",getenv("DATE_RUN_NUMBER"),nevents_physics,nevents_total);
325
326 status=0;
327
328 /* export file to FXS */
329 if (daqDA_FES_storeFile(FILE_RUN, "RUNLevel")) {
330 status=-2;
331 }
332 if (daqDA_FES_storeFile(FILE_TOTAL, "DELAYS")) {
333 status=-2;
334 }
335
336
337 return status;
338}
339