]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Connection to tree in Notify() seems to be safer.
authormorsch <morsch@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 9 Nov 2009 12:39:02 +0000 (12:39 +0000)
committermorsch <morsch@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 9 Nov 2009 12:39:02 +0000 (12:39 +0000)
ANALYSIS/AliMultiEventInputHandler.cxx
ANALYSIS/AliMultiEventInputHandler.h

index d36854dc35c4b16e296b4786b33b67c46f2a89a6..7c6da2987ae57aa51c8a8d77019fa0afef527413 100644 (file)
@@ -100,16 +100,22 @@ Bool_t AliMultiEventInputHandler::Init(TTree* tree, Option_t* /*opt*/)
 {
     // Initialisation necessary for each new tree
     fTree = tree;
-    fTree->GetEntry(0);
-    
     if (!fTree) return kFALSE;
-    // Get pointer to AOD event
-    fEventBuffer[0]->ReadFromTree(fTree, "");
+    for (Int_t i = 0; i < fBufferSize; i++) 
+       fEventBuffer[i]->Clear();
     fIndex     = 0;
     fNBuffered = 1;
     return kTRUE;
 }
 
+
+Bool_t AliMultiEventInputHandler::Notify(const char */*path*/)
+{
+    // Connect to new tree
+    fEventBuffer[0]->ReadFromTree(fTree, "");
+    return (kTRUE);
+}
+
 Bool_t AliMultiEventInputHandler::BeginEvent(Long64_t /*entry*/)
 {
     // Actions before analysis of each event 
@@ -130,7 +136,7 @@ Bool_t AliMultiEventInputHandler::FinishEvent()
     
     fIndex %= fBufferSize;
     AliInfo(Form("Connecting buffer entry %5d", fIndex));
-
+    fEventBuffer[fIndex]->Clear();
     fEventBuffer[fIndex]->ReadFromTree(fTree, "reconnect");
 
     fNBuffered++;
@@ -139,7 +145,6 @@ Bool_t AliMultiEventInputHandler::FinishEvent()
     return (kTRUE);
 }
 
-
 AliVEvent* AliMultiEventInputHandler::GetEvent(Int_t iev) const
 {
     // Get event number iev from buffer
index 76512c35a286597e54263761dd906fb4b42dfcc0..0b492a554c7ac1584d1bfd504d42c837dbcced1a 100644 (file)
@@ -37,7 +37,8 @@ class AliMultiEventInputHandler : public AliInputEventHandler {
     virtual Bool_t Init(TTree* tree, Option_t* /*opt*/);
     virtual Bool_t FinishEvent();
     virtual Bool_t BeginEvent(Long64_t /*entry*/);
-    
+    virtual Bool_t Notify() { return AliInputEventHandler::Notify();}
+    virtual Bool_t Notify(const char */*path*/);
  private:
     AliMultiEventInputHandler(const AliMultiEventInputHandler& handler);             
     AliMultiEventInputHandler& operator=(const AliMultiEventInputHandler& handler);