]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RALICE/AliJet.cxx
Add the possibility to select events with CINT1B/CMUS1B triggers when filling histogr...
[u/mrichter/AliRoot.git] / RALICE / 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 // Class AliJet
20 // Creation and investigation of a jet of particle tracks.
21 // An AliJet can be constructed by adding AliTracks.
22 //
23 // To provide maximal flexibility to the user, two modes of track storage
24 // are provided by means of the memberfunction SetTrackCopy().
25 //
26 // a) SetTrackCopy(0) (which is the default).
27 //    Only the pointers of the 'added' tracks are stored.
28 //    This mode is typically used by making jet studies based on a fixed list
29 //    of tracks which stays under user control or is contained for instance
30 //    in an AliEvent.  
31 //    In this way the AliJet just represents a 'logical structure' for the
32 //    physics analysis which can be embedded in e.g. an AliEvent or AliVertex.
33 //
34 //    Note :
35 //    Modifications made to the original tracks also affect the AliTrack objects
36 //    which are stored in the AliJet. 
37 //
38 // b) SetTrackCopy(1).
39 //    Of every 'added' track a private copy will be made of which the pointer
40 //    will be stored.
41 //    In this way the AliJet represents an entity on its own and modifications
42 //    made to the original tracks do not affect the AliTrack objects which are
43 //    stored in the AliJet. 
44 //    This mode will allow 'adding' many different AliTracks into an AliJet by
45 //    creating only one AliTrack instance in the main programme and using the
46 //    AliTrack::Reset() and AliTrack parameter setting memberfunctions.
47 //
48 // See also the documentation provided for the memberfunction SetOwner(). 
49 //
50 // Coding example to make 2 jets j1 and j2.
51 // ----------------------------------------
52 // j1 contains the AliTracks t1 and t2
53 // j2 contains 10 different AliTracks via tx
54 //
55 // AliTrack t1,t2;
56 //  ...
57 //  ... // code to fill the AliTrack data
58 //  ...
59 // AliJet j1();
60 // j1.AddTrack(t1);
61 // j1.AddTrack(t2);
62 //
63 // AliJet j2();
64 // j2.SetTrackCopy(1);
65 // AliTrack* tx=new AliTrack();
66 // for (Int_t i=0; i<10; i++)
67 // {
68 //  ...
69 //  ... // code to set momentum etc... of the track tx
70 //  ...
71 //  j2.AddTrack(tx);
72 //  tx->Reset();
73 // }
74 //
75 // j1.Data();
76 // j2.Data("sph");
77 //
78 // Float_t e1=j1.GetEnergy();
79 // Float_t pnorm=j1->GetMomentum();
80 // Ali3Vector p=j1->Get3Momentum();
81 // Float_t m=j1.GetInvmass();
82 // Int_t ntk=j1.GetNtracks();
83 // AliTrack* tj=j1.GetTrack(1);
84 //
85 // delete tx;
86 //
87 // Note : By default all quantities are in GeV, GeV/c or GeV/c**2
88 //        but the user can indicate the usage of a different scale
89 //        for the energy-momentum units via the SetEscale() memberfunction.
90 //        The actual energy-momentum unit scale can be obtained via the
91 //        GetEscale() memberfunction.
92 //
93 //--- Author: Nick van Eijndhoven 10-jul-1997 UU-SAP Utrecht
94 //- Modified: NvE $Date$ UU-SAP Utrecht
95 ///////////////////////////////////////////////////////////////////////////
96
97 #include <cstdlib>
98
99 #include "Riostream.h"
100 #include "TObjArray.h"
101
102 #include "AliJet.h"
103 #include "AliPositionObj.h"
104 #include "AliSignal.h"
105  
106 ClassImp(AliJet) // Class implementation to enable ROOT I/O
107  
108 AliJet::AliJet() : TNamed(),Ali4Vector()
109 {
110 // Default constructor
111 // All variables initialised to 0
112 // Initial maximum number of tracks is set to the default value
113  Init();
114  Reset();
115  SetNtinit();
116 }
117 ///////////////////////////////////////////////////////////////////////////
118 void AliJet::Init()
119 {
120 // Initialisation of pointers etc...
121  fTracks=0;
122  fNtinit=0;
123  fTrackCopy=0;
124  fRef=0;
125  fSelected=0;
126  fEscale=1;
127 }
128 ///////////////////////////////////////////////////////////////////////////
129 AliJet::AliJet(Int_t n) : TNamed(),Ali4Vector()
130 {
131 // Create a jet to hold initially a maximum of n tracks
132 // All variables initialised to 0
133  Init();
134  Reset();
135  if (n > 0)
136  {
137   SetNtinit(n);
138  }
139  else
140  {
141   cout << endl;
142   cout << " *AliJet* Initial max. number of tracks entered : " << n << endl;
143   cout << " This is invalid. Default initial maximum will be used." << endl;
144   cout << endl;
145   SetNtinit();
146  }
147 }
148 ///////////////////////////////////////////////////////////////////////////
149 AliJet::~AliJet()
150 {
151 // Default destructor
152  if (fTracks)
153  {
154   delete fTracks;
155   fTracks=0;
156  }
157  if (fRef)
158  {
159   delete fRef;
160   fRef=0;
161  }
162  if (fSelected)
163  {
164   delete fSelected;
165   fSelected=0;
166  }
167 }
168 ///////////////////////////////////////////////////////////////////////////
169 void AliJet::SetOwner(Bool_t own)
170 {
171 // Set ownership of all added objects. 
172 // The default parameter is own=kTRUE.
173 //
174 // Invokation of this memberfunction also sets all the copy modes
175 // (e.g. TrackCopy & co.) according to the value of own.
176 //
177 // This function (with own=kTRUE) is particularly useful when reading data
178 // from a tree/file, since Reset() will then actually remove all the
179 // added objects from memory irrespective of the copy mode settings
180 // during the tree/file creation process. In this way it provides a nice way
181 // of preventing possible memory leaks in the reading/analysis process.
182 //
183 // In addition this memberfunction can also be used as a shortcut to set all
184 // copy modes in one go during a tree/file creation process.
185 // However, in this case the user has to take care to only set/change the
186 // ownership (and copy mode) for empty objects (e.g. newly created objects
187 // or after invokation of the Reset() memberfunction) otherwise it will
188 // very likely result in inconsistent destructor behaviour.
189
190  Int_t mode=1;
191  if (!own) mode=0;
192  if (fTracks) fTracks->SetOwner(own);
193  fTrackCopy=mode;
194 }
195 ///////////////////////////////////////////////////////////////////////////
196 AliJet::AliJet(const AliJet& j) : TNamed(j),Ali4Vector(j)
197 {
198 // Copy constructor
199  fNtinit=j.fNtinit;
200  fNtmax=j.fNtmax;
201  fQ=j.fQ;
202  fEscale=j.fEscale;
203  fNtrk=j.fNtrk;
204  fTrackCopy=j.fTrackCopy;
205  fUserId=j.fUserId;
206  if (j.fRef) fRef=new AliPositionObj(*(j.fRef));
207
208  fSelected=0;
209
210  fTracks=0;
211  if (fNtrk)
212  {
213   fTracks=new TObjArray(fNtmax);
214   if (fTrackCopy) fTracks->SetOwner();
215  }
216
217  for (Int_t i=1; i<=fNtrk; i++)
218  {
219   AliTrack* tx=j.GetTrack(i);
220   if (fTrackCopy)
221   {
222    fTracks->Add(tx->Clone());
223   }
224   else
225   {
226    fTracks->Add(tx);
227   }
228  } 
229 }
230 ///////////////////////////////////////////////////////////////////////////
231 void AliJet::SetNtinit(Int_t n)
232 {
233 // Set the initial maximum number of tracks for this jet
234  fNtinit=n;
235  fNtmax=n;
236
237  if (fTracks)
238  {
239   delete fTracks;
240   fTracks=0;
241  }
242  if (fRef)
243  {
244   delete fRef;
245   fRef=0;
246  }
247 }
248 ///////////////////////////////////////////////////////////////////////////
249 void AliJet::Reset()
250 {
251 // Reset all variables to 0
252 // The max. number of tracks is set to the initial value again
253 // Note : The scale for the energy/momentum units will not be changed.
254  fNtrk=0;
255  fQ=0;
256  fUserId=0;
257  Double_t a[4]={0,0,0,0};
258  SetVector(a,"sph");
259  if (fNtinit > 0) SetNtinit(fNtinit);
260 }
261 ///////////////////////////////////////////////////////////////////////////
262 void AliJet::AddTrack(AliTrack& t)
263 {
264 // Add a track to the jet.
265 // In case the maximum number of tracks has been reached
266 // space will be extended to hold an additional amount of tracks as
267 // was initially reserved.
268 // See SetTrackCopy() to tailor the functionality of the stored structures.
269 //
270 // Notes :
271 // -------
272 // In case a private copy is made, this is performed via the Clone() memberfunction.
273 // All AliTrack and derived classes have the default TObject::Clone() memberfunction.
274 // However, derived classes generally contain an internal data structure which may
275 // include pointers to other objects. Therefore it is recommended to provide
276 // for all derived classes a specific copy constructor and override the default Clone()
277 // memberfunction using this copy constructor.
278 // An example for this may be seen from AliTrack.   
279 //
280 // In case NO private copy is made, a check will be performed if this
281 // specific track is already present in the jet.
282 // If this is the case, no action is performed to prevent multiple
283 // additions of the same track.
284
285
286  AddTrack(t,1);
287 }
288 ///////////////////////////////////////////////////////////////////////////
289 void AliJet::AddTrack(AliTrack& t,Int_t copy)
290 {
291 // Internal memberfunction to actually add a track to the jet.
292 // In case the maximum number of tracks has been reached
293 // space will be extended to hold an additional amount of tracks as
294 // was initially reserved.
295 //
296 // If copy=0 NO copy of the track will be made, irrespective of the setting
297 // of the TrackCopy flag.
298 // This allows a proper treatment of automatically generated connecting
299 // tracks between vertices.
300 //
301 // In case NO copy of the track is made, a check will be performed if this
302 // specific track is already present in the jet.
303 // If this is the case, no action is performed to prevent multiple
304 // additions of the same track.
305 //
306 // Note :
307 // In case a private copy is made, this is performed via the Clone() memberfunction.
308
309  if (!fTracks)
310  {
311   fTracks=new TObjArray(fNtmax);
312   if (fTrackCopy) fTracks->SetOwner();
313  }
314  else if (!fTrackCopy || !copy) // Check if this track is already present
315  {
316   for (Int_t i=0; i<fNtrk; i++)
317   {
318    AliTrack* tx=(AliTrack*)fTracks->At(i);
319    if (tx == &t) return;
320   }
321  }
322
323  if (fNtrk == fNtmax) // Check if maximum track number is reached
324  {
325   fNtmax+=fNtinit;
326   fTracks->Expand(fNtmax);
327  }
328  
329  // Add the track to this jet
330  fNtrk++;
331  if (fTrackCopy && copy)
332  {
333   fTracks->Add(t.Clone());
334  }
335  else
336  {
337   fTracks->Add(&t);
338  }
339
340  fQ+=t.GetCharge();
341
342  Ali4Vector p4=(Ali4Vector)t;
343  Float_t tscale=t.GetEscale();
344  if ((tscale/fEscale > 1.1) || (fEscale/tscale > 1.1)) p4=p4*(tscale/fEscale);
345  (*this)+=p4;
346
347 }
348 ///////////////////////////////////////////////////////////////////////////
349 void AliJet::Data(TString f,TString u)
350 {
351 // Provide jet information within the coordinate frame f
352 //
353 // The string argument "u" allows to choose between different angular units
354 // in case e.g. a spherical frame is selected.
355 // u = "rad" : angles provided in radians
356 //     "deg" : angles provided in degrees
357 //
358 // The defaults are f="car" and u="rad".
359
360  const char* name=GetName();
361  const char* title=GetTitle();
362
363  cout << " *AliJet::Data*";
364  if (strlen(name))  cout << " Name : " << GetName();
365  if (strlen(title)) cout << " Title : " << GetTitle();
366  cout << endl;
367  cout << " Id : " << fUserId << " Invmass : " << GetInvmass() << " Charge : " << fQ
368       << " Momentum : " << GetMomentum() << " Energy scale : " << fEscale << " GeV" << endl;
369
370  ShowTracks(0);
371
372  Ali4Vector::Data(f,u); 
373
374 ///////////////////////////////////////////////////////////////////////////
375 void AliJet::List(TString f,TString u)
376 {
377 // Provide jet and primary track information within the coordinate frame f
378 //
379 // The string argument "u" allows to choose between different angular units
380 // in case e.g. a spherical frame is selected.
381 // u = "rad" : angles provided in radians
382 //     "deg" : angles provided in degrees
383 //
384 // The defaults are f="car" and u="rad".
385
386  Data(f,u); // Information of the current jet
387  if (fRef)   { cout << " Ref-point   :"; fRef->Data(f,u); }
388
389  // The tracks of this jet
390  AliTrack* t; 
391  for (Int_t it=1; it<=fNtrk; it++)
392  {
393   t=GetTrack(it);
394   if (t)
395   {
396    cout << "  ---Track no. " << it << endl;
397    cout << " ";
398    t->Data(f,u); 
399   }
400   else
401   {
402    cout << " *AliJet::List* Error : No track present." << endl; 
403   }
404  }
405
406 ///////////////////////////////////////////////////////////////////////////
407 void AliJet::ListAll(TString f,TString u)
408 {
409 // Provide jet and prim.+sec. track information within the coordinate frame f
410 //
411 // The string argument "u" allows to choose between different angular units
412 // in case e.g. a spherical frame is selected.
413 // u = "rad" : angles provided in radians
414 //     "deg" : angles provided in degrees
415 //
416 // The defaults are f="car" and u="rad".
417
418  Data(f,u); // Information of the current jet
419  if (fRef)   { cout << " Ref-point   :"; fRef->Data(f,u); }
420
421  // The tracks of this jet
422  AliTrack* t; 
423  for (Int_t it=1; it<=fNtrk; it++)
424  {
425   t=GetTrack(it);
426   if (t)
427   {
428    cout << "  ---Track no. " << it << endl;
429    cout << " ";
430    t->ListAll(f,u); 
431   }
432   else
433   {
434    cout << " *AliJet::List* Error : No track present." << endl; 
435   }
436  }
437
438 ///////////////////////////////////////////////////////////////////////////
439 Int_t AliJet::GetNtracks(Int_t idmode,Int_t chmode,Int_t pcode)
440 {
441 // Provide the number of user selected tracks in this jet based on the
442 // idmode, chmode and pcode selections as specified by the user.
443 // For specification of the selection parameters see GetTracks().
444 // The default parameters correspond to no selection, which implies
445 // that invokation of GetNtracks() just returns the total number of
446 // tracks registered in this jet.
447 //
448 // Note : In case certain selections are specified, this function
449 //        invokes GetTracks(idmode,chmode,pcode) to determine the
450 //        number of tracks corresponding to the selections.
451 //        When the jet contains a large number of tracks, invokation
452 //        of GetTracks(idmode,chmode,pcode) and subsequently invoking
453 //        GetEntries() for the resulting TObjArray* might be slightly
454 //        faster.
455
456  Int_t n=0;
457  if (idmode==0 && chmode==2 && pcode==0)
458  {
459   return fNtrk;
460  }
461  else
462  {
463   TObjArray* arr=GetTracks(idmode,chmode,pcode);
464   if (arr) n=arr->GetEntries();
465   return n;
466  }
467 }
468 ///////////////////////////////////////////////////////////////////////////
469 Int_t AliJet::GetNtracks(TString name)
470 {
471 // Provide the number of tracks with the specified name.
472 //
473 // Note :
474 // ------
475 // This facility invokes the corresponding GetTracks memberfunction
476 // and as such may result in overwriting existing track selection
477 // arrays. Please refer to the docs of GetTracks for further details.
478   
479  Int_t n=0;
480
481  TObjArray* arr=GetTracks(name);
482  if (arr) n=arr->GetEntries();
483  return n;
484 }
485 ///////////////////////////////////////////////////////////////////////////
486 Double_t AliJet::GetEnergy(Float_t scale)
487 {
488 // Return the total energy of the jet.
489 // By default the energy is returned in the units as it was stored in the jet
490 // structure. However, the user can select a different energy unit scale by
491 // specification of the scale parameter.
492 // The convention is that scale=1 corresponds to GeV, so specification
493 // of scale=0.001 will provide the energy in MeV.
494 // The error can be obtained by invoking GetResultError() after
495 // invokation of GetEnergy().
496  Double_t E=GetScalar();
497  if (E>0)
498  {
499   if (scale>0)
500   {
501    E*=fEscale/scale;
502    fDresult*=fEscale/scale;
503   }
504   return E;
505  }
506  else
507  {
508   return 0;
509  }
510 }
511 ///////////////////////////////////////////////////////////////////////////
512 Double_t AliJet::GetMomentum(Float_t scale)
513 {
514 // Return the value of the total jet 3-momentum
515 // By default the momentum is returned in the units as it was stored in the jet
516 // structure. However, the user can select a different momentum unit scale by
517 // specification of the scale parameter.
518 // The convention is that scale=1 corresponds to GeV/c, so specification
519 // of scale=0.001 will provide the momentum in MeV/c.
520 // The error can be obtained by invoking GetResultError() after
521 // invokation of GetMomentum().
522
523  Double_t norm=fV.GetNorm();
524  fDresult=fV.GetResultError();
525  if (scale>0)
526  {
527   norm*=fEscale/scale;
528   fDresult*=fEscale/scale;
529  }
530  return norm;
531 }
532 ///////////////////////////////////////////////////////////////////////////
533 Ali3Vector AliJet::Get3Momentum(Float_t scale) const
534 {
535 // Return the the total jet 3-momentum
536 // By default the components of the 3-momentum are returned in the units
537 // as they were stored in the jet structure.
538 // However, the user can select a different momentum unit scale for the
539 // components by specification of the scale parameter.
540 // The convention is that scale=1 corresponds to GeV/c, so specification
541 // of scale=0.001 will provide the 3-momentum in MeV/c.
542
543  Ali3Vector p=Get3Vector();
544  if (scale>0) p*=fEscale/scale;
545  return p;
546 }
547 ///////////////////////////////////////////////////////////////////////////
548 Double_t AliJet::GetInvmass(Float_t scale)
549 {
550 // Return the invariant mass of the jet.
551 // By default the mass is returned in the units as it was stored in the jet
552 // structure. However, the user can select a different mass unit scale by
553 // specification of the scale parameter.
554 // The convention is that scale=1 corresponds to GeV/c**2, so specification
555 // of scale=0.001 will provide the invariant mass in MeV/c**2.
556 // The error can be obtained by invoking GetResultError() after
557 // invokation of GetInvmass().
558
559  Double_t inv=Dot(*this);
560  Double_t dinv=GetResultError();
561  Double_t dm=0;
562  if (inv >= 0)
563  {
564  Double_t m=sqrt(inv);
565  if (m) dm=dinv/(2.*m);
566  if (scale>0)
567  {
568   m*=fEscale/scale;
569   dm*=fEscale/scale;
570  }
571  fDresult=dm;
572  return m;
573  }
574  else
575  {
576   fDresult=dm;
577   return 0;
578  }
579 }
580 ///////////////////////////////////////////////////////////////////////////
581 Float_t AliJet::GetCharge() const
582 {
583 // Return the total charge of the jet
584  return fQ;
585 }
586 ///////////////////////////////////////////////////////////////////////////
587 AliTrack* AliJet::GetTrack(Int_t i) const
588 {
589 // Return the i-th track of this jet
590
591  if (!fTracks) return 0;
592
593  if (i<=0 || i>fNtrk)
594  {
595   cout << " *AliJet*::GetTrack* Invalid argument i : " << i
596        << " Ntrk = " << fNtrk << endl;
597   return 0;
598  }
599  else
600  {
601   return (AliTrack*)fTracks->At(i-1);
602  }
603 }
604 ///////////////////////////////////////////////////////////////////////////
605 AliTrack* AliJet::GetIdTrack(Int_t id) const
606 {
607 // Return the track with user identifier "id" of this jet
608  if (!fTracks) return 0;
609
610  AliTrack* tx=0;
611  for (Int_t i=0; i<fNtrk; i++)
612  {
613   tx=(AliTrack*)fTracks->At(i);
614   if (id == tx->GetId()) return tx;
615  }
616  return 0; // No matching id found
617 }
618 ///////////////////////////////////////////////////////////////////////////
619 TObjArray* AliJet::GetTracks(Int_t idmode,Int_t chmode,Int_t pcode)
620 {
621 // Provide references to user selected tracks based on the idmode, chmode
622 // and pcode selections as specified by the user.
623 //
624 // The following selection combinations are available :
625 // ----------------------------------------------------
626 // idmode = -1 ==> Select tracks with negative user identifier "id"
627 //           0 ==> No selection on user identifier
628 //           1 ==> Select tracks with positive user identifier "id"
629 //
630 // chmode = -1 ==> Select tracks with negative charge
631 //           0 ==> Select neutral tracks
632 //           1 ==> Select tracks with positive charge
633 //           2 ==> No selection on charge
634 //           3 ==> Select all charged tracks
635 //
636 // pcode  =  0 ==> No selection on particle code
637 //           X ==> Select tracks with particle code +X or -X
638 //                 This allows selection of both particles and anti-particles
639 //                 in case of PDG particle codes.
640 //                 Selection of either particles or anti-particles can be
641 //                 obtained in combination with the "chmode" selector.
642 //
643 // Examples :
644 // ----------
645 // idmode=-1 chmode=0 pcode=0   : Selection of all neutral tracks with negative id.
646 // idmode=0  chmode=2 pcode=211 : Selection of all charged pions (PDG convention).
647 // idmode=0  chmode=1 pcode=321 : Selection of all positive kaons (PDG convention).
648 //
649 // The default values are idmode=0 chmode=2 pcode=0 (i.e. no selections applied).
650 //
651 // Notes :
652 // -------
653 // 1) In case the user has labeled simulated tracks with negative id and
654 //    reconstructed tracks with positive id, this memberfunction provides
655 //    easy access to either all simulated or reconstructed tracks.
656 // 2) Subsequent invokations of this memberfunction with e.g. chmode=-1 and chmode=1
657 //    provides a convenient way to investigate particle pairs with opposite charge
658 //    (e.g. for invariant mass analysis).
659 // 3) The selected track pointers are returned via a multi-purpose array,
660 //    which will be overwritten by subsequent selections.
661 //    In case the selected track list is to be used amongst other selections,
662 //    the user is advised to store the selected track pointers in a local
663 //    TObjArray or TRefArray.  
664
665  if (fSelected)
666  {
667   fSelected->Clear();
668  }
669  else
670  {
671   fSelected=new TObjArray();
672  }
673
674  if (!fTracks) return fSelected;
675
676  AliTrack* tx=0;
677  Int_t code=0;
678  Int_t id=0;
679  Float_t q=0;
680  for (Int_t i=0; i<fNtrk; i++)
681  {
682   tx=(AliTrack*)fTracks->At(i);
683   if (!tx) continue;
684
685   code=tx->GetParticleCode();
686   if (pcode && abs(pcode)!=abs(code)) continue;
687
688   id=tx->GetId();
689   if (idmode==-1 && id>=0) continue;
690   if (idmode==1 && id<=0) continue;
691
692   q=tx->GetCharge();
693   if (chmode==-1 && q>=0) continue;
694   if (chmode==0 && fabs(q)>1e-10) continue;
695   if (chmode==1 && q<=0) continue;
696   if (chmode==3 && fabs(q)<1e-10) continue;
697
698   fSelected->Add(tx);
699  }
700
701  return fSelected;
702 }
703 ///////////////////////////////////////////////////////////////////////////
704 TObjArray* AliJet::GetTracks(TString name)
705 {
706 // Provide references to all tracks with the specified name.
707 //
708 // Notes :
709 // -------
710 // 1) In case the user has labeled reconstructed tracks with the name of
711 //    the applied reconstruction algorithm, this memberfunction provides
712 //    easy access to all tracks reconstructed by a certain method.
713 // 2) The selected track pointers are returned via a multi-purpose array,
714 //    which will be overwritten by subsequent selections.
715 //    In case the selected track list is to be used amongst other selections,
716 //    the user is advised to store the selected track pointers in a local
717 //    TObjArray or TRefArray.  
718
719  if (fSelected)
720  {
721   fSelected->Clear();
722  }
723  else
724  {
725   fSelected=new TObjArray();
726  }
727
728  if (!fTracks) return fSelected;
729
730  AliTrack* tx=0;
731  TString s;
732  for (Int_t i=0; i<fNtrk; i++)
733  {
734   tx=(AliTrack*)fTracks->At(i);
735   if (!tx) continue;
736
737   s=tx->GetName();
738   if (s == name) fSelected->Add(tx);
739  }
740
741  return fSelected;
742 }
743 ///////////////////////////////////////////////////////////////////////////
744 void AliJet::RemoveTracks(TString name)
745 {
746 // Remove all tracks with the specified name.
747 // If name="*" all tracks will be removed.
748 //
749 // Note :
750 // ------
751 // In case the user has labeled reconstructed tracks with the name of
752 // the applied reconstruction algorithm, this memberfunction provides
753 // easy removal of all tracks reconstructed by a certain method.
754
755  if (!fTracks) return;
756
757  AliTrack* tx=0;
758  TString s;
759  TObject* obj=0;
760  for (Int_t i=0; i<fNtrk; i++)
761  {
762   tx=(AliTrack*)fTracks->At(i);
763   if (!tx) continue;
764
765   s=tx->GetName();
766   if (s==name || name=="*")
767   {
768    obj=fTracks->Remove(tx);
769    if (obj && fTracks->IsOwner()) delete tx;
770   }
771  }
772  fTracks->Compress();
773  fNtrk=fTracks->GetEntries();
774 }
775 ///////////////////////////////////////////////////////////////////////////
776 void AliJet::RemoveTracks(Int_t idmode,Int_t chmode,Int_t pcode)
777 {
778 // Remove user selected tracks based on the idmode, chmode and pcode
779 // selections as specified by the user.
780 // For defintions of these selections see the corresponding GetTracks()
781 // memberfunction.
782
783  if (!fTracks) return;
784
785  TObjArray* arr=GetTracks(idmode,chmode,pcode);
786  if (!arr) return;
787  
788  Int_t ntk=arr->GetEntries();
789  if (!ntk) return;
790
791  AliTrack* tx=0;
792  TObject* obj=0;
793  for (Int_t i=0; i<ntk; i++)
794  {
795   tx=(AliTrack*)arr->At(i);
796   if (!tx) continue;
797
798   obj=fTracks->Remove(tx);
799   if (obj && fTracks->IsOwner()) delete tx;
800  }
801  fTracks->Compress();
802  fNtrk=fTracks->GetEntries();
803  arr->Clear();
804 }
805 ///////////////////////////////////////////////////////////////////////////
806 void AliJet::ShowTracks(Int_t mode)
807 {
808 // Provide an overview of the available tracks.
809 // The argument mode determines the amount of information as follows :
810 // mode = 0 ==> Only printout of the number of tracks
811 //        1 ==> Provide a listing with 1 line of info for each track
812 //
813 // The default is mode=1.
814 //
815  Int_t ntk=GetNtracks();
816  if (ntk)
817  {
818   if (!mode)
819   {
820    cout << " There are " << ntk << " tracks available." << endl; 
821   }
822   else
823   {
824    cout << " The following " << ntk << " tracks are available :" << endl; 
825    for (Int_t i=1; i<=ntk; i++)
826    {
827     AliTrack* tx=GetTrack(i);
828     if (tx)
829     {
830      const char* name=tx->GetName();
831      const char* title=tx->GetTitle();
832      cout << " Track : " << i;
833      cout << " Id : " << tx->GetId();
834      cout << " Q : " << tx->GetCharge() << " m : " << tx->GetMass() << " p : " << tx->GetMomentum();
835      if (strlen(name)) cout << " Name : " << name;
836      if (strlen(title)) cout << " Title : " << title;
837      cout << endl;
838     }
839    }
840   }
841  }
842  else
843  {
844   cout << " No tracks are present." << endl;
845  }
846 }
847 ///////////////////////////////////////////////////////////////////////////
848 Double_t AliJet::GetPt(Float_t scale)
849 {
850 // Provide the transverse momentum value w.r.t. z-axis.
851 // By default the value is returned in the units as it was stored in the jet
852 // structure. However, the user can select a different momentum unit scale by
853 // specification of the scale parameter.
854 // The convention is that scale=1 corresponds to GeV/c, so specification
855 // of scale=0.001 will provide the transverse momentum in MeV/c.
856 // The error on the value can be obtained by GetResultError()
857 // after invokation of GetPt().
858  Ali3Vector v;
859  v=GetVecTrans();
860  Double_t norm=v.GetNorm();
861  fDresult=v.GetResultError();
862  if (scale>0)
863  {
864   norm*=fEscale/scale;
865   fDresult*=fEscale/scale;
866  }
867
868  return norm;
869 }
870 ///////////////////////////////////////////////////////////////////////////
871 Double_t AliJet::GetPl(Float_t scale)
872 {
873 // Provide the longitudinal momentum value w.r.t. z-axis.
874 // By default the value is returned in the units as it was stored in the jet
875 // structure. However, the user can select a different momentum unit scale by
876 // specification of the scale parameter.
877 // The convention is that scale=1 corresponds to GeV/c, so specification
878 // of scale=0.001 will provide the longitudinal momentum in MeV/c.
879 // Note : the returned value can also be negative.
880 // The error on the value can be obtained by GetResultError()
881 // after invokation of GetPl().
882
883  Ali3Vector v;
884  v=GetVecLong();
885
886  Double_t pl=v.GetNorm();
887  fDresult=v.GetResultError();
888
889  Double_t a[3];
890  v.GetVector(a,"sph");
891  if (cos(a[1])<0) pl=-pl;
892  if (scale>0)
893  {
894   pl*=fEscale/scale;
895   fDresult*=fEscale/scale;
896  }
897
898  return pl;
899 }
900 ///////////////////////////////////////////////////////////////////////////
901 Double_t AliJet::GetEt(Float_t scale)
902 {
903 // Provide transverse energy value w.r.t. z-axis.
904 // By default the value is returned in the units as it was stored in the jet
905 // structure. However, the user can select a different energy unit scale by
906 // specification of the scale parameter.
907 // The convention is that scale=1 corresponds to GeV, so specification
908 // of scale=0.001 will provide the transverse energy in MeV.
909 // The error on the value can be obtained by GetResultError()
910 // after invokation of GetEt().
911
912  Double_t et=GetScaTrans();
913  if (scale>0)
914  {
915   et*=fEscale/scale;
916   fDresult*=fEscale/scale;
917  }
918
919  return et;
920 }
921 ///////////////////////////////////////////////////////////////////////////
922 Double_t AliJet::GetEl(Float_t scale)
923 {
924 // Provide longitudinal energy value w.r.t. z-axis.
925 // By default the value is returned in the units as it was stored in the jet
926 // structure. However, the user can select a different energy unit scale by
927 // specification of the scale parameter.
928 // The convention is that scale=1 corresponds to GeV, so specification
929 // of scale=0.001 will provide the longitudinal energy in MeV.
930 // Note : the returned value can also be negative.
931 // The error on the value can be obtained by GetResultError()
932 // after invokation of GetEl().
933
934  Double_t el=GetScaLong();
935  if (scale>0)
936  {
937   el*=fEscale/scale;
938   fDresult*=fEscale/scale;
939  }
940
941  return el;
942 }
943 ///////////////////////////////////////////////////////////////////////////
944 Double_t AliJet::GetMt(Float_t scale)
945 {
946 // Provide transverse mass value w.r.t. z-axis.
947 // By default the value is returned in the units as it was stored in the jet
948 // structure. However, the user can select a different energy unit scale by
949 // specification of the scale parameter.
950 // The convention is that scale=1 corresponds to GeV, so specification
951 // of scale=0.001 will provide the transverse mass in MeV.
952 // The error on the value can be obtained by GetResultError()
953 // after invokation of GetMt().
954  Double_t pt=GetPt();
955  Double_t dpt=GetResultError();
956  Double_t m=GetInvmass();
957  Double_t dm=GetResultError();
958
959  Double_t mt=sqrt(pt*pt+m*m);
960  Double_t dmt2=0;
961  if (mt) dmt2=(pow((pt*dpt),2)+pow((m*dm),2))/(mt*mt);
962
963  fDresult=sqrt(dmt2);
964  if (scale>0)
965  {
966   mt*=fEscale/scale;
967   fDresult*=fEscale/scale;
968  }
969  return mt;
970 }
971 ///////////////////////////////////////////////////////////////////////////
972 Double_t AliJet::GetRapidity()
973 {
974 // Provide rapidity value w.r.t. z-axis.
975 // The error on the value can be obtained by GetResultError()
976 // after invokation of GetRapidity().
977 // Note : Also GetPseudoRapidity() is available since this class is
978 //        derived from Ali4Vector.
979  Double_t e=GetEnergy();
980  Double_t de=GetResultError();
981  Double_t pl=GetPl();
982  Double_t dpl=GetResultError();
983  Double_t sum=e+pl;
984  Double_t dif=e-pl;
985
986  Double_t y=9999,dy2=0;
987  if (sum && dif) y=0.5*log(sum/dif);
988
989  if (sum*dif) dy2=(1./(sum*dif))*(pow((pl*de),2)+pow((e*dpl),2));
990
991  fDresult=sqrt(dy2);
992  return y;
993 }
994 ///////////////////////////////////////////////////////////////////////////
995 void AliJet::SetTrackCopy(Int_t j)
996 {
997 // (De)activate the creation of private copies of the added tracks.
998 // j=0 ==> No private copies are made; pointers of original tracks are stored.
999 // j=1 ==> Private copies of the tracks are made and these pointers are stored.
1000 //
1001 // Note : Once the storage contains pointer(s) to AliTrack(s) one cannot
1002 //        change the TrackCopy mode anymore.
1003 //        To change the TrackCopy mode for an existing AliJet containing
1004 //        tracks one first has to invoke Reset().
1005  if (!fTracks)
1006  {
1007   if (j==0 || j==1)
1008   {
1009    fTrackCopy=j;
1010   }
1011   else
1012   {
1013    cout << "*AliJet::SetTrackCopy* Invalid argument : " << j << endl;
1014   }
1015  }
1016  else
1017  {
1018   cout << "*AliJet::SetTrackCopy* Storage already contained tracks."
1019        << "  ==> TrackCopy mode not changed." << endl; 
1020  }
1021 }
1022 ///////////////////////////////////////////////////////////////////////////
1023 Int_t AliJet::GetTrackCopy() const
1024 {
1025 // Provide value of the TrackCopy mode.
1026 // 0 ==> No private copies are made; pointers of original tracks are stored.
1027 // 1 ==> Private copies of the tracks are made and these pointers are stored.
1028  return fTrackCopy;
1029 }
1030 ///////////////////////////////////////////////////////////////////////////
1031 void AliJet::SetId(Int_t id)
1032 {
1033 // Set a user defined identifier for this jet.
1034  fUserId=id;
1035 }
1036 ///////////////////////////////////////////////////////////////////////////
1037 Int_t AliJet::GetId() const
1038 {
1039 // Provide the user defined identifier of this jet.
1040  return fUserId;
1041 }
1042 ///////////////////////////////////////////////////////////////////////////
1043 void AliJet::SetReferencePoint(AliPosition& p)
1044 {
1045 // Store the position of the jet reference-point.
1046 // The reference-point of a jet provides a means to define a generic
1047 // space-time location for the jet as a whole.
1048 // This doesn't have to be necessarily the location where all the constituent
1049 // tracks originate (e.g. a bundle of parallel tracks doesn't have such
1050 // a location). As such the meaning of this reference-point is different from
1051 // a normal vertex position and allows to provide complimentary information. 
1052 // This reference point is the preferable point to start e.g. extrapolations
1053 // and investigate coincidences in space and/or time.
1054  if (fRef) delete fRef;
1055  fRef=new AliPositionObj(p);
1056 }
1057 ///////////////////////////////////////////////////////////////////////////
1058 AliPosition* AliJet::GetReferencePoint()
1059 {
1060 // Provide the position of the jet reference-point.
1061 // The reference-point of a jet provides a means to define a generic
1062 // space-time location for the jet as a whole.
1063 // This doesn't have to be necessarily the location where all the constituent
1064 // tracks originate (e.g. a bundle of parallel tracks doesn't have such
1065 // a location). As such the meaning of this reference-point is different from
1066 // a normal vertex position and allows to provide complimentary information. 
1067 // This reference point is the preferable point to start e.g. extrapolations
1068 // and investigate coincidences in space and/or time.
1069  return fRef;
1070 }
1071 ///////////////////////////////////////////////////////////////////////////
1072 TObjArray* AliJet::SortTracks(Int_t mode,TObjArray* tracks)
1073 {
1074 // Order the references to an array of tracks by looping over the input array "tracks"
1075 // and checking the value of a certain observable.
1076 // The ordered array is returned as a TObjArray.
1077 // In case tracks=0 (default), the registered tracks of the current jet are used. 
1078 // Note that the original track array is not modified.
1079 // Via the "mode" argument the user can specify the observable to be checked upon
1080 // and specify whether sorting should be performed in decreasing order (mode<0)
1081 // or in increasing order (mode>0).
1082 //
1083 // The convention for the observable selection is the following :
1084 // mode : 1 ==> Number of signals associated to the track
1085 //        2 ==> Track energy
1086 //        3 ==> Track momentum
1087 //        4 ==> Mass of the track
1088 //        5 ==> Transverse momentum of the track
1089 //        6 ==> Longitudinal momentum of the track
1090 //        7 ==> Transverse energy of the track
1091 //        8 ==> Longitudinal energy of the track
1092 //        9 ==> Transverse mass of the track
1093 //       10 ==> Track rapidity
1094 //       11 ==> Pseudo-rapidity of the track
1095 //       12 ==> Charge of the track
1096 //       13 ==> Probability of the track hypothesis
1097 //
1098 // The default is mode=-1.
1099 //
1100 // Note : This sorting routine uses a common area in memory, which is used
1101 //        by various other sorting facilities as well.
1102 //        This means that the resulting sorted TObjArray may be overwritten
1103 //        when another sorting is invoked.
1104 //        To retain the sorted list of pointers, the user is advised to copy
1105 //        the pointers contained in the returned TObjArray into a private
1106 //        TObjArray instance.
1107
1108  if (fSelected)
1109  {
1110   delete fSelected;
1111   fSelected=0;
1112  }
1113
1114  if (!tracks) tracks=fTracks;
1115  
1116  if (!mode || abs(mode)>13 || !tracks) return fSelected;
1117
1118  Int_t ntracks=tracks->GetEntries();
1119  if (!ntracks)
1120  {
1121   return fSelected;
1122  }
1123  else
1124  {
1125   fSelected=new TObjArray(ntracks);
1126  }
1127
1128  Double_t val1,val2; // Values of the observable to be tested upon
1129  
1130  Int_t nord=0;
1131  for (Int_t i=0; i<ntracks; i++) // Loop over all tracks of the array
1132  {
1133   AliTrack* tx=(AliTrack*)tracks->At(i);
1134
1135   if (!tx) continue;
1136  
1137   if (nord == 0) // store the first track at the first ordered position
1138   {
1139    nord++;
1140    fSelected->AddAt(tx,nord-1);
1141    continue;
1142   }
1143  
1144   for (Int_t j=0; j<=nord; j++) // put track in the right ordered position
1145   {
1146    if (j == nord) // track has smallest (mode<0) or largest (mode>0) observable value seen so far
1147    {
1148     nord++;
1149     fSelected->AddAt(tx,j); // add track at the end
1150     break; // go for next track
1151    }
1152
1153    val1=0;
1154    val2=0;
1155
1156    switch (abs(mode))
1157    {
1158     case 1:
1159      val1=tx->GetNsignals();
1160      val2=((AliTrack*)fSelected->At(j))->GetNsignals();
1161      break;
1162     case 2:
1163      val1=tx->GetEnergy(1);
1164      val2=((AliTrack*)fSelected->At(j))->GetEnergy(1);
1165      break;
1166     case 3:
1167      val1=tx->GetMomentum(1);
1168      val2=((AliTrack*)fSelected->At(j))->GetMomentum(1);
1169      break;
1170     case 4:
1171      val1=tx->GetMass(1);
1172      val2=((AliTrack*)fSelected->At(j))->GetMass(1);
1173      break;
1174     case 5:
1175      val1=tx->GetPt(1);
1176      val2=((AliTrack*)fSelected->At(j))->GetPt(1);
1177      break;
1178     case 6:
1179      val1=tx->GetPl(1);
1180      val2=((AliTrack*)fSelected->At(j))->GetPl(1);
1181      break;
1182     case 7:
1183      val1=tx->GetEt(1);
1184      val2=((AliTrack*)fSelected->At(j))->GetEt(1);
1185      break;
1186     case 8:
1187      val1=tx->GetEl(1);
1188      val2=((AliTrack*)fSelected->At(j))->GetEl(1);
1189      break;
1190     case 9:
1191      val1=tx->GetMt(1);
1192      val2=((AliTrack*)fSelected->At(j))->GetMt(1);
1193      break;
1194     case 10:
1195      val1=tx->GetRapidity();
1196      val2=((AliTrack*)fSelected->At(j))->GetRapidity();
1197      break;
1198     case 11:
1199      val1=tx->GetPseudoRapidity();
1200      val2=((AliTrack*)fSelected->At(j))->GetPseudoRapidity();
1201      break;
1202     case 12:
1203      val1=tx->GetCharge();
1204      val2=((AliTrack*)fSelected->At(j))->GetCharge();
1205      break;
1206     case 13:
1207      val1=tx->GetProb();
1208      val2=((AliTrack*)fSelected->At(j))->GetProb();
1209      break;
1210    }
1211
1212    if (mode<0 && val1 <= val2) continue;
1213    if (mode>0 && val1 >= val2) continue;
1214  
1215    nord++;
1216    for (Int_t k=nord-1; k>j; k--) // create empty position
1217    {
1218     fSelected->AddAt(fSelected->At(k-1),k);
1219    }
1220    fSelected->AddAt(tx,j); // put track at empty position
1221    break; // go for next track
1222   }
1223  }
1224  return fSelected;
1225 }
1226 ///////////////////////////////////////////////////////////////////////////
1227 Double_t AliJet::GetDistance(AliPosition* p,Float_t scale)
1228 {
1229 // Provide distance of the current jet to the position p.
1230 // The error on the result can be obtained as usual by invoking
1231 // GetResultError() afterwards. 
1232 //
1233 // By default the distance will be provided in the metric unit scale of
1234 // the AliPosition p.
1235 // However, the user can select a different metric unit scale by
1236 // specification of the scale parameter.
1237 // The convention is that scale=1 corresponds to meter, so specification
1238 // of scale=0.01 will provide the distance in cm.
1239 // As such it is possible to obtain a correctly computed distance even in case
1240 // the jet parameters have a different unit scale.
1241 // However, it is recommended to work always with one single unit scale.
1242 //
1243 // Note : In case of incomplete information, a distance value of -1 is
1244 //        returned.
1245  
1246  Double_t dist=-1.;
1247  fDresult=0.;
1248
1249  if (!p) return dist;
1250
1251  // Obtain a defined position on this jet
1252  AliPosition* rx=fRef;
1253
1254  if (!rx) return dist;
1255
1256  Ali3Vector pj=Get3Momentum();
1257
1258  if (pj.GetNorm() <= 0.) return dist;
1259
1260  AliTrack tj;
1261  tj.Set3Momentum(pj);
1262  tj.SetReferencePoint(*rx);
1263  dist=tj.GetDistance(p,scale);
1264  fDresult=tj.GetResultError();
1265  return dist;
1266 }
1267 ///////////////////////////////////////////////////////////////////////////
1268 Double_t AliJet::GetDistance(AliTrack* t,Float_t scale)
1269 {
1270 // Provide distance of the current jet to the track t.
1271 // The error on the result can be obtained as usual by invoking
1272 // GetResultError() afterwards. 
1273 //
1274 // By default the distance will be provided in the metric unit scale of
1275 // the current jet.
1276 // However, the user can specify a required metric unit scale by specification
1277 // of the scale parameter.
1278 // The convention is that scale=1 corresponds to meter, so specification
1279 // of scale=0.01 will provide the distance in cm.
1280 // As such it is possible to obtain a correctly computed distance even in case
1281 // the jet and track parameters have a different unit scale.
1282 // However, it is recommended to work always with one single unit scale.
1283 //
1284 // Note : In case of incomplete information, a distance value of -1 is
1285 //        returned.
1286  
1287  Double_t dist=-1.;
1288  fDresult=0.;
1289
1290  if (!t) return dist;
1291
1292  // Obtain a defined position on this jet
1293  AliPosition* rx=fRef;
1294
1295  if (!rx) return dist;
1296
1297  Ali3Vector pj=Get3Momentum();
1298
1299  if (pj.GetNorm() <= 0.) return dist;
1300
1301  AliTrack tj;
1302  tj.Set3Momentum(pj);
1303  tj.SetReferencePoint(*rx);
1304  dist=tj.GetDistance(t,scale);
1305  fDresult=tj.GetResultError();
1306  return dist;
1307 }
1308 ///////////////////////////////////////////////////////////////////////////
1309 Double_t AliJet::GetDistance(AliJet* j,Float_t scale)
1310 {
1311 // Provide distance of the current jet to the jet j.
1312 // The error on the result can be obtained as usual by invoking
1313 // GetResultError() afterwards. 
1314 //
1315 // By default the distance will be provided in the metric unit scale of
1316 // the current jet.
1317 // This implies that the results of j1.GetDistance(j2) and j2.GetDistance(j1)
1318 // may be numerically different in case j1 and j2 have different metric units.
1319 // However, the user can specify a required metric unit scale by specification
1320 // of the scale parameter.
1321 // The convention is that scale=1 corresponds to meter, so specification
1322 // of scale=0.01 will provide the distance in cm.
1323 // As such it is possible to obtain a correctly computed distance even in case
1324 // the jet parameters have a different unit scale.
1325 // However, it is recommended to work always with one single unit scale.
1326 //
1327 // Note : In case of incomplete information, a distance value of -1 is
1328 //        returned.
1329  
1330  Double_t dist=-1.;
1331  fDresult=0.;
1332
1333  if (!j) return dist;
1334
1335  // Obtain a defined position on jet j
1336  AliPosition* rx=j->GetReferencePoint();
1337
1338  if (!rx) return dist;
1339
1340  Ali3Vector pj=j->Get3Momentum();
1341
1342  if (pj.GetNorm() <= 0.) return dist;
1343
1344  AliTrack tj;
1345  tj.Set3Momentum(pj);
1346  tj.SetReferencePoint(*rx);
1347  dist=GetDistance(tj,scale);
1348  return dist;
1349 }
1350 ///////////////////////////////////////////////////////////////////////////
1351 Int_t AliJet::GetNsignals() const
1352 {
1353 // Provide the number of signals associated to the jet tracks.
1354 // Note : Multiple occurrences of the same signal are only counted once. 
1355
1356  if (fNtrk<1) return 0;
1357
1358  TObjArray arr;
1359
1360  Int_t n=0;
1361  AliTrack* tx=0;
1362  Int_t exists=0;
1363  for (Int_t i=1; i<=fNtrk; i++)
1364  {
1365   tx=GetTrack(i);
1366   for (Int_t j=1; j<=tx->GetNsignals(); j++)
1367   {
1368    AliSignal* sx=tx->GetSignal(j);
1369    if (!sx) continue;
1370    exists=0;
1371    for (Int_t k=0; k<arr.GetEntries(); k++)
1372    {
1373     if (sx==(AliSignal*)arr.At(k))
1374     {
1375      exists=1;
1376      break;
1377     } 
1378    }
1379    if (!exists) arr.Add(sx);
1380   } 
1381  }
1382  n=arr.GetEntries();
1383  return n;
1384 }
1385 ///////////////////////////////////////////////////////////////////////////
1386 void AliJet::SetEscale(Float_t scale)
1387 {
1388 // Indicate the energy/momentum scale as used by the user.
1389 // The convention is that scale=1 indicates values in units
1390 // of GeV, GeV/c or GeV/c**2.
1391 // So, in case one decides to store values in units of MeV, MeV/c or MeV/c**2
1392 // the scale indicator should be set to scale=0.001.
1393 //
1394 // By default scale=1 is set in the constructor.
1395
1396  if (scale>0)
1397  {
1398   fEscale=scale;
1399  }
1400  else
1401  {
1402   cout << " *AliJet::SetEscale* Invalid scale value : " << scale << endl;
1403  }
1404 }
1405 ///////////////////////////////////////////////////////////////////////////
1406 Float_t AliJet::GetEscale() const
1407 {
1408 // Provide the energy/momentum scale as used by the user.
1409 // The convention is that scale=1 indicates values in units
1410 // of GeV, GeV/c or GeV/c**2.
1411 // So, a value of scale=0.001 indicates that energy/momentum values are
1412 // stored in units of MeV, MeV/c or MeV/c**2.
1413  return fEscale;
1414 }
1415 ///////////////////////////////////////////////////////////////////////////
1416 TObject* AliJet::Clone(const char* name) const
1417 {
1418 // Make a deep copy of the current object and provide the pointer to the copy.
1419 // This memberfunction enables automatic creation of new objects of the
1420 // correct type depending on the object type, a feature which may be very useful
1421 // for containers when adding objects in case the container owns the objects.
1422 // This feature allows e.g. AliVertex to store either AliJet objects or
1423 // objects derived from AliJet via the AddJet memberfunction, provided
1424 // these derived classes also have a proper Clone memberfunction. 
1425
1426  AliJet* jet=new AliJet(*this);
1427  if (name)
1428  {
1429   if (strlen(name)) jet->SetName(name);
1430  }
1431  return jet;
1432 }
1433 ///////////////////////////////////////////////////////////////////////////