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