]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RAW/AliRawDB.cxx
No need to liknk with lhapdf, pythia6 and microcern libraries
[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 <TBranch.h>
40
41 #include "AliESDEvent.h"
42 #include "AliRawEvent.h"
43 #include "AliRawDataArray.h"
44 #include "AliRawEventHeaderBase.h"
45 #include "AliRawEquipment.h"
46 #include "AliRawEquipmentHeader.h"
47 #include "AliStats.h"
48
49 #include "AliRawDB.h"
50
51
52 ClassImp(AliRawDB)
53
54 const char *AliRawDB::fgkAliRootTag = "$Rev$";
55
56 // Split TPC into 9 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(AliRawEvent *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 AliRawDataArray*[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 AliRawDataArray(nDDLsPerBranch);
86   }
87
88   fDetRawData[AliDAQ::kNDetectors] = new AliRawDataArray*[fgkDetBranches[AliDAQ::kNDetectors]];
89   for (Int_t iBranch = 0; iBranch < fgkDetBranches[AliDAQ::kNDetectors]; iBranch++)
90     fDetRawData[AliDAQ::kNDetectors][iBranch] = new AliRawDataArray(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    // Create raw data TTree
263    MakeTree();
264
265    return kTRUE;
266 }
267
268 static void BranchResetBit(TBranch *b) 
269 {
270   // Reset MapObject on this branch and all the sub-branches
271
272   b->ResetBit( kBranchObject | kBranchAny ); // Or in newer ROOT: b->ResetBit( kMapObject )
273   TIter next( b->GetListOfBranches() );
274   TBranch *sub = 0;
275   while ( (sub = (TBranch*)next() ) ) {
276     BranchResetBit( sub );
277   }
278 }
279
280 //______________________________________________________________________________
281 void AliRawDB::MakeTree()
282 {
283    // Create ROOT Tree object container.
284
285    fTree = new TTree("RAW", Form("ALICE raw-data tree (%s)", GetAliRootTag()));
286    fTree->SetAutoSave(21000000000LL);  // autosave when 21 Gbyte written
287
288    fTree->BranchRef();
289
290    // splitting 29.6 MB/s, no splitting 35.3 MB/s on P4 2GHz 15k SCSI
291    //Int_t split   = 1;
292    Int_t split   = 0;
293    TBranch *b = fTree->Branch("rawevent", "AliRawEvent", &fEvent, fBasketSize, split);
294    BranchResetBit(b);
295
296    // Make brach for each sub-detector
297    for (Int_t iDet = 0; iDet < AliDAQ::kNDetectors; iDet++) {
298      for (Int_t iBranch = 0; iBranch < fgkDetBranches[iDet]; iBranch++) {
299        b = fTree->Branch(Form("%s%d",AliDAQ::DetectorName(iDet),iBranch),"AliRawDataArray",
300                          &fDetRawData[iDet][iBranch],fBasketSize,split);
301        BranchResetBit(b);
302      }
303    }
304    // Make special branch for unrecognized raw-data payloads
305    for (Int_t iBranch = 0; iBranch < fgkDetBranches[AliDAQ::kNDetectors]; iBranch++) {
306      b = fTree->Branch(Form("Common%d",iBranch),"AliRawDataArray",
307                        &fDetRawData[AliDAQ::kNDetectors][iBranch],fBasketSize,split);
308      BranchResetBit(b);
309    }
310
311    // Create tree which will contain the HLT ESD information
312
313    if (fESD) {
314      fESDTree = new TTree("esdTree", Form("ALICE HLT ESD tree (%s)", GetAliRootTag()));
315      fESDTree->SetAutoSave(21000000000LL);  // autosave when 21 Gbyte written
316      split   = 0;
317      fESDTree->Branch("ESD", "AliESDEvent", &fESD, fBasketSize, split);
318    }
319
320 }
321
322 //______________________________________________________________________________
323 Long64_t AliRawDB::Close()
324 {
325    // Close raw DB.
326    if (!fRawDB) return 0;
327
328    if (!fRawDB->IsOpen()) return 0;
329
330    fRawDB->cd();
331
332    // Write the tree.
333    Bool_t error = kFALSE;
334    if (fTree->Write() == 0)
335      error = kTRUE;
336    if (fESDTree)
337      if (fESDTree->Write() == 0)
338        error = kTRUE;
339
340    // Close DB, this also deletes the fTree
341    fRawDB->Close();
342
343    Long64_t filesize = fRawDB->GetEND();
344
345    if (fDeleteFiles) {
346       gSystem->Unlink(fRawDB->GetName());
347       delete fRawDB;
348       fRawDB = 0;
349       if(!error)
350         return filesize;
351       else
352         return -1;
353    }
354
355    delete fRawDB;
356    fRawDB = 0;
357    if(!error)
358      return filesize;
359    else
360      return -1;
361 }
362
363 //______________________________________________________________________________
364 Int_t AliRawDB::Fill()
365 {
366    // Fill the trees and return the number of written bytes
367
368   for (Int_t iDet = 0; iDet < (AliDAQ::kNDetectors + 1); iDet++)
369     for (Int_t iBranch = 0; iBranch < fgkDetBranches[iDet]; iBranch++)
370       fDetRawData[iDet][iBranch]->ClearData();
371
372    // Move the raw-data payloads to the corresponding branches
373   for(Int_t iSubEvent = 0; iSubEvent < fEvent->GetNSubEvents(); iSubEvent++) {
374     AliRawEvent *subEvent = fEvent->GetSubEvent(iSubEvent);
375     for(Int_t iEquipment = 0; iEquipment < subEvent->GetNEquipments(); iEquipment++) {
376       AliRawEquipment *equipment = subEvent->GetEquipment(iEquipment);
377       Int_t iDet = AliDAQ::kNDetectors;
378       Int_t iBranch = 0; // can we split somehow the unrecognized data??? For the moment - no
379       if(equipment->GetEquipmentHeader()->GetEquipmentSize()) {
380         UInt_t eqId = equipment->GetEquipmentHeader()->GetId();
381         Int_t ddlIndex = -1;
382         iDet = AliDAQ::DetectorIDFromDdlID(eqId,ddlIndex);
383         if (iDet < 0 || iDet >= AliDAQ::kNDetectors)
384           iDet = AliDAQ::kNDetectors;
385         else
386           iBranch = (ddlIndex * fgkDetBranches[iDet])/AliDAQ::NumberOfDdls(iDet);
387       }
388       equipment->SetRawDataRef(fDetRawData[iDet][iBranch]);
389     }
390   }
391
392    Double_t bytes = fRawDB->GetBytesWritten();
393    Bool_t error = kFALSE;
394    if (fTree->Fill() == -1)
395      error = kTRUE;
396    if (fESDTree) 
397      if (fESDTree->Fill() == -1)
398        error = kTRUE;
399    if(!error)
400      return Int_t(fRawDB->GetBytesWritten() - bytes);
401    else
402      return -1;
403 }
404
405 //______________________________________________________________________________
406 Long64_t AliRawDB::GetTotalSize()
407 {
408    // Return the total size of the trees
409   Long64_t total = 0;
410
411   {
412     Int_t skey = 0;
413     TDirectory *dir = fTree->GetDirectory();
414     if (dir) {
415       TKey *key = dir->GetKey(fTree->GetName());
416       if (key) skey = key->GetKeylen();
417     }
418     total += (Long64_t)skey + fTree->GetZipBytes();
419   }
420
421   if(fESDTree)
422     {
423       Int_t skey = 0;
424       TDirectory *dir = fESDTree->GetDirectory();
425       if (dir) {
426         TKey *key = dir->GetKey(fESDTree->GetName());
427         if (key) skey = key->GetKeylen();
428       }
429       total += (Long64_t)skey + fESDTree->GetZipBytes();
430     }
431
432   return total;
433 }
434
435 //______________________________________________________________________________
436 Long64_t AliRawDB::AutoSave()
437 {
438   // Auto-save the raw-data and
439   // esd (if any) trees
440
441   Long64_t nbytes = fTree->AutoSave();
442
443   if (fESDTree) nbytes += fESDTree->AutoSave();
444
445   return nbytes;
446 }
447
448 //______________________________________________________________________________
449 void AliRawDB::WriteStats(AliStats* stats)
450 {
451    // Write stats to raw DB, local run DB and global MySQL DB.
452
453    AliRawEventHeaderBase &header = *GetEvent()->GetHeader();
454
455    // Write stats into RawDB
456    TDirectory *ds = gDirectory;
457    GetDB()->cd();
458    stats->SetEvents(GetEvents());
459    stats->SetLastId(header.GetP("Id")[0]);
460    stats->SetFileSize(GetBytesWritten());
461    stats->SetCompressionFactor(GetCompressionFactor());
462    stats->SetEndTime();
463    stats->Write("stats");
464    ds->cd();
465 }
466
467 //______________________________________________________________________________
468 Bool_t AliRawDB::NextFile(const char* fileName)
469 {
470    // Close te current file and open a new one.
471    // Returns kFALSE in case opening failed.
472
473    Close();
474
475    if (!Create(fileName)) return kFALSE;
476    return kTRUE;
477 }
478
479 //______________________________________________________________________________
480 Float_t AliRawDB::GetCompressionFactor() const
481 {
482    // Return compression factor.
483
484    if (fTree->GetZipBytes() == 0.)
485       return 1.0;
486    else
487       return fTree->GetTotBytes()/fTree->GetZipBytes();
488 }
489
490 //______________________________________________________________________________
491 const char *AliRawDB::GetAliRootTag()
492 {
493   // Return the aliroot tag (version)
494   // used to generate the raw data file.
495   // Stored in the raw-data file title.
496
497   static TString version = fgkAliRootTag;
498   version.Remove(TString::kBoth,'$');
499   version.ReplaceAll("Rev","AliRoot version");
500
501   return version.Data();
502 }
503
504 //______________________________________________________________________________
505 Bool_t AliRawDB::WriteGuidFile(TString &guidFileFolder)
506 {
507   // Write the guid file
508   // in the specified folder or
509   // in the folder where the raw data
510   // file is.
511
512    TString guidFileName;
513    if (!guidFileFolder.IsNull()) {
514      guidFileName = guidFileFolder;
515
516      TString pathStr = fRawDB->GetName();
517      TObjArray *pathArr = pathStr.Tokenize('/');
518      guidFileName.Append("/");
519      guidFileName.Append(((TObjString *)pathArr->Last())->String());
520      pathArr->Delete();
521      delete pathArr;
522    }
523    else
524      guidFileName = fRawDB->GetName();
525
526    guidFileName += ".guid";
527
528    ofstream fguid(guidFileName.Data());
529    if (!fguid.is_open()) {
530      Error("WriteGuidFile", "failure to open guid file %s", guidFileName.Data());
531      return kFALSE;
532    }
533    TString guid = fRawDB->GetUUID().AsString();
534    fguid << "guid: \t" << guid.Data();
535    fguid.close();
536
537    return kTRUE;
538 }