]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - STEER/AliAODEvent.cxx
Removing osolete variable. Adding title to trigger configuration obj also
[u/mrichter/AliRoot.git] / STEER / AliAODEvent.cxx
index 2518e7be4366deedecd03d1897bd4c783051047b..1ed88d1cb3160cc28a20d9600d430f2edfb5aa3d 100644 (file)
 //     Author: Markus Oldenburg, CERN
 //-------------------------------------------------------------------------
 
+#include <TROOT.h>
 #include <TTree.h>
+#include <TFolder.h>
+#include <TFriendElement.h>
+#include <TProcessID.h>
+#include <TCollection.h>
 
 #include "AliAODEvent.h"
 #include "AliAODHeader.h"
@@ -45,6 +50,7 @@ ClassImp(AliAODEvent)
 AliAODEvent::AliAODEvent() :
   AliVEvent(),
   fAODObjects(new TList()),
+  fAODFolder(0),
   fHeader(0),
   fTracks(0),
   fVertices(0),
@@ -60,11 +66,85 @@ AliAODEvent::AliAODEvent() :
   // default constructor
 }
 
+//______________________________________________________________________________
+AliAODEvent::AliAODEvent(const AliAODEvent& aod):
+  AliVEvent(aod),
+  fAODObjects(new TList()),
+  fAODFolder(new TFolder()),
+  fHeader(new AliAODHeader(*aod.fHeader)),
+  fTracks(new TClonesArray(*aod.fTracks)),
+  fVertices(new TClonesArray(*aod.fVertices)),
+  fV0s(new TClonesArray(*aod.fV0s)),
+  fTracklets(new AliAODTracklets(*aod.fTracklets)),
+  fJets(new TClonesArray(*aod.fJets)),
+  fEmcalCells(new AliAODCaloCells(*aod.fEmcalCells)),
+  fPhosCells(new AliAODCaloCells(*aod.fPhosCells)),
+  fCaloClusters(new TClonesArray(*aod.fCaloClusters)),
+  fFmdClusters(new TClonesArray(*aod.fFmdClusters)),
+  fPmdClusters(new TClonesArray(*aod.fPmdClusters))
+{
+  // Copy constructor
+  AddObject(fHeader);
+  AddObject(fTracks);
+  AddObject(fVertices);
+  AddObject(fV0s);
+  AddObject(fTracklets);
+  AddObject(fJets);
+  AddObject(fEmcalCells);
+  AddObject(fPhosCells);
+  AddObject(fCaloClusters);
+  AddObject(fFmdClusters);
+  AddObject(fPmdClusters);
+
+  GetStdContent();
+}
+
+//______________________________________________________________________________
+AliAODEvent & AliAODEvent::operator=(const AliAODEvent& aod) {
+
+    // Assignment operator
+
+    if(&aod == this) return *this;
+    AliVEvent::operator=(aod);
+
+    fAODObjects      = new TList();
+    fAODFolder       = new TFolder();
+    fHeader          = new AliAODHeader(*aod.fHeader);
+    fTracks          = new TClonesArray(*aod.fTracks);
+    fVertices        = new TClonesArray(*aod.fVertices);
+    fV0s             = new TClonesArray(*aod.fV0s);
+    fTracklets       = new AliAODTracklets(*aod.fTracklets);
+    fJets            = new TClonesArray(*aod.fJets);
+    fEmcalCells      = new AliAODCaloCells(*aod.fEmcalCells);
+    fPhosCells       = new AliAODCaloCells(*aod.fPhosCells);
+    fCaloClusters    = new TClonesArray(*aod.fCaloClusters);
+    fFmdClusters     = new TClonesArray(*aod.fFmdClusters);
+    fPmdClusters     = new TClonesArray(*aod.fPmdClusters);
+    
+    fAODObjects = new TList();
+    
+    AddObject(fHeader);
+    AddObject(fTracks);
+    AddObject(fVertices);
+    AddObject(fV0s);
+    AddObject(fTracklets);
+    AddObject(fJets);
+    AddObject(fEmcalCells);
+    AddObject(fPhosCells);
+    AddObject(fCaloClusters);
+    AddObject(fFmdClusters);
+    AddObject(fPmdClusters);
+    GetStdContent();
+    return *this;
+}
+
+
 //______________________________________________________________________________
 AliAODEvent::~AliAODEvent() 
 {
 // destructor
     delete fAODObjects;
+    delete fAODFolder;
 }
 
 //______________________________________________________________________________
@@ -115,10 +195,25 @@ void AliAODEvent::CreateStdContent()
 
   // read back pointers
   GetStdContent();
-
+  CreateStdFolders();
   return;
 }
 
+void  AliAODEvent::MakeEntriesReferencable()
+{
+    // Make all entries referencable in a subsequent process
+    //
+    TIter next(fAODObjects);
+    TObject* obj;
+    while ((obj = next()))
+    {
+       if(obj->InheritsFrom("TCollection"))
+           {
+               AssignIDtoCollection((TCollection*)obj);
+           }
+    }
+}
+
 //______________________________________________________________________________
 void AliAODEvent::SetStdNames()
 {
@@ -140,6 +235,25 @@ void AliAODEvent::SetStdNames()
   }
 } 
 
+void AliAODEvent::CreateStdFolders()
+{
+    // Create the standard folder structure
+    fAODFolder = gROOT->GetRootFolder()->AddFolder("AOD", "AOD");
+    if(fAODObjects->GetEntries()==kAODListN){
+       for(int i = 0;i < fAODObjects->GetEntries();i++){
+           TObject *fObj = fAODObjects->At(i);
+           if(fObj->InheritsFrom("TClonesArray")){
+               fAODFolder->AddFolder(fAODListName[i], fAODListName[i], (TCollection*) fObj);
+           } else {
+               fAODFolder->AddFolder(fAODListName[i], fAODListName[i], 0);
+           }
+       }
+    }
+    else{
+       printf("%s:%d CreateStdFolders() Wrong number of Std Entries \n",(char*)__FILE__,__LINE__);
+    }
+} 
+
 //______________________________________________________________________________
 void AliAODEvent::GetStdContent()
 {
@@ -168,7 +282,6 @@ void AliAODEvent::ResetStd(Int_t trkArrSize,
                           Int_t pmdClusSize)
 {
   // deletes content of standard arrays and resets size 
-
   fTracks->Delete();
   if (trkArrSize > fTracks->GetSize()) 
     fTracks->Expand(trkArrSize);
@@ -219,6 +332,47 @@ void AliAODEvent::ClearStd()
   fPmdClusters   ->Clear();
 }
 
+//_________________________________________________________________
+Int_t AliAODEvent::GetPHOSClusters(TRefArray *clusters) const
+{
+  // fills the provided TRefArray with all found phos clusters
+  
+  clusters->Clear();
+  
+  AliAODCaloCluster *cl = 0;
+  for (Int_t i = 0; i < GetNCaloClusters() ; i++) {
+    
+    if ( (cl = GetCaloCluster(i)) ) {
+      if (cl->IsPHOSCluster()){
+       clusters->Add(cl);
+       //AliDebug(1,Form("IsPHOS cluster %d Size: %d \n",i,clusters->GetEntriesFast()));
+      }
+    }
+  }
+  return clusters->GetEntriesFast();
+}
+
+//_________________________________________________________________
+Int_t AliAODEvent::GetEMCALClusters(TRefArray *clusters) const
+{
+  // fills the provided TRefArray with all found emcal clusters
+
+  clusters->Clear();
+
+  AliAODCaloCluster *cl = 0;
+  for (Int_t i = 0; i < GetNCaloClusters(); i++) {
+
+    if ( (cl = GetCaloCluster(i)) ) {
+      if (cl->IsEMCALCluster()){
+       clusters->Add(cl);
+       //AliDebug(1,Form("IsEMCAL cluster %d Size: %d \n",i,clusters->GetEntriesFast()));
+      }
+    }
+  }
+  return clusters->GetEntriesFast();
+}
+
+
 //______________________________________________________________________________
 Int_t AliAODEvent::GetMuonTracks(TRefArray *muonTracks) const
 {
@@ -237,9 +391,9 @@ Int_t AliAODEvent::GetMuonTracks(TRefArray *muonTracks) const
 }
 
 
-void AliAODEvent::ReadFromTree(TTree *tree)
+void AliAODEvent::ReadFromTree(TTree *tree, Option_t* opt /*= ""*/)
 {
-  // connects aod event to tree
+  // Connects aod event to tree
   
   if(!tree){
     Printf("%s %d AliAODEvent::ReadFromTree() Zero Pointer to Tree \n",(char*)__FILE__,__LINE__);
@@ -254,12 +408,12 @@ void AliAODEvent::ReadFromTree(TTree *tree)
   if(aodEvent){
     // Check if already connected to tree
     TList* connectedList = (TList*) (tree->GetUserInfo()->FindObject("AODObjectsConnectedToTree"));
-    if (connectedList) {
-      // If connected use the connected list if objects
-      fAODObjects->Delete();
-      fAODObjects = connectedList;
-      GetStdContent(); 
-      return;
+    if (connectedList && (strcmp(opt, "reconnect"))) {
+       // If connected use the connected list if objects
+       fAODObjects->Delete();
+       fAODObjects = connectedList;
+       GetStdContent(); 
+       return;
     }
     // Connect to tree
     // prevent a memory leak when reading back the TList
@@ -273,28 +427,53 @@ void AliAODEvent::ReadFromTree(TTree *tree)
       printf("%s %d AliAODEvent::ReadFromTree() TList contains less than the standard contents %d < %d \n",
             (char*)__FILE__,__LINE__,fAODObjects->GetEntries(),kAODListN);
     }
-    // set the branch addresses
+    //
+    // Let's find out whether we have friends
+    TList* friendL = tree->GetTree()->GetListOfFriends();
+    if (friendL) 
+    {
+       TIter next(friendL);
+       TFriendElement* fe;
+       while ((fe = (TFriendElement*)next())){
+           aodEvent = (AliAODEvent*)(fe->GetTree()->GetUserInfo()->FindObject("AliAODEvent"));
+           if (!aodEvent) {
+               printf("No UserInfo on tree \n");
+           } else {
+
+               TList* objL = (TList*)(aodEvent->GetList()->Clone());
+               printf("Get list of object from tree %d !!\n", objL->GetEntries());
+               TIter nextobject(objL);
+               TObject* obj =  0;
+               while((obj = nextobject()))
+               {
+                   printf("Adding object from friend %s !\n", obj->GetName());
+                   fAODObjects->Add(obj);
+               } // object "branch" loop
+           } // has userinfo  
+       } // friend loop
+    } // has friends   
+           
+
+// set the branch addresses
     TIter next(fAODObjects);
     TNamed *el;
     while((el=(TNamed*)next())){
       TString bname(el->GetName());
       // check if branch exists under this Name
-      TBranch *br = tree->GetBranch(bname.Data());
+      TBranch *br = tree->GetTree()->GetBranch(bname.Data());
       if(br){
        tree->SetBranchAddress(bname.Data(),fAODObjects->GetObjectRef(el));
-      }
-      else{
-       br = tree->GetBranch(Form("%s.",bname.Data()));
-       if(br){
-         tree->SetBranchAddress(Form("%s.",bname.Data()),fAODObjects->GetObjectRef(el));
-       }
-       else{
-         printf("%s %d AliAODEvent::ReadFromTree() No Branch found with Name %s or %s. \n",
-                (char*)__FILE__,__LINE__,bname.Data(),bname.Data());
-       }       
+      } else {
+         br = tree->GetBranch(Form("%s.",bname.Data()));
+         if(br){
+             tree->SetBranchAddress(Form("%s.",bname.Data()),fAODObjects->GetObjectRef(el));
+         }
+         else{
+             printf("%s %d AliAODEvent::ReadFromTree() No Branch found with Name %s. \n",
+                    (char*)__FILE__,__LINE__,bname.Data());
+         }     
       }
     }
-    
     GetStdContent();
     // when reading back we are not owner of the list 
     // must not delete it
@@ -328,3 +507,15 @@ void AliAODEvent::Print(Option_t *) const
   
   return;
 }
+
+void AliAODEvent::AssignIDtoCollection(TCollection* col)
+{
+    // Static method which assigns a ID to each object in a collection
+    // In this way the objects are marked as referenced and written with 
+    // an ID. This has the advantage that TRefs to this objects can be 
+    // written by a subsequent process.
+    TIter next(col);
+    TObject* obj;
+    while ((obj = next()))
+       TProcessID::AssignID(obj);
+}