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