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