]> git.uio.no Git - u/mrichter/AliRoot.git/blob - JETAN/AliJet.cxx
First commit of new jet reconstruction and analysis code to be used for the
[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 //---------------------------------------------------------------------
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"
27 ClassImp(AliJet)
28  
29  
30 ////////////////////////////////////////////////////////////////////////
31
32 AliJet::AliJet() 
33 {
34   //
35   // Constructor
36   //
37   fJets = new TClonesArray("TLorentzVector",1000);
38   fNInput=0;
39   fNJets=0;
40   fInJet=TArrayI();
41   fMultiplicities=TArrayI();
42
43
44 ////////////////////////////////////////////////////////////////////////
45
46 AliJet::~AliJet()
47 {
48   // destructor
49   if (fJets) {
50     fJets->Delete();
51     delete fJets;
52   }
53 }
54
55 ////////////////////////////////////////////////////////////////////////
56
57 Bool_t AliJet::OutOfRange(Int_t i, const char *s) const
58 {
59   // checks if i is a valid index. s= name of calling method
60   if (i >= fNJets || i < 0) {
61     cout << s << " Index " << i << " out of range" << endl;
62     return kTRUE;
63   }
64   return kFALSE;
65 }
66
67 ////////////////////////////////////////////////////////////////////////
68
69 TLorentzVector* AliJet::GetJet(Int_t i)
70 {
71   // returns i-jet
72   if (OutOfRange(i, "AliJet::GetJet:")) return 0;
73
74   TLorentzVector *lv = (TLorentzVector*) fJets->At(i);
75   return lv; 
76 }
77
78 ////////////////////////////////////////////////////////////////////////
79
80 Int_t AliJet::GetMultiplicity(Int_t i)
81 {
82   // gets multiplicity of i-jet
83   if (OutOfRange(i, "AliJet::GetMultiplicity:")) return 0;
84   return fMultiplicities[i];
85 }
86
87 ////////////////////////////////////////////////////////////////////////
88
89 Double_t AliJet::GetPx(Int_t i)
90 {
91 // Get Px component of jet i
92   if (OutOfRange(i, "AliJet::GetPx:")) return -1e30;
93
94   TLorentzVector *lv = (TLorentzVector*) fJets->At(i);
95   return lv->Px();
96 }
97
98 ////////////////////////////////////////////////////////////////////////
99
100 Double_t AliJet::GetPy(Int_t i)
101 {
102 // Get Py component of jet i
103   if (OutOfRange(i, "AliJet::GetPy:")) return -1e30;
104
105   TLorentzVector *lv = (TLorentzVector*) fJets->At(i);
106   return lv->Py();
107 }
108
109 ////////////////////////////////////////////////////////////////////////
110
111 Double_t AliJet::GetPz(Int_t i)
112 {
113 // Get Pz component of jet i
114   if (OutOfRange(i, "AliJet::GetPz:")) return -1e30;
115
116   TLorentzVector *lv = (TLorentzVector*) fJets->At(i);
117   return lv->Pz();
118 }
119
120 ////////////////////////////////////////////////////////////////////////
121
122 Double_t AliJet::GetP(Int_t i)
123 {
124 // Get momentum of jet i
125   if (OutOfRange(i, "AliJet::GetP:")) return -1e30;
126
127   TLorentzVector *lv = (TLorentzVector*) fJets->At(i);
128   return lv->P();
129 }
130
131 ////////////////////////////////////////////////////////////////////////
132
133 Double_t AliJet::GetE(Int_t i)
134 {
135 // Get energy of jet i
136   if (OutOfRange(i, "AliJet::GetE:")) return -1e30;
137
138   TLorentzVector *lv = (TLorentzVector*) fJets->At(i);
139   return lv->E();
140 }
141
142 ////////////////////////////////////////////////////////////////////////
143
144 Double_t AliJet::GetPt(Int_t i)
145 {
146 // Get transverse momentum of jet i
147   if (OutOfRange(i, "AliJet::GetPt:")) return -1e30;
148
149   TLorentzVector *lv = (TLorentzVector*) fJets->At(i);
150   return lv->Pt();
151 }
152
153 ////////////////////////////////////////////////////////////////////////
154
155 Double_t AliJet::GetEta(Int_t i)
156 {
157 // Get eta of jet i
158   if (OutOfRange(i, "AliJet::GetEta:")) return -1e30;
159
160   TLorentzVector *lv = (TLorentzVector*) fJets->At(i);
161   return lv->Eta();
162 }
163
164 ////////////////////////////////////////////////////////////////////////
165
166 Double_t AliJet::GetPhi(Int_t i)
167 {
168 // Get phi of jet i
169   if (OutOfRange(i, "AliJet::GetPhi:")) return -1e30;
170
171   TLorentzVector *lv = (TLorentzVector*) fJets->At(i);
172   return ( (lv->Phi() < 0) ? (lv->Phi()) + 2. * TMath::Pi() : lv->Phi());
173 }
174
175 ////////////////////////////////////////////////////////////////////////
176
177 Double_t AliJet::GetTheta(Int_t i)
178 {
179 // Get theta of jet i
180   if (OutOfRange(i, "AliJet::GetTheta:")) return -1e30;
181
182   TLorentzVector *lv = (TLorentzVector*) fJets->At(i);
183   return lv->Theta();
184 }
185
186 ////////////////////////////////////////////////////////////////////////
187
188 Double_t AliJet::GetMass(Int_t i)
189 {
190 // Get invariant mass of jet i
191   if (OutOfRange(i, "AliJet::GetMass:")) return -1e30;
192
193   TLorentzVector *lv = (TLorentzVector*) fJets->At(i);
194   return lv->M();
195 }
196
197 ////////////////////////////////////////////////////////////////////////
198  
199
200 void AliJet::AddJet(Double_t px, Double_t py, Double_t pz, Double_t e)
201 {
202 // Add new jet to the list
203   new ((*fJets)[fNJets++]) TLorentzVector(px,py,pz,e);
204 }
205
206 ////////////////////////////////////////////////////////////////////////
207
208 void AliJet::SetInJet(Int_t* j)
209 {
210   // set information of which input object belongs
211   // to each jet. filled in by AliJetFinder
212   if (fNInput>0) fInJet.Set(fNInput, j);
213 }
214
215 ////////////////////////////////////////////////////////////////////////
216
217 void AliJet::SetMultiplicities(Int_t* m)
218 {
219   // set information of jet multiplicities
220   // filled in by AliJetFinder
221   if (fNJets>0) fMultiplicities.Set(fNJets, m);
222 }
223
224 ////////////////////////////////////////////////////////////////////////
225
226 void AliJet::ClearJets(Option_t *option)
227 {
228   // reset all values
229   fJets->Clear(option);
230   fNInput=0;
231   fNJets=0;
232   fMultiplicities.Set(0);
233   fInJet.Set(0); 
234 }
235
236 ////////////////////////////////////////////////////////////////////////
237
238 void AliJet::PrintJets()
239 {
240 // Print jet information
241   if (fNJets == 0) {
242     cout << " AliJet::PrintJets: There are no jets in this event " << endl;
243     return;
244   }
245   cout << " AliJet::PrintJets: There are " << fNJets
246        << " jets in this event" << endl;
247   for(Int_t i=0;i<fNJets;i++) {
248     cout << "   Jet " << i << " (px,py,pz,en)=(" << GetPx(i)
249          << "," << GetPy(i)
250          << "," << GetPz(i)
251          << "," << GetE(i) 
252          << ")" << endl;
253     cout << "         (pt,eta,phi)=(" << GetPt(i)
254          << "," << GetEta(i)
255          << "," << GetPhi(i) << ")" << endl;
256     cout << "         # of tracks =" << GetMultiplicity(i) << endl;
257   }
258 }