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