]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PHOS/PHOSPEDda.cxx
Detector Algorithm for pedestal runs.
[u/mrichter/AliRoot.git] / PHOS / PHOSPEDda.cxx
CommitLineData
9d12ffd1 1/*\r
2contact: Boris.Polishchuk@cern.ch\r
3reference run: /alice/data/2010/LHC10a_PHOS/000112189/raw/10000112189003.10.root\r
4run type: PEDESTAL\r
5DA type: MON \r
6number of events needed: 200\r
7input 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\r
8Output files: PHOS_PED.root\r
9*/\r
10\r
11\r
12#include "event.h"\r
13#include "monitor.h"\r
14\r
15extern "C" {\r
16#include "daqDA.h"\r
17}\r
18\r
19#include <stdio.h>\r
20#include <stdlib.h>\r
21\r
22#include <TSystem.h>\r
23#include <TROOT.h>\r
24#include <TPluginManager.h>\r
25\r
26#include <TStyle.h>\r
27#include <TFile.h>\r
28#include <TH1F.h>\r
29#include <TH2F.h>\r
30#include <TString.h>\r
31#include "AliRawReader.h"\r
32#include "AliCaloRawStreamV3.h"\r
33#include "AliLog.h"\r
34\r
35int main(int argc, char **argv) \r
36{\r
37 // Read raw data, decode it to samples,\r
38 // calculate pedestals from presamples, \r
39 // evaluate the signal amplitude as a maximum sample, \r
40 // and fill histograms with pedestals and amplitudes\r
41\r
42 gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",\r
43 "*",\r
44 "TStreamerInfo",\r
45 "RIO",\r
46 "TStreamerInfo()");\r
47 \r
48 AliLog::SetGlobalDebugLevel(0) ;\r
49 AliLog::SetGlobalLogLevel(AliLog::kFatal);\r
50\r
51 int status;\r
52 \r
53 if (argc!=2) {\r
54 printf("Wrong number of arguments\n");\r
55 return -1;\r
56 }\r
57\r
58 /* Retrieve mapping files from DAQ DB */ \r
59 const char* mapFiles[20] = {\r
60 "Mod0RCU0.data",\r
61 "Mod0RCU1.data",\r
62 "Mod0RCU2.data",\r
63 "Mod0RCU3.data",\r
64 "Mod1RCU0.data",\r
65 "Mod1RCU1.data",\r
66 "Mod1RCU2.data",\r
67 "Mod1RCU3.data",\r
68 "Mod2RCU0.data",\r
69 "Mod2RCU1.data",\r
70 "Mod2RCU2.data",\r
71 "Mod2RCU3.data",\r
72 "Mod3RCU0.data",\r
73 "Mod3RCU1.data",\r
74 "Mod3RCU2.data",\r
75 "Mod3RCU3.data",\r
76 "Mod4RCU0.data",\r
77 "Mod4RCU1.data",\r
78 "Mod4RCU2.data",\r
79 "Mod4RCU3.data"\r
80 };\r
81\r
82 for(Int_t iFile=0; iFile<20; iFile++) {\r
83 int failed = daqDA_DB_getFile(mapFiles[iFile], mapFiles[iFile]);\r
84 if(failed) { \r
85 printf("Cannot retrieve file %s from DAQ DB. Exit.\n",mapFiles[iFile]);\r
86 return -1;\r
87 }\r
88 }\r
89 \r
90 /* Open mapping files */\r
91 AliAltroMapping *mapping[20];\r
92 TString path = "./";\r
93\r
94 path += "Mod";\r
95 TString path2;\r
96 TString path3;\r
97 Int_t iMap = 0;\r
98\r
99 for(Int_t iMod = 0; iMod < 5; iMod++) {\r
100 path2 = path;\r
101 path2 += iMod;\r
102 path2 += "RCU";\r
103\r
104 for(Int_t iRCU=0; iRCU<4; iRCU++) {\r
105 path3 = path2;\r
106 path3 += iRCU;\r
107 path3 += ".data";\r
108 mapping[iMap] = new AliCaloAltroMapping(path3.Data());\r
109 iMap++;\r
110 }\r
111 }\r
112 \r
113 /* define data source : this is argument 1 */ \r
114 status=monitorSetDataSource( argv[1] );\r
115 if (status!=0) {\r
116 printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));\r
117 return -1;\r
118 }\r
119 \r
120 /* declare monitoring program */\r
121 status=monitorDeclareMp( __FILE__ );\r
122 if (status!=0) {\r
123 printf("monitorDeclareMp() failed : %s\n",monitorDecodeError(status));\r
124 return -1;\r
125 }\r
126 \r
127 /* define wait event timeout - 1s max */\r
128 monitorSetNowait();\r
129 monitorSetNoWaitNetworkTimeout(1000);\r
130 \r
131 /* init some counters */\r
132 int nevents_physics=0;\r
133 int nevents_total=0;\r
134\r
135 AliRawReader * reader = NULL;\r
136 AliCaloRawStreamV3* stream = NULL;\r
137 \r
138 TString baseNamePed ="hPed";\r
139 TString baseTitlePed="Ped in cell (";\r
140 const char* sgain[3]={"LG","HG", "TRU"};\r
141\r
142 const Int_t caloFlagMax=3,modMax=5,cellXMax=64,cellZMax=56;\r
143 TH1F *hPed[5][3][64][56] = {};\r
144\r
145 TH1F *hPedHiMean1m2 = new TH1F("hPedHiMean1m2","Mean pedestals in module 2, high gain" ,100,0.,100.);\r
146 TH1F *hPedHiRMS1m2 = new TH1F("hPedHiRMS1m2" ,"RMS pedestals in module 2, high gain" ,100,0.,50.);\r
147 TH1F *hPedLoMean1m2 = new TH1F("hPedLoMean1m2","Mean pedestals in module 2, low gain" ,100,0.,100.);\r
148 TH1F *hPedLoRMS1m2 = new TH1F("hPedLoRMS1m2" ,"RMS pedestals in module 2, low gain" ,100,0.,50.);\r
149 TH1F *hPedTRUMean1m2 = new TH1F("hPedTRUMean1m2","Mean pedestals in module 2, TRU" ,1000,0.,1000.);\r
150 TH1F *hPedTRURMS1m2 = new TH1F("hPedTRURMS1m2" ,"RMS pedestals in module 2, TRU" ,100,0.,50.);\r
151\r
152 TH1F *hPedHiMean1m3 = new TH1F("hPedHiMean1m3","Mean pedestals in module 3, high gain" ,100,0.,100.);\r
153 TH1F *hPedHiRMS1m3 = new TH1F("hPedHiRMS1m3" ,"RMS pedestals in module 3, high gain" ,100,0.,50.);\r
154 TH1F *hPedLoMean1m3 = new TH1F("hPedLoMean1m3","Mean pedestals in module 3, low gain" ,100,0.,100.);\r
155 TH1F *hPedLoRMS1m3 = new TH1F("hPedLoRMS1m3" ,"RMS pedestals in module 3, low gain" ,100,0.,50.);\r
156 TH1F *hPedTRUMean1m3 = new TH1F("hPedTRUMean1m3","Mean pedestals in module 3, TRU" ,1000,0.,1000.);\r
157 TH1F *hPedTRURMS1m3 = new TH1F("hPedTRURMS1m3" ,"RMS pedestals in module 3, TRU" ,100,0.,50.);\r
158\r
159 TH1F *hPedHiMean1m4 = new TH1F("hPedHiMean1m4","Mean pedestals in module 4, high gain" ,100,0.,100.);\r
160 TH1F *hPedHiRMS1m4 = new TH1F("hPedHiRMS1m4" ,"RMS pedestals in module 4, high gain" ,100,0.,50.);\r
161 TH1F *hPedLoMean1m4 = new TH1F("hPedLoMean1m4","Mean pedestals in module 4, low gain" ,100,0.,100.);\r
162 TH1F *hPedLoRMS1m4 = new TH1F("hPedLoRMS1m4" ,"RMS pedestals in module 4, low gain" ,100,0.,50.);\r
163 TH1F *hPedTRUMean1m4 = new TH1F("hPedTRUMean1m4","Mean pedestals in module 4, TRU" ,1000,0.,1000.);\r
164 TH1F *hPedTRURMS1m4 = new TH1F("hPedTRURMS1m4" ,"RMS pedestals in module 4, TRU" ,100,0.,50.);\r
165\r
166 hPedHiMean1m2->Sumw2();\r
167 hPedHiRMS1m2 ->Sumw2();\r
168 hPedLoMean1m2->Sumw2();\r
169 hPedLoRMS1m2 ->Sumw2();\r
170 hPedTRUMean1m2->Sumw2();\r
171 hPedTRURMS1m2 ->Sumw2();\r
172 hPedHiMean1m3->Sumw2();\r
173 hPedHiRMS1m3 ->Sumw2();\r
174 hPedLoMean1m3->Sumw2();\r
175 hPedLoRMS1m3 ->Sumw2();\r
176 hPedTRUMean1m3->Sumw2();\r
177 hPedTRURMS1m3 ->Sumw2();\r
178 hPedHiMean1m4->Sumw2();\r
179 hPedHiRMS1m4 ->Sumw2();\r
180 hPedLoMean1m4->Sumw2();\r
181 hPedLoRMS1m4 ->Sumw2();\r
182 hPedTRUMean1m4->Sumw2();\r
183 hPedTRURMS1m4 ->Sumw2();\r
184\r
185 TH2F *hPedHiMeanm2 = new TH2F("hPedHiMeanm2","Mean pedestals in module 2, high gain",\r
186 cellXMax,0.,cellXMax, cellZMax,0.,cellZMax);\r
187 TH2F *hPedHiRMSm2 = new TH2F("hPedHiRMSm2" ,"R.M.S. of pedestals in module 2, high gain",\r
188 cellXMax,0.,cellXMax, cellZMax,0.,cellZMax);\r
189 TH2F *hPedHiNumm2 = new TH2F("hPedHiNumm2" ,"Number of pedestals in module 2, high gain",\r
190 cellXMax,0.,cellXMax, cellZMax,0.,cellZMax);\r
191 TH2F *hPedLoMeanm2 = new TH2F("hPedLoMeanm2","Mean pedestals in module 2, low gain",\r
192 cellXMax,0.,cellXMax, cellZMax,0.,cellZMax);\r
193 TH2F *hPedLoRMSm2 = new TH2F("hPedLoRMSm2" ,"R.M.S. of pedestals in module 2, low gain",\r
194 cellXMax,0.,cellXMax, cellZMax,0.,cellZMax);\r
195 TH2F *hPedLoNumm2 = new TH2F("hPedLoNumm2" ,"Number of pedestals in module 2, low gain",\r
196 cellXMax,0.,cellXMax, cellZMax,0.,cellZMax);\r
197\r
198 TH2F *hPedHiMeanm3 = new TH2F("hPedHiMeanm3","Mean pedestals in module 3, high gain",\r
199 cellXMax,0.,cellXMax, cellZMax,0.,cellZMax);\r
200 TH2F *hPedHiRMSm3 = new TH2F("hPedHiRMSm3" ,"R.M.S. of pedestals in module 3, high gain",\r
201 cellXMax,0.,cellXMax, cellZMax,0.,cellZMax);\r
202 TH2F *hPedHiNumm3 = new TH2F("hPedHiNumm3" ,"Number of pedestals in module 3, high gain",\r
203 cellXMax,0.,cellXMax, cellZMax,0.,cellZMax);\r
204 TH2F *hPedLoMeanm3 = new TH2F("hPedLoMeanm3","Mean pedestals in module 3, low gain",\r
205 cellXMax,0.,cellXMax, cellZMax,0.,cellZMax);\r
206 TH2F *hPedLoRMSm3 = new TH2F("hPedLoRMSm3" ,"R.M.S. of pedestals in module 3, low gain",\r
207 cellXMax,0.,cellXMax, cellZMax,0.,cellZMax);\r
208 TH2F *hPedLoNumm3 = new TH2F("hPedLoNumm3" ,"Number of pedestals in module 3, low gain",\r
209 cellXMax,0.,cellXMax, cellZMax,0.,cellZMax);\r
210\r
211 TH2F *hPedHiMeanm4 = new TH2F("hPedHiMeanm4","Mean pedestals in module 4, high gain",\r
212 cellXMax,0.,cellXMax, cellZMax,0.,cellZMax);\r
213 TH2F *hPedHiRMSm4 = new TH2F("hPedHiRMSm4" ,"R.M.S. of pedestals in module 4, high gain",\r
214 cellXMax,0.,cellXMax, cellZMax,0.,cellZMax);\r
215 TH2F *hPedHiNumm4 = new TH2F("hPedHiNumm4" ,"Number of pedestals in module 4, high gain",\r
216 cellXMax,0.,cellXMax, cellZMax,0.,cellZMax);\r
217 TH2F *hPedLoMeanm4 = new TH2F("hPedLoMeanm4","Mean pedestals in module 4, low gain",\r
218 cellXMax,0.,cellXMax, cellZMax,0.,cellZMax);\r
219 TH2F *hPedLoRMSm4 = new TH2F("hPedLoRMSm4" ,"R.M.S. of pedestals in module 4, low gain",\r
220 cellXMax,0.,cellXMax, cellZMax,0.,cellZMax);\r
221 TH2F *hPedLoNumm4 = new TH2F("hPedLoNumm4" ,"Number of pedestals in module 4, low gain",\r
222 cellXMax,0.,cellXMax, cellZMax,0.,cellZMax);\r
223\r
224 TH1I *hNBunches = new TH1I("hNBunches","Number of bunches",10,0,10);\r
225 TH2I *hHWaddr = new TH2I("hHWaddr","DDL is vs HW addr",216,0,216,4096,0,4096);\r
226 TH1I *hModule = new TH1I("hModule" ,"Module number", 5,0.,5);\r
227\r
228 Int_t runNum=0;\r
229 Int_t module,cellX,cellZ,caloFlag;\r
230\r
231 /* main loop (infinite) */\r
232 for(;;) {\r
233 struct eventHeaderStruct *event;\r
234 eventTypeType eventT;\r
235 \r
236 /* check shutdown condition */\r
237 if (daqDA_checkShutdown()) {break;}\r
238 \r
239 /* get next event (blocking call until timeout) */\r
240 status=monitorGetEventDynamic((void **)&event);\r
241 if (status==MON_ERR_EOF) {\r
242 printf ("End of File detected\n");\r
243 break; /* end of monitoring file has been reached */\r
244 }\r
245 \r
246 if (status!=0) {\r
247 printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));\r
248 break;\r
249 }\r
250 \r
251 /* retry if got no event */\r
252 if (event==NULL) {\r
253 continue;\r
254 }\r
255\r
256 /* use event - here, just write event id to result file */\r
257 eventT=event->eventType;\r
258 \r
259 if (eventT==PHYSICS_EVENT) {\r
260\r
261 reader = new AliRawReaderDate((void*)event);\r
262 stream = new AliCaloRawStreamV3(reader,"PHOS",mapping);\r
263\r
264 while (reader->NextEvent()) {\r
265 runNum = reader->GetRunNumber();\r
266 while (stream->NextDDL()) {\r
267 while (stream->NextChannel()) {\r
268 module = stream->GetModule();\r
269 cellX = stream->GetCellX();\r
270 cellZ = stream->GetCellZ();\r
271 caloFlag = stream->GetCaloFlag();\r
272 if (caloFlag!=0 && caloFlag!=1) continue;\r
273\r
274 hHWaddr->Fill(stream->GetDDLNumber(),stream->GetHWAddress());\r
275 hModule->Fill(module);\r
276 if (!hPed[module][caloFlag][cellX][cellZ]) {\r
277 TString name = baseNamePed;\r
278 TString title = baseTitlePed;\r
279 name +="_g"; name +=caloFlag;\r
280 name +="_m"; name +=module;\r
281 name +="_x"; name +=cellX;\r
282 name +="_z"; name +=cellZ;\r
283\r
284 title +=module; title +=",";\r
285 title +=cellX; title +=",";\r
286 title +=cellZ; title +="), ";\r
287 title +=sgain[caloFlag];\r
288 \r
289 Int_t nx,xmin,xmax;\r
290 if (caloFlag==0 || caloFlag==1) {\r
291 nx=100;\r
292 xmin=0.;\r
293 xmax=100.;\r
294 }\r
295 else {\r
296 nx=1000;\r
297 xmin=0.;\r
298 xmax=1000.;\r
299 }\r
300 hPed[module][caloFlag][cellX][cellZ] = new TH1F(name,title,100,0.,100.);\r
301 hPed[module][caloFlag][cellX][cellZ]->Sumw2();\r
302 hPed[module][caloFlag][cellX][cellZ]->SetMarkerStyle(20);\r
303 hPed[module][caloFlag][cellX][cellZ]->SetOption("eph");\r
304 }\r
305\r
306 Int_t nBunches = 0;\r
307 while (stream->NextBunch()) {\r
308 nBunches++;\r
309 const UShort_t *sig = stream->GetSignals();\r
310 Int_t sigLength = stream->GetBunchLength();\r
311 for (Int_t i = 0; i < sigLength; i++) {\r
312 hPed[module][caloFlag][cellX][cellZ]->Fill(sig[i]);\r
313 }\r
314 }\r
315 hNBunches->Fill(nBunches);\r
316 } // end of NextChannel()\r
317\r
318 } // end of NextDDL()\r
319 } // end of nextEvent()\r
320 \r
321 // Fill 2-dim histograms for mean, rms and n pedestals\r
322 \r
323 for (Int_t mod=2; mod<=4; mod++) {\r
324 for (Int_t caloFlag=0; caloFlag<2; caloFlag++) {\r
325 for (Int_t cellX=0; cellX<cellXMax; cellX++) {\r
326 for (Int_t cellZ=0; cellZ<cellZMax; cellZ++) {\r
327 if (hPed[mod][caloFlag][cellX][cellZ] != 0) {\r
328 if (caloFlag == 0) {\r
329 if (mod==2) {\r
330 hPedLoMean1m2->Fill( hPed[mod][caloFlag][cellX][cellZ]->GetMean());\r
331 hPedLoRMS1m2 ->Fill( hPed[mod][caloFlag][cellX][cellZ]->GetRMS());\r
332 hPedLoMeanm2 ->Fill( cellX, cellZ, hPed[mod][caloFlag][cellX][cellZ]->GetMean());\r
333 hPedLoRMSm2 ->Fill( cellX, cellZ, hPed[mod][caloFlag][cellX][cellZ]->GetRMS());\r
334 hPedLoNumm2 ->Fill( cellX, cellZ, hPed[mod][caloFlag][cellX][cellZ]->GetEntries());\r
335 }\r
336 else if (mod==3) {\r
337 hPedLoMean1m3->Fill( hPed[mod][caloFlag][cellX][cellZ]->GetMean());\r
338 hPedLoRMS1m3 ->Fill( hPed[mod][caloFlag][cellX][cellZ]->GetRMS());\r
339 hPedLoMeanm3 ->Fill( cellX, cellZ, hPed[mod][caloFlag][cellX][cellZ]->GetMean());\r
340 hPedLoRMSm3 ->Fill( cellX, cellZ, hPed[mod][caloFlag][cellX][cellZ]->GetRMS());\r
341 hPedLoNumm3 ->Fill( cellX, cellZ, hPed[mod][caloFlag][cellX][cellZ]->GetEntries());\r
342 }\r
343 else if (mod==4) {\r
344 hPedLoMean1m4->Fill( hPed[mod][caloFlag][cellX][cellZ]->GetMean());\r
345 hPedLoRMS1m4 ->Fill( hPed[mod][caloFlag][cellX][cellZ]->GetRMS());\r
346 hPedLoMeanm4 ->Fill( cellX, cellZ, hPed[mod][caloFlag][cellX][cellZ]->GetMean());\r
347 hPedLoRMSm4 ->Fill( cellX, cellZ, hPed[mod][caloFlag][cellX][cellZ]->GetRMS());\r
348 hPedLoNumm4 ->Fill( cellX, cellZ, hPed[mod][caloFlag][cellX][cellZ]->GetEntries());\r
349 }\r
350 }\r
351 else if (caloFlag == 1) {\r
352 if (mod==2) {\r
353 hPedHiMean1m2->Fill( hPed[mod][caloFlag][cellX][cellZ]->GetMean());\r
354 hPedHiRMS1m2 ->Fill( hPed[mod][caloFlag][cellX][cellZ]->GetRMS() );\r
355 hPedHiMeanm2 ->Fill( cellX, cellZ, hPed[mod][caloFlag][cellX][cellZ]->GetMean());\r
356 hPedHiRMSm2 ->Fill( cellX, cellZ, hPed[mod][caloFlag][cellX][cellZ]->GetRMS());\r
357 hPedHiNumm2 ->Fill( cellX, cellZ, hPed[mod][caloFlag][cellX][cellZ]->GetEntries());\r
358 }\r
359 if (mod==3) {\r
360 hPedHiMean1m3->Fill( hPed[mod][caloFlag][cellX][cellZ]->GetMean());\r
361 hPedHiRMS1m3 ->Fill( hPed[mod][caloFlag][cellX][cellZ]->GetRMS() );\r
362 hPedHiMeanm3 ->Fill( cellX, cellZ, hPed[mod][caloFlag][cellX][cellZ]->GetMean());\r
363 hPedHiRMSm3 ->Fill( cellX, cellZ, hPed[mod][caloFlag][cellX][cellZ]->GetRMS());\r
364 hPedHiNumm3 ->Fill( cellX, cellZ, hPed[mod][caloFlag][cellX][cellZ]->GetEntries());\r
365 }\r
366 if (mod==4) {\r
367 hPedHiMean1m4->Fill( hPed[mod][caloFlag][cellX][cellZ]->GetMean());\r
368 hPedHiRMS1m4 ->Fill( hPed[mod][caloFlag][cellX][cellZ]->GetRMS());\r
369 hPedHiMeanm4 ->Fill( cellX, cellZ, hPed[mod][caloFlag][cellX][cellZ]->GetMean());\r
370 hPedHiRMSm4 ->Fill( cellX, cellZ, hPed[mod][caloFlag][cellX][cellZ]->GetRMS());\r
371 hPedHiNumm4 ->Fill( cellX, cellZ, hPed[mod][caloFlag][cellX][cellZ]->GetEntries());\r
372 }\r
373 }\r
374 else if (caloFlag == 2) {\r
375 if (mod==2) {\r
376 hPedTRUMean1m2->Fill( hPed[mod][caloFlag][cellX][cellZ]->GetMean());\r
377 hPedTRURMS1m2 ->Fill( hPed[mod][caloFlag][cellX][cellZ]->GetRMS() );\r
378 }\r
379 if (mod==3) {\r
380 hPedTRUMean1m3->Fill( hPed[mod][caloFlag][cellX][cellZ]->GetMean());\r
381 hPedTRURMS1m3 ->Fill( hPed[mod][caloFlag][cellX][cellZ]->GetRMS() );\r
382 }\r
383 if (mod==4) {\r
384 hPedTRUMean1m4->Fill( hPed[mod][caloFlag][cellX][cellZ]->GetMean());\r
385 hPedTRURMS1m4 ->Fill( hPed[mod][caloFlag][cellX][cellZ]->GetRMS() );\r
386 }\r
387 }\r
388 }\r
389 }\r
390 }\r
391 }\r
392 }\r
393 \r
394 delete rawReader; \r
395 delete stream;\r
396 nevents_physics++;\r
397 } // end of if (eventT==PHYSICS_EVENT)\r
398 \r
399 nevents_total++;\r
400 \r
401 /* free resources */\r
402 free(event);\r
403 \r
404 /* exit when last event received, no need to wait for TERM signal */\r
405 if (eventT==END_OF_RUN) {\r
406 printf("EOR event detected\n");\r
407 break;\r
408 }\r
409\r
410 } // end of inf. loop over events\r
411 \r
412 for(Int_t i = 0; i < 20; i++) delete mapping[i];\r
413 \r
414 // Write existing histograms to a root file\r
415 \r
416 TString fileName = "PHOS_PED.root";\r
417 TFile *file = new TFile(fileName,"RECREATE");\r
418 \r
419 for (Int_t mod=2; mod<=3; mod++) {\r
420 for (Int_t caloFlag=0; caloFlag<caloFlagMax; caloFlag++) {\r
421 for (Int_t mod=0; mod<modMax; mod++) {\r
422 for (Int_t cellX=0; cellX<cellXMax; cellX++) {\r
423 for (Int_t cellZ=0; cellZ<cellZMax; cellZ++) {\r
424 if (hPed[mod][caloFlag][cellX][cellZ] != 0)\r
425 hPed[mod][caloFlag][cellX][cellZ]->Write();\r
426 }\r
427 }\r
428 }\r
429 }\r
430 }\r
431 \r
432 hPedHiMean1m2->Write();\r
433 hPedHiRMS1m2 ->Write();\r
434 hPedLoMean1m2->Write();\r
435 hPedLoRMS1m2 ->Write();\r
436 hPedHiMeanm2 ->Write();\r
437 hPedHiRMSm2 ->Write();\r
438 hPedHiNumm2 ->Write();\r
439 hPedLoMeanm2 ->Write();\r
440 hPedLoRMSm2 ->Write();\r
441 hPedLoNumm2 ->Write();\r
442 hPedTRUMean1m2->Write();\r
443 hPedTRURMS1m2 ->Write();\r
444 \r
445 hPedHiMean1m3->Write();\r
446 hPedHiRMS1m3 ->Write();\r
447 hPedLoMean1m3->Write();\r
448 hPedLoRMS1m3 ->Write();\r
449 hPedHiMeanm3 ->Write();\r
450 hPedHiRMSm3 ->Write();\r
451 hPedHiNumm3 ->Write();\r
452 hPedLoMeanm3 ->Write();\r
453 hPedLoRMSm3 ->Write();\r
454 hPedLoNumm3 ->Write();\r
455 hPedTRUMean1m3->Write();\r
456 hPedTRURMS1m3 ->Write();\r
457\r
458 hPedHiMean1m4->Write();\r
459 hPedHiRMS1m4 ->Write();\r
460 hPedLoMean1m4->Write();\r
461 hPedLoRMS1m4 ->Write();\r
462 hPedHiMeanm4 ->Write();\r
463 hPedHiRMSm4 ->Write();\r
464 hPedHiNumm4 ->Write();\r
465 hPedLoMeanm4 ->Write();\r
466 hPedLoRMSm4 ->Write();\r
467 hPedLoNumm4 ->Write();\r
468 hPedTRUMean1m4->Write();\r
469 hPedTRURMS1m4 ->Write();\r
470\r
471 hNBunches ->Write();\r
472 hHWaddr ->Write();\r
473 hModule ->Write();\r
474 \r
475 file->Close();\r
476 \r
477 /* Store output files to the File Exchange Server */\r
478 daqDA_FES_storeFile(fileName.Data(),"PED");\r
479 \r
480 printf("%d physics events of %d total processed.\n",nevents_physics,nevents_total);\r
481 return status;\r
482}\r