]> git.uio.no Git - u/mrichter/AliRoot.git/blob - JETAN/AliAnalysisTaskJets.cxx
Moved MakePrimary to EMCRecPoint to rpduce correct order of primaries
[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 /* $Id$ */
17  
18 #include <TROOT.h>
19 #include <TSystem.h>
20 #include <TInterpreter.h>
21 #include <TChain.h>
22 #include <TFile.h>
23 #include <TList.h>
24
25 #include "AliAnalysisTaskJets.h"
26 #include "AliAnalysisManager.h"
27 #include "AliJetFinder.h"
28 #include "AliJetHistos.h"
29 #include "AliESDEvent.h"
30 #include "AliESD.h"
31 #include "AliAODEvent.h"
32 #include "AliAODHandler.h"
33 #include "AliMCEventHandler.h"
34 #include "AliMCEvent.h"
35 #include "AliStack.h"
36
37
38 ClassImp(AliAnalysisTaskJets)
39
40 ////////////////////////////////////////////////////////////////////////
41
42 AliAnalysisTaskJets::AliAnalysisTaskJets():
43     fDebug(0),
44     fJetFinder(0x0),
45     fChain(0x0),
46     fESD(0x0),
47     fAOD(0x0),
48     fTreeA(0x0),
49     fHistos(0x0),
50     fListOfHistos(0x0)
51 {
52   // Default constructor
53 }
54
55 AliAnalysisTaskJets::AliAnalysisTaskJets(const char* name):
56     AliAnalysisTask(name, "AnalysisTaskJets"),
57     fDebug(0),
58     fJetFinder(0x0),
59     fChain(0x0),
60     fESD(0x0),
61     fAOD(0x0),
62     fTreeA(0x0),
63     fHistos(0x0),
64     fListOfHistos(0x0)
65 {
66   // Default constructor
67     DefineInput (0, TChain::Class());
68     DefineOutput(0, TTree::Class());
69     DefineOutput(1, TList::Class());
70 }
71
72 void AliAnalysisTaskJets::CreateOutputObjects()
73 {
74 // Create the output container
75 //
76 //  Default AOD
77     if (fDebug > 1) printf("AnalysisTaskJets::CreateOutPutData() \n");
78     AliAODHandler* handler = (AliAODHandler*) ((AliAnalysisManager::GetAnalysisManager())->GetOutputEventHandler());
79     
80     fAOD   = handler->GetAOD();
81     fTreeA = handler->GetTree();
82     fJetFinder->ConnectAOD(fAOD);
83 //
84 //  Histograms
85     OpenFile(1);
86     fListOfHistos = new TList();
87     fHistos       = new AliJetHistos();
88     fHistos->AddHistosToList(fListOfHistos);
89 }
90
91 void AliAnalysisTaskJets::Init()
92 {
93     // Initialization
94     if (fDebug > 1) printf("AnalysisTaskJets::Init() \n");
95
96     // Call configuration file
97     gROOT->LoadMacro("ConfigJetAnalysis.C");
98     fJetFinder = (AliJetFinder*) gInterpreter->ProcessLine("ConfigJetAnalysis()");
99     // Initialise Jet Analysis
100     fJetFinder->Init();
101     // Write header information to local file
102     fJetFinder->WriteHeaders();
103 }
104
105 void AliAnalysisTaskJets::ConnectInputData(Option_t */*option*/)
106 {
107 // Connect the input data
108     if (fDebug > 1) printf("AnalysisTaskJets::ConnectInputData() \n");
109     fChain = (TChain*)GetInputData(0);
110     fESD = new AliESDEvent();
111     fESD->ReadFromTree(fChain);
112
113     AliMCEventHandler*    mcTruth = (AliMCEventHandler*) 
114         ((AliAnalysisManager::GetAnalysisManager())->GetMCtruthEventHandler());
115
116     fJetFinder->GetReader()->SetInputEvent(fESD, fAOD, mcTruth);
117 }
118
119 void AliAnalysisTaskJets::Exec(Option_t */*option*/)
120 {
121 // Execute analysis for current event
122 //
123     AliMCEventHandler*    mctruth = (AliMCEventHandler*) 
124         ((AliAnalysisManager::GetAnalysisManager())->GetMCtruthEventHandler());
125     if (mctruth) {
126         AliStack* stack = mctruth->MCEvent()->Stack();
127         printf("AliAnalysisTaskJets: Number of tracks on stack %5d\n", stack->GetNtrack());
128     }
129     
130     AliESD* old = fESD->GetAliESDOld();
131     if (old) {
132         fESD->CopyFromOldESD();
133         old->Reset();
134     }
135     
136     Long64_t ientry = fChain->GetReadEntry();
137     if (fDebug > 1) printf("Analysing event # %5d\n", (Int_t) ientry);
138     fJetFinder->ProcessEvent(ientry);
139
140     // Fill control histos
141     fHistos->FillHistos(fAOD->GetJets());
142     
143     // Post the data
144     PostData(0, fTreeA);
145     PostData(1, fListOfHistos);
146 }
147
148 void AliAnalysisTaskJets::Terminate(Option_t */*option*/)
149 {
150 // Terminate analysis
151 //
152     if (fDebug > 1) printf("AnalysisJets: Terminate() \n");
153     // if (fJetFinder) fJetFinder->FinishRun();
154 }
155