]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/DA/TOFphysicsda.cxx
corrected TRD/TOF MV position
[u/mrichter/AliRoot.git] / TOF / DA / TOFphysicsda.cxx
1 /*
2  
3  TOF DA for online calibration
4  
5  Version: 1.0
6  Contact: Roberto.Preghenella@bo.infn.it
7  Reference Run: 115401
8  Run Type: PHYSICS
9  DA Type: MON
10  Number of events needed:
11  Input Files: no input
12  Output Files: TOFdaHits.root TOFdaReadout.root
13  Event types used: PHYSICS_EVENT
14  
15  */
16
17 #define FILE_HITS "TOFdaHits.root"
18 #define FILE_READOUT "TOFdaReadout.root"
19
20 #define EVENT_DETECTOR_TOF 5
21 #define READOUT_INFO_HISTO 1
22
23 // DATE
24 #include "event.h"
25 #include "monitor.h"
26 #include "daqDA.h"
27
28 #include <stdio.h>
29 #include <stdlib.h>
30 #include <errno.h>
31
32 //AliRoot
33 #include "TROOT.h"
34 #include "AliTOFRawStream.h"
35 #include "AliRawReaderDate.h"
36 #include "AliRawReader.h"
37 #include "AliDAQ.h"
38 #include "AliTOFHitData.h"
39 #include "AliTOFHitDataBuffer.h"
40 #include "AliTOFDaConfigHandler.h"
41 #include "AliTOFHitField.h"
42 #include "AliLog.h"
43 #include "AliTOFGeometry.h"
44 #include "AliTOFDecoderV2.h"
45 #include "AliTOFDecoderSummaryData.h"
46 #include "AliTOFDRMSummaryData.h"
47 #include "AliTOFTRMSummaryData.h"
48 #include "AliTOFChainSummaryData.h"
49 #include "AliTOFTDCHitBuffer.h"
50 #include "AliTOFTDCHit.h"
51 #include "AliTOFTDCErrorBuffer.h"
52 #include "AliTOFTDCError.h"
53
54 //ROOT
55 #include "TFile.h"
56 #include "TKey.h"
57 #include "TH2S.h"
58 #include "TObject.h"
59 #include "TMath.h"
60 #include "TSystem.h"
61 #include "TROOT.h"
62 #include "TPluginManager.h"
63 #include "TSAXParser.h"
64 #include "TTree.h"
65
66 /* Main routine
67  Arguments:
68  1- monitoring data source
69  */
70 int
71 main(int argc, char **argv)
72 {
73   
74   printf("in main\n"); // remove me
75   fflush(stdout); // remove me
76     
77     /* magic line from Rene */
78     gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
79                                           "*",
80                                           "TStreamerInfo",
81                                           "RIO",
82                                           "TStreamerInfo()");
83     
84     
85     /* log start of process */
86     printf("TOF DA started\n");
87   fflush(stdout); // remove me
88     
89     /* check that we got some arguments = list of files */
90     if (argc!=2) {
91         printf("Wrong number of arguments\n");
92         return -1;
93     }
94     
95     /*
96      * CONFIG
97      */
98     
99     printf("retrieve config file\n"); // remove me
100   fflush(stdout); // remove me
101
102
103     /* retrieve config file */
104     int getConfigFile = daqDA_DB_getFile("TOFPhysicsConfig.xml","TOFPhysicsConfig.xml");
105     if (getConfigFile != 0){
106         printf("Failed to retrieve config file from DB! returning...\n");
107         return -1;
108     }
109
110     printf("parse config file\n"); // remove me
111   fflush(stdout); // remove me
112
113     /* parse config file */
114     AliTOFDaConfigHandler* tofHandler = new AliTOFDaConfigHandler();
115     TSAXParser *parser = new TSAXParser();
116     parser->ConnectToHandler("AliTOFDaConfigHandler", tofHandler);
117     if (parser->ParseFile("./TOFPhysicsConfig.xml") != 0) {
118         printf("Failed parsing config file! retunring... \n");
119         return -1;
120     }
121     /* setup config params */
122     Int_t meanMultiplicity = tofHandler->GetMeanMultiplicity(); /* average expected TOF multiplicity */
123     Int_t maxHits = tofHandler->GetMaxHits(); /* max number of hits to be collected */
124     printf("current settings:\n");
125     printf(" - meanMultiplicity = %d\n", meanMultiplicity);
126     printf(" - maxHits          = %d\n", maxHits);
127   fflush(stdout); // remove me
128     /* constants */
129     const Int_t nChannels = 157248;
130     Int_t noiseCheckTrigger = 10; /* first noise check after 10 events */
131     Float_t meanChannelRate = (Float_t)meanMultiplicity / (Float_t)nChannels; /* average expected channel rate (hits/event) */
132     Float_t noiseThreshold = 10. * meanChannelRate; /* noise threshold (hits/event) */
133     Int_t minNoiseHits = 10; /* min number of channel hits to check noise */
134     /* counters and flags */
135     Int_t nPhysicsEvents, nCollectedPhysicsEvents, totHits;
136     Int_t nChHits[nChannels];
137     Bool_t inhibitCollection;
138     Bool_t noiseFlag[nChannels];
139     /* variables */
140     Int_t ddl, slot, trm, chain, tdc, channel, index, timebin, totbin, deltaBC, l0l1latency, det[5], dummy;
141     Float_t noiseHitThreshold;
142     
143     /*
144      * INIT
145      */
146     
147     printf("init\n"); // remove me
148   fflush(stdout); // remove me
149
150     /* init counters and flags */
151     nPhysicsEvents = 0;
152     nCollectedPhysicsEvents = 0;
153     totHits = 0;
154     inhibitCollection = kFALSE;
155     for (Int_t ich = 0; ich < nChannels; ich++) {
156         nChHits[ich] = 0;
157         noiseFlag[ich] = kFALSE;
158     }
159     
160     /* TOF raw data handling */
161     AliTOFRawStream *rawStream = new AliTOFRawStream();
162     AliTOFDecoderV2 *decoder = rawStream->GetDecoderV2();
163     AliTOFDecoderSummaryData *decodersd;
164     AliTOFDRMSummaryData *drmsd;
165     //  AliTOFLTMSummaryData *ltmsd;
166     AliTOFTRMSummaryData *trmsd;
167     AliTOFChainSummaryData *chainsd;
168     AliTOFTDCHitBuffer *hitBuffer;
169     AliTOFTDCHit *hit;
170     AliTOFTDCErrorBuffer *errorBuffer;
171     AliTOFTDCError *error;
172     UShort_t errorFlags;
173     UChar_t *data = 0x0;
174     Int_t dataSize;
175     Int_t dataWords;
176     Int_t currentDDL;
177     const AliRawDataHeader *currentCDH;
178     const AliRawDataHeaderV3 *currentCDHV3;
179     Int_t currentMiniEventID;
180     Int_t currentEventID1;
181     Int_t currentL0BCID ;
182     Int_t currentBunchID;
183     Bool_t skipTRM, skipChain;
184     Int_t trmIndex, chainIndex, tdcIndex;
185     Double_t chainEff;
186     
187     printf("open files and create histos\n"); // remove me
188   fflush(stdout); // remove me
189
190     /* open HITS output file */
191     TFile *fileOutHits = new TFile(FILE_HITS, "RECREATE");
192     /* create hit field data structure */
193     AliTOFHitField *hitField = new AliTOFHitField();
194     /* create temporary tree */
195     TTree *tempTree = new TTree("tempTree", "temporary tree");
196     tempTree->Branch("hit", "AliTOFHitField", &hitField);
197     /* create output tree */
198     TTree *outTree = new TTree("hitTree", "hit tree");
199     outTree->Branch("hit", "AliTOFHitField", &hitField);
200     
201     /* open READOUT output file */
202     TFile *fileOutReadout = new TFile(FILE_READOUT, "RECREATE");
203     /* create chain readout efficiency histo */
204     TH1F *hChainEfficiency = new TH1F("hChainEfficiency", "Chain efficiency;chain;efficiency", 1440, 0., 1440.);
205     
206 #if READOUT_INFO_HISTO
207     /* create TRM data histo */
208     TH1F *hTRMData = new TH1F("hTRMData", "TRM data;TRM;frequency", 720, 0., 720.);
209     /* create TRM empty event frequency histo */
210     TH1F *hTRMEmptyEvent = new TH1F("hTRMEmptyEvent", "TRM empty event error;TRM;frequency", 720, 0., 720.);
211     /* create TRM bad event counter frequency histo */
212     TH1F *hTRMBadEventCounter = new TH1F("hTRMBadEventCounter", "TRM bad event counter;TRM;frequency", 720, 0., 720.);
213     /* create TRM bad CRC frequency histo */
214     TH1F *hTRMBadCRC = new TH1F("hTRMBadCRC", "TRM bad CRC;TRM;frequency", 720, 0., 720.);
215     
216     /* create chain data histo */
217     TH1F *hChainData = new TH1F("hChainData", "Chain data;chain;frequency", 1440, 0., 1440.);
218     /* create chain bad status frequency histo */
219     TH1F *hChainBadStatus = new TH1F("hChainBadStatus", "Chain bad status;chain;frequency", 1440, 0., 1440.);
220     /* create chain bad event counter frequency histo */
221     TH1F *hChainBadEventCounter = new TH1F("hChainBadEventCounter", "Chain bad event counter;chain;status;frequency", 1440, 0., 1440.);
222     
223     /* create TDC error frequency histo */
224     TH1F *hTDCError = new TH1F("hTDCError", "TDC error;TDC;frequency", 21600, 0., 21600.);
225     /* create TDC error flags frequency histo */
226     TH2F *hTDCErrorFlags = new TH2F("hTDCErrorFlags", "TDC error flags;TDC;error flag;frequency", 21600, 0., 21600., 15, 0., 15);
227 #endif
228     
229     
230     /*
231      * ONLINE MONITOR
232      */
233     
234     printf("monitor declare table\n"); // remove me
235   fflush(stdout); // remove me
236
237     AliLog::SetGlobalLogLevel(AliLog::kFatal);
238     struct eventHeaderStruct *event;
239     int ret;
240     /* define monitoring table */
241     char *monTable[5] = {
242       (char *)"ALL", (char *)"no",
243       (char *)"PHY", (char *)"yes",
244       (char *)NULL
245     };
246     ret = monitorDeclareTable(monTable);
247     if (ret != 0) {
248         printf("monitorDeclareTable() failed: %s\n", monitorDecodeError(ret));
249         return -1;
250     }
251
252     printf("monitor set data source\n"); // remove me 
253   fflush(stdout); // remove me
254
255     /* define data source : this is argument 1 */
256     ret = monitorSetDataSource(argv[1]);
257     if (ret != 0) {
258         printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(ret));
259         return -1;
260     }
261
262     printf("monitor declare mp\n"); // remove me
263   fflush(stdout); // remove me
264
265     /* declare monitoring program */
266     ret = monitorDeclareMp("TOFdaPhysics");
267     if (ret != 0) {
268         printf("monitorDeclareMp() failed : %s\n",monitorDecodeError(ret));
269         return -1;
270     }
271
272     printf("monitor set no wait\n"); // remove me
273   fflush(stdout); // remove me
274
275     /* define wait event timeout - 1s max */
276     monitorSetNowait();
277
278     printf("monitor set no wait network timeout\n"); // remove me
279   fflush(stdout); // remove me
280
281     monitorSetNoWaitNetworkTimeout(1000);
282     
283     /* variables */
284
285     printf("endless loop\n"); // remove me
286   fflush(stdout); // remove me
287
288     /* loop over events */
289     while (1) {
290   fflush(stdout); // remove me
291         
292         /* check shutdown condition */
293         if (daqDA_checkShutdown()) break;
294         
295         /*
296          * NOISE CHECK
297          */
298         
299         /* check inhibit collection */
300         if (!inhibitCollection) {
301             /* check number of events and check noise */
302             if (nCollectedPhysicsEvents >= noiseCheckTrigger || totHits >= maxHits) {
303                 noiseHitThreshold = noiseThreshold * nCollectedPhysicsEvents;
304                 printf("noise check triggered after %d events: threshold is %f hits\n", nCollectedPhysicsEvents, noiseHitThreshold);
305                 /* loop over all channels */
306                 for (Int_t ich = 0; ich < nChannels; ich++) {
307                     /* check */
308                     if (nChHits[ich] < minNoiseHits || noiseFlag[ich] || nChHits[ich] < noiseHitThreshold) continue;
309                     printf("channel %d tagged as noisy (%d hits): disabled\n", ich, nChHits[ich]);
310                     noiseFlag[ich] = kTRUE;
311                     totHits -= nChHits[ich];
312                 } /* end of loop over all channels */
313                 /* set new noise check trigger value */
314                 noiseCheckTrigger *= 10;
315             } /* end of noise check */
316   fflush(stdout); // remove me
317             
318             /* inhibit hit collection when maximum number of hits exceeded */
319             if (totHits >= maxHits) {
320                 printf("maximum number of hits exceeded (%d): inhibit hit collection\n", maxHits);
321                 inhibitCollection = kTRUE;
322             }
323         }
324         
325         /*
326          * GET EVENT
327          */
328         
329         /* get next event (blocking call until timeout) */
330         ret = monitorGetEventDynamic((void **)&event);
331         if (ret == MON_ERR_EOF) {
332             printf ("End of File detected\n");
333             break; /* end of monitoring file has been reached */
334         }
335         if (ret != 0) {
336             printf("monitorGetEventDynamic() failed (ret=%d errno=%d): %s\n", ret, errno, monitorDecodeError(ret));
337             break;
338         }
339         /* retry if got no event */
340         if (event==NULL) continue;
341         /* check TOF in partecipating detectors */
342         if (!TEST_DETECTOR_IN_PATTERN(event->eventDetectorPattern, EVENT_DETECTOR_TOF)) {
343             free(event);
344             continue;
345         }
346         /* check event type */
347         if (event->eventType != PHYSICS_EVENT) {
348             printf("not a physics event: %d\n", event->eventType);
349             free(event);
350             continue;
351         }
352         /* increment number of physics events */
353         nPhysicsEvents++;
354         if (!inhibitCollection) nCollectedPhysicsEvents++;
355         
356         /*
357          * DECODE EVENT
358          */
359         
360         /* create and setup raw reader */
361         AliRawReader *rawReader = new AliRawReaderDate((void *)event);
362         rawReader->Reset();
363         rawReader->Select("TOF", 0, AliDAQ::NumberOfDdls("TOF") - 1);
364         /* setup raw stream */
365         rawStream->SetRawReader(rawReader);
366         
367         /* loop over DDLs - rawReader->ReadHeader() */
368         while (rawReader->ReadHeader()) {
369             
370             /* read equipment data */
371             dataSize = rawReader->GetDataSize();
372             data = new UChar_t[dataSize];
373             if (!rawReader->ReadNext(data, dataSize)){
374                 delete [] data;
375                 continue;
376             }
377             
378             /* decode data */
379             dataWords = dataSize / 4;
380             decoder->Decode((UInt_t *)data, dataWords);
381             delete [] data;
382             
383             /* read equipment info */
384             currentDDL = rawReader->GetDDLID();
385             currentCDH = rawReader->GetDataHeader();
386             currentCDHV3 = rawReader->GetDataHeaderV3();
387             currentMiniEventID = currentCDH ? currentCDH->GetMiniEventID() : -1;
388             currentMiniEventID = currentCDHV3 ? currentCDHV3->GetMiniEventID(): currentMiniEventID;
389             currentEventID1 = currentCDH ? currentCDH->GetEventID1() : -1;
390             currentEventID1 = currentCDHV3? currentCDHV3->GetEventID1() : currentEventID1;
391             
392             /* read decoder summary data */
393             decodersd = decoder->GetDecoderSummaryData();
394             
395             /* check DRM header/trailer */
396             drmsd = decodersd->GetDRMSummaryData();
397             if (!drmsd->GetHeader() || !drmsd->GetTrailer()) continue;
398             
399             /* get DRM data */
400             currentL0BCID = drmsd->GetL0BCID();
401             
402             /* loop over TRM to get hits */
403             for (Int_t itrm = 0; itrm < 10; itrm++) {
404                 trmsd = drmsd->GetTRMSummaryData(itrm);
405                 trmIndex = itrm + 10 * currentDDL;
406                 skipTRM = kFALSE;
407                 
408                 /* check header/trailer */
409                 if (!trmsd->GetHeader() || !trmsd->GetTrailer()) continue;
410                 
411 #if READOUT_INFO_HISTO
412                 /* fill TRM data */
413                 hTRMData->Fill(trmIndex);
414                 /* fill TRM empty event */
415                 if (trmsd->GetEBit() != 0) {
416                     hTRMEmptyEvent->Fill(trmIndex);
417                     skipTRM = kTRUE;
418                 }
419                 /* fill TRM bad event counter */
420                 if (trmsd->GetEventCounter() != drmsd->GetLocalEventCounter()) {
421                     hTRMBadEventCounter->Fill(trmIndex);
422                     skipTRM = kTRUE;
423                 }
424                 /* fill TRM bad CRC */
425                 if (trmsd->GetEventCRC() != trmsd->GetDecoderCRC()) {
426                     hTRMBadCRC->Fill(trmIndex);
427                     skipTRM = kTRUE;
428                 }
429 #else
430                 /* check bad condition and skip TRM */
431                 if ( trmsd->GetEBit() != 0 ||
432                     trmsd->GetEventCounter() != drmsd->GetLocalEventCounter() ||
433                     trmsd->GetEventCRC() != trmsd->GetDecoderCRC() ) continue;
434 #endif
435                 
436                 /* loop over chains */
437                 for (Int_t ichain = 0; ichain < 2; ichain++) {
438                     chainsd = trmsd->GetChainSummaryData(ichain);
439                     chainIndex = ichain + 2 * itrm + 20 * currentDDL;
440                     skipChain = kFALSE;
441                     
442                     /* check header/trailer */
443                     if (!chainsd->GetHeader() || !chainsd->GetTrailer()) continue;
444                     currentBunchID = chainsd->GetBunchID();
445                     hitBuffer = chainsd->GetTDCPackedHitBuffer();
446                     errorBuffer = chainsd->GetTDCErrorBuffer();
447                     
448 #if READOUT_INFO_HISTO
449                     /* fill chain data */
450                     hChainData->Fill(chainIndex);
451                     /* check chain bad status */
452                     if (chainsd->GetStatus() != 0) {
453                         hChainBadStatus->Fill(chainIndex);
454                         skipChain = kTRUE;
455                     }
456                     /* check chain bad event counter */
457                     if (chainsd->GetEventCounter() != drmsd->GetLocalEventCounter()) {
458                         hChainBadEventCounter->Fill(chainIndex);
459                         skipChain = kTRUE;
460                     }
461                     /* fill TDC error frequency histo */
462                     for (Int_t ierr = 0; ierr < errorBuffer->GetEntries(); ierr++) {
463                         error = errorBuffer->GetError(ierr);
464                         tdc = error->GetTDCID();
465                         tdcIndex = tdc + 15 * ichain + 30 * itrm + 300 * currentDDL;
466                         hTDCError->Fill(tdcIndex);
467                         errorFlags = error->GetErrorFlags();
468                         for (Int_t ierflg = 0; ierflg < 15; ierflg++)
469                             if (errorFlags & (1 << ierflg))
470                                 hTDCErrorFlags->Fill(tdcIndex, ierflg);
471                     }
472 #else
473                     /* check bad condition and skip chain */
474                     if ( chainsd->GetStatus() != 0 ||
475                         chainsd->GetEventCounter() != drmsd->GetLocalEventCounter() ) continue;
476 #endif
477                     
478                     /*
479                      * CHAIN READOUT EFFICIENCY
480                      */
481                     
482                     /* compute number of available channels removing TDCs in error */
483                     chainEff = (120. - 8. * errorBuffer->GetEntries()) / 120.;
484                     /* fill chain readout efficiency histo */
485                     if (!skipTRM && !skipChain)
486                         hChainEfficiency->Fill(chainIndex, chainEff);
487                     
488                     /*
489                      * HIT MANIPULATION
490                      */
491                     
492                     /* check inhibit collection */
493                     if (inhibitCollection) continue;
494                     
495                     /* loop over hits in buffer */
496                     for (Int_t ihit = 0; ihit < hitBuffer->GetEntries(); ihit++) {
497                         
498                         /* get hit */
499                         hit = hitBuffer->GetHit(ihit);
500                         
501                         /* get channel info */
502                         ddl = currentDDL;
503                         slot = trmsd->GetSlotID();
504                         trm = slot - 3;
505                         chain = chainsd->GetChain();
506                         tdc = hit->GetTDCID();
507                         channel = hit->GetChan();
508                         /* get index */
509                         rawStream->EquipmentId2VolumeId(ddl, slot, chain, tdc, channel, det);
510                         dummy = det[4];
511                         det[4] = det[3];
512                         det[3] = dummy;
513                         /* check valid index */
514                         if (det[0] < 0 || det[0] > 17 ||
515                             det[1] < 0 || det[1] > 5 ||
516                             det[2] < 0 || det[2] > 18 ||
517                             det[3] < 0 || det[3] > 1 ||
518                             det[4] < 0 || det[4] > 47) continue;
519                         index = AliTOFGeometry::GetIndex(det);
520                         
521                         /* check noise flag */
522                         if (noiseFlag[index]) continue;
523                         /* increment number of channel hits and total hits */
524                         nChHits[index]++;
525                         totHits++;
526                         /* get signal info */
527                         timebin = hit->GetHitTime();
528                         totbin = hit->GetTOTWidth();
529                         deltaBC = currentBunchID - currentEventID1;
530                         l0l1latency = currentMiniEventID - currentL0BCID;
531                         /* set hit field data */
532                         hitField->SetIndex(index);
533                         hitField->SetTimeBin(timebin);
534                         hitField->SetTOTBin(totbin);
535                         hitField->SetDeltaBC(deltaBC);
536                         hitField->SetL0L1Latency(l0l1latency);
537                         /* fill temp tree */
538                         tempTree->Fill();
539                         
540                     } /* end of loop over hits in buffer */
541                 } /* end of loop over chains */
542             } /* end of loop over TRMs */
543         } /* end of loop over DDLs - rawReader->ReadHeader() */
544         
545         /* delete raw reader */
546         delete rawReader;
547         /* free event */
548         free(event);
549         
550     } /* end of loop over events */
551     
552     /* final noise check */
553     noiseHitThreshold = noiseThreshold * nCollectedPhysicsEvents;
554     printf("final noise check after collectiong %d events: threshold is %f hits\n", nCollectedPhysicsEvents, noiseHitThreshold);
555     /* loop over all channels */
556     for (Int_t ich = 0; ich < nChannels; ich++) {
557         /* check */
558         if (nChHits[ich] < minNoiseHits || noiseFlag[ich] || nChHits[ich] < noiseHitThreshold) continue;
559         printf("channel %d tagged as noisy (%d hits): disabled\n", ich, nChHits[ich]);
560         noiseFlag[ich] = kTRUE;
561         totHits -= nChHits[ich];
562     } /* end of loop over all channels */
563     
564     /* copy hits into output tree from temp tree */
565     printf("copy hits from temporary tree into output tree\n");
566     printf("temporary tree contains %d hits\n", (Int_t)tempTree->GetEntries());
567     for (Int_t ihit = 0; ihit < tempTree->GetEntries(); ihit++) {
568         /* get entry */
569         tempTree->GetEntry(ihit);
570         /* check noise flag */
571         if (noiseFlag[hitField->GetIndex()]) continue;
572         /* fill output tree */
573         outTree->Fill();
574     } /* end of copy hits into output tree from temp tree */
575     printf("output tree contains %d hits\n", (Int_t)outTree->GetEntries());
576     
577     /* write output tree on HITS file */
578     fileOutHits->cd();
579     outTree->Write();
580     fileOutHits->Close();
581     /* export file to FXS */
582     if (daqDA_FES_storeFile(FILE_HITS, "HITS"))
583         return -2;
584     
585 #if READOUT_INFO_HISTO
586     hTRMData->Sumw2();
587     hTRMEmptyEvent->Sumw2();
588     hTRMBadEventCounter->Sumw2();
589     hTRMBadCRC->Sumw2();
590     
591     hChainData->Sumw2();
592     hChainBadStatus->Sumw2();
593     hChainBadEventCounter->Sumw2();
594     
595     hTDCError->Sumw2();
596     hTDCErrorFlags->Sumw2();
597     
598     /* divide histos */
599     if (nPhysicsEvents > 0) {
600         hTRMEmptyEvent->Divide(hTRMData);
601         hTRMBadEventCounter->Divide(hTRMData);
602         hTRMBadCRC->Divide(hTRMData);
603         hTRMData->Scale(1. / nPhysicsEvents);
604         
605         hChainBadStatus->Divide(hChainData);
606         hChainBadEventCounter->Divide(hChainData);
607         hChainData->Scale(1. / nPhysicsEvents);
608         
609         hTDCError->Scale(1. / nPhysicsEvents);
610         hTDCErrorFlags->Scale(1. / nPhysicsEvents);
611     }
612 #endif
613     
614     /* scale chain efficiency by number of physics events */
615     printf("found %d physics events\n", nPhysicsEvents);
616     hChainEfficiency->Sumw2();
617     if (nPhysicsEvents > 0) hChainEfficiency->Scale(1. / (nPhysicsEvents));
618     
619     /* write efficiency histo on READOUT file */
620     fileOutReadout->cd();
621 #if READOUT_INFO_HISTO
622     hTRMData->Write();
623     hTRMEmptyEvent->Write();
624     hTRMBadEventCounter->Write();
625     hTRMBadCRC->Write();
626     
627     hChainData->Write();
628     hChainBadStatus->Write();
629     hChainBadEventCounter->Write();
630     
631     hTDCError->Write();
632     hTDCErrorFlags->Write();
633 #endif
634     hChainEfficiency->Write();
635     fileOutReadout->Close();
636     /* export file to FXS */
637     if (daqDA_FES_storeFile(FILE_READOUT, "READOUT"))
638         return -2;
639     
640     return 0;
641 }