]> git.uio.no Git - u/mrichter/AliRoot.git/blob - JETAN/AliAnalysisTaskJets.cxx
Updates needed for the Analysis framework.
[u/mrichter/AliRoot.git] / JETAN / AliAnalysisTaskJets.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15  
16 #include <TROOT.h>
17 #include <TSystem.h>
18 #include <TInterpreter.h>
19 #include <TChain.h>
20 #include <TFile.h>
21 #include <TH1.h>
22
23 #include "AliAnalysisTaskJets.h"
24 #include "AliJetFinder.h"
25 #include "AliESD.h"
26
27 ClassImp(AliAnalysisTaskJets)
28
29 ////////////////////////////////////////////////////////////////////////
30
31 AliAnalysisTaskJets::AliAnalysisTaskJets():
32     fDebug(0),
33     fJetFinder(0x0),
34     fChain(0x0),
35     fESD(0x0),
36     fTreeJ(0x0)
37 {
38   // Default constructor
39 }
40
41 AliAnalysisTaskJets::AliAnalysisTaskJets(const char* name):
42     AliAnalysisTask(name, "AnalysisTaskJets"),
43     fDebug(0),
44     fJetFinder(0x0),
45     fChain(0x0),
46     fESD(0x0),
47     fTreeJ(0x0)
48 {
49   // Default constructor
50     DefineInput (0, TChain::Class());
51     DefineOutput(0, TTree::Class());
52 }
53
54 void AliAnalysisTaskJets::CreateOutputObjects()
55 {
56 // Create the output container
57     OpenFile(0);
58     fTreeJ = fJetFinder->MakeTreeJ("TreeJ");
59 }
60
61 void AliAnalysisTaskJets::Init()
62 {
63     // Initialization
64     if (fDebug > 1) printf("AnalysisTaskJets::Init() \n");
65
66     // Call configuration file
67     gROOT->LoadMacro("ConfigJetAnalysis.C");
68     fJetFinder = (AliJetFinder*) gInterpreter->ProcessLine("ConfigJetAnalysis()");
69     // Initialise Jet Analysis
70     fJetFinder->Init();
71     // Write header information to local file
72     fJetFinder->WriteHeaders();
73 }
74
75 void AliAnalysisTaskJets::ConnectInputData(Option_t */*option*/)
76 {
77 // Connect the input data
78 //
79     if (fDebug > 1) printf("AnalysisTaskJets::ConnectInputData() \n");
80     fChain = (TChain*)GetInputData(0);
81
82     char ** address = (char **)GetBranchAddress(0, "ESD");
83     if (address)     {
84         
85 // Branch has been already connected
86         fESD = (AliESD*)(*address);
87     }
88     else     {
89 // First task taking the branch enables it
90         fESD = new AliESD();
91         SetBranchAddress(0, "ESD", &fESD);
92     }
93     
94     fJetFinder->ConnectTree(fChain, fESD);
95 }
96
97 void AliAnalysisTaskJets::Exec(Option_t */*option*/)
98 {
99 // Execute analysis for current event
100 //
101     Long64_t ientry = fChain->GetReadEntry();
102     if (fDebug > 1) printf("Analysing event # %5d\n", (Int_t) ientry);
103     fJetFinder->ProcessEvent(ientry);
104     PostData(0, fTreeJ);
105 }
106
107 void AliAnalysisTaskJets::Terminate(Option_t */*option*/)
108 {
109 // Terminate analysis
110 //
111     if (fDebug > 1) printf("AnalysisJets: Terminate() \n");
112     // if (fJetFinder) fJetFinder->FinishRun();
113 }
114