]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/TOFda.cxx
Fixed memory leaks for #86360: High memory consumption in 2.76TeV p+p RAW reco jobs
[u/mrichter/AliRoot.git] / TOF / TOFda.cxx
1 /*
2
3 TOF DA for online calibration
4
5 Contact: Chiara.Zampolli@bo.infn.it
6          Roberto.Preghenella@bo.infn.it
7
8 Run Type: PHYSICS
9 DA Type: MON
10 Number of events needed:
11 Input Files: no input
12 Output Files: TOFdaHits.root
13 Event types used: PHYSICS_EVENT
14
15 */
16
17 #define FILE_HITS "TOFdaHits.root"
18 #define FILE_CALIB "TOFdaCalib.root"
19
20 // DATE
21 #include "event.h"
22 #include "monitor.h"
23 #include "daqDA.h"
24
25 #include <stdio.h>
26 #include <stdlib.h>
27 #include <errno.h>
28
29 //AliRoot
30 #include "TROOT.h"
31 #include "AliTOFRawStream.h"
32 #include "AliRawReaderDate.h"
33 #include "AliRawReader.h"
34 #include "AliDAQ.h"
35 #include "AliTOFHitData.h"
36 #include "AliTOFHitDataBuffer.h"
37 #include "AliTOFDaConfigHandler.h"
38 #include "AliTOFHitField.h"
39 #include "AliLog.h"
40 #include "AliTOFGeometry.h"
41
42 //ROOT
43 #include "TFile.h"
44 #include "TKey.h"
45 #include "TH2S.h"
46 #include "TObject.h"
47 #include "TMath.h"
48 #include "TSystem.h"
49 #include "TROOT.h"
50 #include "TPluginManager.h"
51 #include "TSAXParser.h"
52 #include "TTree.h"
53
54 /* Main routine
55       Arguments: 
56       1- monitoring data source
57 */
58 int 
59 main(int argc, char **argv) 
60 {
61   
62   /* magic line from Rene */
63   gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
64                                         "*",
65                                         "TStreamerInfo",
66                                         "RIO",
67                                         "TStreamerInfo()");
68   
69
70   /* log start of process */
71   printf("TOF DA started\n");  
72   
73   /* check that we got some arguments = list of files */
74   if (argc!=2) {
75     printf("Wrong number of arguments\n");
76     return -1;
77   }
78
79   /*
80    * CONFIG
81    */
82   
83   /* retrieve config file */
84   int getConfigFile = daqDA_DB_getFile("TOFPhysicsConfig.xml","TOFPhysicsConfig.xml");
85   if (getConfigFile != 0){
86     printf("Failed to retrieve config file from DB! returning...\n");
87     return -1;
88   }
89   /* parse config file */
90   AliTOFDaConfigHandler* tofHandler = new AliTOFDaConfigHandler();
91   TSAXParser *parser = new TSAXParser();
92   parser->ConnectToHandler("AliTOFDaConfigHandler", tofHandler);  
93   if (parser->ParseFile("./TOFPhysicsConfig.xml") != 0) {
94     printf("Failed parsing config file! retunring... \n");
95     return -1;
96   }
97   /* setup config params */
98   Int_t meanMultiplicity = tofHandler->GetMeanMultiplicity(); /* average expected TOF multiplicity */
99   Int_t maxHits = tofHandler->GetMaxHits(); /* max number of hits to be collected */
100   printf("current settings:\n");
101   printf(" - meanMultiplicity = %d\n", meanMultiplicity);
102   printf(" - maxHits          = %d\n", maxHits);
103   /* constants */
104   const Int_t nChannels = 157248;
105   Int_t noiseCheckTrigger = 10; /* first noise check after 10 events */
106   Float_t meanChannelRate = (Float_t)meanMultiplicity / (Float_t)nChannels; /* average expected channel rate (hits/event) */
107   Float_t noiseThreshold = 10. * meanChannelRate; /* noise threshold (hits/event) */
108   Int_t minNoiseHits = 10; /* min number of channel hits to check noise */
109   /* counters and flags */
110   Int_t nPhysicsEvents, nCalibEvents, totHits;
111   Int_t nChHits[nChannels];
112   Bool_t inhibitCollection;
113   Bool_t noiseFlag[nChannels];
114   /* variables */
115   Int_t nhits, ddl, slot, trm, chain, tdc, channel, index, timebin, totbin, deltaBC, l0l1latency, det[5], dummy;
116   Float_t noiseHitThreshold;
117
118   /*
119    * INIT 
120    */
121
122   /* init counters and flags */
123   nPhysicsEvents = 0;
124   nCalibEvents = 0;
125   totHits = 0;
126   inhibitCollection = kFALSE;
127   for (Int_t ich = 0; ich < nChannels; ich++) {
128     nChHits[ich] = 0;
129     noiseFlag[ich] = kFALSE;
130   }
131
132   /* TOF raw data handling */
133   AliTOFRawStream *rawStream = new AliTOFRawStream();
134   AliTOFHitDataBuffer *pdb = NULL;
135   AliTOFHitData *hit = NULL;
136   
137   /* open HITS output file */
138   TFile *fileOutHits = new TFile(FILE_HITS, "RECREATE"); 
139   /* create hit field data structure */
140   AliTOFHitField *hitField = new AliTOFHitField();
141   /* create temporary tree */
142   TTree *tempTree = new TTree("tempTree", "temporary tree");
143   tempTree->Branch("hit", "AliTOFHitField", &hitField);
144   /* create output tree */
145   TTree *outTree = new TTree("hitTree", "hit tree");
146   outTree->Branch("hit", "AliTOFHitField", &hitField);
147
148   /* open CALIB output file */
149   TFile *fileOutCalib = new TFile(FILE_CALIB, "RECREATE"); 
150   /* create calib hit histo */
151   TH1F *hCalibHit = new TH1F("hCalibHit", "Calibration events;index;N_{hits}/N_{events}", nChannels, 0., nChannels);
152
153   /*
154    * ONLINE MONITOR
155    */
156
157   AliLog::SetGlobalLogLevel(AliLog::kFatal);
158   struct eventHeaderStruct *event;
159   int ret;
160   /* define data source : this is argument 1 */  
161   ret = monitorSetDataSource(argv[1]);
162   if (ret != 0) {
163     printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(ret));
164     return -1;
165   }
166   /* declare monitoring program */
167   ret = monitorDeclareMp("tofDA");
168   if (ret != 0) {
169     printf("monitorDeclareMp() failed : %s\n",monitorDecodeError(ret));
170     return -1;
171   }
172   /* define wait event timeout - 1s max */
173   monitorSetNowait();
174   monitorSetNoWaitNetworkTimeout(1000);
175
176   /* loop over events */
177   while (1) {
178     
179     /* check shutdown condition */
180     if (daqDA_checkShutdown()) break;
181     
182     /*
183      * NOISE CHECK
184      */
185
186     /* check inhibit collection */
187     if (!inhibitCollection) {
188       /* check number of events and check noise */
189       if (nPhysicsEvents >= noiseCheckTrigger || totHits >= maxHits) {
190         noiseHitThreshold = noiseThreshold * nPhysicsEvents;
191         printf("noise check triggered after %d events: threshold is %f hits\n", nPhysicsEvents, noiseHitThreshold);
192         /* loop over all channels */
193         for (Int_t ich = 0; ich < nChannels; ich++) {
194           /* check */
195           if (nChHits[ich] < minNoiseHits || noiseFlag[ich] || nChHits[ich] < noiseHitThreshold) continue;
196           printf("channel %d tagged as noisy (%d hits): disabled\n", ich, nChHits[ich]);
197           noiseFlag[ich] = kTRUE;
198           totHits -= nChHits[ich];
199         } /* end of loop over all channels */
200         /* set new noise check trigger value */
201         noiseCheckTrigger *= 10;
202       } /* end of noise check */    
203       
204       /* inhibit hit collection when maximum number of hits exceeded */
205       if (totHits >= maxHits) {
206         printf("maximum number of hits exceeded (%d): inhibit hit collection\n", maxHits);
207         inhibitCollection = kTRUE;
208       }
209     }
210     
211     /*
212      * GET EVENT
213      */
214
215     /* get next event (blocking call until timeout) */
216     ret = monitorGetEventDynamic((void **)&event);
217     if (ret == MON_ERR_EOF) {
218       printf ("End of File detected\n");
219       break; /* end of monitoring file has been reached */
220     }
221     if (ret != 0) {
222       printf("monitorGetEventDynamic() failed (ret=%d errno=%d): %s\n", ret, errno, monitorDecodeError(ret));
223       break;
224     }
225     /* retry if got no event */
226     if (event==NULL) continue;
227     /* check event type */
228     if (event->eventType != PHYSICS_EVENT && event->eventType != CALIBRATION_EVENT) {
229       free(event);
230       continue;
231     }
232     /* check inhibit collection */
233     if (event->eventType == PHYSICS_EVENT && inhibitCollection) {
234       free(event);
235       continue;
236     }
237     /* increment number of physics events */
238     if (event->eventType == PHYSICS_EVENT) nPhysicsEvents++;
239     /* increment number of calib events */
240     if (event->eventType == CALIBRATION_EVENT) nCalibEvents++;
241     
242     /*
243      * DECODE EVENT
244      */
245
246     /* create raw reader */
247     AliRawReader *rawReader = new AliRawReaderDate((void *)event);
248     /* setup raw stream */
249     rawStream->SetRawReader(rawReader);
250     /* reset buffers */
251     rawStream->ResetBuffers();
252     /* decode */
253     rawStream->DecodeDDL(0, AliDAQ::NumberOfDdls("TOF") - 1, 0);
254
255     /*
256      * HIT MANIPULATION
257      */
258
259     /* loop over DDLs */
260     for (Int_t iddl = 0; iddl < AliDAQ::NumberOfDdls("TOF"); iddl++) {
261       /* get packed-data buffer */
262       pdb = rawStream->GetPackedDataBuffer(iddl);
263       nhits = pdb->GetEntries();
264       /* loop over hits in buffer */
265       for (Int_t ihit = 0; ihit < nhits; ihit++) {
266         /* get hit */
267         hit = pdb->GetHit(ihit);
268         /* get channel info */
269         ddl = iddl;
270         slot = hit->GetSlotID();
271         trm = slot - 3;
272         chain = hit->GetChain();
273         tdc = hit->GetTDC();
274         channel = hit->GetChan();
275         /* get index */
276         rawStream->EquipmentId2VolumeId(ddl, slot, chain, tdc, channel, det);
277         dummy = det[4];
278         det[4] = det[3];
279         det[3] = dummy;
280         /* check valid index */
281         if (det[0] < 0 || det[0] > 17 ||
282             det[1] < 0 || det[1] > 5 ||
283             det[2] < 0 || det[2] > 18 ||
284             det[3] < 0 || det[3] > 1 ||
285             det[4] < 0 || det[4] > 47) continue;
286         index = AliTOFGeometry::GetIndex(det);
287
288         /* switch event type */
289         switch (event->eventType) {
290
291           /*
292            * PHYSICS EVENT
293            */
294
295         case PHYSICS_EVENT:
296           /* check noise flag */
297           if (noiseFlag[index]) continue;
298           /* increment number of channel hits and total hits */
299           nChHits[index]++;
300           totHits++;
301           /* get signal info */
302           timebin = hit->GetTimeBin();
303           totbin = hit->GetTOTBin();
304           deltaBC = hit->GetDeltaBunchID();
305           l0l1latency = hit->GetL0L1Latency();
306           /* set hit field data */
307           hitField->SetIndex(index);
308           hitField->SetTimeBin(timebin);
309           hitField->SetTOTBin(totbin);
310           hitField->SetDeltaBC(deltaBC);
311           hitField->SetL0L1Latency(l0l1latency);
312           /* fill temp tree */
313           tempTree->Fill();
314           break;
315
316           /*
317            * CALIBRATION EVENT
318            */
319
320         case CALIBRATION_EVENT:
321           /* fill calib hit histo */
322           hCalibHit->Fill(index);
323           break;
324           
325         } /* end of switch event type */
326       } /* end of loop over hits in buffer */
327     } /* end of loop over DDLs */
328     
329     /* delete raw reader */
330     delete rawReader;
331     /* free event */
332     free(event);
333     
334   } /* end of loop over events */
335
336   /* final noise check */
337   noiseHitThreshold = noiseThreshold * nPhysicsEvents;
338   printf("final noise check after %d events: threshold is %f hits\n", nPhysicsEvents, noiseHitThreshold);
339   /* loop over all channels */
340   for (Int_t ich = 0; ich < nChannels; ich++) {
341     /* check */
342     if (nChHits[ich] < minNoiseHits || noiseFlag[ich] || nChHits[ich] < noiseHitThreshold) continue;
343     printf("channel %d tagged as noisy (%d hits): disabled\n", ich, nChHits[ich]);
344     noiseFlag[ich] = kTRUE;
345     totHits -= nChHits[ich];
346   } /* end of loop over all channels */
347   
348   /* copy hits into output tree from temp tree */
349   printf("copy hits from temporary tree into output tree\n");
350   printf("temporary tree contains %d hits\n", (Int_t)tempTree->GetEntries());
351   for (Int_t ihit = 0; ihit < tempTree->GetEntries(); ihit++) {
352     /* get entry */
353     tempTree->GetEntry(ihit);
354     /* check noise flag */
355     if (noiseFlag[hitField->GetIndex()]) continue;
356     /* fill output tree */
357     outTree->Fill();
358   } /* end of copy hits into output tree from temp tree */
359   printf("output tree contains %d hits\n", (Int_t)outTree->GetEntries());
360
361   /* write output tree on HITS file */
362   fileOutHits->cd();
363   outTree->Write();
364   fileOutHits->Close();
365   /* export file to FXS */
366   if (daqDA_FES_storeFile(FILE_HITS, "HITS"))
367     return -2;
368
369   /* scale calib hit histo by number of calib events */
370   printf("found %d calibration events\n", nCalibEvents);
371   hCalibHit->Sumw2();
372   if (nCalibEvents > 0)
373     hCalibHit->Scale(1. / nCalibEvents);
374   
375   /* write calib hit histo on CALIB file */
376   fileOutCalib->cd();
377   hCalibHit->Write();
378   fileOutCalib->Close();
379   /* export file to FXS */
380   if (daqDA_FES_storeFile(FILE_CALIB, "CALIB"))
381     return -2;
382
383   return 0;
384 }