]> git.uio.no Git - u/mrichter/AliRoot.git/blob - JETAN/AliAnalysisTaskJets.cxx
pid of TRefArray equal to that of the contained objects.
[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 #include <Riostream.h> 
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 "AliJetHeader.h"
29 #include "AliJetReader.h"
30 #include "AliJetReaderHeader.h"
31 #include "AliJetHistos.h"
32 #include "AliESDEvent.h"
33 #include "AliESD.h"
34 #include "AliAODEvent.h"
35 #include "AliAODHandler.h"
36 #include "AliMCEventHandler.h"
37 #include "AliESDInputHandler.h"
38 #include "AliMCEvent.h"
39 #include "AliStack.h"
40
41
42 ClassImp(AliAnalysisTaskJets)
43
44 ////////////////////////////////////////////////////////////////////////
45
46 AliAnalysisTaskJets::AliAnalysisTaskJets():
47   AliAnalysisTaskSE(),
48   fConfigFile("ConfigJetAnalysis.C"),
49   fNonStdBranch(""), 
50   fJetFinder(0x0),
51   fHistos(0x0),
52   fListOfHistos(0x0),
53   fChain(0x0),
54   fOpt(0)
55 {
56   // Default constructor
57 }
58
59 AliAnalysisTaskJets::AliAnalysisTaskJets(const char* name):
60   AliAnalysisTaskSE(name),
61   fConfigFile("ConfigJetAnalysis.C"),
62   fNonStdBranch(""),
63   fJetFinder(0x0),
64   fHistos(0x0),
65   fListOfHistos(0x0),
66   fChain(0x0),
67   fOpt(0)
68 {
69   // Default constructor
70   DefineOutput(1, TList::Class());
71 }
72
73 AliAnalysisTaskJets::AliAnalysisTaskJets(const char* name, TChain* chain):
74   AliAnalysisTaskSE(name),
75   fConfigFile("ConfigJetAnalysis.C"),
76   fNonStdBranch(""),
77   fJetFinder(0x0),
78   fHistos(0x0),
79   fListOfHistos(0x0),
80   fChain(chain),
81   fOpt(0)
82 {
83   // Default constructor
84   DefineOutput(1, TList::Class());
85 }
86
87 //----------------------------------------------------------------
88 void AliAnalysisTaskJets::UserCreateOutputObjects()
89 {
90   // Create the output container
91   //
92   if (fDebug > 1) printf("AnalysisTaskJets::CreateOutPutData() \n");
93   
94   if(fNonStdBranch.Length()==0)
95     {
96       //  Connec default AOD to jet finder
97       fJetFinder->ConnectAOD(AODEvent());
98     }
99   else
100     {
101       // Create a new branch for jets...
102       // how is this reset? -> cleared in the UserExec....
103       // Can this be handled by the framework?
104       TClonesArray *tca = new TClonesArray("AliAODJet", 0);
105       tca->SetName(fNonStdBranch);
106       AddAODBranch("TClonesArray",&tca);
107       fJetFinder->ConnectAODNonStd(AODEvent(), fNonStdBranch.Data()); 
108     }
109   AliAODJetEventBackground* evBkg = fJetFinder->GetEventBackground();
110   if (evBkg) AddAODBranch("AliAODJetEventBackground",&evBkg);
111
112   // Histograms
113   OpenFile(1);
114   fListOfHistos = new TList();
115   fHistos       = new AliJetHistos();
116   fHistos->AddHistosToList(fListOfHistos);
117   
118   // Add the JetFinderInformation to the Outputlist
119   AliJetHeader *fH = fJetFinder->GetHeader();
120   
121   // Compose a characteristic output name
122   // with the name of the output branch
123   if(fH) {
124     if(fNonStdBranch.Length()==0) {
125       fH->SetName("AliJetHeader_jets");
126     }
127     else {
128       fH->SetName(Form("AliJetHeader_%s",fNonStdBranch.Data()));
129     }
130   }
131   OutputTree()->GetUserInfo()->Add(fH);
132 }
133
134 //----------------------------------------------------------------
135 void AliAnalysisTaskJets::Init()
136 {
137   // Initialization
138   if (fDebug > 1) printf("AnalysisTaskJets::Init() \n");
139   
140   // Call configuration file
141   if (fConfigFile.Length()) {
142     gROOT->LoadMacro(fConfigFile);
143     fJetFinder = (AliJetFinder*) gInterpreter->ProcessLine("ConfigJetAnalysis()");
144   }
145   AliJetReaderHeader *header = (AliJetReaderHeader*)fJetFinder->GetReader()->GetReaderHeader();
146   fOpt = header->GetDetector();
147
148   // Initialise Jet Analysis
149   if(fOpt == 0) fJetFinder->Init();
150   else fJetFinder->InitTask(fChain); // V2
151   
152   // Write header information to local file
153   fJetFinder->WriteHeaders();
154 }
155
156 //----------------------------------------------------------------
157 void AliAnalysisTaskJets::UserExec(Option_t */*option*/)
158 {
159   // Execute analysis for current event
160   //
161   // Fill control histos
162   TClonesArray* jarray = 0;
163
164   if(fNonStdBranch.Length()==0) {
165     jarray = AODEvent()->GetJets();
166   }
167   else {
168     jarray = dynamic_cast<TClonesArray*>(AODEvent()->FindListObject(fNonStdBranch.Data()));
169     jarray->Delete();    // this is our responsibility, clear before filling again
170   }
171   if (dynamic_cast<AliAODEvent*>(InputEvent()) !=  0) {
172       fJetFinder->GetReader()->SetInputEvent(InputEvent(), InputEvent(), MCEvent());
173   } else {
174       fJetFinder->GetReader()->SetInputEvent(InputEvent(), AODEvent(), MCEvent());
175   }
176   
177   
178
179   if(fOpt==0) fJetFinder->ProcessEvent();
180   else fJetFinder->ProcessEvent2();    // V2
181  
182   // Fill control histos
183   fHistos->FillHistos(jarray);
184
185   // Post the data
186   PostData(1, fListOfHistos);
187 }
188
189
190 //*************************************************************
191
192 void AliAnalysisTaskJets::Terminate(Option_t */*option*/)
193 {
194 // Terminate analysis
195 //
196     if (fDebug > 1) printf("AnalysisJets: Terminate() \n");
197 //    if (fJetFinder) fJetFinder->FinishRun();
198 }
199