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