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