]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Method to get an approximate output file size is added. AliMDC ProcessEvent method...
authorcvetan <cvetan@f7af4fe6-9843-0410-8265-dc069ae4e863>
Wed, 16 Mar 2005 11:44:00 +0000 (11:44 +0000)
committercvetan <cvetan@f7af4fe6-9843-0410-8265-dc069ae4e863>
Wed, 16 Mar 2005 11:44:00 +0000 (11:44 +0000)
RAW/AliMDC.cxx
RAW/AliMDC.h
RAW/AliRawDB.cxx
RAW/AliRawDB.h
RAW/mdc.cxx
RAW/mdc.h

index 9da6b27cc6fe11ffb667502141e6bb2f915fc592..ed929a6c4456470c38156101e5215208151e780f 100644 (file)
@@ -227,6 +227,18 @@ Int_t AliMDC::ProcessEvent(void* event, Bool_t isIovecArray)
 // or, if isIovecArray is kTRUE, a pointer to an array of iovecs with one
 // iovec per subevent (used by the event builder).
 // The return value is the number of written bytes or an error code
+  const UInt_t kFileSizeWarningLevel = 1800000000;
+  const UInt_t kFileSizeErrorLevel   = 1900000000;
+
+  UInt_t currentFileSize = GetTotalSize();
+  if(currentFileSize > kFileSizeErrorLevel) {
+    Error("ProcessEvent", "file size (%u) exceeds the limit "
+         , currentFileSize);
+    return kErrFileSize;
+  }
+  if(currentFileSize > kFileSizeWarningLevel)
+    Warning("ProcessEvent", "file size (%u) is close to the limit "
+           , currentFileSize);
 
   Int_t status;
   char* data = (char*) event;
@@ -373,6 +385,16 @@ Int_t AliMDC::ProcessEvent(void* event, Bool_t isIovecArray)
   return nBytes;
 }
 
+//______________________________________________________________________________
+Int_t AliMDC::GetTotalSize()
+{
+// return the total current raw DB file size
+
+  if (!fRawDB) return -1;
+
+  return fRawDB->GetTotalSize();
+}
+
 //______________________________________________________________________________
 Int_t AliMDC::Close()
 {
index ee7815bced7a764a38a02cc1c9de91222dc28447..2e238e010709ab5168b2cf39f60e2ec10ff9694b 100644 (file)
@@ -47,7 +47,8 @@ public:
                     kErrSubHeader = -4, 
                     kErrDataSize = -5, 
                     kErrEquipmentHeader = -6, 
-                    kErrEquipment = -7 };
+                    kErrEquipment = -7,
+                     kErrFileSize = -8 };
 
    AliMDC(Int_t compress, Bool_t deleteFiles, 
          EFilterMode filterMode = kFilterTransparent, 
@@ -58,6 +59,7 @@ public:
 
    Int_t      Open(EWriteMode mode, const char* fileName);
    Int_t      ProcessEvent(void* event, Bool_t isIovecArray = kFALSE);
+   Int_t      GetTotalSize();
    Int_t      Close();
 
    Int_t      Run(const char* inputFile, Bool_t loop,
index 5259a34447903bf96d9b31238461d9c82a60d587..34bcbfa48e15830d1baccb52bdc008ae205e1af2 100644 (file)
@@ -25,6 +25,7 @@
 #include <errno.h>
 
 #include <TSystem.h>
+#include <TKey.h>
 
 #ifdef ALI_DATE
 #include "event.h"
@@ -315,6 +316,44 @@ Int_t AliRawDB::Fill()
    return Int_t(fRawDB->GetBytesWritten() - bytes);
 }
 
+//______________________________________________________________________________
+Int_t AliRawDB::GetTotalSize()
+{
+   // Return the total size of the trees
+  Int_t total = 0;
+
+  {
+    Int_t skey = 0;
+    TDirectory *dir = fTree->GetDirectory();
+    if (dir) {
+      TKey *key = dir->GetKey(fTree->GetName());
+      if (key) skey = key->GetKeylen();
+    }
+    total += skey;
+    if (fTree->GetZipBytes() > 0) total += fTree->GetTotBytes();
+    TBuffer b(TBuffer::kWrite,10000);
+    TTree::Class()->WriteBuffer(b,fTree);
+    total += b.Length();
+  }
+
+  if(fESDTree)
+    {
+      Int_t skey = 0;
+      TDirectory *dir = fESDTree->GetDirectory();
+      if (dir) {
+       TKey *key = dir->GetKey(fESDTree->GetName());
+       if (key) skey = key->GetKeylen();
+      }
+      total += skey;
+      if (fESDTree->GetZipBytes() > 0) total += fESDTree->GetTotBytes();
+      TBuffer b(TBuffer::kWrite,10000);
+      TTree::Class()->WriteBuffer(b,fESDTree);
+      total += b.Length();
+    }
+
+  return total;
+}
+
 //______________________________________________________________________________
 void AliRawDB::WriteStats(AliStats* stats)
 {
index d223f848f574bef7d63811ce26c427a7d3088271..bd178fe36604dc8d17a405e1797fdc11c09317e2 100644 (file)
@@ -49,6 +49,7 @@ public:
    virtual Bool_t      Create(const char* fileName = NULL);
    virtual void        Close();
    Int_t               Fill();
+   Int_t               GetTotalSize();
 
    void         WriteStats(AliStats* stats);
 
index a3797e9cfecd37f118d102e7e56ea51a5244de47..d142c454524027b9a18f8a1b05fed3f4fc675ef8 100644 (file)
@@ -52,6 +52,13 @@ int alimdcProcessEvent(void* alimdc, void* event, int isIovecArray)
   return ((AliMDC*)alimdc)->ProcessEvent(event, isIovecArray);
 }
 
+int alimdcGetTotalFileSize(void* alimdc)
+{
+// return the total current file size
+
+  return ((AliMDC*)alimdc)->GetTotalSize();
+}
+
 int alimdcClose(void* alimdc)
 {
 // close the raw DB
index 46860d36a4f809c849c7ba60b00a1b684db38ec0..38851f8cb3f7b27724dae392233ab238763b3afc 100644 (file)
--- a/RAW/mdc.h
+++ b/RAW/mdc.h
@@ -21,6 +21,7 @@ void* alimdcCreate(int compress, int filterMode,
                   double maxSizeTagDB, const char* fileNameTagDB);
 int   alimdcOpen(void* alimdc, int mode, const char* fileName);
 int   alimdcProcessEvent(void* alimdc, void* event, int isIovecArray);
+int   alimdcGetTotalFileSize(void* alimdc);
 int   alimdcClose(void* alimdc);
 void  alimdcDelete(void* alimdc);
 void  alimdcEnableDebug();