]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/PHOSLEDda.cxx
Mapping error 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: /castor/cern.ch/alice/phos/2007/10/04/18/07000008249001.1000.root
5 run type: STANDALONE
6 DA type: MON 
7 number of events needed: 1000
8 input files: RCU0.data  RCU1.data  RCU2.data  RCU3.data
9 Output files: PHOS_Module2_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
34
35 /* Main routine
36       Arguments: 
37       1- monitoring data source
38 */
39 int main(int argc, char **argv) {
40
41   gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
42                                         "*",
43                                         "TStreamerInfo",
44                                         "RIO",
45                                         "TStreamerInfo()");
46
47   int status;
48   
49   if (argc!=2) {
50     printf("Wrong number of arguments\n");
51     return -1;
52   }
53
54   /* Retrieve ZS parameters from DAQ DB */
55   const char* zsfile = "zs.txt";
56   int failZS = daqDA_DB_getFile(zsfile, zsfile);
57   
58   Int_t offset,threshold;
59   
60   if(!failZS) {
61     FILE *f = fopen(zsfile,"r");
62     int scan = fscanf(f,"%d %d",&offset,&threshold);
63   }
64
65   /* Retrieve mapping files from DAQ DB */ 
66   const char* mapFiles[4] = {"RCU0.data","RCU1.data","RCU2.data","RCU3.data"};
67
68   for(Int_t iFile=0; iFile<4; iFile++) {
69     int failed = daqDA_DB_getFile(mapFiles[iFile], mapFiles[iFile]);
70     if(failed) { 
71       printf("Cannot retrieve file %s from DAQ DB. Exit.\n",mapFiles[iFile]);
72       return -1;
73     }
74   }
75   
76   /* Open mapping files */
77   AliAltroMapping *mapping[4];
78   TString path = "./";
79   path += "RCU";
80   TString path2;
81   for(Int_t i = 0; i < 4; i++) {
82     path2 = path;
83     path2 += i;
84     path2 += ".data";
85     mapping[i] = new AliCaloAltroMapping(path2.Data());
86   }
87   
88
89   /* define data source : this is argument 1 */  
90   status=monitorSetDataSource( argv[1] );
91   if (status!=0) {
92     printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
93     return -1;
94   }
95
96
97   /* declare monitoring program */
98   status=monitorDeclareMp( __FILE__ );
99   if (status!=0) {
100     printf("monitorDeclareMp() failed : %s\n",monitorDecodeError(status));
101     return -1;
102   }
103
104
105   /* define wait event timeout - 1s max */
106   monitorSetNowait();
107   monitorSetNoWaitNetworkTimeout(1000);
108   
109    /* init some counters */
110   int nevents_physics=0;
111   int nevents_total=0;
112   
113   AliRawReader *rawReader = NULL;
114   AliPHOSRcuDA1* dAs[5];
115   
116   for(Int_t iMod=0; iMod<5; iMod++) {
117     dAs[iMod] = 0;
118   }
119
120   Float_t e[64][56][2];
121   Float_t t[64][56][2];
122   
123   for(Int_t iX=0; iX<64; iX++) {
124     for(Int_t iZ=0; iZ<56; iZ++) {
125       for(Int_t iGain=0; iGain<2; iGain++) {
126         e[iX][iZ][iGain] = 0.;
127         t[iX][iZ][iGain] = 0.;
128       }
129     }
130   }
131   
132   Int_t cellX    = -1;
133   Int_t cellZ    = -1;
134   Int_t nBunches =  0;
135   Int_t nFired   = -1;
136   Int_t sigStart, sigLength;
137   Int_t caloFlag;
138   
139
140   TH1I fFiredCells("fFiredCells","Number of fired cells per event",100,0,1000);
141
142   /* main loop (infinite) */
143   for(;;) {
144     struct eventHeaderStruct *event;
145     eventTypeType eventT;
146   
147     /* check shutdown condition */
148     if (daqDA_checkShutdown()) {break;}
149     
150     /* get next event (blocking call until timeout) */
151     status=monitorGetEventDynamic((void **)&event);
152     if (status==MON_ERR_EOF) {
153       printf ("End of File detected\n");
154       break; /* end of monitoring file has been reached */
155     }
156     
157     if (status!=0) {
158       printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
159       break;
160     }
161
162     /* retry if got no event */
163     if (event==NULL) {
164       continue;
165     }
166
167
168     /* use event - here, just write event id to result file */
169     eventT=event->eventType;
170     
171     if (eventT==PHYSICS_EVENT) {
172       
173       nFired = 0;
174       
175       rawReader = new AliRawReaderDate((void*)event);
176       AliCaloRawStreamV3 stream(rawReader,"PHOS",mapping);
177       AliPHOSRawFitterv0 fitter;
178       fitter.SubtractPedestals(kTRUE); // assume that data is non-ZS
179
180       if(!failZS) {
181         fitter.SubtractPedestals(kFALSE);
182         fitter.SetAmpOffset(offset);
183         fitter.SetAmpThreshold(threshold);
184       }
185       
186       while (stream.NextDDL()) {
187         while (stream.NextChannel()) {
188           
189           cellX    = stream.GetCellX();
190           cellZ    = stream.GetCellZ();
191           caloFlag = stream.GetCaloFlag();  // 0=LG, 1=HG, 2=TRU
192           
193           if(caloFlag!=0 && caloFlag!=1) continue; //TRU data!
194           
195           // In case of oscillating signals with ZS, a channel can have several bunches
196           nBunches = 0;
197           while (stream.NextBunch()) {
198             nBunches++;
199             if (nBunches > 1) continue;
200             sigStart  = stream.GetStartTimeBin();
201             sigLength = stream.GetBunchLength();
202             fitter.SetChannelGeo(stream.GetModule(),cellX,cellZ,caloFlag);
203             fitter.Eval(stream.GetSignals(),sigStart,sigLength);
204           } // End of NextBunch()
205           
206           if(nBunches>1) continue;
207           
208           e[cellX][cellZ][caloFlag] = fitter.GetEnergy();
209           t[cellX][cellZ][caloFlag] = fitter.GetTime();
210           
211           if(caloFlag==1 && fitter.GetEnergy()>40)
212             nFired++;
213         }
214         
215         if(stream.GetModule()<0 || stream.GetModule()>4) continue;
216         
217         if(dAs[stream.GetModule()])
218           dAs[stream.GetModule()]->FillHistograms(e,t);
219         else {
220           dAs[stream.GetModule()] = new AliPHOSRcuDA1(stream.GetModule(),-1,0);
221           dAs[stream.GetModule()]->FillHistograms(e,t);
222         } 
223         
224         for(Int_t iX=0; iX<64; iX++) {
225           for(Int_t iZ=0; iZ<56; iZ++) {
226             for(Int_t iGain=0; iGain<2; iGain++) {
227               e[iX][iZ][iGain] = 0.;
228               t[iX][iZ][iGain] = 0.;
229             }
230           }
231         }
232         
233       }
234       
235       fFiredCells.Fill(nFired);
236       
237       delete rawReader;     
238       nevents_physics++;
239     }
240     
241     nevents_total++;
242     
243     /* free resources */
244     free(event);
245     
246     /* exit when last event received, no need to wait for TERM signal */
247     if (eventT==END_OF_RUN) {
248       printf("EOR event detected\n");
249       break;
250     }
251   }
252   
253   for(Int_t i = 0; i < 4; i++) delete mapping[i];  
254   
255   /* Be sure that all histograms are saved */
256   
257   const TH2F* h2=0;
258   const TH1F* h1=0;
259   char localfile[128];
260
261   Int_t nGood=0;    // >10 entries in peak
262   Int_t nMax=-111;  // max. number of entries in peak
263   Int_t iXmax=-1;
264   Int_t iZmax=-1;
265   Int_t iModMax=-1;
266   
267   for(Int_t iMod=0; iMod<5; iMod++) {
268     if(!dAs[iMod]) continue;
269   
270     printf("DA1 for module %d detected.\n",iMod);
271     sprintf(localfile,"PHOS_Module%d_LED.root",iMod);
272     TFile* f = new TFile(localfile,"recreate");
273     
274     for(Int_t iX=0; iX<64; iX++) {
275       for(Int_t iZ=0; iZ<56; iZ++) {
276         
277         h1 = dAs[iMod]->GetHgLgRatioHistogram(iX,iZ); // High Gain/Low Gain ratio
278         if(h1) {
279           if(h1->GetMaximum()>10.) nGood++;
280           if(h1->GetMaximum()>nMax) {
281             nMax = (Int_t)h1->GetMaximum(); iXmax=iX; iZmax=iZ; iModMax=iMod;
282           }
283           h1->Write(); 
284         }
285       
286         for(Int_t iGain=0; iGain<2; iGain++) {
287           h2 = dAs[iMod]->GetTimeEnergyHistogram(iX,iZ,iGain); // Time vs Energy
288           if(h2) h2->Write();
289         }
290       
291       }
292     }
293     
294     fFiredCells.Write();
295     f->Close();
296     
297     /* Store output files to the File Exchange Server */
298     daqDA_FES_storeFile(localfile,"LED");
299   }
300   
301   printf("%d physics events of %d total processed.\n",nevents_physics,nevents_total);
302   printf("%d histograms has >10 entries in maximum, max. is %d entries ",nGood,nMax);
303   printf("(module,iX,iZ)=(%d,%d,%d)",iModMax,iXmax,iZmax);
304   printf("\n");
305
306   return status;
307 }