]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Use pointers for output TTrees and use AliESDInputHandler (Mihaela)
authordainese <dainese@f7af4fe6-9843-0410-8265-dc069ae4e863>
Tue, 1 Apr 2008 12:14:11 +0000 (12:14 +0000)
committerdainese <dainese@f7af4fe6-9843-0410-8265-dc069ae4e863>
Tue, 1 Apr 2008 12:14:11 +0000 (12:14 +0000)
PWG3/AliAnalysisTaskVertexingHF.cxx
PWG3/AliAnalysisTaskVertexingHF.h
PWG3/AliAnalysisTaskVertexingHFTest.C

index 1160d7c8fe667fd92815aaeb39c84602feade54f..671f33999a0e02fcf0a39b512a95ab757c594462 100644 (file)
@@ -29,6 +29,7 @@
 
 #include "AliAnalysisTask.h"
 #include "AliESDEvent.h"
+#include "AliESDInputHandler.h"
 #include "AliAnalysisVertexingHF.h"
 #include "AliAnalysisTaskVertexingHF.h"
 
@@ -54,16 +55,36 @@ fTrees(0)
   DefineOutput(3, TTree::Class());
 }
 //________________________________________________________________________
+AliAnalysisTaskVertexingHF::~AliAnalysisTaskVertexingHF()
+{
+  // Destructor
+  if (fTrees) delete [] fTrees;
+}  
+//________________________________________________________________________
 void AliAnalysisTaskVertexingHF::ConnectInputData(Option_t *)
 {
   //Implementation of AliAnalysisTask::ConnectInputData
 
   Info("ConnectInputData","ConnectInputData of task %s\n",GetName());
   fChain = (TChain*)GetInputData(0);
-  fESD = new AliESDEvent();
-  fESD->ReadFromTree(fChain);
 
+  if (!fChain) {
+    printf("ERROR: Could not read chain from input slot 0\n");
+    return;
+  } 
+
+  //fESD = new AliESDEvent();
+  //fESD->ReadFromTree(fChain); 
+  
+  AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
+  
+  if (!esdH) {
+    printf("ERROR: Could not get ESDInputHandler\n");
+  } else {
+    fESD = esdH->GetEvent();
+    if (!fESD) printf("ERROR: Could not get ESD object\n"); 
+ }
   return;
 }
 //________________________________________________________________________
@@ -72,22 +93,18 @@ void AliAnalysisTaskVertexingHF::CreateOutputObjects()
   //Implementation of AliAnalysisTask::CreateOutputObjects
   //4 output trees (D0 in 2-prongs, J/Psi to e+e-, 3-prongs (D+, Ds, /\c), D0 in 4-prongs)
 
-  fTrees = new TTree[4];
+  fTrees = new TTree*[4];
   AliAODRecoDecayHF2Prong *rd2=0;
   AliAODRecoDecayHF3Prong *rd3=0;
   AliAODRecoDecayHF4Prong *rd4=0;
-  fTrees[0].SetName("NameD0toKpi");
-  fTrees[0].SetTitle("TitleD0toKpi");
-  fTrees[1].SetName("NameJPSItoEle");
-  fTrees[1].SetTitle("TitleJPSItoEle");
-  fTrees[2].SetName("NameCharmto3Prong");
-  fTrees[2].SetTitle("TitleCharmto3Prong");
-  fTrees[3].SetName("NameD0to4Prong");
-  fTrees[3].SetTitle("TitleD0to4Prong");
-  fTrees[0].Branch("D0toKpi","AliAODRecoDecayHF2Prong",&rd2);
-  fTrees[1].Branch("JPSItoEle","AliAODRecoDecayHF2Prong",&rd2);
-  fTrees[2].Branch("Charmto3Prong","AliAODRecoDecayHF3Prong",&rd3);
-  fTrees[3].Branch("D0to4Prong","AliAODRecoDecayHF4Prong",&rd4);
+  fTrees[0] = new TTree("NameD0toKpi","TitleD0toKpi");
+  fTrees[0]->Branch("D0toKpi","AliAODRecoDecayHF2Prong",&rd2);
+  fTrees[1] = new TTree("NameJPSItoEle", "TitleJPSItoEle");
+  fTrees[1]->Branch("JPSItoEle","AliAODRecoDecayHF2Prong",&rd2);
+  fTrees[2] = new TTree("NameCharmto3Prong", "TitleCharmto3Prong");
+  fTrees[2]->Branch("Charmto3Prong","AliAODRecoDecayHF3Prong",&rd3);
+  fTrees[3] = new TTree("NameD0to4Prong","TitleD0to4Prong");
+  fTrees[3]->Branch("D0to4Prong","AliAODRecoDecayHF4Prong",&rd4);
 
   return;
 }
@@ -107,19 +124,15 @@ void AliAnalysisTaskVertexingHF::LocalInit()
 void AliAnalysisTaskVertexingHF::Exec(Option_t *)
 {
   //Performs heavy flavor vertexing
-
-  //Transition to new ESD format (from AliRoot v4-06 ?) :
-  fESD->GetAliESDOld();
-  if (fESD->GetAliESDOld()) fESD->CopyFromOldESD();
   
   //Heavy flavor vertexing :
   fVHF->FindCandidates(fESD,fTrees);
   
   //Post final data. It will be written to a file with option "RECREATE"
-  PostData(0, &fTrees[0]);
-  PostData(1, &fTrees[1]);
-  PostData(2, &fTrees[2]);
-  PostData(3, &fTrees[3]);
+  PostData(0, fTrees[0]);
+  PostData(1, fTrees[1]);
+  PostData(2, fTrees[2]);
+  PostData(3, fTrees[3]);
 
   return;
 }      
index 48e926d74ba3854e6edf5eda99298edf6fb12738..552896bb5da6c3f94a2e95da05c83332dfbc5749 100644 (file)
@@ -26,6 +26,7 @@ class AliAnalysisTaskVertexingHF : public AliAnalysisTask
 
   AliAnalysisTaskVertexingHF() : AliAnalysisTask(), fESD(0), fChain(0), fVHF(0), fTrees(0) {}
   AliAnalysisTaskVertexingHF(const char *name);
+  virtual ~AliAnalysisTaskVertexingHF();
 
   AliAnalysisTaskVertexingHF(const AliAnalysisTaskVertexingHF &source);
   AliAnalysisTaskVertexingHF& operator=(const AliAnalysisTaskVertexingHF& source); 
@@ -38,12 +39,12 @@ class AliAnalysisTaskVertexingHF : public AliAnalysisTask
   
  private:
 
-  AliESDEvent            *fESD;   //ESD
-  TChain                 *fChain; //Chain
+  AliESDEvent            *fESD;   //!ESD
+  TChain                 *fChain; //!Chain
   AliAnalysisVertexingHF *fVHF;    //Vertexer heavy flavour
-  TTree                  *fTrees; //Output trees (D0 in 2-prongs, J/Psi to e+e-, 3-prongs (D+, Ds, Lc), D0 in 4-prongs)
+  TTree                  **fTrees; //Output trees (D0 in 2-prongs, J/Psi to e+e-, 3-prongs (D+, Ds, Lc), D0 in 4-prongs)
   
-  ClassDef(AliAnalysisTaskVertexingHF,1); //AliAnalysisTask for the reconstruction of heavy-flavour decay candidates
+  ClassDef(AliAnalysisTaskVertexingHF,2); //AliAnalysisTask for the reconstruction of heavy-flavour decay candidates
 };
 
 #endif
index 8854cccc8bad08701e5cda91120204112c3e6201..a7053196b73f35aeb456ededcdae56c174fbd29c 100644 (file)
@@ -38,13 +38,14 @@ void AliAnalysisTaskVertexingHFTest() {
   TChain* chain = tagAna->QueryTags(runCuts,lhcCuts,detCuts,eventCuts);
   */
  
-  //Temporary solution to avoid memory leaks : 
-  chain->SetBranchStatus("*FMD*",0);
-  chain->SetBranchStatus("*CaloClusters*",0);
-
   //Create tasks
   AliAnalysisManager *analManager = new AliAnalysisManager("myAnalysisManager");
   //analManager->SetDebugLevel(10);
+  // ESD input handler
+  AliESDInputHandler *esdHandler = new AliESDInputHandler();
+  esdHandler->SetInactiveBranches("FMD CaloCluster");
+  analManager->SetInputEventHandler(esdHandler);
+
   AliAnalysisTaskVertexingHF *task1 = new AliAnalysisTaskVertexingHF("myTask");
   analManager->AddTask(task1);
   //Create containers for input/output
@@ -58,7 +59,6 @@ void AliAnalysisTaskVertexingHFTest() {
   analManager->ConnectOutput(task1,1,coutput2);
   analManager->ConnectOutput(task1,2,coutput3);
   analManager->ConnectOutput(task1,3,coutput4);
-  cinput1->SetData(chain);
  
   printf("CHAIN HAS %d ENTRIES\n",(Int_t)chain->GetEntries());
   if (analManager->InitAnalysis()) {