]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Added functionality for ESD merging and copying of TClonesArrays, not yet enabled
authorrichterm <richterm@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 8 Sep 2008 08:27:48 +0000 (08:27 +0000)
committerrichterm <richterm@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 8 Sep 2008 08:27:48 +0000 (08:27 +0000)
HLT/configure.ac
HLT/rec/AliHLTEsdManagerImplementation.cxx
HLT/rec/AliHLTEsdManagerImplementation.h

index 57ec0486ee2cb002e32023287561df3276c4ffe2..72042a972b07d9444fb45966310b468f0ecb9810 100644 (file)
@@ -252,9 +252,32 @@ if test ! "x$have_aliroot" = "xno" ; then
     ALIROOT_LIBS=
   fi
 
+  dnl
+  dnl ESD supports non-std content
+  dnl
+  have_esd_nonstd=no
+  AC_RUN_IFELSE([AC_LANG_PROGRAM([#include <AliESDEvent.h>
+                                  #include <AliExternalTrackParam.h>
+                                  #include <TTree.h>
+                                  #include <TClonesArray.h>],
+                                 [AliESDEvent esd;
+                                  esd.CreateStdContent();
+                                  TTree* tree=new TTree("esdTree", "Tree with HLT ESD objects");
+                                  esd.WriteToTree(tree);
+                                  TClonesArray* a=new TClonesArray("AliExternalTrackParam");
+                                  a->SetName("SomeObject");
+                                  esd.AddObject(a);
+                                  if (!tree->FindBranch("SomeObject")) return 1;
+                                  return 0;])],
+                 [have_esd_nonstd=yes], 
+                [AC_DEFINE(HAVE_NOT_ESD_NONSTD)])
+  AC_MSG_CHECKING(whether ESD supports non standard content)
+  AC_MSG_RESULT([$have_esd_nonstd])
+
   dnl
   dnl ESD copy function added May 9 2008 rev 25667
   dnl
+  if test "x$have_esd_nonstd" != "xyes"; then
   have_esd_copy=no
   AC_LINK_IFELSE([AC_LANG_PROGRAM([#include <AliESDEvent.h>],
                                  [AliESDEvent esd1;
@@ -264,6 +287,7 @@ if test ! "x$have_aliroot" = "xno" ; then
                 [AC_DEFINE(HAVE_NOT_ESD_COPY)])
   AC_MSG_CHECKING(for ESD assignment operator)
   AC_MSG_RESULT([$have_esd_copy])
+  fi
 
   dnl
   dnl AliRawReaderMemory support for multiple buffers added
index 936ecb15e914c31c30fc4a379473e3825da73472..c56493195452b6dffb05134bcd2736484e1b40e9 100644 (file)
@@ -27,6 +27,7 @@
 #include "AliESDEvent.h"
 #include "AliHLTMessage.h"
 #include "AliESDEvent.h"
+#include "AliESDtrack.h"
 #include "TFile.h"
 #include "TTree.h"
 #include "TClass.h"
@@ -117,7 +118,15 @@ int AliHLTEsdManagerImplementation::WriteESD(const AliHLTUInt8_t* pBuffer, AliHL
          fESDs.push_back(newEntry);
        }
        if (tgtesd) {
-#ifndef HAVE_NOT_ESD_COPY
+#if 0 // later extension !defined(HAVE_NOT_ESD_NONSTD)
+         entry=Find(dt);
+         if (entry) {
+           entry->CopyNonEmptyObjects(tgtesd, pESD);
+         } else {
+           HLTError("internal mismatch, can not create list entry");
+           iResult=-ENOMEM;
+         }
+#elif !defined(HAVE_NOT_ESD_COPY)
          *tgtesd=*pESD;
 #else //HAVE_NOT_ESD_COPY
          static bool warningPrinted=false;
@@ -217,7 +226,8 @@ AliHLTEsdManagerImplementation::AliHLTEsdListEntry::AliHLTEsdListEntry(AliHLTCom
   fDt(dt),
   fpFile(NULL),
   fpTree(NULL),
-  fpEsd(NULL)
+  fpEsd(NULL),
+  fPrefix()
 {
   // see header file for class documentation
 }
@@ -251,15 +261,11 @@ int AliHLTEsdManagerImplementation::AliHLTEsdListEntry::WriteESD(AliESDEvent* pS
 #ifndef HAVE_NOT_ESD_COPY
   if (fName.IsNull()) {
     // this is the first event, create the file name
-    TString origin;
-    origin.Insert(0, fDt.fOrigin, kAliHLTComponentDataTypefOriginSize);
-    origin.Remove(TString::kTrailing, ' ');
-    origin.ToUpper();
     fName="";
     if (!fDirectory.IsNull()) {
       fName+=fDirectory; fName+="/";
     }
-    fName+="AliHLT"; fName+=origin;
+    fName+="Ali"; fName+=GetPrefix();
     if (fDt!=kAliHLTDataTypeESDObject &&
        fDt!=kAliHLTDataTypeESDTree) {
 
@@ -586,3 +592,99 @@ const char* AliHLTEsdManagerImplementation::AliHLTEsdListEntry::GetFileName() co
   // see header file for class documentation
   return fName.Data();
 }
+
+const char* AliHLTEsdManagerImplementation::AliHLTEsdListEntry::GetPrefix()
+{
+  // see header file for class documentation
+  if (fPrefix.IsNull()) {
+    fPrefix.Insert(0, fDt.fOrigin, kAliHLTComponentDataTypefOriginSize);
+    fPrefix.Remove(TString::kTrailing, ' ');
+    fPrefix.ToUpper();
+    if (!fPrefix.Contains("HLT")) {
+      fPrefix.Insert(0, "HLT");
+    }
+  }
+  return fPrefix.Data();
+}
+
+int AliHLTEsdManagerImplementation::AliHLTEsdListEntry::CopyNonEmptyObjects(AliESDEvent* pTgt, AliESDEvent* pSrc)
+{
+  // see header file for class documentation
+  int iResult=0;
+  if (!pTgt || !pSrc) return -EINVAL;
+
+  const char* defaultPrefix="HLT";
+  TIter next(pSrc->GetList());
+  TObject* pSrcObject=NULL;
+  TString name;
+  while ((pSrcObject=next())) {
+    if(!pSrcObject->InheritsFrom("TCollection")){
+      // simple objects
+    } else if(pSrcObject->InheritsFrom("TClonesArray")){
+      TClonesArray* pTClA=dynamic_cast<TClonesArray*>(pSrcObject);
+      if (pTClA!=NULL && pTClA->GetEntriesFast()>0) {
+       bool bMakeNewName=true;
+       name=defaultPrefix;
+       name+=pTClA->GetName();
+       TObject* pTgtObject=pTgt->GetList()->FindObject(name);
+       TClonesArray* pTgtArray=NULL;
+       if (bMakeNewName=(pTgtObject!=NULL) && pTgtObject->InheritsFrom("TClonesArray")){
+         pTgtArray=dynamic_cast<TClonesArray*>(pTgtObject);
+         if (pTgtArray) {
+           TString classType=pTClA->Class()->GetName();
+           if (classType.CompareTo(pTgtArray->Class()->GetName())==0) {
+             if (pTgtArray->GetEntries()==0) {
+               bMakeNewName=false;
+             } else {
+               HLTWarning("TClonesArray \"%s\"  in target ESD %p is already filled with %d entries",
+                          name.Data(), pTgt, pTgtArray->GetEntries());
+             }
+           } else {
+             HLTWarning("TClonesArray \"%s\" exists in target ESD %p, but describes incompatible class type %s instead of %s",
+                        name.Data(), pTgt, pTgtArray->GetClass()->GetName(), pTClA->GetClass()->GetName());
+           }
+         } else {
+           HLTError("internal error: dynamic cast failed for object %s %p", pTgtObject->GetName(), pTgtObject);
+         }
+       } else if (pTgtObject) {
+         HLTWarning("object \"%s\" does already exist in target ESD %p, but is %s rather than TClonesArray",
+                    name.Data(), pTgt, pTgtObject->Class()->GetName());
+         // TODO: temporary solution, think about a general naming scheme and add
+         // the name as property of AliHLTEsdListEntry
+       }
+
+       if (bMakeNewName) {
+         pTgtArray=NULL;
+         int count=1;
+         while (pTgt->GetList()->FindObject(name)) {
+           name.Form("%sESD_%s", GetPrefix(), pTClA->GetName());
+           if (count++>1) {
+             name+=Form("%d", count);
+           }
+         }
+         HLTWarning("adding new TClonesArray \"%s\" because of conflicts", name.Data());
+       }
+
+       if (pTgtArray) {
+         pTgtArray->ExpandCreate(pTClA->GetEntries());
+         HLTInfo("Expanding TClonesArray \"%s\" to %d elements", pTgtArray->GetClass()->GetName(), pTClA->GetEntries());
+       } else {
+         pTgtArray=new TClonesArray(pTClA->GetClass(), pTClA->GetEntries());
+         pTgtArray->ExpandCreate(pTClA->GetEntries());
+         pTgtArray->SetName(name);
+         pTgt->AddObject(pTgtArray);
+         HLTInfo("Adding TClonesArray \"%s\" with %d elements to ESD %p", name.Data(), pTClA->GetEntries(), pTgt);
+       }
+
+       if (pTgtArray) {
+         for(int i=0; i<pTClA->GetEntriesFast(); ++i){
+           (*pTClA)[i]->Copy(*((*pTgtArray)[i]));
+         }
+       } else {
+         iResult=-ENOMEM;
+       }
+      }
+    }
+  }
+  return iResult;
+}
index 830d9680420f17756bc11024fad489b75284b139..a8c553048de8857eb9626fe3bb860fd64602a08a 100644 (file)
@@ -107,6 +107,13 @@ class AliHLTEsdManagerImplementation : public AliHLTEsdManager {
      */
     int WriteESD(AliESDEvent* pESD, int eventno=-1);
 
+    /**
+     * Copy non-empty objects of the source ESD to the target ESD.
+     * The first implementation just copies the entries of type TClonesArray which
+     * are not empty.
+     */
+    int CopyNonEmptyObjects(AliESDEvent* pTgt, AliESDEvent* pSrc);
+
     /**
      * Set the target directory for the ESD file.
      */
@@ -122,6 +129,13 @@ class AliHLTEsdManagerImplementation : public AliHLTEsdManager {
      */
     const char* GetFileName() const;
 
+    /**
+     * Get the object name prefix generated from the data origin
+     * The prefix is added to the names of the ESD objects when copied to the
+     * master ESD.
+     */
+    const char* GetPrefix();
+
     bool operator==(AliHLTComponentDataType dt) const;
 
   private:
@@ -151,6 +165,8 @@ class AliHLTEsdManagerImplementation : public AliHLTEsdManager {
     TTree* fpTree; //!transient
     /** the esd to fill into the tree */
     AliESDEvent* fpEsd; //!transient
+    /** Prefix for generated ESD objects in the master ESD */
+    TString fPrefix; //!transient
   };
 
   typedef vector<AliHLTEsdListEntry*> AliHLTEsdPList;