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