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