1 // @(#)alimdc:$Name$:$Id$
2 // Author: Fons Rademakers 26/11/99
4 /**************************************************************************
5 * Copyright(c) 1998-2003, ALICE Experiment at CERN, All rights reserved. *
7 * Author: The ALICE Off-line Project. *
8 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
19 //////////////////////////////////////////////////////////////////////////
23 //////////////////////////////////////////////////////////////////////////
31 #include "AliRawEvent.h"
32 #include "AliRawEventHeaderBase.h"
41 //______________________________________________________________________________
42 AliRawDB::AliRawDB(AliRawEvent *event,
45 const char* fileName) :
58 // Create a new raw DB
61 if (!Create(fileName))
66 //______________________________________________________________________________
67 AliRawDB::AliRawDB(const AliRawDB& rawDB): TObject(rawDB)
71 Fatal("AliRawDB", "copy constructor not implemented");
74 //______________________________________________________________________________
75 AliRawDB& AliRawDB::operator = (const AliRawDB& /*rawDB*/)
77 // assignment operator
79 Fatal("operator =", "assignment operator not implemented");
83 //______________________________________________________________________________
84 Bool_t AliRawDB::FSHasSpace(const char *fs) const
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.
89 Long_t id, bsize, blocks, bfree;
91 if (gSystem->GetFsInfo(fs, &id, &bsize, &blocks, &bfree) == 1) {
92 Error("FSHasSpace", "could not stat file system %s", fs);
96 // Leave 5 percent of diskspace free
97 Double_t avail = Double_t(bfree) * 0.95;
98 if (avail*bsize > fMaxSize)
101 Warning("FSHasSpace", "no space on file system %s", fs);
105 //______________________________________________________________________________
106 const char *AliRawDB::GetFileName() const
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.
113 static TString fname;
114 static Bool_t fstoggle = kFALSE;
116 TString fs = fstoggle ? fFS2 : fFS1;
119 TString hostname = gSystem->HostName();
121 if ((pos = hostname.Index(".")) != kNPOS)
122 hostname.Remove(pos);
124 if (!FSHasSpace(fs)) {
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
135 fname = fs + "/" + hostname + "_";
136 fname += dt.GetDate();
138 fname += dt.GetTime();
141 fstoggle = !fstoggle;
146 //______________________________________________________________________________
147 void AliRawDB::SetFS(const char* fs1, const char* fs2)
149 // set the file system location
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);
162 gSystem->ResetErrno();
163 gSystem->MakeDirectory(fs2);
164 if (gSystem->GetErrno() && gSystem->GetErrno() != EEXIST) {
165 SysError("SetFS", "mkdir %s", fs2);
170 //______________________________________________________________________________
171 Bool_t AliRawDB::Create(const char* fileName)
173 // Create a new raw DB.
175 const Int_t kMaxRetry = 1;
176 const Int_t kMaxSleep = 1; // seconds
177 const Int_t kMaxSleepLong = 10; // seconds
181 if (fStop) return kFALSE;
183 const char *fname = fileName;
184 if (!fname) fname = GetFileName();
186 Error("Create", "error getting raw DB file name");
192 fRawDB = TFile::Open(fname, GetOpenOption(),
193 Form("ALICE MDC%d raw DB", kMDC), fCompress,
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);
202 Error("Create", "failure to open file %s after %d tries", fname, kMaxRetry);
206 Warning("Create", "succeeded to open file after %d retries", retry);
208 if (fRawDB->IsZombie()) {
209 if (fRawDB->GetErrno() == ENOSPC ||
210 fRawDB->GetErrno() == 1018 || // SECOMERR
211 fRawDB->GetErrno() == 1027) { // SESYSERR
212 fRawDB->ResetErrno();
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
219 Error("Create", "file %s is zombie", fname);
220 fRawDB->ResetErrno();
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);
229 Error("Create", "failure to open file %s after %d tries", fname, kMaxRetry);
233 // Create raw data TTree
239 //______________________________________________________________________________
240 void AliRawDB::MakeTree()
242 // Create ROOT Tree object container.
244 fTree = new TTree("RAW", Form("ALICE MDC%d raw data tree", kMDC));
245 fTree->SetAutoSave(2000000000); // autosave when 2 Gbyte written
247 Int_t bufsize = 256000;
248 // splitting 29.6 MB/s, no splitting 35.3 MB/s on P4 2GHz 15k SCSI
251 fTree->Branch("rawevent", "AliRawEvent", &fEvent, bufsize, split);
253 // Create tree which will contain the HLT ESD information
256 fESDTree = new TTree("esdTree", Form("ALICE MDC%d HLT ESD tree", kMDC));
257 fESDTree->SetAutoSave(2000000000); // autosave when 2 Gbyte written
259 fESDTree->Branch("ESD", "AliESD", &fESD, bufsize, split);
264 //______________________________________________________________________________
265 Int_t AliRawDB::Close()
268 if (!fRawDB) return 0;
270 if (!fRawDB->IsOpen()) return 0;
275 Bool_t error = kFALSE;
276 if (fTree->Write() == 0)
279 if (fESDTree->Write() == 0)
282 // Close DB, this also deletes the fTree
285 Int_t filesize = fRawDB->GetEND();
288 gSystem->Unlink(fRawDB->GetName());
297 // Create semaphore to say this file is finished
298 Int_t tfd = ::creat(Form("%s.done", fRawDB->GetName()), 0644);
309 //______________________________________________________________________________
310 Int_t AliRawDB::Fill()
312 // Fill the trees and return the number of written bytes
314 Double_t bytes = fRawDB->GetBytesWritten();
315 Bool_t error = kFALSE;
316 if (fTree->Fill() == -1)
319 if (fESDTree->Fill() == -1)
322 return Int_t(fRawDB->GetBytesWritten() - bytes);
327 //______________________________________________________________________________
328 Int_t AliRawDB::GetTotalSize()
330 // Return the total size of the trees
335 TDirectory *dir = fTree->GetDirectory();
337 TKey *key = dir->GetKey(fTree->GetName());
338 if (key) skey = key->GetKeylen();
341 if (fTree->GetZipBytes() > 0) total += fTree->GetTotBytes();
342 TBuffer b(TBuffer::kWrite,10000);
343 TTree::Class()->WriteBuffer(b,fTree);
350 TDirectory *dir = fESDTree->GetDirectory();
352 TKey *key = dir->GetKey(fESDTree->GetName());
353 if (key) skey = key->GetKeylen();
356 if (fESDTree->GetZipBytes() > 0) total += fESDTree->GetTotBytes();
357 TBuffer b(TBuffer::kWrite,10000);
358 TTree::Class()->WriteBuffer(b,fESDTree);
365 //______________________________________________________________________________
366 void AliRawDB::WriteStats(AliStats* stats)
368 // Write stats to raw DB, local run DB and global MySQL DB.
370 AliRawEventHeaderBase &header = *GetEvent()->GetHeader();
372 // Write stats into RawDB
373 TDirectory *ds = gDirectory;
375 stats->SetEvents(GetEvents());
376 stats->SetLastId(header.Get("RunNb"), header.GetP("Id")[0]);
377 stats->SetFileSize(GetBytesWritten());
378 stats->SetCompressionFactor(GetCompressionFactor());
380 stats->Write("stats");
384 //______________________________________________________________________________
385 Bool_t AliRawDB::NextFile(const char* fileName)
387 // Close te current file and open a new one.
388 // Returns kFALSE in case opening failed.
392 if (!Create(fileName)) return kFALSE;
396 //______________________________________________________________________________
397 Float_t AliRawDB::GetCompressionFactor() const
399 // Return compression factor.
401 if (fTree->GetZipBytes() == 0.)
404 return fTree->GetTotBytes()/fTree->GetZipBytes();