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