]> git.uio.no Git - u/mrichter/AliRoot.git/blob - JETAN/AliAnalysisTaskJets.cxx
limit printouts
[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 "AliAODJetEventBackground.h"
36 #include "AliAODHandler.h"
37 #include "AliMCEventHandler.h"
38 #include "AliESDInputHandler.h"
39 #include "AliMCEvent.h"
40 #include "AliStack.h"
41
42
43 ClassImp(AliAnalysisTaskJets)
44
45 ////////////////////////////////////////////////////////////////////////
46
47 AliAnalysisTaskJets::AliAnalysisTaskJets():
48   AliAnalysisTaskSE(),
49   fConfigFile("ConfigJetAnalysis.C"),
50   fNonStdBranch(""), 
51   fNonStdFile(""),
52   fJetFinder(0x0),
53   fHistos(0x0),
54   fAODExtension(0x0),
55   fListOfHistos(0x0),
56   fChain(0x0),
57   fOpt(0),
58   fReadAODFromOutput(0)
59 {
60   // Default constructor
61 }
62
63 AliAnalysisTaskJets::AliAnalysisTaskJets(const char* name):
64   AliAnalysisTaskSE(name),
65   fConfigFile("ConfigJetAnalysis.C"),
66   fNonStdBranch(""),
67   fNonStdFile(""),
68   fJetFinder(0x0),
69   fHistos(0x0),
70   fAODExtension(0x0),
71   fListOfHistos(0x0),
72   fChain(0x0),
73   fOpt(0),
74   fReadAODFromOutput(0)
75 {
76   // Default constructor
77   DefineOutput(1, TList::Class());
78 }
79
80 AliAnalysisTaskJets::AliAnalysisTaskJets(const char* name, TChain* chain):
81   AliAnalysisTaskSE(name),
82   fConfigFile("ConfigJetAnalysis.C"),
83   fNonStdBranch(""),
84   fNonStdFile(""),
85   fJetFinder(0x0),
86   fHistos(0x0),
87   fAODExtension(0x0),
88   fListOfHistos(0x0),
89   fChain(chain),
90   fOpt(0),
91   fReadAODFromOutput(0)
92 {
93   // Default constructor
94   DefineOutput(1, TList::Class());
95 }
96
97 //----------------------------------------------------------------
98 void AliAnalysisTaskJets::UserCreateOutputObjects()
99 {
100   // Create the output container
101   //
102   if (fDebug > 1) printf("AnalysisTaskJets::CreateOutPutData() \n");
103
104   if(fNonStdBranch.Length()==0)
105     {
106       //  Connec default AOD to jet finder
107       // create a new branch for the background
108       if(!AODEvent()->FindListObject(AliAODJetEventBackground::StdBranchName())){
109         AliAODJetEventBackground* evBkg = new AliAODJetEventBackground();
110         evBkg->SetName(AliAODJetEventBackground::StdBranchName());
111         AddAODBranch("AliAODJetEventBackground",&evBkg);
112       }
113       fJetFinder->ConnectAOD(AODEvent());
114     }
115   else
116     {
117       // Create a new branch for jets...
118       // how is this reset? -> cleared in the UserExec....
119       // Can this be handled by the framework?
120       // here we can also have the case that the brnaches are written to a separate file
121
122       TClonesArray *tca = new TClonesArray("AliAODJet", 0);
123       tca->SetName(fNonStdBranch.Data());
124       AddAODBranch("TClonesArray",&tca,fNonStdFile.Data());
125       if(!AODEvent()->FindListObject(Form("%s_%s",AliAODJetEventBackground::StdBranchName(),fNonStdBranch.Data()))){
126         AliAODJetEventBackground* evBkg = new AliAODJetEventBackground();
127         evBkg->SetName(Form("%s_%s",AliAODJetEventBackground::StdBranchName(),fNonStdBranch.Data()));
128         AddAODBranch("AliAODJetEventBackground",&evBkg,fNonStdFile.Data());
129       }
130
131       if(fNonStdFile.Length()!=0){
132         // 
133         // case that we have an AOD extension we need to fetch the jets from the extended output
134         // we identifay the extension aod event by looking for the branchname
135         AliAODHandler *aodH = dynamic_cast<AliAODHandler*>(AliAnalysisManager::GetAnalysisManager()->GetOutputEventHandler());
136         TObjArray* extArray = aodH->GetExtensions();
137         if (extArray) {
138           TIter next(extArray);
139           while ((fAODExtension=(AliAODExtension*)next())){
140             TObject *obj = fAODExtension->GetAOD()->FindListObject(fNonStdBranch.Data());
141             if(fDebug>10){
142               Printf("%s:%d Dumping..",(char*)__FILE__,__LINE__);
143               fAODExtension->GetAOD()->Dump();
144             }
145             if(obj){
146               if(fDebug>1)Printf("AODExtension found for %s",fNonStdBranch.Data());
147               break;
148             }
149             fAODExtension = 0;
150           }
151         }
152         if(fAODExtension)fJetFinder->ConnectAODNonStd(fAODExtension->GetAOD(), fNonStdBranch.Data()); 
153       }
154       else{
155         fJetFinder->ConnectAODNonStd(AODEvent(), fNonStdBranch.Data()); 
156       }
157     }
158
159
160   // Histograms
161   OpenFile(1);
162   fListOfHistos = new TList();
163   fHistos       = new AliJetHistos();
164   fHistos->AddHistosToList(fListOfHistos);
165   
166   // Add the JetFinderInformation to the Outputlist
167   AliJetHeader *fH = fJetFinder->GetHeader();
168   
169   // Compose a characteristic output name
170   // with the name of the output branch
171   if(fH) {
172     if(fNonStdBranch.Length()==0) {
173       fH->SetName("AliJetHeader_jets");
174     }
175     else {
176       fH->SetName(Form("AliJetHeader_%s",fNonStdBranch.Data()));
177     }
178   }
179   if(!fAODExtension)OutputTree()->GetUserInfo()->Add(fH);
180   else fAODExtension->GetTree()->GetUserInfo()->Add(fH);
181 }
182
183 //----------------------------------------------------------------
184 void AliAnalysisTaskJets::Init()
185 {
186   // Initialization
187   if (fDebug > 1) printf("AnalysisTaskJets::Init() \n");
188   
189   // Call configuration file
190   if (fConfigFile.Length()) {
191     gROOT->LoadMacro(fConfigFile);
192     fJetFinder = (AliJetFinder*) gInterpreter->ProcessLine("ConfigJetAnalysis()");
193   }
194   AliJetReaderHeader *header = (AliJetReaderHeader*)fJetFinder->GetReader()->GetReaderHeader();
195   fOpt = header->GetDetector();
196
197   // Initialise Jet Analysis
198   if(fOpt == 0) fJetFinder->Init();
199   else fJetFinder->InitTask(fChain); // V2
200 }
201
202
203 //----------------------------------------------------------------
204 void AliAnalysisTaskJets::UserExec(Option_t */*option*/)
205 {
206   // Execute analysis for current event
207   //
208   // Fill control histos
209   TClonesArray* jarray = 0;
210   AliAODJetEventBackground*  evBkg;
211
212   if(fNonStdBranch.Length()==0) {
213     jarray = AODEvent()->GetJets();
214     evBkg = (AliAODJetEventBackground*)(AODEvent()->FindListObject(AliAODJetEventBackground::StdBranchName()));
215     evBkg->Reset();
216   }
217   else {
218     jarray = (TClonesArray*)(AODEvent()->FindListObject(fNonStdBranch.Data()));
219     if(!jarray)jarray = (TClonesArray*)(fAODExtension->GetAOD()->FindListObject(fNonStdBranch.Data()));
220     jarray->Delete();    // this is our responsibility, clear before filling again
221     evBkg = (AliAODJetEventBackground*)(AODEvent()->FindListObject(Form("%s_%s",AliAODJetEventBackground::StdBranchName(),fNonStdBranch.Data())));
222     if(!evBkg)  evBkg = (AliAODJetEventBackground*)(fAODExtension->GetAOD()->FindListObject(Form("%s_%s",AliAODJetEventBackground::StdBranchName(),fNonStdBranch.Data())));
223     evBkg->Reset();
224   }
225
226   if (dynamic_cast<AliAODEvent*>(InputEvent()) !=  0 && !fReadAODFromOutput) {
227 // AOD is input event..........................................V                                       
228       fJetFinder->GetReader()->SetInputEvent(InputEvent(), InputEvent(), MCEvent());
229   } else {
230 // AOD is read from output ....................................V      
231       fJetFinder->GetReader()->SetInputEvent(InputEvent(), AODEvent(), MCEvent());
232   }
233   
234   
235
236   if(fOpt==0) fJetFinder->ProcessEvent();
237   else fJetFinder->ProcessEvent2();    // V2
238  
239   // Fill control histos
240   fHistos->FillHistos(jarray);
241
242   // Post the data
243   PostData(1, fListOfHistos);
244   return;
245 }
246
247
248 //*************************************************************
249
250 void AliAnalysisTaskJets::Terminate(Option_t */*option*/)
251 {
252 // Terminate analysis
253 //
254     if (fDebug > 1) printf("AnalysisJets: Terminate() \n");
255 //    if (fJetFinder) fJetFinder->FinishRun();
256 }
257