]> git.uio.no Git - u/mrichter/AliRoot.git/blame - JETAN/AliJetFinder.cxx
Fast EMCAL simulation option added.
[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 **************************************************************************/
15
16
17//---------------------------------------------------------------------
18// Jet finder base class
19// manages the search for jets
20// Author: jgcn@mda.cinvestav.mx
21//---------------------------------------------------------------------
22
23#include <Riostream.h>
24#include <TFile.h>
25#include "AliGenPythiaEventHeader.h"
26#include "AliJetFinder.h"
27#include "AliJet.h"
28#include "AliJetReader.h"
29#include "AliJetReaderHeader.h"
30#include "AliJetControlPlots.h"
31#include "AliLeading.h"
32#include "AliHeader.h"
33
34
35ClassImp(AliJetFinder)
36
37////////////////////////////////////////////////////////////////////////
38
39AliJetFinder::AliJetFinder()
40{
41 //
42 // Constructor
43 //
655dbb2b 44 fOut = 0x0;
45 fJets = new AliJet();
99e5fe42 46 fGenJets = new AliJet();
47 fLeading = new AliLeading();
655dbb2b 48 fReader = 0x0;
49 fPlots = 0x0;
99e5fe42 50 SetPlotMode(kFALSE);
51}
52
53
54////////////////////////////////////////////////////////////////////////
55
56AliJetFinder::~AliJetFinder()
57{
58 //
59 // destructor
60 //
61
62 // here reset and delete jets
63 fJets->ClearJets();
64 delete fJets;
65 fGenJets->ClearJets();
66 delete fGenJets;
67 // close file
68 if (fOut) {
69 fOut->Close();
70 fOut->Delete();
71 }
72 delete fOut;
73 // reset and delete control plots
74 if (fPlots) delete fPlots;
ae8bbe84 75 // delete fLeading;
99e5fe42 76}
77
78
79////////////////////////////////////////////////////////////////////////
80
81void AliJetFinder::SetOutputFile(const char *name)
82{
83 // opens output file
84 fOut = new TFile(name,"recreate");
85}
86
87
88////////////////////////////////////////////////////////////////////////
89
90void AliJetFinder::PrintJets()
91{
92//
93// Print jet information
94 cout << " Jets found with jet algorithm:" << endl;
95 fJets->PrintJets();
96 cout << " Jets found by pythia:" << endl;
97 fGenJets->PrintJets();
98}
99
100
101////////////////////////////////////////////////////////////////////////
102
103void AliJetFinder::SetPlotMode(Bool_t b)
104{
105// Sets the plotting mode
106 fPlotMode=b;
107 if (b && !fPlots) fPlots = new AliJetControlPlots();
108}
109
110////////////////////////////////////////////////////////////////////////
111
112void AliJetFinder::WriteJetsToFile(Int_t i)
113{
114// Writes the jets to file
115 fOut->cd();
116 char hname[30];
117 sprintf(hname,"TreeJ%d",i);
118 TTree* jetT = new TTree(hname,"AliJet");
119 jetT->Branch("FoundJet",&fJets,1000);
120 jetT->Branch("GenJet",&fGenJets,1000);
121 jetT->Branch("LeadingPart",&fLeading,1000);
122 jetT->Fill();
123 jetT->Write(hname);
124 delete jetT;
125}
126
127////////////////////////////////////////////////////////////////////////
128
129void AliJetFinder::WriteRHeaderToFile()
130{
131 // write reader header
132 fOut->cd();
133 AliJetReaderHeader *rh = fReader->GetReaderHeader();
134 rh->Write();
135}
136
137
138////////////////////////////////////////////////////////////////////////
139
140void AliJetFinder::GetGenJets()
141{
142// Get the generated jet information from mc header
99e5fe42 143 AliHeader* alih = fReader->GetAliHeader();
144 if (alih == 0) return;
145 AliGenEventHeader * genh = alih->GenEventHeader();
146 if (genh == 0) return;
147 Int_t nj =((AliGenPythiaEventHeader*)genh)->NTriggerJets();
148 Int_t* m = new Int_t[nj];
149 Int_t* k = new Int_t[nj];
150 for (Int_t i=0; i< nj; i++) {
151 Float_t p[4];
152 ((AliGenPythiaEventHeader*)genh)->TriggerJet(i,p);
153 fGenJets->AddJet(p[0],p[1],p[2],p[3]);
154 m[i]=1;
155 k[i]=i;
156 }
655dbb2b 157 fGenJets->SetNinput(nj);
99e5fe42 158 fGenJets->SetMultiplicities(m);
159 fGenJets->SetInJet(k);
160}
161
162////////////////////////////////////////////////////////////////////////
163
164void AliJetFinder::Run()
165{
166 // do some initialization
167 Init();
168
169 // connect files
170 fReader->OpenInputFiles();
171
172 // write headers
173 if (fOut) {
174 fOut->cd();
175 WriteRHeaderToFile();
176 WriteJHeaderToFile();
177 }
178
179 // loop over events
180 Int_t nFirst,nLast;
181 nFirst = fReader->GetReaderHeader()->GetFirstEvent();
182 nLast = fReader->GetReaderHeader()->GetLastEvent();
183 // loop over events
184 for (Int_t i=nFirst;i<nLast;i++) {
185 fReader->FillMomentumArray(i);
186 fLeading->FindLeading(fReader);
187 GetGenJets();
188 FindJets();
189 if (fOut) WriteJetsToFile(i);
190 if (fPlots) fPlots->FillHistos(fJets,fReader);
191 fLeading->Reset();
192 fGenJets->ClearJets();
193 Reset();
194 }
195 // write out
196 if (fPlots) {
197 fPlots->Normalize();
198 fPlots->PlotHistos();
199 }
200 if (fOut) {
201 fOut->cd();
202 fPlots->Write();
203 fOut->Close();
204 }
205}
206