]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RAW/AliRawDB.cxx
One run number entry in the raw data file statistics
[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
27 #include <TSystem.h>
28 #include <TKey.h>
29
30 #include "AliESD.h"
31 #include "AliRawEvent.h"
32 #include "AliRawEventHeaderBase.h"
33 #include "AliStats.h"
34
35 #include "AliRawDB.h"
36
37
38 ClassImp(AliRawDB)
39
40 const char *AliRawDB::fgkAliRootTag = "$Name$";
41
42 //______________________________________________________________________________
43 AliRawDB::AliRawDB(AliRawEvent *event,
44                    AliESD *esd, 
45                    Int_t compress,
46                    const char* fileName) :
47   fRawDB(NULL),
48   fTree(NULL),
49   fEvent(event),
50   fESDTree(NULL),
51   fESD(esd),
52   fCompress(compress),
53   fMaxSize(-1),
54   fFS1(""),
55   fFS2(""),
56   fDeleteFiles(kFALSE),
57   fStop(kFALSE)
58 {
59    // Create a new raw DB
60
61    if (fileName) {
62       if (!Create(fileName))
63          MakeZombie();
64    }
65 }
66
67 //______________________________________________________________________________
68 Bool_t AliRawDB::FSHasSpace(const char *fs) const
69 {
70    // Check for at least fMaxSize bytes of free space on the file system.
71    // If the space is not available return kFALSE, kTRUE otherwise.
72
73    Long_t id, bsize, blocks, bfree;
74
75    if (gSystem->GetFsInfo(fs, &id, &bsize, &blocks, &bfree) == 1) {
76       Error("FSHasSpace", "could not stat file system %s", fs);
77       return kFALSE;
78    }
79
80    // Leave 5 percent of diskspace free
81    Double_t avail = Double_t(bfree) * 0.95;
82    if (avail*bsize > fMaxSize)
83       return kTRUE;
84
85    Warning("FSHasSpace", "no space on file system %s", fs);
86    return kFALSE;
87 }
88
89 //______________________________________________________________________________
90 const char *AliRawDB::GetFileName() const
91 {
92    // Return filename based on hostname and date and time. This will make
93    // each file unique. Also makes sure (via FSHasSpace()) that there is
94    // enough space on the file system to store the file. Returns 0 in
95    // case of error or interrupt signal.
96
97    static TString fname;
98    static Bool_t  fstoggle = kFALSE;
99
100    TString fs = fstoggle ? fFS2 : fFS1;
101    TDatime dt;
102
103    TString hostname = gSystem->HostName();
104    Int_t pos;
105    if ((pos = hostname.Index(".")) != kNPOS)
106       hostname.Remove(pos);
107
108    if (!FSHasSpace(fs)) {
109       while (1) {
110          fstoggle = !fstoggle;
111          fs = fstoggle ? fFS2 : fFS1;
112          if (FSHasSpace(fs)) break;
113          Info("GetFileName", "sleeping 30 seconds before retrying...");
114          gSystem->Sleep(30000);   // sleep for 30 seconds
115          if (fStop) return 0;
116       }
117    }
118
119    fname = fs + "/" + hostname + "_";
120    fname += dt.GetDate();
121    fname += "_";
122    fname += dt.GetTime();
123    fname += ".root";
124
125    fstoggle = !fstoggle;
126
127    return fname;
128 }
129
130 //______________________________________________________________________________
131 void AliRawDB::SetFS(const char* fs1, const char* fs2)
132 {
133 // set the file system location
134
135   fFS1 = fs1;
136   if (fs1 && !fFS1.Contains(":")) {
137     gSystem->ResetErrno();
138     gSystem->MakeDirectory(fs1);
139     if (gSystem->GetErrno() && gSystem->GetErrno() != EEXIST) {
140       SysError("SetFS", "mkdir %s", fs1);
141     }
142   }
143
144   fFS2 = fs2;
145   if (fs2) {
146     gSystem->ResetErrno();
147     gSystem->MakeDirectory(fs2);
148     if (gSystem->GetErrno() && gSystem->GetErrno() != EEXIST) {
149       SysError("SetFS", "mkdir %s", fs2);
150     }
151   }
152 }
153
154 //______________________________________________________________________________
155 Bool_t AliRawDB::Create(const char* fileName)
156 {
157    // Create a new raw DB.
158
159    const Int_t kMaxRetry = 1;
160    const Int_t kMaxSleep = 1;      // seconds
161    const Int_t kMaxSleepLong = 10; // seconds
162    Int_t retry = 0;
163
164 again:
165    if (fStop) return kFALSE;
166
167    const char *fname = fileName;
168    if (!fname) fname = GetFileName();
169    if (!fname) {
170       Error("Create", "error getting raw DB file name");
171       return kFALSE;
172    }
173
174    retry++;
175
176    fRawDB = TFile::Open(fname, GetOpenOption(),
177                         Form("ALICE raw-data file (%s)", GetAliRootTag()), fCompress,
178                         GetNetopt());
179    if (!fRawDB) {
180       if (retry < kMaxRetry) {
181          Warning("Create", "failure to open file, sleeping %d %s before retrying...",
182                  kMaxSleep, kMaxSleep==1 ? "second" : "seconds");
183          gSystem->Sleep(kMaxSleep*1000);
184          goto again;
185       }
186       Error("Create", "failure to open file %s after %d tries", fname, kMaxRetry);
187       return kFALSE;
188    }
189    if (retry > 1)
190       Warning("Create", "succeeded to open file after %d retries", retry);
191
192    if (fRawDB->IsZombie()) {
193       if (fRawDB->GetErrno() == ENOSPC ||
194           fRawDB->GetErrno() == 1018   ||   // SECOMERR
195           fRawDB->GetErrno() == 1027) {     // SESYSERR
196          fRawDB->ResetErrno();
197          delete fRawDB;
198          Warning("Create", "file is a zombie (no space), sleeping %d %s before retrying...",
199                  kMaxSleepLong, kMaxSleepLong==1 ? "second" : "seconds");
200          gSystem->Sleep(kMaxSleepLong*1000);   // sleep 10 seconds before retrying
201          goto again;
202       }
203       Error("Create", "file %s is zombie", fname);
204       fRawDB->ResetErrno();
205       delete fRawDB;
206       fRawDB = 0;
207       if (retry < kMaxRetry) {
208          Warning("Create", "file is a zombie, sleeping %d %s before retrying...",
209                  kMaxSleep, kMaxSleep==1 ? "second" : "seconds");
210          gSystem->Sleep(kMaxSleep*1000);
211          goto again;
212       }
213       Error("Create", "failure to open file %s after %d tries", fname, kMaxRetry);
214       return kFALSE;
215    }
216
217    // Create raw data TTree
218    MakeTree();
219
220    return kTRUE;
221 }
222
223 //______________________________________________________________________________
224 void AliRawDB::MakeTree()
225 {
226    // Create ROOT Tree object container.
227
228    fTree = new TTree("RAW", Form("ALICE raw-data tree (%s)", GetAliRootTag()));
229    fTree->SetAutoSave(2000000000);  // autosave when 2 Gbyte written
230
231    Int_t bufsize = 256000;
232    // splitting 29.6 MB/s, no splitting 35.3 MB/s on P4 2GHz 15k SCSI
233    //Int_t split   = 1;
234    Int_t split   = 0;
235    fTree->Branch("rawevent", "AliRawEvent", &fEvent, bufsize, split);
236
237    // Create tree which will contain the HLT ESD information
238
239    if (fESD) {
240      fESDTree = new TTree("esdTree", Form("ALICE HLT ESD tree (%s)", GetAliRootTag()));
241      fESDTree->SetAutoSave(2000000000);  // autosave when 2 Gbyte written
242      split   = 0;
243      fESDTree->Branch("ESD", "AliESD", &fESD, bufsize, split);
244    }
245
246 }
247
248 //______________________________________________________________________________
249 Int_t AliRawDB::Close()
250 {
251    // Close raw DB.
252    if (!fRawDB) return 0;
253
254    if (!fRawDB->IsOpen()) return 0;
255
256    fRawDB->cd();
257
258    // Write the tree.
259    Bool_t error = kFALSE;
260    if (fTree->Write() == 0)
261      error = kTRUE;
262    if (fESDTree)
263      if (fESDTree->Write() == 0)
264        error = kTRUE;
265
266    // Close DB, this also deletes the fTree
267    fRawDB->Close();
268
269    Int_t filesize = fRawDB->GetEND();
270
271    if (fDeleteFiles) {
272       gSystem->Unlink(fRawDB->GetName());
273       delete fRawDB;
274       fRawDB = 0;
275       if(!error)
276         return filesize;
277       else
278         return -1;
279    }
280
281    delete fRawDB;
282    fRawDB = 0;
283    if(!error)
284      return filesize;
285    else
286      return -1;
287 }
288
289 //______________________________________________________________________________
290 Int_t AliRawDB::Fill()
291 {
292    // Fill the trees and return the number of written bytes
293
294    Double_t bytes = fRawDB->GetBytesWritten();
295    Bool_t error = kFALSE;
296    if (fTree->Fill() == -1)
297      error = kTRUE;
298    if (fESDTree) 
299      if (fESDTree->Fill() == -1)
300        error = kTRUE;
301    if(!error)
302      return Int_t(fRawDB->GetBytesWritten() - bytes);
303    else
304      return -1;
305 }
306
307 //______________________________________________________________________________
308 Int_t AliRawDB::GetTotalSize()
309 {
310    // Return the total size of the trees
311   Int_t total = 0;
312
313   {
314     Int_t skey = 0;
315     TDirectory *dir = fTree->GetDirectory();
316     if (dir) {
317       TKey *key = dir->GetKey(fTree->GetName());
318       if (key) skey = key->GetKeylen();
319     }
320     total += skey;
321     if (fTree->GetZipBytes() > 0) total += fTree->GetTotBytes();
322     TBuffer b(TBuffer::kWrite,10000);
323     TTree::Class()->WriteBuffer(b,fTree);
324     total += b.Length();
325   }
326
327   if(fESDTree)
328     {
329       Int_t skey = 0;
330       TDirectory *dir = fESDTree->GetDirectory();
331       if (dir) {
332         TKey *key = dir->GetKey(fESDTree->GetName());
333         if (key) skey = key->GetKeylen();
334       }
335       total += skey;
336       if (fESDTree->GetZipBytes() > 0) total += fESDTree->GetTotBytes();
337       TBuffer b(TBuffer::kWrite,10000);
338       TTree::Class()->WriteBuffer(b,fESDTree);
339       total += b.Length();
340     }
341
342   return total;
343 }
344
345 //______________________________________________________________________________
346 void AliRawDB::WriteStats(AliStats* stats)
347 {
348    // Write stats to raw DB, local run DB and global MySQL DB.
349
350    AliRawEventHeaderBase &header = *GetEvent()->GetHeader();
351
352    // Write stats into RawDB
353    TDirectory *ds = gDirectory;
354    GetDB()->cd();
355    stats->SetEvents(GetEvents());
356    stats->SetLastId(header.GetP("Id")[0]);
357    stats->SetFileSize(GetBytesWritten());
358    stats->SetCompressionFactor(GetCompressionFactor());
359    stats->SetEndTime();
360    stats->Write("stats");
361    ds->cd();
362 }
363
364 //______________________________________________________________________________
365 Bool_t AliRawDB::NextFile(const char* fileName)
366 {
367    // Close te current file and open a new one.
368    // Returns kFALSE in case opening failed.
369
370    Close();
371
372    if (!Create(fileName)) return kFALSE;
373    return kTRUE;
374 }
375
376 //______________________________________________________________________________
377 Float_t AliRawDB::GetCompressionFactor() const
378 {
379    // Return compression factor.
380
381    if (fTree->GetZipBytes() == 0.)
382       return 1.0;
383    else
384       return fTree->GetTotBytes()/fTree->GetZipBytes();
385 }
386
387 //______________________________________________________________________________
388 const char *AliRawDB::GetAliRootTag()
389 {
390   // Return the aliroot tag (version)
391   // used to generate the raw data file.
392   // Stored in the raw-data file title.
393
394   TString version = fgkAliRootTag;
395   version.Remove(TString::kBoth,'$');
396   version.ReplaceAll("Name","AliRoot version");
397
398   return version.Data();
399 }