]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/PHOSLEDda.cxx
bugfix: external interface was calling AliHLTComponent::Init twice since r27483
[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 "AliPHOSRawDecoder.h"
31 #include "AliCaloAltroMapping.h"
32
33
34 /* Main routine
35       Arguments: 
36       1- monitoring data source
37 */
38 int main(int argc, char **argv) {
39
40   gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
41                                         "*",
42                                         "TStreamerInfo",
43                                         "RIO",
44                                         "TStreamerInfo()");
45
46   int status;
47   
48   if (argc!=2) {
49     printf("Wrong number of arguments\n");
50     return -1;
51   }
52
53   /* Retrieve mapping files from DAQ DB */ 
54   const char* mapFiles[4] = {"RCU0.data","RCU1.data","RCU2.data","RCU3.data"};
55
56   for(Int_t iFile=0; iFile<4; iFile++) {
57     int failed = daqDA_DB_getFile(mapFiles[iFile], mapFiles[iFile]);
58     if(failed) { 
59       printf("Cannot retrieve file %s from DAQ DB. Exit.\n",mapFiles[iFile]);
60       return -1;
61     }
62   }
63   
64   /* Open mapping files */
65   AliAltroMapping *mapping[4];
66   TString path = "./";
67   path += "RCU";
68   TString path2;
69   for(Int_t i = 0; i < 4; i++) {
70     path2 = path;
71     path2 += i;
72     path2 += ".data";
73     mapping[i] = new AliCaloAltroMapping(path2.Data());
74   }
75   
76
77   /* define data source : this is argument 1 */  
78   status=monitorSetDataSource( argv[1] );
79   if (status!=0) {
80     printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
81     return -1;
82   }
83
84
85   /* declare monitoring program */
86   status=monitorDeclareMp( __FILE__ );
87   if (status!=0) {
88     printf("monitorDeclareMp() failed : %s\n",monitorDecodeError(status));
89     return -1;
90   }
91
92
93   /* define wait event timeout - 1s max */
94   monitorSetNowait();
95   monitorSetNoWaitNetworkTimeout(1000);
96   
97    /* init some counters */
98   int nevents_physics=0;
99   int nevents_total=0;
100
101   Int_t fMod = 2; // module 2!
102   AliRawReader *rawReader = NULL;
103
104   AliPHOSRcuDA1 da1(fMod,-1,0); // DA1 for module2, no input/output file
105   
106   Float_t e[64][56][2];
107   Float_t t[64][56][2];
108
109   Int_t gain = -1;
110   Int_t X = -1;
111   Int_t Z = -1;
112
113   /* main loop (infinite) */
114   for(;;) {
115     struct eventHeaderStruct *event;
116     eventTypeType eventT;
117   
118     /* check shutdown condition */
119     if (daqDA_checkShutdown()) {break;}
120     
121     /* get next event (blocking call until timeout) */
122     status=monitorGetEventDynamic((void **)&event);
123     if (status==MON_ERR_EOF) {
124       printf ("End of File detected\n");
125       break; /* end of monitoring file has been reached */
126     }
127     
128     if (status!=0) {
129       printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
130       break;
131     }
132
133     /* retry if got no event */
134     if (event==NULL) {
135       continue;
136     }
137
138
139     /* use event - here, just write event id to result file */
140     eventT=event->eventType;
141     
142     if (eventT==PHYSICS_EVENT) {
143       
144       for(Int_t iX=0; iX<64; iX++) {
145         for(Int_t iZ=0; iZ<56; iZ++) {
146           for(Int_t iGain=0; iGain<2; iGain++) {
147             e[iX][iZ][iGain] = 0.;
148             t[iX][iZ][iGain] = 0.;
149           }
150         }
151       }
152
153       rawReader = new AliRawReaderDate((void*)event);
154 //       AliPHOSRawDecoderv1 dc(rawReader,mapping);
155       AliPHOSRawDecoder dc(rawReader,mapping);
156       dc.SubtractPedestals(kTRUE);
157       
158       while(dc.NextDigit()) {
159
160         X = dc.GetRow() - 1;
161         Z = dc.GetColumn() - 1;
162
163         if(dc.IsLowGain()) gain = 0;
164         else
165           gain = 1;
166         
167         e[X][Z][gain] = dc.GetEnergy();
168         t[X][Z][gain] = dc.GetTime();
169         
170       }
171       
172       da1.FillHistograms(e,t);
173       
174       delete rawReader;     
175       nevents_physics++;
176     }
177     
178     nevents_total++;
179     
180     /* free resources */
181     free(event);
182     
183     /* exit when last event received, no need to wait for TERM signal */
184     if (eventT==END_OF_RUN) {
185       printf("EOR event detected\n");
186       break;
187     }
188   }
189   
190   for(Int_t i = 0; i < 4; i++) delete mapping[i];  
191   
192   /* Be sure that all histograms are saved */
193
194   char localfile[128];
195   sprintf(localfile,"PHOS_Module%d_LED.root",fMod);
196   TFile* f = new TFile(localfile,"recreate");
197   
198   const TH2F* h2=0;
199   const TH1F* h1=0;
200
201   Int_t nGood=0;    // >10 entries in peak
202   Int_t nMax=-111;  // max. number of entries in peak
203   Int_t iXmax=-1;
204   Int_t iZmax=-1;
205
206   for(Int_t iX=0; iX<64; iX++) {
207     for(Int_t iZ=0; iZ<56; iZ++) {
208       
209       h1 = da1.GetHgLgRatioHistogram(iX,iZ); // High Gain/Low Gain ratio
210       if(h1) {
211         if(h1->GetMaximum()>10.) nGood++;
212         if(h1->GetMaximum()>nMax) {nMax = (Int_t)h1->GetMaximum(); iXmax=iX; iZmax=iZ;}
213         h1->Write(); 
214       }
215       
216       for(Int_t iGain=0; iGain<2; iGain++) {
217         h2 = da1.GetTimeEnergyHistogram(iX,iZ,iGain); // Time vs Energy
218         if(h2) h2->Write();
219       }
220       
221     }
222   }
223   
224   f->Close();
225   
226   /* Store output files to the File Exchange Server */
227   char fesfileid[128];
228   
229   for(Int_t iMod=0; iMod<5; iMod++) {
230     sprintf(localfile,"PHOS_Module%d_LED.root",iMod);
231     sprintf(fesfileid,"PHOS_Module%d_LED",iMod);
232     daqDA_FES_storeFile(localfile,fesfileid);
233   }
234   
235   printf("%d physics events of %d total processed.\n",nevents_physics,nevents_total);
236   printf("%d histograms has >10 entries in maximum, max. is %d entries ",nGood,nMax);
237   printf("(module,iX,iZ)=(%d,%d,%d)",fMod,iXmax,iZmax);
238   printf("\n");
239
240   return status;
241 }