]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Correct memory leaks.
authorgconesab <gconesab@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 11 Mar 2010 14:20:41 +0000 (14:20 +0000)
committergconesab <gconesab@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 11 Mar 2010 14:20:41 +0000 (14:20 +0000)
Correct wagon.

13 files changed:
PWG4/PartCorrBase/AliAnaPartCorrBaseClass.cxx
PWG4/PartCorrBase/AliAnaPartCorrBaseClass.h
PWG4/PartCorrBase/AliAnaPartCorrMaker.cxx
PWG4/PartCorrBase/AliAnalysisTaskParticleCorrelation.cxx
PWG4/PartCorrDep/AliAnaCalorimeterQA.cxx
PWG4/PartCorrDep/AliAnaElectron.cxx
PWG4/PartCorrDep/AliAnaExample.cxx
PWG4/PartCorrDep/AliAnaOmegaToPi0Gamma.cxx
PWG4/PartCorrDep/AliAnaOmegaToPi0Gamma.h
PWG4/PartCorrDep/AliAnaParticleIsolation.cxx
PWG4/PartCorrDep/AliAnaPhoton.cxx
PWG4/PartCorrDep/AliAnaPi0EbE.cxx
PWG4/macros/AddTaskPartCorr.C

index b248b03cb30e6575b55da6895f43e1b3e063e1bb..24f9cd8d460370a228c6f38f48a216f5b79e11f5 100755 (executable)
@@ -50,6 +50,7 @@ ClassImp(AliAnaPartCorrBaseClass)
     fAODObjArrayName(""), fAddToHistogramsName(""),
     fAODCaloCells(0x0),//fAODCaloClusters(0x0),  
     fCaloPID(0x0), fFidCut(0x0), fIC(0x0),fMCUtils(0x0), fNMS(0x0),
+       //fAnaOutContainer(0x0),
     fHistoPtBins(0),   fHistoPtMax(0.),   fHistoPtMin(0.),
     fHistoPhiBins(0),  fHistoPhiMax(0.),  fHistoPhiMin(0.),
     fHistoEtaBins(0),  fHistoEtaMax(0.),  fHistoEtaMin(0.),
@@ -57,13 +58,7 @@ ClassImp(AliAnaPartCorrBaseClass)
        fHistoAsymBins(0), fHistoAsymMax(0.), fHistoAsymMin(0.)
 {
   //Default Ctor
-  
-  fReader  = new AliCaloTrackReader();
-  fCaloPID = new AliCaloPID();
-  fFidCut  = new AliFiducialCut();
-  fIC      = new AliIsolationCut();
-  fMCUtils = new AliMCAnalysisUtils();
-  
+    
   //Initialize parameters
   InitParameters();
 }
@@ -82,6 +77,7 @@ AliAnaPartCorrBaseClass::AliAnaPartCorrBaseClass(const AliAnaPartCorrBaseClass &
   //fAODCaloClusters(new TClonesArray(*abc.fAODCaloClusters)),
   fAODCaloCells(new AliAODCaloCells(*abc.fAODCaloCells)),
   fCaloPID(abc.fCaloPID), fFidCut(abc.fFidCut), fIC(abc.fIC),fMCUtils(abc.fMCUtils), fNMS(abc.fNMS),
+  //fAnaOutContainer(abc.fAnaOutContainer),
   fHistoPtBins(abc.fHistoPtBins),     fHistoPtMax(abc.fHistoPtMax),     fHistoPtMin(abc.fHistoPtMin),
   fHistoPhiBins(abc.fHistoPhiBins),   fHistoPhiMax(abc.fHistoPhiMax),   fHistoPhiMin(abc.fHistoPhiMin),
   fHistoEtaBins(abc.fHistoEtaBins),   fHistoEtaMax(abc.fHistoEtaMax),   fHistoEtaMin(abc.fHistoEtaMin),
@@ -117,7 +113,9 @@ AliAnaPartCorrBaseClass & AliAnaPartCorrBaseClass::operator = (const AliAnaPartC
   fIC      = abc.fIC;
   fMCUtils = abc.fMCUtils;
   fNMS     = abc.fNMS;
-  
+       
+  //fAnaOutContainer     = abc.fAnaOutContainer;
+       
   fInputAODBranch      = new TClonesArray(*abc.fInputAODBranch) ;
   fInputAODName        = abc.fInputAODName;
   fOutputAODBranch     = new TClonesArray(*abc.fOutputAODBranch) ;
@@ -160,6 +158,11 @@ AliAnaPartCorrBaseClass::~AliAnaPartCorrBaseClass()
     delete fAODCaloCells ;
   }
   
+//  if(fAnaOutContainer){
+//     fAnaOutContainer->Clear() ; 
+//     delete fAnaOutContainer ;
+//  }
+               
   if(fReader)  delete fReader ;
   if(fCaloPID) delete fCaloPID ;
   if(fFidCut)  delete fFidCut ;
@@ -471,16 +474,21 @@ void AliAnaPartCorrBaseClass::InitParameters()
   fMinPt = 0.1  ; //Min pt in particle analysis
   fMaxPt = 300. ; //Max pt in particle analysis
 
-  fCaloPID = new AliCaloPID ;  
-  fFidCut = new AliFiducialCut;
-  fIC = new AliIsolationCut;
-  fNMS = new AliNeutralMesonSelection;
-  fNewAOD = kFALSE ;
-  fOutputAODName = "PartCorr";
-  fOutputAODClassName = "AliAODPWG4Particle";
-  fInputAODName = "PartCorr";
+  fReader  = new AliCaloTrackReader();
+  fCaloPID = new AliCaloPID();
+  fFidCut  = new AliFiducialCut();
+  fIC      = new AliIsolationCut();
+  fMCUtils = new AliMCAnalysisUtils(); 
+  fNMS     = new AliNeutralMesonSelection;
+  
+  //fAnaOutContainer = new TList();
+       
+  fNewAOD              = kFALSE ;
+  fOutputAODName       = "PartCorr";
+  fOutputAODClassName  = "AliAODPWG4Particle";
+  fInputAODName        = "PartCorr";
   fAddToHistogramsName = "";
-  fAODObjArrayName="Ref";
+  fAODObjArrayName     = "Ref";
          
   //Histogram settings
   fHistoPtBins    = 240 ;
index ce88f8103af32c7d9ea3c7d45d7087efd2fb6909..8c4322104c1a0fb0bf3f6cfa44fb70563dfdb0c4 100755 (executable)
@@ -47,8 +47,9 @@ public:
   virtual void ConnectAODEMCALCells();
   virtual void ConnectInputOutputAODBranches();
   
-  virtual TList * GetCreateOutputObjects() { return (new TList) ;}
-  
+  virtual TList * GetCreateOutputObjects()      { return (new TList) ;}
+  //virtual TList * GetAnalysisOutputContainer()  { return fAnaOutContainer ;} 
+       
   virtual void AddToHistogramsName(TString add) { fAddToHistogramsName = add; }  
   virtual TString GetAddedHistogramsStringToName() {return fAddToHistogramsName ;}
   
@@ -232,6 +233,8 @@ public:
   AliMCAnalysisUtils       * fMCUtils; // MonteCarlo Analysis utils 
   AliNeutralMesonSelection * fNMS;     // Neutral Meson Selection
   
+  //TList * fAnaOutContainer;  // Temporal histogram output container, contents to be added to the main container passed to the main analysis frame
+
   //Histograms binning and range    
   Int_t   fHistoPtBins   ;  // Number of bins in pt axis
   Float_t fHistoPtMax    ;  // Maximum value of pt histogram range
index 010afe001aed71104e631bf4a96196d31d586d50..57913712b2e53e51f871b6754004e9edf368378d 100755 (executable)
@@ -237,8 +237,11 @@ void AliAnaPartCorrMaker::ProcessEvent(const Int_t iEntry, const char * currentF
   }
   //Each event needs an empty branch
   Int_t nAODBranches = fAODBranchList->GetEntries();
-  for(Int_t iaod = 0; iaod < nAODBranches; iaod++)
-               fAODBranchList->At(iaod)->Clear();
+  for(Int_t iaod = 0; iaod < nAODBranches; iaod++){
+         //fAODBranchList->At(iaod)->Clear();
+         TClonesArray *tca = dynamic_cast<TClonesArray*> (fAODBranchList->At(iaod));
+         if(tca) tca->Delete();
+  }
 
   //Tell the reader to fill the data in the 3 detector lists
   Bool_t ok = fReader->FillInputEvent(iEntry, currentFileName);
index 76afbdfdc02d92fc0e65dede61d19a971e35aad1..3ac473f035d95daa2a7555cb25ef83f90ea1e5ea 100755 (executable)
@@ -92,7 +92,7 @@ void AliAnalysisTaskParticleCorrelation::UserCreateOutputObjects()
   //Histograms container
   OpenFile(1);
   fOutputContainer = fAna->GetOutputContainer();
-  fOutputContainer->SetOwner();
+  fOutputContainer->SetOwner(kTRUE);
   if (DebugLevel() > 1) printf("AliAnalysisTaskParticleCorrelation::UserCreateOutputObjects() - End\n");
  
   PostData(1,fOutputContainer);
index d40994134c54bcfe319ece665a618548eba56b3f..97dfffa18f85f663f3644c48028e2e053f844304 100755 (executable)
@@ -1862,6 +1862,7 @@ void  AliAnaCalorimeterQA::MakeAnalysisFillHistograms()
                        printf("AliAnaCalorimeterQA::MakeAnalysisFillHistograms() - module %d calo %s clusters %d\n", imod, fCalorimeter.Data(), nClustersInModule[imod]); 
                fhNClustersMod[imod]->Fill(nClustersInModule[imod]);
        }
+       delete [] nClustersInModule;
        
        //----------------------------------------------------------
        // CALOCELLS
@@ -2038,7 +2039,8 @@ void  AliAnaCalorimeterQA::MakeAnalysisFillHistograms()
                        printf("AliAnaCalorimeterQA::MakeAnalysisFillHistograms() - module %d calo %s cells %d\n", imod, fCalorimeter.Data(), nCellsInModule[imod]); 
                fhNCellsMod[imod]->Fill(nCellsInModule[imod]) ;
        }
-
+       delete [] nCellsInModule;
+       
        if(GetDebug() > 0)
                printf("AliAnaCalorimeterQA::MakeAnalysisFillHistograms() - End \n");
 }
@@ -2521,11 +2523,11 @@ void AliAnaCalorimeterQA::CorrelateCalorimeters(TRefArray* refArray){
        if(GetReader()->GetDataType()==AliCaloTrackReader::kESD) {
                if(fCalorimeter == "EMCAL"){ 
                        ((AliESDEvent*)GetReader()->GetInputEvent())->GetPHOSClusters(caloClustersPHOS);
-                       caloClustersEMCAL = refArray;
+                       caloClustersEMCAL = new TRefArray(*refArray);
                }
                else if(fCalorimeter == "PHOS") { 
                        ((AliESDEvent*)GetReader()->GetInputEvent())->GetEMCALClusters (caloClustersEMCAL);
-                       caloClustersEMCAL = refArray;
+                       caloClustersPHOS = new TRefArray(*refArray);
                }
                
                //Fill histograms with clusters
@@ -2607,6 +2609,9 @@ void AliAnaCalorimeterQA::CorrelateCalorimeters(TRefArray* refArray){
                }
        }//AOD  
        
+       delete caloClustersEMCAL;
+       delete caloClustersPHOS;
+       
 }
 
 //______________________________________________________________________________
index 86aea75fe78429b37e5ffd53c2fe90b5033e5544..fa42e5f85b0a4ade7a730032ef0e548ba7f04164 100755 (executable)
@@ -92,7 +92,7 @@ AliAnaElectron::AliAnaElectron()
   //reco electrons from various sources\r
   fhPhiConversion(0),fhEtaConversion(0),\r
   //for comparisons with tracking detectors\r
-  fhPtTrack(0),
+  fhPtTrack(0),\r
   fhPtNPEBHadron(0),\r
   //for computing efficiency of B-jet tags\r
   fhBJetPt1x4(0),fhBJetPt2x3(0),fhBJetPt3x2(0),\r
@@ -662,8 +662,6 @@ void  AliAnaElectron::MakeAnalysisFillAOD()
   // Also fill some QA histograms\r
   //\r
 \r
-  TObjArray *cl = new TObjArray();\r
-\r
   Double_t bfield = 0.;\r
   if(GetReader()->GetDataType() != AliCaloTrackReader::kMC) bfield = GetReader()->GetBField();\r
 \r
@@ -672,7 +670,8 @@ void  AliAnaElectron::MakeAnalysisFillAOD()
     printf("This class not yet implemented for PHOS\n");\r
     abort();\r
   }\r
-  cl = GetAODEMCAL();\r
+  \r
+  TObjArray *cl = GetAODEMCAL();\r
   \r
   ////////////////////////////////////////////////\r
   //Start from tracks and get associated clusters \r
index 56f4f4b206001d51b36c8fd360dfad651102ef80..f25e5f24707939cd1994a224271dcc159875a7cd 100755 (executable)
@@ -206,7 +206,7 @@ void  AliAnaExample::MakeAnalysisFillAOD()
   }
   
   //Get List with tracks or clusters  
-  TObjArray * partList = new TObjArray();
+  TObjArray * partList = 0x0;
   if(fDetector == "CTS") partList = GetAODCTS();
   else if(fDetector == "EMCAL") partList = GetAODEMCAL();
   else if(fDetector == "PHOS") partList = GetAODPHOS();
index 771502b994586336aeb8cdcca72e19f34ef050e4..f72b5433447e773d3e2e604841375b790f92a2a3 100644 (file)
@@ -44,10 +44,9 @@ ClassImp(AliAnaOmegaToPi0Gamma)
 AliAnaOmegaToPi0Gamma::AliAnaOmegaToPi0Gamma() : AliAnaPartCorrBaseClass(),\r
 fInputAODGamma(0),fInputAODPi0(0), fInputAODGammaName(""),\r
 fEventsList(0), \r
-fVtxZCut(0), fCent(0), fRp(0), fBadChDist(0),\r
 fNVtxZBin(0), fNCentBin(0), fNRpBin(0), fNBadChDistBin(0), fNpid(0),\r
-fNmaxMixEv(0), fPi0Mass(0),\r
-fPi0MassWindow(0),fPi0OverOmegaPtCut(0),fGammaOverOmegaPtCut(0),\r
+fNmaxMixEv(0),fVtxZCut(0), fCent(0), fRp(0), fBadChDist(0),\r
+fPi0Mass(0), fPi0MassWindow(0),fPi0OverOmegaPtCut(0),fGammaOverOmegaPtCut(0),\r
 fhEtalon(0),\r
 fRealOmega(0), fMixAOmega(0),\r
 fMixBOmega(0), fMixCOmega(0),\r
@@ -67,11 +66,10 @@ fInputAODGamma(new TClonesArray (*ex.fInputAODGamma)),
 fInputAODPi0(new TClonesArray (*ex.fInputAODPi0)),\r
 fInputAODGammaName(ex.fInputAODGammaName),\r
 fEventsList(ex.fEventsList), \r
-fVtxZCut(ex.fVtxZCut), fCent(ex.fCent), fRp(ex.fRp), fBadChDist(ex.fBadChDist),\r
 fNVtxZBin(ex.fNVtxZBin), fNCentBin(ex.fNCentBin), fNRpBin(ex.fNRpBin),\r
-fNBadChDistBin(ex.fNBadChDistBin),fNpid(ex.fNpid),\r
-fNmaxMixEv(ex.fNmaxMixEv), fPi0Mass(ex.fPi0Mass),\r
-fPi0MassWindow(ex.fPi0MassWindow),\r
+fNBadChDistBin(ex.fNBadChDistBin),fNpid(ex.fNpid),fNmaxMixEv(ex.fNmaxMixEv), \r
+fVtxZCut(ex.fVtxZCut), fCent(ex.fCent), fRp(ex.fRp), fBadChDist(ex.fBadChDist),\r
+fPi0Mass(ex.fPi0Mass),fPi0MassWindow(ex.fPi0MassWindow),\r
 fPi0OverOmegaPtCut(ex.fPi0OverOmegaPtCut),fGammaOverOmegaPtCut(ex.fGammaOverOmegaPtCut),\r
 fhEtalon(ex.fhEtalon),\r
 fRealOmega(ex.fRealOmega), fMixAOmega(ex.fMixAOmega),\r
index 6c88adfdfe9f05e7d654b9a8a2fa747abaec2166..2194718a3b0f8c86686118e1bd7c57e3639a42c1 100644 (file)
@@ -54,19 +54,18 @@ class AliAnaOmegaToPi0Gamma : public AliAnaPartCorrBaseClass {
   TString fInputAODGammaName;    //Input AOD gamma name
   TList ** fEventsList;          //event list for mixing 
 
-  Double_t *fVtxZCut;            //[fNVtxZBin] vtertx z cut
-  Double_t *fCent;               //[fNCentBin] centrality cut
-  Double_t *fRp;                 //[fNRpBin] reaction plane cut
-  Int_t *fBadChDist;             //[fNBadChDistBin] bad channel dist
-
   Int_t fNVtxZBin;               //Number of vertex z cut
   Int_t fNCentBin;               //Number of centrality cut
   Int_t fNRpBin;                 //Number of reaction plane cut
   Int_t fNBadChDistBin;          //Number of bad channel dist cut
   Int_t fNpid;                   //Number of PID cut
-
   Int_t fNmaxMixEv;              //buffer size events to be mixed
+
+  Double_t *fVtxZCut;            //[fNVtxZBin] vtertx z cut
+  Double_t *fCent;               //[fNCentBin] centrality cut
+  Double_t *fRp;                 //[fNRpBin] reaction plane cut
+  Int_t *fBadChDist;             //[fNBadChDistBin] bad channel dist
+
   Double_t fPi0Mass;             //nominal pi0 mass
   Double_t fPi0MassWindow;       //pi0 mass windows
   Double_t fPi0OverOmegaPtCut;   //pi0 Pt over omega pt cut
index 11ad48c52ec9b2e40c06a93971a51c32a6b5d6d0..f736bc9605a9f1c7a6896499959930a138079a6e 100755 (executable)
@@ -695,7 +695,7 @@ void  AliAnaParticleIsolation::MakeAnalysisFillAOD()
   Int_t n = 0, nfrac = 0;
   Bool_t isolated = kFALSE ; 
   Float_t coneptsum = 0 ;
-  TObjArray * pl = new TObjArray()
+  TObjArray * pl = 0x0; 
   
   //Select the calorimeter for candidate isolation with neutral particles
   if(fCalorimeter == "PHOS")
index 5707597c8cc341fbc0fe54c69d5e13566a837151..53c72d3556cdbd32de716688ad898fd8a7d6e127 100755 (executable)
@@ -165,7 +165,8 @@ TList *  AliAnaPhoton::GetCreateOutputObjects()
   // store them in outputContainer
   TList * outputContainer = new TList() ; 
   outputContainer->SetName("PhotonHistos") ; 
-  
+  outputContainer->SetOwner(kFALSE);
+       
   Int_t nptbins  = GetHistoPtBins();
   Int_t nphibins = GetHistoPhiBins();
   Int_t netabins = GetHistoEtaBins();
@@ -434,7 +435,7 @@ void  AliAnaPhoton::MakeAnalysisFillAOD()
 {
   //Do analysis and fill aods
   //Search for photons in fCalorimeter 
-  TObjArray * pl = new TObjArray()
+  TObjArray * pl = 0x0
   
   //Get vertex for photon momentum calculation
   Double_t vertex[]  = {0,0,0} ; //vertex 
@@ -617,6 +618,8 @@ void  AliAnaPhoton::MakeAnalysisFillAOD()
     
   }//loop
   
+  delete [] indexConverted;
+       
   if(GetDebug() > 1) printf("AliAnaPhoton::MakeAnalysisFillAOD()  End fill AODs \n");  
   
 }
index 99b4ef84cd0d43784d877774960fc664422c8668..58de786335121c5a8a59b33e791b8de7bab3ebb3 100755 (executable)
@@ -446,7 +446,7 @@ void  AliAnaPi0EbE::MakeShowerShapeIdentification()
 {
   //Search for pi0 in fCalorimeter with shower shape analysis 
   
-  TObjArray * pl = new TObjArray
+  TObjArray * pl = 0x0
   
   //Get vertex for photon momentum calculation
   Double_t vertex[]  = {0,0,0} ; //vertex 
index cd9193e03931adf07ebebc70b2796d302c5e4239..c1e6d89938bf7736f949b0dd220eff08930d30a0 100644 (file)
@@ -38,7 +38,7 @@ AliAnalysisTaskParticleCorrelation *AddTaskPartCorr(TString data, TString calori
    AliCaloTrackReader * reader = 0x0;
    if(data.Contains("AOD")) reader = new AliCaloTrackAODReader();
    else if(data=="ESD") reader = new AliCaloTrackESDReader();
-   else if(data=="MC" && dataType == "ESD") reader = new AliCaloTrackMCReader();
+   else if(data=="MC" && inputDataType == "ESD") reader = new AliCaloTrackMCReader();
    reader->SetDebug(-1);//10 for lots of messages
    reader->SwitchOnCTS();
    //reader->SetDeltaAODFileName("");
@@ -157,7 +157,7 @@ AliAnalysisTaskParticleCorrelation *AddTaskPartCorr(TString data, TString calori
     anaphoton1->SwitchOnFiducialCut();
     anaphoton1->SetFiducialCut(fidCut1stYear);
   }
-       
+  anaphoton1->SwitchOnTrackMatchRejection();
   if(!data.Contains("delta")) anaphoton1->SetOutputAODName(Form("PhotonsForIM%s",calorimeter.Data()));
   else                        anaphoton1->SetInputAODName (Form("PhotonsForIM%s",calorimeter.Data()));
   //Set Histograms bins and ranges
@@ -463,9 +463,9 @@ AliAnalysisTaskParticleCorrelation *AddTaskPartCorr(TString data, TString calori
   anaomegaToPi0Gamma->SetNVtxZ(2);
   anaomegaToPi0Gamma->SetNBadChDist(3);
   anaomegaToPi0Gamma->SetNEventsMixed(4);
-  if(calorimeter.Data()=="PHOS")
+  if(calorimeter=="PHOS")
            anaomegaToPi0Gamma->SetPi0MassPeakWidthCut(0.008); // PHOS
-  else if(calorimeter.Data()=="EMCAL")
+  else if(calorimeter=="EMCAL")
            anaomegaToPi0Gamma->SetPi0MassPeakWidthCut(0.012); // EMCAL 
   anaomegaToPi0Gamma->SetHistoPtRangeAndNBins(0, 20, 200) ;
   anaomegaToPi0Gamma->SetHistoMassRangeAndNBins(0, 1, 200) ;
@@ -481,17 +481,17 @@ AliAnalysisTaskParticleCorrelation *AddTaskPartCorr(TString data, TString calori
   maker->SetReader(reader);//pointer to reader
   //if(!data.Contains("delta")) maker->AddAnalysis(qa,0);
   maker->AddAnalysis(anaphoton1,0);
-  maker->AddAnalysis(anapi0,1);
-  maker->AddAnalysis(anaphoton2,2);
-  maker->AddAnalysis(anaisol,3);
-  maker->AddAnalysis(anacorrjet,4);
-  maker->AddAnalysis(anacorrhadron,5);
-  maker->AddAnalysis(anacorrisohadron,6);
-  maker->AddAnalysis(anapi0ebe,7);
-  maker->AddAnalysis(anaisolpi0,8);
-  maker->AddAnalysis(anacorrhadronpi0,9);
-  maker->AddAnalysis(anacorrhadronisopi0,10);
-  maker->AddAnalysis(anaomegaToPi0Gamma,11);   
+//  maker->AddAnalysis(anapi0,1);
+//  maker->AddAnalysis(anaphoton2,2);
+//  maker->AddAnalysis(anaisol,3);
+//  maker->AddAnalysis(anacorrjet,4);
+//  maker->AddAnalysis(anacorrhadron,5);
+//  maker->AddAnalysis(anacorrisohadron,6);
+//  maker->AddAnalysis(anapi0ebe,7);
+//  maker->AddAnalysis(anaisolpi0,8);
+//  maker->AddAnalysis(anacorrhadronpi0,9);
+//  maker->AddAnalysis(anacorrhadronisopi0,10);
+//  maker->AddAnalysis(anaomegaToPi0Gamma,11);   
   maker->SetAnaDebug(-1)  ;
   maker->SwitchOnHistogramsMaker()  ;
   if(data.Contains("delta")) maker->SwitchOffAODsMaker()  ;