]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Updates needed for the Analysis framework.
authormorsch <morsch@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 3 May 2007 14:21:05 +0000 (14:21 +0000)
committermorsch <morsch@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 3 May 2007 14:21:05 +0000 (14:21 +0000)
JETAN/AliAnalysisTaskJets.cxx
JETAN/AliAnalysisTaskJets.h
JETAN/AliJetESDReader.cxx
JETAN/AliJetFinder.cxx
JETAN/AliUA1JetFinderV1.cxx

index 9abfd8804ca380aa6ed471bed51512b5044d44a6..8b1b75d9df1364c73e45bb4f1433e843d21a27c5 100644 (file)
@@ -17,6 +17,7 @@
 #include <TSystem.h>
 #include <TInterpreter.h>
 #include <TChain.h>
+#include <TFile.h>
 #include <TH1.h>
 
 #include "AliAnalysisTaskJets.h"
@@ -53,6 +54,7 @@ AliAnalysisTaskJets::AliAnalysisTaskJets(const char* name):
 void AliAnalysisTaskJets::CreateOutputObjects()
 {
 // Create the output container
+    OpenFile(0);
     fTreeJ = fJetFinder->MakeTreeJ("TreeJ");
 }
 
@@ -66,6 +68,8 @@ void AliAnalysisTaskJets::Init()
     fJetFinder = (AliJetFinder*) gInterpreter->ProcessLine("ConfigJetAnalysis()");
     // Initialise Jet Analysis
     fJetFinder->Init();
+    // Write header information to local file
+    fJetFinder->WriteHeaders();
 }
 
 void AliAnalysisTaskJets::ConnectInputData(Option_t */*option*/)
@@ -73,18 +77,21 @@ void AliAnalysisTaskJets::ConnectInputData(Option_t */*option*/)
 // Connect the input data
 //
     if (fDebug > 1) printf("AnalysisTaskJets::ConnectInputData() \n");
-    char ** address = (char **)GetBranchAddress(0, "ESD");
+    fChain = (TChain*)GetInputData(0);
 
+    char ** address = (char **)GetBranchAddress(0, "ESD");
     if (address)     {
+       
+// Branch has been already connected
        fESD = (AliESD*)(*address);
     }
     else     {
+// First task taking the branch enables it
        fESD = new AliESD();
-       SetBranchAddress(0, "ESD", &fESD); // first task taking the branch enables it
+       SetBranchAddress(0, "ESD", &fESD);
     }
-    fChain = (TChain*)GetInputData(0);
+    
     fJetFinder->ConnectTree(fChain, fESD);
-    fJetFinder->WriteHeaders();
 }
 
 void AliAnalysisTaskJets::Exec(Option_t */*option*/)
@@ -92,7 +99,7 @@ void AliAnalysisTaskJets::Exec(Option_t */*option*/)
 // Execute analysis for current event
 //
     Long64_t ientry = fChain->GetReadEntry();
-    if (fDebug > 1) printf("Analysing event # %5d \n", (Int_t) ientry);
+    if (fDebug > 1) printf("Analysing event # %5d\n", (Int_t) ientry);
     fJetFinder->ProcessEvent(ientry);
     PostData(0, fTreeJ);
 }
@@ -102,7 +109,6 @@ void AliAnalysisTaskJets::Terminate(Option_t */*option*/)
 // Terminate analysis
 //
     if (fDebug > 1) printf("AnalysisJets: Terminate() \n");
-   
-    if (fJetFinder) fJetFinder->FinishRun();
+    // if (fJetFinder) fJetFinder->FinishRun();
 }
 
index 02dba0470bf8b97fea796c821b66d3bbda800b3f..2b48078c5c78ff4fc7c12be1020cb7c2a00f5919 100644 (file)
@@ -20,6 +20,7 @@ class AliAnalysisTaskJets : public AliAnalysisTask
     virtual void ConnectInputData(Option_t *option = "");
     virtual void CreateOutputObjects();
     virtual void Init();
+    virtual void LocalInit() {Init();}
     virtual void Exec(Option_t *option);
     virtual void Terminate(Option_t *option);
     virtual void SetDebugLevel(Int_t level) {fDebug = level;}
index 5ccba901a75b16555a0422d11c015549d243d85e..af7e1677344067ed6c3531fc468d98bed1def5e3 100755 (executable)
@@ -148,7 +148,7 @@ Bool_t AliJetESDReader::FillMomentumArray(Int_t event)
   ClearArray();
   fDebug = fReaderHeader->GetDebug();
   // get event from chain
-  fChain->GetTree()->GetEntry(event);
+  // fChain->GetTree()->GetEntry(event);
   
   if (!fESD) {
       return kFALSE;
index 97855f190ea22d7f6aa2adda2422e41ab379cbfc..67d76bec640512ea3ef851d1e3baae34c6c08a19 100644 (file)
@@ -103,8 +103,6 @@ void AliJetFinder::SetPlotMode(Bool_t b)
 TTree* AliJetFinder::MakeTreeJ(char* name)
 {
     // Create the tree for reconstructed jets
-    fOut = new TFile("jets.root","recreate");
-    fOut->cd();
     fTreeJ = new TTree(name, "AliJet");
     fTreeJ->Branch("FoundJet",   &fJets,   1000);
     fTreeJ->Branch("GenJet",     &fGenJets,1000);
@@ -117,9 +115,8 @@ TTree* AliJetFinder::MakeTreeJ(char* name)
 void AliJetFinder::WriteRHeaderToFile()
 {
   // write reader header
-  fOut->cd();
-  AliJetReaderHeader *rh = fReader->GetReaderHeader();
-  rh->Write();
+    AliJetReaderHeader *rh = fReader->GetReaderHeader();
+    rh->Write();
 }
 
 ////////////////////////////////////////////////////////////////////////
@@ -194,11 +191,10 @@ void AliJetFinder::ConnectTree(TTree* tree, TObject* data)
 void AliJetFinder::WriteHeaders()
 {
     // Write the Headers
-    if (fOut) {
-       fOut->cd();
-       WriteRHeaderToFile();
-       WriteJHeaderToFile();
-    }
+    TFile* f = new TFile("jets_local.root", "recreate");
+    WriteRHeaderToFile();
+    WriteJHeaderToFile();
+    f->Close();
 }
 
 
@@ -209,15 +205,17 @@ Bool_t AliJetFinder::ProcessEvent(Long64_t entry)
 //
     Int_t debug  = fReader->GetReaderHeader()->GetDebug();
     if (debug > 0) printf("<<<<< Processing Event %5d >>>>> \n", (Int_t) entry);
+    
     Bool_t ok = fReader->FillMomentumArray(entry);
     if (!ok) return kFALSE;
+
+    // Leading particles
     fLeading->FindLeading(fReader);
+    // Jets
     FindJets();
-    if (fOut) {
-       fOut->cd();
-       fTreeJ->Fill();
-    }
-    
+    // Fill the tree
+    fTreeJ->Fill();
+
     if (fPlots) fPlots->FillHistos(fJets);
     fLeading->Reset();
     fGenJets->ClearJets();
index 2fa3f27da37e0c86bd584d82037ec0d9f948fe37..e1e959346ebe7a4b7551b6efefdb3eb6f28772b9 100644 (file)
@@ -161,6 +161,8 @@ void AliUA1JetFinderV1::FindJets()
   // add jets to list
   Int_t* idxjets = new Int_t[nj];
   Int_t nselectj = 0;
+  printf("Found %d jets \n", nj);
+  
   for(Int_t kj=0; kj<nj; kj++){
      if ((etaJet[kj] > (header->GetJetEtaMax())) ||
           (etaJet[kj] < (header->GetJetEtaMin())) ||
@@ -802,7 +804,6 @@ void AliUA1JetFinderV1::Reset()
 void AliUA1JetFinderV1::WriteJHeaderToFile()
 {
   AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
-  fOut->cd();
   header->Write();
 }