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