]> git.uio.no Git - u/mrichter/AliRoot.git/blame - JETAN/AliJet.cxx
Updated version of the mutiplicity reconstruction with tracklets (Massimo)
[u/mrichter/AliRoot.git] / JETAN / AliJet.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// Jet class
18// Stores the output of a jet algorithm
19// Author: jgcn@mda.cinvestav.mx
20//---------------------------------------------------------------------
21
22#include <Riostream.h>
23#include <TClonesArray.h>
24#include <TLorentzVector.h>
25
26#include "AliJet.h"
27ClassImp(AliJet)
1b7d5d7e 28
29AliJet::AliJet():
30 fNInput(0),
31 fNJets(0),
32 fEtAvg(0),
33 fInJet(0),
34 fMultiplicities(0),
35 fNCells(0),
36 fPtFromSignal(0),
37 fJets(0),
38 fEtaIn(0),
39 fPhiIn(0),
40 fPtIn(0)
99e5fe42 41{
1b7d5d7e 42 // Default constructor
99e5fe42 43 fJets = new TClonesArray("TLorentzVector",1000);
1b7d5d7e 44 fInJet = TArrayI();
45 fPtIn = TArrayF();
46 fEtaIn = TArrayF();
47 fPhiIn = TArrayF();
48 fPtFromSignal = TArrayF();
49 fMultiplicities = TArrayI();
50 fNCells = TArrayI();
99e5fe42 51}
52
53////////////////////////////////////////////////////////////////////////
54
55AliJet::~AliJet()
56{
57 // destructor
58 if (fJets) {
59 fJets->Delete();
60 delete fJets;
61 }
62}
63
64////////////////////////////////////////////////////////////////////////
65
66Bool_t AliJet::OutOfRange(Int_t i, const char *s) const
67{
83a444b1 68 // checks if i is a valid index. s = name of calling method
99e5fe42 69 if (i >= fNJets || i < 0) {
70 cout << s << " Index " << i << " out of range" << endl;
71 return kTRUE;
72 }
73 return kFALSE;
74}
75
76////////////////////////////////////////////////////////////////////////
77
78TLorentzVector* AliJet::GetJet(Int_t i)
79{
80 // returns i-jet
81 if (OutOfRange(i, "AliJet::GetJet:")) return 0;
99e5fe42 82 TLorentzVector *lv = (TLorentzVector*) fJets->At(i);
83 return lv;
84}
85
86////////////////////////////////////////////////////////////////////////
87
83a444b1 88Int_t AliJet::GetMultiplicity(Int_t i) const
99e5fe42 89{
90 // gets multiplicity of i-jet
91 if (OutOfRange(i, "AliJet::GetMultiplicity:")) return 0;
92 return fMultiplicities[i];
93}
94
95////////////////////////////////////////////////////////////////////////
96
83a444b1 97Int_t AliJet::GetNCell(Int_t i) const
98{
99 // gets number of cell of i-jet
100 if (OutOfRange(i, "AliJet::GetNCell:")) return 0;
101 return fNCells[i];
102}
103
104////////////////////////////////////////////////////////////////////////
105
99e5fe42 106Double_t AliJet::GetPx(Int_t i)
107{
108// Get Px component of jet i
109 if (OutOfRange(i, "AliJet::GetPx:")) return -1e30;
99e5fe42 110 TLorentzVector *lv = (TLorentzVector*) fJets->At(i);
111 return lv->Px();
112}
113
114////////////////////////////////////////////////////////////////////////
115
116Double_t AliJet::GetPy(Int_t i)
117{
118// Get Py component of jet i
119 if (OutOfRange(i, "AliJet::GetPy:")) return -1e30;
99e5fe42 120 TLorentzVector *lv = (TLorentzVector*) fJets->At(i);
121 return lv->Py();
122}
123
124////////////////////////////////////////////////////////////////////////
125
126Double_t AliJet::GetPz(Int_t i)
127{
128// Get Pz component of jet i
129 if (OutOfRange(i, "AliJet::GetPz:")) return -1e30;
99e5fe42 130 TLorentzVector *lv = (TLorentzVector*) fJets->At(i);
131 return lv->Pz();
132}
133
134////////////////////////////////////////////////////////////////////////
135
136Double_t AliJet::GetP(Int_t i)
137{
138// Get momentum of jet i
139 if (OutOfRange(i, "AliJet::GetP:")) return -1e30;
99e5fe42 140 TLorentzVector *lv = (TLorentzVector*) fJets->At(i);
141 return lv->P();
142}
143
144////////////////////////////////////////////////////////////////////////
145
146Double_t AliJet::GetE(Int_t i)
147{
148// Get energy of jet i
149 if (OutOfRange(i, "AliJet::GetE:")) return -1e30;
99e5fe42 150 TLorentzVector *lv = (TLorentzVector*) fJets->At(i);
151 return lv->E();
152}
153
154////////////////////////////////////////////////////////////////////////
155
156Double_t AliJet::GetPt(Int_t i)
157{
158// Get transverse momentum of jet i
159 if (OutOfRange(i, "AliJet::GetPt:")) return -1e30;
99e5fe42 160 TLorentzVector *lv = (TLorentzVector*) fJets->At(i);
161 return lv->Pt();
162}
163
164////////////////////////////////////////////////////////////////////////
165
166Double_t AliJet::GetEta(Int_t i)
167{
168// Get eta of jet i
169 if (OutOfRange(i, "AliJet::GetEta:")) return -1e30;
99e5fe42 170 TLorentzVector *lv = (TLorentzVector*) fJets->At(i);
171 return lv->Eta();
172}
173
174////////////////////////////////////////////////////////////////////////
175
176Double_t AliJet::GetPhi(Int_t i)
177{
178// Get phi of jet i
179 if (OutOfRange(i, "AliJet::GetPhi:")) return -1e30;
99e5fe42 180 TLorentzVector *lv = (TLorentzVector*) fJets->At(i);
181 return ( (lv->Phi() < 0) ? (lv->Phi()) + 2. * TMath::Pi() : lv->Phi());
182}
183
184////////////////////////////////////////////////////////////////////////
185
186Double_t AliJet::GetTheta(Int_t i)
187{
188// Get theta of jet i
189 if (OutOfRange(i, "AliJet::GetTheta:")) return -1e30;
99e5fe42 190 TLorentzVector *lv = (TLorentzVector*) fJets->At(i);
191 return lv->Theta();
192}
193
194////////////////////////////////////////////////////////////////////////
195
196Double_t AliJet::GetMass(Int_t i)
197{
198// Get invariant mass of jet i
199 if (OutOfRange(i, "AliJet::GetMass:")) return -1e30;
99e5fe42 200 TLorentzVector *lv = (TLorentzVector*) fJets->At(i);
201 return lv->M();
202}
203
204////////////////////////////////////////////////////////////////////////
205
206
207void AliJet::AddJet(Double_t px, Double_t py, Double_t pz, Double_t e)
208{
209// Add new jet to the list
210 new ((*fJets)[fNJets++]) TLorentzVector(px,py,pz,e);
211}
212
213////////////////////////////////////////////////////////////////////////
214
215void AliJet::SetInJet(Int_t* j)
216{
217 // set information of which input object belongs
83a444b1 218 // to each jet. If zero, object was not assigned to
219 // a jet, if n,positive, it was assiged to jet n
220 // if n, negative, it is within cone of jet n, but
221 // it did not passed the user cuts. filled in by AliJetFinder
99e5fe42 222 if (fNInput>0) fInJet.Set(fNInput, j);
223}
224
225////////////////////////////////////////////////////////////////////////
226
83a444b1 227void AliJet::SetEtaIn(Float_t* r)
228{
229 if (fNInput>0) fEtaIn.Set(fNInput, r);
230}
231
232////////////////////////////////////////////////////////////////////////
233
234void AliJet::SetPtIn(Float_t* pt)
235{
236 if (fNInput>0) fPtIn.Set(fNInput, pt);
237}
238
239////////////////////////////////////////////////////////////////////////
240
241void AliJet::SetPhiIn(Float_t* x)
242{
243 if (fNInput>0) fPhiIn.Set(fNInput, x);
244}
245
246////////////////////////////////////////////////////////////////////////
247
365b766d 248void AliJet::SetPtFromSignal(Float_t* p)
249{
250 // set information of percentage of pt of jets
251 // coming from signal (ie Pythia)
252 if (fNJets>0) fPtFromSignal.Set(fNJets, p);
253}
254
255////////////////////////////////////////////////////////////////////////
256
99e5fe42 257void AliJet::SetMultiplicities(Int_t* m)
258{
259 // set information of jet multiplicities
260 // filled in by AliJetFinder
261 if (fNJets>0) fMultiplicities.Set(fNJets, m);
262}
263
264////////////////////////////////////////////////////////////////////////
265
83a444b1 266void AliJet::SetNCells(Int_t* n)
267{
268 if (fNJets>0) fNCells.Set(fNJets, n);
269}
270
271////////////////////////////////////////////////////////////////////////
272
99e5fe42 273void AliJet::ClearJets(Option_t *option)
274{
275 // reset all values
276 fJets->Clear(option);
277 fNInput=0;
278 fNJets=0;
279 fMultiplicities.Set(0);
280 fInJet.Set(0);
83a444b1 281 fPtFromSignal.Set(0);
282 fPhiIn.Set(0);
283 fEtaIn.Set(0);
284 fPtIn.Set(0);
285 fNCells.Set(0);
99e5fe42 286}
287
288////////////////////////////////////////////////////////////////////////
289
290void AliJet::PrintJets()
291{
292// Print jet information
293 if (fNJets == 0) {
294 cout << " AliJet::PrintJets: There are no jets in this event " << endl;
295 return;
296 }
297 cout << " AliJet::PrintJets: There are " << fNJets
298 << " jets in this event" << endl;
299 for(Int_t i=0;i<fNJets;i++) {
300 cout << " Jet " << i << " (px,py,pz,en)=(" << GetPx(i)
301 << "," << GetPy(i)
302 << "," << GetPz(i)
303 << "," << GetE(i)
304 << ")" << endl;
305 cout << " (pt,eta,phi)=(" << GetPt(i)
306 << "," << GetEta(i)
307 << "," << GetPhi(i) << ")" << endl;
308 cout << " # of tracks =" << GetMultiplicity(i) << endl;
309 }
310}