]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Added possibility to not create aod tclones branch
authorslindal <slindal@f7af4fe6-9843-0410-8265-dc069ae4e863>
Sun, 19 Dec 2010 15:37:46 +0000 (15:37 +0000)
committerslindal <slindal@f7af4fe6-9843-0410-8265-dc069ae4e863>
Sun, 19 Dec 2010 15:37:46 +0000 (15:37 +0000)
Removed delete[] from conversion particle destructor

PWG4/GammaConv/AliAODConversionParticle.cxx
PWG4/GammaConv/AliAnalysisTaskGammaConversion.cxx
PWG4/GammaConv/AliAnalysisTaskGammaConversion.h
PWG4/macros/ConfigGammaConversion.C

index 75f20f09d88dda855b9f623c7a5249f9f30f8d9d..4746b6a67eafbfa37c96b142eb392867cf361941 100644 (file)
@@ -64,7 +64,6 @@ AliAODConversionParticle::AliAODConversionParticle(TLorentzVector & momentum) :
 AliAODConversionParticle::~AliAODConversionParticle() {
   
   ///destructor
-  delete[] fLabel;
 
 }
 
index 3f01c05f14ee9b3eb56a273082e3e75376f397e8..58a822be570f3f06c6f8acbaec4d1e33f8c91cb0 100644 (file)
@@ -129,6 +129,7 @@ AliAnalysisTaskSE(),
   fAODOmega(NULL),
   fAODBranchName("GammaConv"),
   fOutputAODClassName("AliAODConversionParticle"),
+  fKFCreateAOD(kTRUE),
   fKFForceAOD(kFALSE),
   fKFDeltaAODFileName(""),
   fDoNeutralMesonV0MCCheck(kFALSE),
@@ -226,6 +227,7 @@ AliAnalysisTaskGammaConversion::AliAnalysisTaskGammaConversion(const char* name)
   fAODOmega(NULL),
   fAODBranchName("GammaConv"),
   fOutputAODClassName("AliAODConversionParticle"),
+  fKFCreateAOD(kTRUE),
   fKFForceAOD(kFALSE),
   fKFDeltaAODFileName(""),
   fDoNeutralMesonV0MCCheck(kFALSE),
@@ -620,8 +622,11 @@ void AliAnalysisTaskGammaConversion::UserExec(Option_t */*option*/)
   
 
   //Fill Gamma AOD
-  FillAODWithConversionGammas() ; 
-       
+  if(fKFCreateAOD) {
+    FillAODWithConversionGammas() ; 
+  }
+
+
   // Process reconstructed gammas
   if(fDoNeutralMeson == kTRUE){
     ProcessGammasForNeutralMesonAnalysis();
@@ -2390,7 +2395,7 @@ void AliAnalysisTaskGammaConversion::ProcessGammasForNeutralMesonAnalysis(){
 
            Double_t lowMassPi0=0.1;
            Double_t highMassPi0=0.15;
-           if (massTwoGammaCandidate > lowMassPi0 && massTwoGammaCandidate < highMassPi0 ){
+           if ( fKFCreateAOD && (massTwoGammaCandidate > lowMassPi0) && (massTwoGammaCandidate < highMassPi0) ){
              new((*fKFReconstructedPi0sTClone)[fKFReconstructedPi0sTClone->GetEntriesFast()])  AliKFParticle(*twoGammaCandidate);
              fGammav1.push_back(firstGammaIndex);
              fGammav2.push_back(secondGammaIndex);
@@ -3451,29 +3456,33 @@ void AliAnalysisTaskGammaConversion::Terminate(Option_t */*option*/)
 
 void AliAnalysisTaskGammaConversion::UserCreateOutputObjects()
 {
-  //AOD
-  if(!fAODGamma) fAODGamma = new TClonesArray(fOutputAODClassName.Data(), 0);
-  else fAODGamma->Delete();
-  fAODGamma->SetName(Form("%s_gamma", fAODBranchName.Data()));
   
-  if(!fAODPi0) fAODPi0 = new TClonesArray(fOutputAODClassName.Data(), 0);
-  else fAODPi0->Delete();
-  fAODPi0->SetName(Form("%s_Pi0", fAODBranchName.Data()));
-
-  if(!fAODOmega) fAODOmega = new TClonesArray(fOutputAODClassName.Data(), 0);
-  else fAODOmega->Delete();
-  fAODOmega->SetName(Form("%s_Omega", fAODBranchName.Data()));
-
-  //If delta AOD file name set, add in separate file. Else add in standard aod file. 
-  if(fKFDeltaAODFileName.Length() > 0) {
-    AddAODBranch("TClonesArray", &fAODGamma, fKFDeltaAODFileName.Data());
-    AddAODBranch("TClonesArray", &fAODPi0, fKFDeltaAODFileName.Data());
-    AddAODBranch("TClonesArray", &fAODOmega, fKFDeltaAODFileName.Data());
-    AliAnalysisManager::GetAnalysisManager()->RegisterExtraFile(fKFDeltaAODFileName.Data());
-  } else  {
-    AddAODBranch("TClonesArray", &fAODGamma);
-    AddAODBranch("TClonesArray", &fAODPi0);
-    AddAODBranch("TClonesArray", &fAODOmega);
+  if(fKFCreateAOD) {
+
+    //AOD
+    if(!fAODGamma) fAODGamma = new TClonesArray(fOutputAODClassName.Data(), 0);
+    else fAODGamma->Delete();
+    fAODGamma->SetName(Form("%s_gamma", fAODBranchName.Data()));
+    
+    if(!fAODPi0) fAODPi0 = new TClonesArray(fOutputAODClassName.Data(), 0);
+    else fAODPi0->Delete();
+    fAODPi0->SetName(Form("%s_Pi0", fAODBranchName.Data()));
+    
+    if(!fAODOmega) fAODOmega = new TClonesArray(fOutputAODClassName.Data(), 0);
+    else fAODOmega->Delete();
+    fAODOmega->SetName(Form("%s_Omega", fAODBranchName.Data()));
+    
+    //If delta AOD file name set, add in separate file. Else add in standard aod file. 
+    if(GetDeltaAODFileName().Length() > 0) {
+      AddAODBranch("TClonesArray", &fAODGamma, GetDeltaAODFileName().Data());
+      AddAODBranch("TClonesArray", &fAODPi0, GetDeltaAODFileName().Data());
+      AddAODBranch("TClonesArray", &fAODOmega, GetDeltaAODFileName().Data());
+      AliAnalysisManager::GetAnalysisManager()->RegisterExtraFile(GetDeltaAODFileName().Data());
+    } else  {
+      AddAODBranch("TClonesArray", &fAODGamma);
+      AddAODBranch("TClonesArray", &fAODPi0);
+      AddAODBranch("TClonesArray", &fAODOmega);
+    }
   }
 
   // Create the output container
index a4d77ef16c183190b0937f7bf5c5ac05660dad68..bb5b32ba876ccc44a76e449e6cba7c880d79e698 100644 (file)
@@ -131,7 +131,8 @@ class AliAnalysisTaskGammaConversion : public AliAnalysisTaskSE
   Double_t GetMCOpeningAngle(const TParticle* const daughter0,const TParticle* const daughter1) const;
   void CheckV0Efficiency();
   void SetDeltaAODFileName(TString fn) { fKFDeltaAODFileName = fn; };
-               
+  void SetCreateAOD(Bool_t doAod) { fKFCreateAOD = doAod; };
+  TString GetDeltaAODFileName() const { return fKFDeltaAODFileName; };
   //////////////////Chi_c Analysis////////////////////////////
   void GetPID(const AliESDtrack *track, Stat_t &pid, Stat_t &weight);  
   double GetSigmaToVertex(const AliESDtrack* t);
@@ -300,6 +301,7 @@ class AliAnalysisTaskGammaConversion : public AliAnalysisTaskSE
   TClonesArray * fAODOmega; //TTClonesArray for omegas to put in AOD
   TString fAODBranchName; // New AOD branch name
   TString fOutputAODClassName; //Class to use for the AOD
+  Bool_t fKFCreateAOD; //Create the AOD tclones? (regardless if storing or not)
   
   Bool_t fKFForceAOD;  //Set the Analysis Manager FillAOD variable to true every event
   TString fKFDeltaAODFileName; //! File name for delta AOD (if any)
@@ -320,7 +322,7 @@ class AliAnalysisTaskGammaConversion : public AliAnalysisTaskSE
   Int_t fUseMultiplicityBin;
   Int_t fUseCentrality;
   Int_t fUseCentralityBin;
-  ClassDef(AliAnalysisTaskGammaConversion, 16); // Analysis task for gamma conversions
+  ClassDef(AliAnalysisTaskGammaConversion, 17); // Analysis task for gamma conversions
 };
 
 #endif //ALIANALYSISTASKGAMMA_H
index f7cf37f22c759b7907301f59978c7acafeb758e8..e447886530ff86e3a7318ed2d7b7c6e799a843c3 100644 (file)
@@ -1586,13 +1586,13 @@ AliAnalysisTaskGammaConversion* ConfigGammaConversion(TString arguments, AliAnal
     //gammaconversion->SetOutputAODClassName("AliGammaConversionAODObject");
 
     if( kGCrunOnTrain ) {
-
+      
       AliAODHandler * aodHandler = dynamic_cast<AliAODHandler*>(mgr->GetOutputEventHandler());
       if(!aodHandler) {
        ::Error("This task requires an AOD handler");
-       return -1;
+       return NULL;
       }
-
+      
       gammaconversion->SetDeltaAODFileName(kGCDeltaAODFilename);
       
       if(kGCDeltaAODFilename.Length() > 0) {
@@ -1602,23 +1602,25 @@ AliAnalysisTaskGammaConversion* ConfigGammaConversion(TString arguments, AliAnal
     } else {
       if(kGCDeltaAODFilename.Length() == 0 ) {
        cout << "Error:: Need a file name for the AOD"<<endl;
-       return;
+       return NULL;
       }
       AliAODHandler* aodHandler = new AliAODHandler();
       aodHandler->SetOutputFileName(kGCDeltaAODFilename);
       aodHandler->SetCreateNonStandardAOD();
       mgr->SetOutputEventHandler(aodHandler);  
     }
+  }  else {
+    gammaconversion->SetCreateAOD(kFALSE);
   }
 
   // Connect I/O to the task
   mgr->ConnectInput (gammaconversion, 0, cinput1);
   if(mgr->GetCommonOutputContainer())
     mgr->ConnectOutput(gammaconversion, 0, mgr->GetCommonOutputContainer());
-
+  
   mgr->ConnectOutput(gammaconversion, 1, coutput2);
   mgr->ConnectOutput(gammaconversion, 2, coutput3);
-
+  
   if(kGCRunGammaJetTask) {
     AliAnalysisTaskGammaJet * gammaJetTask = new AliAnalysisTaskGammaJet("GammaJetTask");
     if(kGCrunOnTrain) {