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