]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/PHOSLEDda.cxx
Stuff from aldaqdqm09
[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   Int_t cellX    = -1;
124   Int_t cellZ    = -1;
125   Int_t nBunches =  0;
126   Int_t nFired   = -1;
127   Int_t sigStart, sigLength;
128   Int_t caloFlag;
129   
130
131   TH1I fFiredCells("fFiredCells","Number of fired cells per event",100,0,1000);
132
133   /* main loop (infinite) */
134   for(;;) {
135     struct eventHeaderStruct *event;
136     eventTypeType eventT;
137   
138     /* check shutdown condition */
139     if (daqDA_checkShutdown()) {break;}
140     
141     /* get next event (blocking call until timeout) */
142     status=monitorGetEventDynamic((void **)&event);
143     if (status==MON_ERR_EOF) {
144       printf ("End of File detected\n");
145       break; /* end of monitoring file has been reached */
146     }
147     
148     if (status!=0) {
149       printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
150       break;
151     }
152
153     /* retry if got no event */
154     if (event==NULL) {
155       continue;
156     }
157
158
159     /* use event - here, just write event id to result file */
160     eventT=event->eventType;
161     
162     if (eventT==PHYSICS_EVENT) {
163       
164       nFired = 0;
165       
166       rawReader = new AliRawReaderDate((void*)event);
167       AliCaloRawStreamV3 stream(rawReader,"PHOS",mapping);
168       AliPHOSRawFitterv0 fitter;
169       fitter.SubtractPedestals(kTRUE); // assume that data is non-ZS
170
171       if(!failZS) {
172         fitter.SubtractPedestals(kFALSE);
173         fitter.SetAmpOffset(offset);
174         fitter.SetAmpThreshold(threshold);
175       }
176       
177       while (stream.NextDDL()) {
178         while (stream.NextChannel()) {
179           
180           cellX    = stream.GetCellX();
181           cellZ    = stream.GetCellZ();
182           caloFlag = stream.GetCaloFlag();  // 0=LG, 1=HG, 2=TRU
183           
184           // In case of oscillating signals with ZS, a channel can have several bunches
185           nBunches = 0;
186           while (stream.NextBunch()) {
187             nBunches++;
188             if (nBunches > 1) continue;
189             sigStart  = stream.GetStartTimeBin();
190             sigLength = stream.GetBunchLength();
191             fitter.SetChannelGeo(stream.GetModule(),cellX,cellZ,caloFlag);
192             fitter.Eval(stream.GetSignals(),sigStart,sigLength);
193           } // End of NextBunch()
194           
195           if(nBunches>1) continue;
196 //        if (nBunches>1 || caloFlag!=0 || caloFlag!=1 || fitter.GetSignalQuality()>1) continue;
197           
198           e[cellX][cellZ][caloFlag] = fitter.GetEnergy();
199           t[cellX][cellZ][caloFlag] = fitter.GetTime();
200           
201           if(caloFlag==1 && fitter.GetEnergy()>40)
202             nFired++;
203         }
204         
205         if(dAs[stream.GetModule()])
206           dAs[stream.GetModule()]->FillHistograms(e,t);
207         else
208           dAs[stream.GetModule()] = new AliPHOSRcuDA1(stream.GetModule(),-1,0);
209         
210         for(Int_t iX=0; iX<64; iX++) {
211           for(Int_t iZ=0; iZ<56; iZ++) {
212             for(Int_t iGain=0; iGain<2; iGain++) {
213               e[iX][iZ][iGain] = 0.;
214               t[iX][iZ][iGain] = 0.;
215             }
216           }
217         }
218         
219       }
220       
221       fFiredCells.Fill(nFired);
222       
223       delete rawReader;     
224       nevents_physics++;
225     }
226     
227     nevents_total++;
228     
229     /* free resources */
230     free(event);
231     
232     /* exit when last event received, no need to wait for TERM signal */
233     if (eventT==END_OF_RUN) {
234       printf("EOR event detected\n");
235       break;
236     }
237   }
238   
239   for(Int_t i = 0; i < 4; i++) delete mapping[i];  
240   
241   /* Be sure that all histograms are saved */
242   
243   const TH2F* h2=0;
244   const TH1F* h1=0;
245   char localfile[128];
246
247   Int_t nGood=0;    // >10 entries in peak
248   Int_t nMax=-111;  // max. number of entries in peak
249   Int_t iXmax=-1;
250   Int_t iZmax=-1;
251   Int_t iModMax=-1;
252   
253   for(Int_t iMod=0; iMod<5; iMod++) {
254     if(!dAs[iMod]) continue;
255   
256     printf("DA1 for module %d detected.\n",iMod);
257     sprintf(localfile,"PHOS_Module%d_LED.root",iMod);
258     TFile* f = new TFile(localfile,"recreate");
259     
260     for(Int_t iX=0; iX<64; iX++) {
261       for(Int_t iZ=0; iZ<56; iZ++) {
262         
263         h1 = dAs[iMod]->GetHgLgRatioHistogram(iX,iZ); // High Gain/Low Gain ratio
264         if(h1) {
265           if(h1->GetMaximum()>10.) nGood++;
266           if(h1->GetMaximum()>nMax) {
267             nMax = (Int_t)h1->GetMaximum(); iXmax=iX; iZmax=iZ; iModMax=iMod;
268           }
269           h1->Write(); 
270         }
271       
272         for(Int_t iGain=0; iGain<2; iGain++) {
273           h2 = dAs[iMod]->GetTimeEnergyHistogram(iX,iZ,iGain); // Time vs Energy
274           if(h2) h2->Write();
275         }
276       
277       }
278     }
279     
280     fFiredCells.Write();
281     f->Close();
282     
283     /* Store output files to the File Exchange Server */
284     daqDA_FES_storeFile(localfile,"LED");
285   }
286   
287   printf("%d physics events of %d total processed.\n",nevents_physics,nevents_total);
288   printf("%d histograms has >10 entries in maximum, max. is %d entries ",nGood,nMax);
289   printf("(module,iX,iZ)=(%d,%d,%d)",iModMax,iXmax,iZmax);
290   printf("\n");
291
292   return status;
293 }