]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
memory leak fixed in track reconstruction and obsolete code removed
authorcussonno <cussonno@f7af4fe6-9843-0410-8265-dc069ae4e863>
Fri, 22 Jul 2005 07:59:56 +0000 (07:59 +0000)
committercussonno <cussonno@f7af4fe6-9843-0410-8265-dc069ae4e863>
Fri, 22 Jul 2005 07:59:56 +0000 (07:59 +0000)
MUON/AliMUONData.cxx
MUON/AliMUONEventReconstructor.cxx
MUON/AliMUONEventReconstructor.h
MUON/AliMUONTrack.cxx

index c44cab5fd4b441c267ea512032c08c0fa70d256e..9bd25ad095796548e827331f8187778a83cfb036 100644 (file)
@@ -845,14 +845,14 @@ void AliMUONData::ResetRecTracks()
 {
   // Reset tracks information
   fNrectracks = 0;
-  if (fRecTracks) fRecTracks->Clear();
+  if (fRecTracks) fRecTracks->Delete(); // necessary to delete in case of memory allocation
 }
 //____________________________________________________________________________
 void AliMUONData::ResetRecTriggerTracks()
 {
   // Reset tracks information
   fNrectriggertracks = 0;
-  if (fRecTriggerTracks) fRecTriggerTracks->Clear();
+  if (fRecTriggerTracks) fRecTriggerTracks->Delete(); // necessary to delete in case of memory allocation
 }
 //_____________________________________________________________________________
 void AliMUONData::SetTreeAddress(Option_t* option)
index b94f13ffda69fce51ec0878959e38d42ae1217b3..a3eeafce7379e2e5d6b5604db53493030b802da4 100644 (file)
@@ -47,7 +47,6 @@
 #include "AliMUONRawCluster.h"
 #include "AliMUONLocalTrigger.h"
 #include "AliMUONGlobalTrigger.h"
-#include "AliMUONRecoEvent.h"
 #include "AliMUONSegment.h"
 #include "AliMUONTrack.h"
 #include "AliMUONTrackHit.h"
@@ -103,9 +102,6 @@ AliMUONEventReconstructor::AliMUONEventReconstructor(AliLoader* loader)
   // Is 10 the right size ????
   fRecTracksPtr = new TClonesArray("AliMUONTrack", 10);
   fNRecTracks = 0; // really needed or GetEntriesFast sufficient ????
-// trigger tracks
-  fRecTriggerTracksPtr = new TClonesArray("AliMUONTriggerTrack", 10);
-  fNRecTriggerTracks = 0; // really needed or GetEntriesFast sufficient ????
   // Memory allocation for the TClonesArray of hits on reconstructed tracks
   // Is 100 the right size ????
   //AZ fRecTrackHitsPtr = new TClonesArray("AliMUONTrack", 100);
@@ -133,12 +129,7 @@ AliMUONEventReconstructor::AliMUONEventReconstructor(AliLoader* loader)
   fBkgTrackRefParticles = 0;
   fBkgTrackRefTTR = 0;
   fBkgTrackRefEventNumber = -1;
-  
-  // Initialize to 0 pointers to RecoEvent, tree and tree file
-  fRecoEvent = 0;
-  fEventTree = 0;
-  fTreeFile  = 0;
+   
   // initialize loader's
   fLoader = loader;
 
@@ -172,12 +163,6 @@ AliMUONEventReconstructor::operator=(const AliMUONEventReconstructor& rhs)
 AliMUONEventReconstructor::~AliMUONEventReconstructor(void)
 {
   // Destructor for class AliMUONEventReconstructor
-  if (fTreeFile) {
-     fTreeFile->Close();
-     delete fTreeFile;
-  }
-//  if (fEventTree) delete fEventTree;
-  if (fRecoEvent) delete fRecoEvent;
   delete fHitsForRecPtr; // Correct destruction of everything ???? or delete [] ????
   for (Int_t st = 0; st < AliMUONConstants::NTrackingCh()/2; st++)
     delete fSegmentsPtr[st]; // Correct destruction of everything ????
@@ -424,16 +409,6 @@ void AliMUONEventReconstructor::ResetTracks(void)
   return;
 }
 
-  //__________________________________________________________________________
-void AliMUONEventReconstructor::ResetTriggerTracks(void)
-{
-  // To reset the TClonesArray of reconstructed trigger tracks
-    if (fRecTriggerTracksPtr) fRecTriggerTracksPtr->Delete();
-  // Delete in order that the Track destructors are called,
-  // hence the space for the TClonesArray of pointers to TrackHit's is freed
-    fNRecTriggerTracks = 0;
-  return;
-}
 
   //__________________________________________________________________________
 void AliMUONEventReconstructor::ResetTrackHits(void)
@@ -961,7 +936,6 @@ void AliMUONEventReconstructor::ValidateTracksWithTrigger(void)
   AliMUONTrack *track;
   TClonesArray *recTriggerTracks;
   
-  fMUONData->ResetRecTriggerTracks();
   fMUONData->SetTreeAddress("RL");
   fMUONData->GetRecTriggerTracks();
   recTriggerTracks = fMUONData->RecTriggerTracks();
@@ -979,7 +953,6 @@ Bool_t AliMUONEventReconstructor::MakeTriggerTracks(void)
 {
     // To make the trigger tracks from Local Trigger
   AliDebug(1, "Enter MakeTriggerTracks");
-    //    ResetTriggerTracks();
     
     Int_t nTRentries;
     Long_t gloTrigPat;
@@ -1051,9 +1024,11 @@ Bool_t AliMUONEventReconstructor::MakeTriggerTracks(void)
       Float_t thetay = TMath::ATan2( (y21-y11) , (z21-z11) );
 
       recTriggerTrack = new AliMUONTriggerTrack(x11,y11,thetax,thetay,gloTrigPat);
+      
       // since static statement does not work, set gloTrigPat for each track
 
       fMUONData->AddRecTriggerTrack(*recTriggerTrack);
+      delete recTriggerTrack;
     } // end of loop on Local Trigger
     return kTRUE;    
 }
@@ -1617,36 +1592,6 @@ void AliMUONEventReconstructor::EventDumpTrigger(void)
   } 
 }
 
-//__________________________________________________________________________
-void AliMUONEventReconstructor::FillEvent()
-{
-// Create a new AliMUONRecoEvent, fill its track list, then add it as a
-// leaf in the Event branch of TreeRecoEvent tree
-   cout << "Enter FillEvent() ...\n";
-
-   if (!fRecoEvent) {
-      fRecoEvent = new AliMUONRecoEvent();
-   } else {
-      fRecoEvent->Clear();
-   }
-   //save current directory
-   TDirectory *current =  gDirectory;
-   if (!fTreeFile)  fTreeFile  = new TFile("tree_reco.root", "RECREATE");
-   if (!fEventTree) fEventTree = new TTree("TreeRecoEvent", "MUON reconstructed events");
-   //AZif (fRecoEvent->MakeDumpTracks(fRecTracksPtr)) {
-   if (fRecoEvent->MakeDumpTracks(fMuons, fRecTracksPtr, this)) { //AZ
-     if (AliLog::GetGlobalDebugLevel() > 1) fRecoEvent->EventInfo();
-      TBranch *branch = fEventTree->GetBranch("Event");
-      if (!branch) branch = fEventTree->Branch("Event", "AliMUONRecoEvent", &fRecoEvent, 64000);
-      branch->SetAutoDelete();
-      fTreeFile->cd();
-      fEventTree->Fill();
-      fTreeFile->Write();
-   }
-   // restore directory
-   current->cd();
-}
-
 //__________________________________________________________________________
 void AliMUONEventReconstructor::MakeTrackCandidatesK(void)
 {
index 5ba9a30e27d19ca929165c27ea66f7e12e5a8897..90ef83b53a3d887de0a1dcf608fd77d87bfd2cd2 100644 (file)
@@ -18,7 +18,6 @@ class AliMUONSegment;
 class TClonesArray;
 class TFile;
 class TTree;
-class AliMUONRecoEvent;
 class AliMUONData;
 class AliRunLoader;
 class AliLoader;
@@ -71,11 +70,6 @@ class AliMUONEventReconstructor : public TObject {
   void SetNRecTracks(Int_t NRecTracks) {fNRecTracks = NRecTracks;}
   TClonesArray* GetRecTracksPtr(void) const {return fRecTracksPtr;} // Array
  
- // Reconstructed trigger tracks
-   Int_t GetNRecTriggerTracks() const {return fNRecTriggerTracks;} // Number
-   void SetNRecTriggerTracks(Int_t NRecTriggerTracks) {fNRecTriggerTracks = NRecTriggerTracks;}
-   TClonesArray* GetRecTriggerTracksPtr(void) const {return fRecTriggerTracksPtr;} // Array
-
   // Hits on reconstructed tracks
   Int_t GetNRecTrackHits() const {return fNRecTrackHits;} // Number
   void SetNRecTrackHits(Int_t NRecTrackHits) {fNRecTrackHits = NRecTrackHits;}
@@ -172,19 +166,10 @@ class AliMUONEventReconstructor : public TObject {
   TClonesArray *fRecTracksPtr; // pointer to array of reconstructed tracks
   Int_t fNRecTracks; // number of reconstructed tracks
 
-  // Reconstructed trigger tracks
-  TClonesArray *fRecTriggerTracksPtr; // pointer to array of reconstructed trigger tracks
-  Int_t fNRecTriggerTracks; // number of reconstructed trigger tracks
-
   // Track hits on reconstructed tracks
   TClonesArray *fRecTrackHitsPtr; // pointer to array of hits on reconstructed tracks
   Int_t fNRecTrackHits; // number of hits on reconstructed tracks
 
-  // Objects needed for tree writing
-  AliMUONRecoEvent *fRecoEvent; // the reconstructed event
-  TTree            *fEventTree; // tree of reconstructed events
-  TFile            *fTreeFile;  // file where the tree is outputed
-
   // data container
   AliMUONData* fMUONData; // Data container for MUON subsystem 
 
@@ -209,7 +194,6 @@ class AliMUONEventReconstructor : public TObject {
   Bool_t MakeTriggerTracks(void);
   void ResetTrackHits(void);
   void ResetTracks(void);
-  void ResetTriggerTracks(void);
   Int_t MakeTrackCandidatesWithTwoSegments(AliMUONSegment *BegSegment);
   Int_t MakeTrackCandidatesWithOneSegmentAndOnePoint(AliMUONSegment *BegSegment);
   void MakeTrackCandidates(void);
index d490e53512558d0d64f8dd22e989f7f2c62c1b63..bed3e822d1926c0686a3a2f5b6e6f754653ccc16 100644 (file)
@@ -33,6 +33,8 @@
 #include <TObjArray.h>
 #include <TVirtualFitter.h>
 
+#include "AliLog.h"
+
 #include "AliMUONEventReconstructor.h" 
 #include "AliMUONHitForRec.h" 
 #include "AliMUONSegment.h" 
@@ -123,21 +125,19 @@ AliMUONTrack::~AliMUONTrack()
 {
   // Destructor
   if (fTrackHitsPtr) {
-    fTrackHitsPtr->Clear();
-    delete fTrackHitsPtr; // delete the TObjArray of pointers to TrackHit's
+    // delete the TObjArray of pointers to TrackHit's   
+    delete fTrackHitsPtr; 
     fTrackHitsPtr = NULL;
   }
   
   if (fTrackParamAtHit) {
     // delete the TClonesArray of pointers to TrackParam
-    fTrackParamAtHit->Clear();
     delete fTrackParamAtHit;
     fTrackParamAtHit = NULL;
   }
 
   if (fHitForRecAtHit) {
     // delete the TClonesArray of pointers to HitForRec
-    fHitForRecAtHit->Clear();
     delete fHitForRecAtHit;
     fHitForRecAtHit = NULL;
   }
@@ -160,6 +160,7 @@ AliMUONTrack::AliMUONTrack (const AliMUONTrack& theMUONTrack)
     AliMUONTrackHit *trackHit = new AliMUONTrackHit(*(AliMUONTrackHit*)(theMUONTrack.fTrackHitsPtr)->At(index));
     fTrackHitsPtr->Add(trackHit);
   }
+  fTrackHitsPtr->SetOwner(); // nedeed for deleting TClonesArray
 
   // necessary to make a copy of the objects and not only the pointers in TClonesArray.
   fTrackParamAtHit  =  new TClonesArray("AliMUONTrackParam",10);
@@ -196,8 +197,9 @@ AliMUONTrack & AliMUONTrack::operator=(const AliMUONTrack& theMUONTrack)
   // base class assignement
   TObject::operator=(theMUONTrack);
 
-  // fEventReconstructor =  new AliMUONEventReconstructor(*MUONTrack.fEventReconstructor); // is it right ?
-                               // is it right ? NO because it would use dummy copy constructor
+  // fEventReconstructor =  new AliMUONEventReconstructor(*MUONTrack.fEventReconstructor); 
+  // is it right ?
+  // is it right ? NO because it would use dummy copy constructor
   fEventReconstructor =  theMUONTrack.fEventReconstructor;
   fTrackParamAtVertex =  theMUONTrack.fTrackParamAtVertex;
 
@@ -207,6 +209,7 @@ AliMUONTrack & AliMUONTrack::operator=(const AliMUONTrack& theMUONTrack)
     AliMUONTrackHit *trackHit = new AliMUONTrackHit(*(AliMUONTrackHit*)(theMUONTrack.fTrackHitsPtr)->At(index));
     fTrackHitsPtr->Add(trackHit);
   }
+  fTrackHitsPtr->SetOwner();  // nedeed for deleting TClonesArray
 
   // necessary to make a copy of the objects and not only the pointers in TClonesArray.
   fTrackParamAtHit  =  new TClonesArray("AliMUONTrackParam",10);
@@ -545,6 +548,7 @@ void AliMUONTrack::AddHitForRec(AliMUONHitForRec* HitForRec)
     new ((*recTrackHitsPtr)[eventTrackHits]) AliMUONTrackHit(HitForRec);
   this->fEventReconstructor->SetNRecTrackHits(eventTrackHits + 1);
   // track
+  if (fTrackHitsPtr->IsOwner()) AliFatal("fTrackHitsPtr is owner");
   fTrackHitsPtr->Add(trackHit);
   fNTrackHits++;
 }