- Warnings corrected
authormorsch <morsch@f7af4fe6-9843-0410-8265-dc069ae4e863>
Wed, 5 Mar 2008 09:09:08 +0000 (09:09 +0000)
committermorsch <morsch@f7af4fe6-9843-0410-8265-dc069ae4e863>
Wed, 5 Mar 2008 09:09:08 +0000 (09:09 +0000)
- AliAODEvent with standard folder structure

STEER/AliAODCaloCells.cxx
STEER/AliAODCaloCells.h
STEER/AliAODEvent.cxx
STEER/AliAODEvent.h
STEER/AliAODTracklets.cxx
STEER/AliAODTracklets.h

index c2f1539..67764e1 100644 (file)
@@ -34,6 +34,38 @@ AliAODCaloCells::AliAODCaloCells(const char* name, const char* title, AODCells_t
   // TNamed constructor
 }
 
+AliAODCaloCells::AliAODCaloCells(const AliAODCaloCells& cells) :
+    TNamed(cells),
+    fNCells(cells.fNCells),
+    fCellNumber(0),
+    fAmplitude(0),
+    fIsSorted(cells.fIsSorted),
+    fType(cells.fType)
+{
+// Copy constructor
+
+    fCellNumber = new Short_t[fNCells];
+    fAmplitude =  new Double32_t[fNCells];
+
+    for (Int_t i = 0; i < fNCells; i++) {
+       fCellNumber[i]    = cells.fCellNumber[i];
+       fAmplitude[i]     = cells.fAmplitude[i];
+    }
+}
+
+AliAODCaloCells& AliAODCaloCells::operator=(const AliAODCaloCells& cells)
+{
+// Assignment operator
+    if(&cells == this) return *this;
+    TNamed::operator=(cells);
+    fNCells = cells.fNCells;
+    for (Int_t i = 0; i < fNCells; i++) {
+       fCellNumber[i]    = cells.fCellNumber[i];
+       fAmplitude[i]     = cells.fAmplitude[i];
+    }
+    return *this;
+}
+
 AliAODCaloCells::~AliAODCaloCells()
 {
   // destructor
index 091bb5d..89631fc 100644 (file)
@@ -23,6 +23,8 @@ class AliAODCaloCells : public TNamed
 
   AliAODCaloCells();
   AliAODCaloCells(const char* name, const char* title, AODCells_t ttype=kUndef);
+  AliAODCaloCells(const AliAODCaloCells& cells); 
+  AliAODCaloCells& operator=(const AliAODCaloCells& cells);
   
   virtual ~AliAODCaloCells();
   
@@ -48,9 +50,6 @@ class AliAODCaloCells : public TNamed
   Bool_t      fIsSorted;     //! true if cell arrays are sorted by index
   Char_t      fType;         // Cell type
   
- private:
-  AliAODCaloCells(const AliAODCaloCells& tow); 
-  AliAODCaloCells& operator=(const AliAODCaloCells& tow);
   
   ClassDef(AliAODCaloCells, 1);
 };
index ebc422d..348b42c 100644 (file)
@@ -20,7 +20,9 @@
 //     Author: Markus Oldenburg, CERN
 //-------------------------------------------------------------------------
 
+#include <TROOT.h>
 #include <TTree.h>
+#include <TFolder.h>
 
 #include "AliAODEvent.h"
 #include "AliAODHeader.h"
@@ -45,6 +47,7 @@ ClassImp(AliAODEvent)
 AliAODEvent::AliAODEvent() :
   AliVEvent(),
   fAODObjects(new TList()),
+  fAODFolder(0),
   fHeader(0),
   fTracks(0),
   fVertices(0),
@@ -60,11 +63,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,7 +192,7 @@ void AliAODEvent::CreateStdContent()
 
   // read back pointers
   GetStdContent();
-
+  CreateStdFolders();
   return;
 }
 
@@ -140,6 +217,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()
 {
index 268ccc0..d34a1c9 100644 (file)
@@ -30,6 +30,7 @@
 #include "AliAODFmdCluster.h"
 
 class TTree;
+class TFolder;
 
 class AliAODEvent : public AliVEvent {
 
@@ -51,8 +52,8 @@ class AliAODEvent : public AliVEvent {
   AliAODEvent();
   virtual ~AliAODEvent();
 
-  //AliAODEvent(const AliAODEvent& aodevent);  // not implemented
-  //AliAODEvent& operator=(const AliAODEvent& aodevent);  // not implemented
+  AliAODEvent(const AliAODEvent& aodevent);             
+  AliAODEvent& operator=(const AliAODEvent& aodevent);  
 
   void          AddObject(TObject *obj);
   void          RemoveObject(TObject *obj);
@@ -159,6 +160,7 @@ class AliAODEvent : public AliVEvent {
   void    CreateStdContent();
   void    SetStdNames();
   void    GetStdContent();
+  void    CreateStdFolders();
   void    ResetStd(Int_t trkArrSize = 0, 
                   Int_t vtxArrSize = 0, 
                   Int_t v0ArrSize = 0, 
@@ -174,8 +176,9 @@ class AliAODEvent : public AliVEvent {
 
  private :
 
-  TList *fAODObjects; //  list of AODObjects
-  
+  TList   *fAODObjects; //  list of AODObjects
+  TFolder *fAODFolder;  //  folder structure of branches
   // standard content
   AliAODHeader    *fHeader;       //! event information
   TClonesArray    *fTracks;       //! charged tracks
@@ -191,7 +194,7 @@ class AliAODEvent : public AliVEvent {
 
   static const char* fAODListName[kAODListN]; //!
 
-  ClassDef(AliAODEvent,3);
+  ClassDef(AliAODEvent, 4);
 };
 
 #endif
index 82892d1..412141d 100644 (file)
@@ -35,6 +35,41 @@ AliAODTracklets::AliAODTracklets(const char* name, const char* title) : TNamed(n
   // TNamed constructor
 }
 
+AliAODTracklets::AliAODTracklets(const AliAODTracklets& tracklet) :
+    TNamed(tracklet),
+    fNTracks(tracklet.fNTracks),
+    fTheta(0),
+    fPhi(0),
+    fDeltaPhi(0),
+    fLabels(0)
+{
+// Copy constructor
+    fTheta = new Double32_t[fNTracks];
+    fPhi = new Double32_t[fNTracks];
+    fDeltaPhi = new Double32_t[fNTracks];
+    fLabels = new Int_t[fNTracks];
+    for (Int_t i = 0; i < fNTracks; i++) {
+       fTheta[i]    = tracklet.fTheta[i];
+       fPhi[i]      = tracklet.fPhi[i];
+       fDeltaPhi[i] = tracklet.fDeltaPhi[i];
+       fLabels[i]   = tracklet.fLabels[i];
+    }
+}
+
+AliAODTracklets& AliAODTracklets::operator=(const AliAODTracklets& tracklet)
+{
+// Assignment operator
+    if(&tracklet == this) return *this;
+    TNamed::operator=(tracklet);
+    fNTracks = tracklet.fNTracks;
+        for (Int_t i = 0; i < fNTracks; i++) {
+       fTheta[i]    = tracklet.fTheta[i];
+       fPhi[i]      = tracklet.fPhi[i];
+       fDeltaPhi[i] = tracklet.fDeltaPhi[i];
+       fLabels[i]   = tracklet.fLabels[i];
+    }
+}
+
 void AliAODTracklets::CreateContainer(Int_t nTracks)
 {
   // function that creates container to store tracklets
@@ -54,6 +89,7 @@ void AliAODTracklets::CreateContainer(Int_t nTracks)
   fLabels = new Int_t[fNTracks];
 }
 
+
 AliAODTracklets::~AliAODTracklets()
 {
   // destructor
index 0e0bd3b..db98812 100644 (file)
@@ -19,6 +19,8 @@ class AliAODTracklets : public TNamed
  public:
   AliAODTracklets();
   AliAODTracklets(const char* name, const char* title);
+  AliAODTracklets(const AliAODTracklets& evt); 
+  AliAODTracklets& operator=(const AliAODTracklets& evt);
 
   virtual ~AliAODTracklets();
 
@@ -40,9 +42,6 @@ class AliAODTracklets : public TNamed
   Double32_t *fDeltaPhi;     //[fNTracks] array with delta phi values
   Int_t      *fLabels;       //[fNTracks] array with labels of tracklets
 
- private:
-  AliAODTracklets(const AliAODTracklets& evt); 
-  AliAODTracklets& operator=(const AliAODTracklets& evt);
 
   ClassDef(AliAODTracklets, 2);
 };