]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RAW/AliRawDB.cxx
-updates (ShinIchi)
[u/mrichter/AliRoot.git] / RAW / AliRawDB.cxx
1 // @(#)alimdc:$Name:  $:$Id$
2 // Author: Fons Rademakers  26/11/99
3
4 /**************************************************************************
5  * Copyright(c) 1998-2003, ALICE Experiment at CERN, All rights reserved. *
6  *                                                                        *
7  * Author: The ALICE Off-line Project.                                    *
8  * Contributors are mentioned in the code where appropriate.              *
9  *                                                                        *
10  * Permission to use, copy, modify and distribute this software and its   *
11  * documentation strictly for non-commercial purposes is hereby granted   *
12  * without fee, provided that the above copyright notice appears in all   *
13  * copies and that both the copyright notice and this permission notice   *
14  * appear in the supporting documentation. The authors make no claims     *
15  * about the suitability of this software for any purpose. It is          *
16  * provided "as is" without express or implied warranty.                  *
17  **************************************************************************/
18
19 //////////////////////////////////////////////////////////////////////////
20 //                                                                      //
21 // AliRawDB                                                             //
22 //                                                                      //
23 //////////////////////////////////////////////////////////////////////////
24
25 #include <errno.h>
26 #include <Riostream.h>
27 #include <RVersion.h>
28
29 #if ROOT_VERSION_CODE >= ROOT_VERSION(5,15,0)
30 #include <TBufferFile.h>
31 #else
32 #include <TBuffer.h>
33 #endif
34
35 #include <TSystem.h>
36 #include <TKey.h>
37
38 #include <TObjString.h>
39
40 #include <TBranch.h>
41
42 #include "AliESDEvent.h"
43 #include "AliRawEventV2.h"
44 #include "AliRawDataArrayV2.h"
45 #include "AliRawEventHeaderBase.h"
46 #include "AliRawEquipmentHeader.h"
47
48 #include "AliRawDB.h"
49
50 using std::ofstream;
51
52 ClassImp(AliRawDB)
53
54 const char *AliRawDB::fgkAliRootTag = "$Rev$";
55
56 // Split TPC into 18 branches in order to avoid problems with big memory
57 // consumption in case of TPC events w/o zero-suppression
58 Int_t AliRawDB::fgkDetBranches[AliDAQ::kNDetectors+1] = {1,1,1,18,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,10,1};
59
60 //______________________________________________________________________________
61 AliRawDB::AliRawDB(AliRawEventV2 *event,
62                    AliESDEvent *esd, 
63                    Int_t compress,
64                    const char* fileName,
65                    Int_t basketsize) :
66   fRawDB(NULL),
67   fTree(NULL),
68   fEvent(event),
69   fESDTree(NULL),
70   fESD(esd),
71   fCompress(compress),
72   fBasketSize(basketsize),
73   fMaxSize(-1),
74   fFS1(""),
75   fFS2(""),
76   fDeleteFiles(kFALSE),
77   fStop(kFALSE)
78 {
79    // Create a new raw DB
80
81   for (Int_t iDet = 0; iDet < AliDAQ::kNDetectors; iDet++) {
82     fDetRawData[iDet] = new AliRawDataArrayV2*[fgkDetBranches[iDet]];
83     Int_t nDDLsPerBranch = AliDAQ::NumberOfDdls(iDet)/fgkDetBranches[iDet];
84     for (Int_t iBranch = 0; iBranch < fgkDetBranches[iDet]; iBranch++)
85       fDetRawData[iDet][iBranch] = new AliRawDataArrayV2(nDDLsPerBranch);
86   }
87
88   fDetRawData[AliDAQ::kNDetectors] = new AliRawDataArrayV2*[fgkDetBranches[AliDAQ::kNDetectors]];
89   for (Int_t iBranch = 0; iBranch < fgkDetBranches[AliDAQ::kNDetectors]; iBranch++)
90     fDetRawData[AliDAQ::kNDetectors][iBranch] = new AliRawDataArrayV2(100);
91
92    if (fileName) {
93       if (!Create(fileName))
94          MakeZombie();
95    }
96 }
97
98
99 //______________________________________________________________________________
100 AliRawDB::~AliRawDB() {
101   // Destructor
102
103   if(Close()==-1) Error("~AliRawDB", "cannot close output file!");
104
105   for (Int_t iDet = 0; iDet < (AliDAQ::kNDetectors + 1); iDet++) {
106     for (Int_t iBranch = 0; iBranch < fgkDetBranches[iDet]; iBranch++)
107       delete fDetRawData[iDet][iBranch];
108     delete [] fDetRawData[iDet];
109   }
110 }
111
112 //______________________________________________________________________________
113 Bool_t AliRawDB::FSHasSpace(const char *fs) const
114 {
115    // Check for at least fMaxSize bytes of free space on the file system.
116    // If the space is not available return kFALSE, kTRUE otherwise.
117
118    Long_t id, bsize, blocks, bfree;
119
120    if (gSystem->GetFsInfo(fs, &id, &bsize, &blocks, &bfree) == 1) {
121       Error("FSHasSpace", "could not stat file system %s", fs);
122       return kFALSE;
123    }
124
125    // Leave 5 percent of diskspace free
126    Double_t avail = Double_t(bfree) * 0.95;
127    if (avail*bsize > fMaxSize)
128       return kTRUE;
129
130    Warning("FSHasSpace", "no space on file system %s", fs);
131    return kFALSE;
132 }
133
134 //______________________________________________________________________________
135 const char *AliRawDB::GetFileName() const
136 {
137    // Return filename based on hostname and date and time. This will make
138    // each file unique. Also makes sure (via FSHasSpace()) that there is
139    // enough space on the file system to store the file. Returns 0 in
140    // case of error or interrupt signal.
141
142    static TString fname;
143    static Bool_t  fstoggle = kFALSE;
144
145    TString fs = fstoggle ? fFS2 : fFS1;
146    TDatime dt;
147
148    TString hostname = gSystem->HostName();
149    Int_t pos;
150    if ((pos = hostname.Index(".")) != kNPOS)
151       hostname.Remove(pos);
152
153    if (!FSHasSpace(fs)) {
154       while (1) {
155          fstoggle = !fstoggle;
156          fs = fstoggle ? fFS2 : fFS1;
157          if (FSHasSpace(fs)) break;
158          Info("GetFileName", "sleeping 30 seconds before retrying...");
159          gSystem->Sleep(30000);   // sleep for 30 seconds
160          if (fStop) return 0;
161       }
162    }
163
164    fname = fs + "/" + hostname + "_";
165    fname += dt.GetDate();
166    fname += "_";
167    fname += dt.GetTime();
168    fname += ".root";
169
170    fstoggle = !fstoggle;
171
172    return fname;
173 }
174
175 //______________________________________________________________________________
176 void AliRawDB::SetFS(const char* fs1, const char* fs2)
177 {
178 // set the file system location
179
180   fFS1 = fs1;
181   if (fs1 && !fFS1.Contains(":")) {
182     gSystem->ResetErrno();
183     gSystem->MakeDirectory(fs1);
184     if (gSystem->GetErrno() && gSystem->GetErrno() != EEXIST) {
185       SysError("SetFS", "mkdir %s", fs1);
186     }
187   }
188
189   fFS2 = fs2;
190   if (fs2) {
191     gSystem->ResetErrno();
192     gSystem->MakeDirectory(fs2);
193     if (gSystem->GetErrno() && gSystem->GetErrno() != EEXIST) {
194       SysError("SetFS", "mkdir %s", fs2);
195     }
196   }
197 }
198
199 //______________________________________________________________________________
200 Bool_t AliRawDB::Create(const char* fileName)
201 {
202    // Create a new raw DB.
203
204    const Int_t kMaxRetry = 1;
205    const Int_t kMaxSleep = 1;      // seconds
206    const Int_t kMaxSleepLong = 10; // seconds
207    Int_t retry = 0;
208
209 again:
210    if (fStop) return kFALSE;
211
212    const char *fname = fileName;
213    if (!fname) fname = GetFileName();
214    if (!fname) {
215       Error("Create", "error getting raw DB file name");
216       return kFALSE;
217    }
218
219    retry++;
220
221    fRawDB = TFile::Open(fname, GetOpenOption(),
222                         Form("ALICE raw-data file (%s)", GetAliRootTag()), fCompress,
223                         GetNetopt());
224    if (!fRawDB) {
225       if (retry < kMaxRetry) {
226          Warning("Create", "failure to open file, sleeping %d %s before retrying...",
227                  kMaxSleep, kMaxSleep==1 ? "second" : "seconds");
228          gSystem->Sleep(kMaxSleep*1000);
229          goto again;
230       }
231       Error("Create", "failure to open file %s after %d tries", fname, kMaxRetry);
232       return kFALSE;
233    }
234    if (retry > 1)
235       Warning("Create", "succeeded to open file after %d retries", retry);
236
237    if (fRawDB->IsZombie()) {
238       if (fRawDB->GetErrno() == ENOSPC ||
239           fRawDB->GetErrno() == 1018   ||   // SECOMERR
240           fRawDB->GetErrno() == 1027) {     // SESYSERR
241          fRawDB->ResetErrno();
242          delete fRawDB;
243          Warning("Create", "file is a zombie (no space), sleeping %d %s before retrying...",
244                  kMaxSleepLong, kMaxSleepLong==1 ? "second" : "seconds");
245          gSystem->Sleep(kMaxSleepLong*1000);   // sleep 10 seconds before retrying
246          goto again;
247       }
248       Error("Create", "file %s is zombie", fname);
249       fRawDB->ResetErrno();
250       delete fRawDB;
251       fRawDB = 0;
252       if (retry < kMaxRetry) {
253          Warning("Create", "file is a zombie, sleeping %d %s before retrying...",
254                  kMaxSleep, kMaxSleep==1 ? "second" : "seconds");
255          gSystem->Sleep(kMaxSleep*1000);
256          goto again;
257       }
258       Error("Create", "failure to open file %s after %d tries", fname, kMaxRetry);
259       return kFALSE;
260    }
261
262    return kTRUE;
263 }
264
265 static void BranchResetBit(TBranch *b) 
266 {
267   // Reset MapObject on this branch and all the sub-branches
268
269   b->ResetBit( kBranchObject | kBranchAny ); // Or in newer ROOT: b->ResetBit( kMapObject )
270   TIter next( b->GetListOfBranches() );
271   TBranch *sub = 0;
272   while ( (sub = (TBranch*)next() ) ) {
273     BranchResetBit( sub );
274   }
275 }
276
277 //______________________________________________________________________________
278 void AliRawDB::MakeTree()
279 {
280    // Create ROOT Tree object container.
281
282    fTree = new TTree("RAW", Form("ALICE raw-data tree (%s)", GetAliRootTag()));
283    fTree->SetAutoSave(21000000000LL);  // autosave when 21 Gbyte written
284
285    fTree->BranchRef();
286
287    Int_t split   = 99;
288    TBranch *b = fTree->Branch("rawevent", "AliRawEventV2", &fEvent, fBasketSize, split);
289    BranchResetBit(b);
290
291    // Make brach for each sub-detector
292    for (Int_t iDet = 0; iDet < AliDAQ::kNDetectors; iDet++) {
293      for (Int_t iBranch = 0; iBranch < fgkDetBranches[iDet]; iBranch++) {
294        b = fTree->Branch(Form("%s%d",AliDAQ::DetectorName(iDet),iBranch),"AliRawDataArrayV2",
295                          &fDetRawData[iDet][iBranch],fBasketSize,split);
296        BranchResetBit(b);
297      }
298    }
299    // Make special branch for unrecognized raw-data payloads
300    for (Int_t iBranch = 0; iBranch < fgkDetBranches[AliDAQ::kNDetectors]; iBranch++) {
301      b = fTree->Branch(Form("Common%d",iBranch),"AliRawDataArrayV2",
302                        &fDetRawData[AliDAQ::kNDetectors][iBranch],fBasketSize,split);
303      BranchResetBit(b);
304    }
305
306    // Create tree which will contain the HLT ESD information
307
308    if (fESD) {
309      fESDTree = new TTree("esdTree", Form("ALICE HLT ESD tree (%s)", GetAliRootTag()));
310      fESDTree->SetAutoSave(21000000000LL);  // autosave when 21 Gbyte written
311      split   = 0;
312      fESDTree->Branch("ESD", "AliESDEvent", &fESD, fBasketSize, split);
313    }
314
315 }
316
317 //______________________________________________________________________________
318 Long64_t AliRawDB::Close()
319 {
320    // Close raw DB.
321    if (!fRawDB) return 0;
322
323    if (!fRawDB->IsOpen()) return 0;
324
325    fRawDB->cd();
326
327    // Write the tree.
328    Bool_t error = kFALSE;
329    if (fTree)
330      if (fTree->Write() == 0)
331        error = kTRUE;
332    if (fESDTree)
333      if (fESDTree->Write() == 0)
334        error = kTRUE;
335
336    // Close DB, this also deletes the fTree
337    fRawDB->Close();
338
339    fTree = NULL;
340
341    Long64_t filesize = fRawDB->GetEND();
342
343    if (fDeleteFiles) {
344       gSystem->Unlink(fRawDB->GetName());
345       delete fRawDB;
346       fRawDB = 0;
347       if(!error)
348         return filesize;
349       else
350         return -1;
351    }
352
353    delete fRawDB;
354    fRawDB = 0;
355    if(!error)
356      return filesize;
357    else
358      return -1;
359 }
360
361 //______________________________________________________________________________
362 Int_t AliRawDB::Fill()
363 {
364   // Fill the trees and return the number of written bytes
365
366   // Create raw data TTree if it not yet done
367   if (!fTree) MakeTree();
368
369    Double_t bytes = fRawDB->GetBytesWritten();
370    Bool_t error = kFALSE;
371    if (fTree->Fill() == -1)
372      error = kTRUE;
373    if (fESDTree) 
374      if (fESDTree->Fill() == -1)
375        error = kTRUE;
376    if(!error)
377      return Int_t(fRawDB->GetBytesWritten() - bytes);
378    else
379      return -1;
380 }
381
382 //______________________________________________________________________________
383 Long64_t AliRawDB::GetTotalSize()
384 {
385    // Return the total size of the trees
386   Long64_t total = 0;
387
388   if (fTree) {
389     Int_t skey = 0;
390     TDirectory *dir = fTree->GetDirectory();
391     if (dir) {
392       TKey *key = dir->GetKey(fTree->GetName());
393       if (key) skey = key->GetKeylen();
394     }
395     total += (Long64_t)skey + fTree->GetZipBytes();
396   }
397
398   if(fESDTree)
399     {
400       Int_t skey = 0;
401       TDirectory *dir = fESDTree->GetDirectory();
402       if (dir) {
403         TKey *key = dir->GetKey(fESDTree->GetName());
404         if (key) skey = key->GetKeylen();
405       }
406       total += (Long64_t)skey + fESDTree->GetZipBytes();
407     }
408
409   return total;
410 }
411
412 //______________________________________________________________________________
413 Long64_t AliRawDB::AutoSave()
414 {
415   // Auto-save the raw-data and
416   // esd (if any) trees
417
418   Long64_t nbytes = fTree->AutoSave();
419
420   if (fESDTree) nbytes += fESDTree->AutoSave();
421
422   return nbytes;
423 }
424
425 //______________________________________________________________________________
426 Bool_t AliRawDB::NextFile(const char* fileName)
427 {
428    // Close te current file and open a new one.
429    // Returns kFALSE in case opening failed.
430
431    Close();
432
433    if (!Create(fileName)) return kFALSE;
434    return kTRUE;
435 }
436
437 //______________________________________________________________________________
438 Float_t AliRawDB::GetCompressionFactor() const
439 {
440    // Return compression factor.
441
442    if (fTree->GetZipBytes() == 0.)
443       return 1.0;
444    else
445       return fTree->GetTotBytes()/fTree->GetZipBytes();
446 }
447
448 //______________________________________________________________________________
449 const char *AliRawDB::GetAliRootTag()
450 {
451   // Return the aliroot tag (version)
452   // used to generate the raw data file.
453   // Stored in the raw-data file title.
454
455   static TString version = fgkAliRootTag;
456   version.Remove(TString::kBoth,'$');
457   version.ReplaceAll("Rev","AliRoot version");
458
459   return version.Data();
460 }
461
462 //______________________________________________________________________________
463 Bool_t AliRawDB::WriteGuidFile(TString &guidFileFolder)
464 {
465   // Write the guid file
466   // in the specified folder or
467   // in the folder where the raw data
468   // file is.
469
470    TString guidFileName;
471    if (!guidFileFolder.IsNull()) {
472      guidFileName = guidFileFolder;
473
474      TString pathStr = fRawDB->GetName();
475      TObjArray *pathArr = pathStr.Tokenize('/');
476      guidFileName.Append("/");
477      guidFileName.Append(((TObjString *)pathArr->Last())->String());
478      pathArr->Delete();
479      delete pathArr;
480    }
481    else
482      guidFileName = fRawDB->GetName();
483
484    guidFileName += ".guid";
485
486    ofstream fguid(guidFileName.Data());
487    if (!fguid.is_open()) {
488      Error("WriteGuidFile", "failure to open guid file %s", guidFileName.Data());
489      return kFALSE;
490    }
491    TString guid = fRawDB->GetUUID().AsString();
492    fguid << "guid: \t" << guid.Data();
493    fguid.close();
494
495    return kTRUE;
496 }
497
498
499 //______________________________________________________________________________
500 void AliRawDB::Reset()
501 {
502   // Clear the raw-data arrays
503   // Should be done before processing the raw-data event
504
505   for (Int_t iDet = 0; iDet < (AliDAQ::kNDetectors + 1); iDet++)
506     for (Int_t iBranch = 0; iBranch < fgkDetBranches[iDet]; iBranch++)
507       fDetRawData[iDet][iBranch]->ClearData();
508 }
509
510 //______________________________________________________________________________
511 AliRawDataArrayV2 *AliRawDB::GetRawDataArray(UInt_t eqSize, UInt_t eqId) const
512 {
513   // Return the corresponding raw-datra array (branch)
514   // depending on the equipment ID
515
516   Int_t iDet = AliDAQ::kNDetectors;
517   Int_t iBranch = 0; // can we split somehow the unrecognized data??? For the moment - no
518   if(eqSize) {
519     Int_t ddlIndex = -1;
520     iDet = AliDAQ::DetectorIDFromDdlID(eqId,ddlIndex);
521     if (iDet < 0 || iDet >= AliDAQ::kNDetectors)
522       iDet = AliDAQ::kNDetectors;
523     else
524       iBranch = (ddlIndex * fgkDetBranches[iDet])/AliDAQ::NumberOfDdls(iDet);
525   }
526
527   return fDetRawData[iDet][iBranch];
528 }
529