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