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