]> git.uio.no Git - u/mrichter/AliRoot.git/blame - JETAN/AliJetFinder.cxx
First functional implementation.
[u/mrichter/AliRoot.git] / JETAN / AliJetFinder.cxx
CommitLineData
99e5fe42 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 **************************************************************************/
83a444b1 15
52ff852a 16/* $Id$ */
17
99e5fe42 18//---------------------------------------------------------------------
19// Jet finder base class
20// manages the search for jets
7d0f353c 21// Authors: jgcn@mda.cinvestav.mx
22// andreas.morsch@cern.ch
99e5fe42 23//---------------------------------------------------------------------
24
25#include <Riostream.h>
26#include <TFile.h>
99e5fe42 27#include "AliJetFinder.h"
28#include "AliJet.h"
29#include "AliJetReader.h"
30#include "AliJetReaderHeader.h"
31#include "AliJetControlPlots.h"
32#include "AliLeading.h"
99e5fe42 33
99e5fe42 34ClassImp(AliJetFinder)
35
1b7d5d7e 36AliJetFinder::AliJetFinder():
eaabc21f 37 fTreeJ(0),
7d0f353c 38 fPlotMode(kFALSE),
39 fJets(0),
40 fGenJets(0),
41 fLeading(0),
42 fReader(0x0),
19e6695b 43 fHeader(0x0),
7d0f353c 44 fPlots(0x0),
45 fOut(0x0)
46
99e5fe42 47{
99e5fe42 48 // Constructor
655dbb2b 49 fJets = new AliJet();
99e5fe42 50 fGenJets = new AliJet();
51 fLeading = new AliLeading();
99e5fe42 52}
53
99e5fe42 54////////////////////////////////////////////////////////////////////////
55
56AliJetFinder::~AliJetFinder()
57{
99e5fe42 58 // destructor
99e5fe42 59 // here reset and delete jets
60 fJets->ClearJets();
61 delete fJets;
62 fGenJets->ClearJets();
63 delete fGenJets;
64 // close file
65 if (fOut) {
66 fOut->Close();
67 fOut->Delete();
68 }
69 delete fOut;
70 // reset and delete control plots
71 if (fPlots) delete fPlots;
99e5fe42 72}
73
99e5fe42 74////////////////////////////////////////////////////////////////////////
75
c457fc2e 76void AliJetFinder::SetOutputFile(const char* name)
99e5fe42 77{
eaabc21f 78 // opens output file
c457fc2e 79 fOut = new TFile(name,"recreate");
99e5fe42 80}
81
99e5fe42 82////////////////////////////////////////////////////////////////////////
83
84void AliJetFinder::PrintJets()
85{
83a444b1 86 // Print jet information
99e5fe42 87 cout << " Jets found with jet algorithm:" << endl;
88 fJets->PrintJets();
89 cout << " Jets found by pythia:" << endl;
90 fGenJets->PrintJets();
91}
92
99e5fe42 93////////////////////////////////////////////////////////////////////////
94
95void AliJetFinder::SetPlotMode(Bool_t b)
96{
83a444b1 97 // Sets the plotting mode
99e5fe42 98 fPlotMode=b;
99 if (b && !fPlots) fPlots = new AliJetControlPlots();
100}
101
102////////////////////////////////////////////////////////////////////////
eaabc21f 103TTree* AliJetFinder::MakeTreeJ(char* name)
99e5fe42 104{
eaabc21f 105 // Create the tree for reconstructed jets
eaabc21f 106 fTreeJ = new TTree(name, "AliJet");
107 fTreeJ->Branch("FoundJet", &fJets, 1000);
108 fTreeJ->Branch("GenJet", &fGenJets,1000);
109 fTreeJ->Branch("LeadingPart",&fLeading,1000);
110 return fTreeJ;
99e5fe42 111}
112
113////////////////////////////////////////////////////////////////////////
114
115void AliJetFinder::WriteRHeaderToFile()
116{
117 // write reader header
76c48857 118 AliJetReaderHeader *rh = fReader->GetReaderHeader();
119 rh->Write();
99e5fe42 120}
121
99e5fe42 122////////////////////////////////////////////////////////////////////////
123
99e5fe42 124
125void AliJetFinder::Run()
126{
7d0f353c 127 // Do some initialization
99e5fe42 128 Init();
99e5fe42 129 // connect files
130 fReader->OpenInputFiles();
131
132 // write headers
7d0f353c 133 WriteHeaders();
99e5fe42 134 // loop over events
b45b0c92 135 Int_t nFirst, nLast, option, debug, arrayInitialised;
99e5fe42 136 nFirst = fReader->GetReaderHeader()->GetFirstEvent();
b45b0c92 137 nLast = fReader->GetReaderHeader()->GetLastEvent();
138
139 option = fReader->GetReaderHeader()->GetDetector();
140 debug = fReader->GetReaderHeader()->GetDebug();
141 arrayInitialised = fReader->GetArrayInitialised();
142
99e5fe42 143 // loop over events
144 for (Int_t i=nFirst;i<nLast;i++) {
145 fReader->FillMomentumArray(i);
146 fLeading->FindLeading(fReader);
7d0f353c 147 fReader->GetGenJets(fGenJets);
b45b0c92 148
149 if (option == 0) { // TPC with fMomentumArray
150 if(debug > 1)
151 printf("In FindJetsTPC() routine: find jets with fMomentumArray !!!\n");
152 FindJetsTPC();
153 } else {
c457fc2e 154 if(debug > 1) printf("In FindJets() routine: find jets with fUnitArray !!!\n");
155 FindJets();
b45b0c92 156 }
eaabc21f 157 if (fOut) {
158 fOut->cd();
159 fTreeJ->Fill();
160 }
161
83a444b1 162 if (fPlots) fPlots->FillHistos(fJets);
99e5fe42 163 fLeading->Reset();
164 fGenJets->ClearJets();
165 Reset();
166 }
b45b0c92 167
99e5fe42 168 // write out
169 if (fPlots) {
170 fPlots->Normalize();
171 fPlots->PlotHistos();
172 }
173 if (fOut) {
174 fOut->cd();
175 fPlots->Write();
176 fOut->Close();
177 }
178}
7d0f353c 179
180
181//
182// The following methods have been added to allow for event steering from the outside
183//
184
f3f3617d 185void AliJetFinder::ConnectTree(TTree* tree, TObject* data)
7d0f353c 186{
187 // Connect the input file
f3f3617d 188 fReader->ConnectTree(tree, data);
7d0f353c 189}
190
191void AliJetFinder::WriteHeaders()
192{
193 // Write the Headers
76c48857 194 TFile* f = new TFile("jets_local.root", "recreate");
195 WriteRHeaderToFile();
196 WriteJHeaderToFile();
197 f->Close();
7d0f353c 198}
199
200
201Bool_t AliJetFinder::ProcessEvent(Long64_t entry)
202{
203//
204// Process one event
205//
eaabc21f 206 Int_t debug = fReader->GetReaderHeader()->GetDebug();
207 if (debug > 0) printf("<<<<< Processing Event %5d >>>>> \n", (Int_t) entry);
76c48857 208
7d0f353c 209 Bool_t ok = fReader->FillMomentumArray(entry);
210 if (!ok) return kFALSE;
76c48857 211
212 // Leading particles
7d0f353c 213 fLeading->FindLeading(fReader);
76c48857 214 // Jets
7d0f353c 215 FindJets();
76c48857 216 // Fill the tree
217 fTreeJ->Fill();
218
7d0f353c 219 if (fPlots) fPlots->FillHistos(fJets);
220 fLeading->Reset();
221 fGenJets->ClearJets();
222 Reset();
223 return kTRUE;
224}
225
226void AliJetFinder::FinishRun()
227{
228 // Finish a run
eaabc21f 229 if (fPlots) {
230 fPlots->Normalize();
231 fPlots->PlotHistos();
232 }
233
234 if (fOut) {
235 fOut->cd();
236 fTreeJ->Write();
237 if (fPlots) {
238 fPlots->Write();
239 }
240 fOut->Close();
241 }
7d0f353c 242}
243