Memory leaks fixed. (M. Bondila)
authormorsch <morsch@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 12 Nov 2001 14:31:00 +0000 (14:31 +0000)
committermorsch <morsch@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 12 Nov 2001 14:31:00 +0000 (14:31 +0000)
EVGEN/AliGenExtFile.cxx
EVGEN/AliGenReaderCwn.cxx
EVGEN/AliGenReaderCwn.h
EVGEN/AliGenReaderTreeK.cxx
EVGEN/AliGenReaderTreeK.h

index eae2d34..e9a83b1 100644 (file)
@@ -15,6 +15,9 @@
 
 /*
 $Log$
+Revision 1.17  2001/11/09 09:12:58  morsch
+Generalization by using AliGenReader object to read particles from file.
+
 Revision 1.16  2001/07/27 17:09:35  morsch
 Use local SetTrack, KeepTrack and SetHighWaterMark methods
 to delegate either to local stack or to stack owned by AliRun.
@@ -163,6 +166,7 @@ void AliGenExtFile::Generate()
         y     < fYMin     || y     > fYMax        )
       {
          printf("\n Not selected %d %f %f %f %f %f", i, theta, phi, pmom, pt, y);
+         delete iparticle;
          continue;
       }
       p[0] = iparticle->Px();
@@ -178,7 +182,8 @@ void AliGenExtFile::Generate()
          }
       }
       SetTrack(fTrackIt,-1,idpart,p,origin,polar,0,kPPrimary,nt);
-  }
+      delete iparticle;
+   }
   TFile *pFile=0;
 // Get AliRun object or create it 
   if (!gAlice) {
index 48223be..2c164ac 100644 (file)
@@ -16,6 +16,9 @@
 
 /*
 $Log$
+Revision 1.1  2001/11/09 09:10:46  morsch
+Realisation of AliGenReader that reads the old cwn event format.
+
 */
 
 // Read the old ALICE event format based on CW-ntuples
@@ -40,6 +43,11 @@ AliGenReaderCwn::AliGenReaderCwn()
     fTreeNtuple = 0;
 }
 
+AliGenReaderCwn::~AliGenReaderCwn()
+{
+    delete fTreeNtuple;
+}
+
 void AliGenReaderCwn::Init() 
 {
 //
@@ -131,6 +139,11 @@ TParticle* AliGenReaderCwn::NextParticle()
 
 
 
+AliGenReaderCwn& AliGenReaderCwn::operator=(const  AliGenReaderCwn& rhs)
+{
+// Assignment operator
+    return *this;
+}
 
 
 
index 9fb4ebb..50ebd7d 100644 (file)
@@ -14,13 +14,13 @@ class AliGenReaderCwn : public AliGenReader
     AliGenReaderCwn();
     
     AliGenReaderCwn(const AliGenReaderCwn &reader){;}
-    virtual ~AliGenReaderCwn(){;}
-    // Initialise 
+    virtual ~AliGenReaderCwn();
+        // Initialise 
     virtual void Init();
     // Read
     virtual Int_t NextEvent();
     virtual TParticle*  NextParticle();
-    AliGenReaderCwn & operator=(const AliGenReader & rhs){return *this;}
+    AliGenReaderCwn & operator=(const AliGenReaderCwn & rhs);
  protected:
     Int_t             fNcurrent;      // points to the next entry
     Int_t             fNparticle;     // particle number in event
index 880c742..13d5944 100644 (file)
@@ -15,6 +15,9 @@
 
 /*
 $Log$
+Revision 1.1  2001/11/09 09:11:24  morsch
+Realisation of AliGenReader that reads the kine tree (TreeK).
+
 */
 #include <TFile.h>
 #include <TTree.h>
@@ -38,6 +41,13 @@ AliGenReaderTreeK::AliGenReaderTreeK():AliGenReader()
     fNparticle      = 0;
     fFile           = 0;
     fBaseFile       = 0;
+    fTreeE          = 0;
+}
+
+AliGenReaderTreeK::~AliGenReaderTreeK() 
+{
+// Destructor
+    delete fTreeE;
 }
 
 void AliGenReaderTreeK::Init() 
@@ -63,13 +73,17 @@ Int_t AliGenReaderTreeK::NextEvent()
     fStack = new AliStack(1000);
     fFile->cd();
 //  Connect treeE
-    TTree* treeE = (TTree*)gDirectory->Get("TE");
-    treeE->ls();
+    if (fTreeE) {
+       delete fTreeE;
+    } else {
+       fTreeE = (TTree*)gDirectory->Get("TE");
+    }
+
     if (fHeader) delete fHeader;
     fHeader = 0;
-    treeE->SetBranchAddress("Header", &fHeader);
+    fTreeE->SetBranchAddress("Header", &fHeader);
 //  Get next event
-    treeE->GetEntry(fNcurrent);
+    fTreeE->GetEntry(fNcurrent);
     fStack = fHeader->Stack();
     fStack->GetEvent(fNcurrent);
     
@@ -94,6 +108,11 @@ TParticle* AliGenReaderTreeK::NextParticle()
 
 
 
+AliGenReaderTreeK& AliGenReaderTreeK::operator=(const  AliGenReaderTreeK& rhs)
+{
+// Assignment operator
+    return *this;
+}
 
 
 
index 29f3df0..6c30527 100644 (file)
@@ -15,15 +15,15 @@ class AliGenReaderTreeK : public AliGenReader
 {
  public:
     AliGenReaderTreeK();
-    
     AliGenReaderTreeK(const AliGenReaderTreeK &reader){;}
-    virtual ~AliGenReaderTreeK(){;}
+    virtual ~AliGenReaderTreeK();
     // Initialise 
     virtual void Init();
     // Read
     virtual Int_t NextEvent();
     virtual TParticle*  NextParticle();
-    AliGenReaderTreeK & operator=(const AliGenReader & rhs){return *this;}
+    AliGenReaderTreeK & operator=(const AliGenReader & rhs);
+    
  protected:
     Int_t             fNcurrent;          // points to the next entry
     Int_t             fNparticle;         // Next particle in list
@@ -32,6 +32,7 @@ class AliGenReaderTreeK : public AliGenReader
     TFile            *fBaseFile;          //! pointer to base file
     AliStack         *fStack;             //! Particle stack
     AliHeader        *fHeader;            //! Pointer to event header
+    TTree            *fTreeE;             //! Pointer to header tree
     ClassDef(AliGenReaderTreeK,1) // Read particles from TreeK
 };
 #endif