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