]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/PHOSLEDda.cxx
Bad design which caused memory leaks in rotation matrices, is corrected.
[u/mrichter/AliRoot.git] / PHOS / PHOSLEDda.cxx
1 /*
2 contact: Boris.Polishchuk@cern.ch
3 link: see comments in the $ALICE_ROOT/PHOS/AliPHOSRcuDA1.cxx
4 reference run: /alice/data/2009/LHC09b_PHOS/000075883/raw/09000075883017.20.root
5 run type: LED
6 DA type: MON 
7 number of events needed: 1000
8 input files: Mod0RCU0.data Mod0RCU1.data Mod0RCU2.data Mod0RCU3.data Mod1RCU0.data Mod1RCU1.data Mod1RCU2.data Mod1RCU3.data Mod2RCU0.data Mod2RCU1.data Mod2RCU2.data Mod2RCU3.data Mod3RCU0.data Mod3RCU1.data Mod3RCU2.data Mod3RCU3.data Mod4RCU0.data Mod4RCU1.data Mod4RCU2.data Mod4RCU3.data zs.txt
9 Output files: PHOS_LED.root.
10 Trigger types used: CALIBRATION_EVENT
11 */
12
13
14 #include "event.h"
15 #include "monitor.h"
16 extern "C" {
17 #include "daqDA.h"
18 }
19
20 #include <stdio.h>
21 #include <stdlib.h>
22
23 #include <TSystem.h>
24 #include <TROOT.h>
25 #include <TPluginManager.h>
26
27 #include "AliRawReader.h"
28 #include "AliRawReaderDate.h"
29 #include "AliPHOSRcuDA1.h"
30 #include "AliPHOSRawFitterv0.h"
31 #include "AliCaloAltroMapping.h"
32 #include "AliCaloRawStreamV3.h"
33 #include "AliLog.h"
34
35
36 /* Main routine
37       Arguments: 
38       1- monitoring data source
39 */
40 int main(int argc, char **argv) {
41
42   gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
43                                         "*",
44                                         "TStreamerInfo",
45                                         "RIO",
46                                         "TStreamerInfo()");
47
48   AliLog::SetGlobalDebugLevel(0) ;
49   AliLog::SetGlobalLogLevel(AliLog::kFatal);
50   
51   int status;
52   
53   if (argc!=2) {
54     printf("Wrong number of arguments\n");
55     return -1;
56   }
57
58   /* Retrieve ZS parameters from DAQ DB */
59   const char* zsfile = "zs.txt";
60   int failZS = daqDA_DB_getFile(zsfile, zsfile);
61   
62   short offset,threshold;
63   
64   if(!failZS) {
65     FILE *f = fopen(zsfile,"r");
66     int scan = fscanf(f,"%d %d",&offset,&threshold);
67   }
68   
69   /* Retrieve mapping files from DAQ DB */ 
70   const char* mapFiles[20] = {
71     "Mod0RCU0.data",
72     "Mod0RCU1.data",
73     "Mod0RCU2.data",
74     "Mod0RCU3.data",
75     "Mod1RCU0.data",
76     "Mod1RCU1.data",
77     "Mod1RCU2.data",
78     "Mod1RCU3.data",
79     "Mod2RCU0.data",
80     "Mod2RCU1.data",
81     "Mod2RCU2.data",
82     "Mod2RCU3.data",
83     "Mod3RCU0.data",
84     "Mod3RCU1.data",
85     "Mod3RCU2.data",
86     "Mod3RCU3.data",
87     "Mod4RCU0.data",
88     "Mod4RCU1.data",
89     "Mod4RCU2.data",
90     "Mod4RCU3.data"
91   };
92
93   for(Int_t iFile=0; iFile<20; iFile++) {
94     int failed = daqDA_DB_getFile(mapFiles[iFile], mapFiles[iFile]);
95     if(failed) { 
96       printf("Cannot retrieve file %s from DAQ DB. Exit.\n",mapFiles[iFile]);
97       return -1;
98     }
99   }
100   
101   /* Open mapping files */
102   AliAltroMapping *mapping[20];
103   TString path = "./";
104
105   path += "Mod";
106   TString path2;
107   TString path3;
108   Int_t iMap = 0;
109
110   for(Int_t iMod = 0; iMod < 5; iMod++) {
111     path2 = path;
112     path2 += iMod;
113     path2 += "RCU";
114
115     for(Int_t iRCU=0; iRCU<4; iRCU++) {
116       path3 = path2;
117       path3 += iRCU;
118       path3 += ".data";
119       mapping[iMap] = new AliCaloAltroMapping(path3.Data());
120       iMap++;
121     }
122   }
123   
124   /* define data source : this is argument 1 */  
125   status=monitorSetDataSource( argv[1] );
126   if (status!=0) {
127     printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
128     return -1;
129   }
130
131
132   /* declare monitoring program */
133   status=monitorDeclareMp( __FILE__ );
134   if (status!=0) {
135     printf("monitorDeclareMp() failed : %s\n",monitorDecodeError(status));
136     return -1;
137   }
138
139
140   /* define wait event timeout - 1s max */
141   monitorSetNowait();
142   monitorSetNoWaitNetworkTimeout(1000);
143   
144    /* init some counters */
145   int nevents_physics=0;
146   int nevents_total=0;
147   
148   AliRawReader *rawReader = NULL;
149   AliPHOSRcuDA1* dAs[5];
150   
151   for(Int_t iMod=0; iMod<5; iMod++) {
152     dAs[iMod] = 0;
153   }
154
155   Float_t e[64][56][2];
156   Float_t t[64][56][2];
157   
158   for(Int_t iX=0; iX<64; iX++) {
159     for(Int_t iZ=0; iZ<56; iZ++) {
160       for(Int_t iGain=0; iGain<2; iGain++) {
161         e[iX][iZ][iGain] = 0.;
162         t[iX][iZ][iGain] = 0.;
163       }
164     }
165   }
166   
167   Int_t cellX    = -1;
168   Int_t cellZ    = -1;
169   Int_t nBunches =  0;
170   Int_t nFired   = -1;
171   Int_t sigStart, sigLength;
172   Int_t caloFlag;
173   
174
175   TH1I fFiredCells("fFiredCells","Number of fired cells per event",100,0,1000);
176
177   /* main loop (infinite) */
178   for(;;) {
179     struct eventHeaderStruct *event;
180     eventTypeType eventT;
181   
182     /* check shutdown condition */
183     if (daqDA_checkShutdown()) {break;}
184     
185     /* get next event (blocking call until timeout) */
186     status=monitorGetEventDynamic((void **)&event);
187     if (status==MON_ERR_EOF) {
188       printf ("End of File detected\n");
189       break; /* end of monitoring file has been reached */
190     }
191     
192     if (status!=0) {
193       printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
194       break;
195     }
196
197     /* retry if got no event */
198     if (event==NULL) {
199       continue;
200     }
201
202
203     /* use event - here, just write event id to result file */
204     eventT=event->eventType;
205     
206     if (eventT==PHYSICS_EVENT) {
207       
208       nFired = 0;
209       
210       rawReader = new AliRawReaderDate((void*)event);
211       AliCaloRawStreamV3 stream(rawReader,"PHOS",mapping);
212       AliPHOSRawFitterv0 fitter;
213       fitter.SubtractPedestals(kTRUE); // assume that data is non-ZS
214
215       if(!failZS) {
216         fitter.SubtractPedestals(kFALSE);
217         fitter.SetAmpOffset(offset);
218         fitter.SetAmpThreshold(threshold);
219       }
220       
221       while (stream.NextDDL()) {
222         while (stream.NextChannel()) {
223           
224           /* Retrieve ZS parameters from data*/
225           if(failZS) {
226             short value = stream.GetAltroCFG1();
227             bool ZeroSuppressionEnabled = (value >> 15) & 0x1;
228             bool AutomaticBaselineSubtraction = (value >> 14) & 0x1;
229             if(ZeroSuppressionEnabled) {
230               offset = (value >> 10) & 0xf;
231               threshold = value & 0x3ff;
232               fitter.SubtractPedestals(kFALSE);
233               fitter.SetAmpOffset(offset);
234               fitter.SetAmpThreshold(threshold);
235             }
236           }
237           
238           cellX    = stream.GetCellX();
239           cellZ    = stream.GetCellZ();
240           caloFlag = stream.GetCaloFlag();  // 0=LG, 1=HG, 2=TRU
241           
242           if(caloFlag!=0 && caloFlag!=1) continue; //TRU data!
243           
244           // In case of oscillating signals with ZS, a channel can have several bunches
245           nBunches = 0;
246           while (stream.NextBunch()) {
247             nBunches++;
248             if (nBunches > 1) continue;
249             sigStart  = stream.GetStartTimeBin();
250             sigLength = stream.GetBunchLength();
251             fitter.SetChannelGeo(stream.GetModule(),cellX,cellZ,caloFlag);
252             fitter.Eval(stream.GetSignals(),sigStart,sigLength);
253           } // End of NextBunch()
254           
255           if(nBunches != 1) continue;
256           
257           e[cellX][cellZ][caloFlag] = fitter.GetEnergy();
258           t[cellX][cellZ][caloFlag] = fitter.GetTime();
259           
260           if(caloFlag==1 && fitter.GetEnergy()>40)
261             nFired++;
262         }
263         
264         if(stream.GetModule()<0 || stream.GetModule()>4) continue;
265         
266         if(dAs[stream.GetModule()])
267           dAs[stream.GetModule()]->FillHistograms(e,t);
268         else {
269           dAs[stream.GetModule()] = new AliPHOSRcuDA1(stream.GetModule(),-1,0);
270           dAs[stream.GetModule()]->FillHistograms(e,t);
271         } 
272         
273         for(Int_t iX=0; iX<64; iX++) {
274           for(Int_t iZ=0; iZ<56; iZ++) {
275             for(Int_t iGain=0; iGain<2; iGain++) {
276               e[iX][iZ][iGain] = 0.;
277               t[iX][iZ][iGain] = 0.;
278             }
279           }
280         }
281         
282       }
283       
284       fFiredCells.Fill(nFired);
285       
286       delete rawReader;     
287       nevents_physics++;
288     }
289     
290     nevents_total++;
291     
292     /* free resources */
293     free(event);
294     
295     /* exit when last event received, no need to wait for TERM signal */
296     if (eventT==END_OF_RUN) {
297       printf("EOR event detected\n");
298       break;
299     }
300   }
301   
302   for(Int_t i = 0; i < 20; i++) delete mapping[i];  
303   
304   /* Be sure that all histograms are saved */
305   
306   const TH2F* h2=0;
307   const TH1F* h1=0;
308   char localfile[128];
309
310   Int_t nGood=0;    // >10 entries in peak
311   Int_t nMax=-111;  // max. number of entries in peak
312   Int_t iXmax=-1;
313   Int_t iZmax=-1;
314   Int_t iModMax=-1;
315   
316   sprintf(localfile,"PHOS_LED.root");
317   TFile* f = new TFile(localfile,"recreate");
318
319   for(Int_t iMod=0; iMod<5; iMod++) {
320     if(!dAs[iMod]) continue;
321   
322     printf("DA1 for module %d detected.\n",iMod);
323     
324     for(Int_t iX=0; iX<64; iX++) {
325       for(Int_t iZ=0; iZ<56; iZ++) {
326         
327         h1 = dAs[iMod]->GetHgLgRatioHistogram(iX,iZ); // High Gain/Low Gain ratio
328         if(h1) {
329           if(h1->GetMaximum()>10.) nGood++;
330           if(h1->GetMaximum()>nMax) {
331             nMax = (Int_t)h1->GetMaximum(); iXmax=iX; iZmax=iZ; iModMax=iMod;
332           }
333           h1->Write(); 
334         }
335       
336         for(Int_t iGain=0; iGain<2; iGain++) {
337           h2 = dAs[iMod]->GetTimeEnergyHistogram(iX,iZ,iGain); // Time vs Energy
338           if(h2) h2->Write();
339         }
340         
341       }
342     }
343     
344     fFiredCells.Write();    
345   }
346   
347   f->Close();
348   
349   /* Store output files to the File Exchange Server */
350   daqDA_FES_storeFile(localfile,"LED");
351   
352   printf("%d physics events of %d total processed.\n",nevents_physics,nevents_total);
353   printf("%d histograms has >10 entries in maximum, max. is %d entries ",nGood,nMax);
354   printf("(module,iX,iZ)=(%d,%d,%d)",iModMax,iXmax,iZmax);
355   printf("\n");
356
357   return status;
358 }