1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 ///////////////////////////////////////////////////////////////////////////
20 // Creation and investigation of a jet of particle tracks.
21 // An AliJet can be constructed by adding AliTracks.
23 // To provide maximal flexibility to the user, two modes of track storage
24 // are provided by means of the memberfunction SetTrackCopy().
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
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.
35 // Modifications made to the original tracks also affect the AliTrack objects
36 // which are stored in the AliJet.
38 // b) SetTrackCopy(1).
39 // Of every 'added' track a private copy will be made of which the pointer
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.
48 // See also the documentation provided for the memberfunction SetOwner().
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
57 // ... // code to fill the AliTrack data
64 // j2.SetTrackCopy(1);
65 // AliTrack* tx=new AliTrack();
66 // for (Int_t i=0; i<10; i++)
69 // ... // code to set momentum etc... of the track tx
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);
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.
93 //--- Author: Nick van Eijndhoven 10-jul-1997 UU-SAP Utrecht
94 //- Modified: NvE $Date$ UU-SAP Utrecht
95 ///////////////////////////////////////////////////////////////////////////
98 #include "Riostream.h"
100 ClassImp(AliJet) // Class implementation to enable ROOT I/O
102 AliJet::AliJet() : TNamed(),Ali4Vector()
104 // Default constructor
105 // All variables initialised to 0
106 // Initial maximum number of tracks is set to the default value
111 ///////////////////////////////////////////////////////////////////////////
114 // Initialisation of pointers etc...
122 ///////////////////////////////////////////////////////////////////////////
123 AliJet::AliJet(Int_t n) : TNamed(),Ali4Vector()
125 // Create a jet to hold initially a maximum of n tracks
126 // All variables initialised to 0
136 cout << " *AliJet* Initial max. number of tracks entered : " << n << endl;
137 cout << " This is invalid. Default initial maximum will be used." << endl;
142 ///////////////////////////////////////////////////////////////////////////
145 // Default destructor
162 ///////////////////////////////////////////////////////////////////////////
163 void AliJet::SetOwner(Bool_t own)
165 // Set ownership of all added objects.
166 // The default parameter is own=kTRUE.
168 // Invokation of this memberfunction also sets all the copy modes
169 // (e.g. TrackCopy & co.) according to the value of own.
171 // This function (with own=kTRUE) is particularly useful when reading data
172 // from a tree/file, since Reset() will then actually remove all the
173 // added objects from memory irrespective of the copy mode settings
174 // during the tree/file creation process. In this way it provides a nice way
175 // of preventing possible memory leaks in the reading/analysis process.
177 // In addition this memberfunction can also be used as a shortcut to set all
178 // copy modes in one go during a tree/file creation process.
179 // However, in this case the user has to take care to only set/change the
180 // ownership (and copy mode) for empty objects (e.g. newly created objects
181 // or after invokation of the Reset() memberfunction) otherwise it will
182 // very likely result in inconsistent destructor behaviour.
186 if (fTracks) fTracks->SetOwner(own);
189 ///////////////////////////////////////////////////////////////////////////
190 AliJet::AliJet(const AliJet& j) : TNamed(j),Ali4Vector(j)
198 fTrackCopy=j.fTrackCopy;
200 if (j.fRef) fRef=new AliPositionObj(*(j.fRef));
207 fTracks=new TObjArray(fNtmax);
208 if (fTrackCopy) fTracks->SetOwner();
211 for (Int_t i=1; i<=fNtrk; i++)
213 AliTrack* tx=j.GetTrack(i);
216 fTracks->Add(tx->Clone());
224 ///////////////////////////////////////////////////////////////////////////
225 void AliJet::SetNtinit(Int_t n)
227 // Set the initial maximum number of tracks for this jet
242 ///////////////////////////////////////////////////////////////////////////
245 // Reset all variables to 0
246 // The max. number of tracks is set to the initial value again
247 // Note : The scale for the energy/momentum units will not be changed.
251 Double_t a[4]={0,0,0,0};
253 if (fNtinit > 0) SetNtinit(fNtinit);
255 ///////////////////////////////////////////////////////////////////////////
256 void AliJet::AddTrack(AliTrack& t)
258 // Add a track to the jet.
259 // In case the maximum number of tracks has been reached
260 // space will be extended to hold an additional amount of tracks as
261 // was initially reserved.
262 // See SetTrackCopy() to tailor the functionality of the stored structures.
266 // In case a private copy is made, this is performed via the Clone() memberfunction.
267 // All AliTrack and derived classes have the default TObject::Clone() memberfunction.
268 // However, derived classes generally contain an internal data structure which may
269 // include pointers to other objects. Therefore it is recommended to provide
270 // for all derived classes a specific copy constructor and override the default Clone()
271 // memberfunction using this copy constructor.
272 // An example for this may be seen from AliTrack.
274 // In case NO private copy is made, a check will be performed if this
275 // specific track is already present in the jet.
276 // If this is the case, no action is performed to prevent multiple
277 // additions of the same track.
282 ///////////////////////////////////////////////////////////////////////////
283 void AliJet::AddTrack(AliTrack& t,Int_t copy)
285 // Internal memberfunction to actually add a track to the jet.
286 // In case the maximum number of tracks has been reached
287 // space will be extended to hold an additional amount of tracks as
288 // was initially reserved.
290 // If copy=0 NO copy of the track will be made, irrespective of the setting
291 // of the TrackCopy flag.
292 // This allows a proper treatment of automatically generated connecting
293 // tracks between vertices.
295 // In case NO copy of the track is made, a check will be performed if this
296 // specific track is already present in the jet.
297 // If this is the case, no action is performed to prevent multiple
298 // additions of the same track.
301 // In case a private copy is made, this is performed via the Clone() memberfunction.
305 fTracks=new TObjArray(fNtmax);
306 if (fTrackCopy) fTracks->SetOwner();
308 else if (!fTrackCopy || !copy) // Check if this track is already present
310 for (Int_t i=0; i<fNtrk; i++)
312 AliTrack* tx=(AliTrack*)fTracks->At(i);
313 if (tx == &t) return;
317 if (fNtrk == fNtmax) // Check if maximum track number is reached
320 fTracks->Expand(fNtmax);
323 // Add the track to this jet
325 if (fTrackCopy && copy)
327 fTracks->Add(t.Clone());
336 Ali4Vector p4=(Ali4Vector)t;
337 Float_t tscale=t.GetEscale();
338 if ((tscale/fEscale > 1.1) || (fEscale/tscale > 1.1)) p4=p4*(tscale/fEscale);
342 ///////////////////////////////////////////////////////////////////////////
343 void AliJet::Data(TString f,TString u)
345 // Provide jet information within the coordinate frame f
347 // The string argument "u" allows to choose between different angular units
348 // in case e.g. a spherical frame is selected.
349 // u = "rad" : angles provided in radians
350 // "deg" : angles provided in degrees
352 // The defaults are f="car" and u="rad".
354 const char* name=GetName();
355 const char* title=GetTitle();
357 cout << " *AliJet::Data*";
358 if (strlen(name)) cout << " Name : " << GetName();
359 if (strlen(title)) cout << " Title : " << GetTitle();
361 cout << " Id : " << fUserId << " Invmass : " << GetInvmass() << " Charge : " << fQ
362 << " Momentum : " << GetMomentum() << " Energy scale : " << fEscale << " GeV" << endl;
366 Ali4Vector::Data(f,u);
368 ///////////////////////////////////////////////////////////////////////////
369 void AliJet::List(TString f,TString u)
371 // Provide jet and primary track information within the coordinate frame f
373 // The string argument "u" allows to choose between different angular units
374 // in case e.g. a spherical frame is selected.
375 // u = "rad" : angles provided in radians
376 // "deg" : angles provided in degrees
378 // The defaults are f="car" and u="rad".
380 Data(f,u); // Information of the current jet
381 if (fRef) { cout << " Ref-point :"; fRef->Data(f,u); }
383 // The tracks of this jet
385 for (Int_t it=1; it<=fNtrk; it++)
390 cout << " ---Track no. " << it << endl;
396 cout << " *AliJet::List* Error : No track present." << endl;
400 ///////////////////////////////////////////////////////////////////////////
401 void AliJet::ListAll(TString f,TString u)
403 // Provide jet and prim.+sec. track information within the coordinate frame f
405 // The string argument "u" allows to choose between different angular units
406 // in case e.g. a spherical frame is selected.
407 // u = "rad" : angles provided in radians
408 // "deg" : angles provided in degrees
410 // The defaults are f="car" and u="rad".
412 Data(f,u); // Information of the current jet
413 if (fRef) { cout << " Ref-point :"; fRef->Data(f,u); }
415 // The tracks of this jet
417 for (Int_t it=1; it<=fNtrk; it++)
422 cout << " ---Track no. " << it << endl;
428 cout << " *AliJet::List* Error : No track present." << endl;
432 ///////////////////////////////////////////////////////////////////////////
433 Int_t AliJet::GetNtracks(Int_t idmode,Int_t chmode,Int_t pcode)
435 // Provide the number of user selected tracks in this jet based on the
436 // idmode, chmode and pcode selections as specified by the user.
437 // For specification of the selection parameters see GetTracks().
438 // The default parameters correspond to no selection, which implies
439 // that invokation of GetNtracks() just returns the total number of
440 // tracks registered in this jet.
442 // Note : In case certain selections are specified, this function
443 // invokes GetTracks(idmode,chmode,pcode) to determine the
444 // number of tracks corresponding to the selections.
445 // When the jet contains a large number of tracks, invokation
446 // of GetTracks(idmode,chmode,pcode) and subsequently invoking
447 // GetEntries() for the resulting TObjArray* might be slightly
451 if (idmode==0 && chmode==2 && pcode==0)
457 TObjArray* arr=GetTracks(idmode,chmode,pcode);
462 ///////////////////////////////////////////////////////////////////////////
463 Double_t AliJet::GetEnergy(Float_t scale)
465 // Return the total energy of the jet.
466 // By default the energy is returned in the units as it was stored in the jet
467 // structure. However, the user can select a different energy unit scale by
468 // specification of the scale parameter.
469 // The convention is that scale=1 corresponds to GeV, so specification
470 // of scale=0.001 will provide the energy in MeV.
471 // The error can be obtained by invoking GetResultError() after
472 // invokation of GetEnergy().
473 Double_t E=GetScalar();
479 fDresult*=fEscale/scale;
488 ///////////////////////////////////////////////////////////////////////////
489 Double_t AliJet::GetMomentum(Float_t scale)
491 // Return the value of the total jet 3-momentum
492 // By default the momentum is returned in the units as it was stored in the jet
493 // structure. However, the user can select a different momentum unit scale by
494 // specification of the scale parameter.
495 // The convention is that scale=1 corresponds to GeV/c, so specification
496 // of scale=0.001 will provide the momentum in MeV/c.
497 // The error can be obtained by invoking GetResultError() after
498 // invokation of GetMomentum().
500 Double_t norm=fV.GetNorm();
501 fDresult=fV.GetResultError();
505 fDresult*=fEscale/scale;
509 ///////////////////////////////////////////////////////////////////////////
510 Ali3Vector AliJet::Get3Momentum(Float_t scale) const
512 // Return the the total jet 3-momentum
513 // By default the components of the 3-momentum are returned in the units
514 // as they were stored in the jet structure.
515 // However, the user can select a different momentum unit scale for the
516 // components by specification of the scale parameter.
517 // The convention is that scale=1 corresponds to GeV/c, so specification
518 // of scale=0.001 will provide the 3-momentum in MeV/c.
520 Ali3Vector p=Get3Vector();
521 if (scale>0) p*=fEscale/scale;
524 ///////////////////////////////////////////////////////////////////////////
525 Double_t AliJet::GetInvmass(Float_t scale)
527 // Return the invariant mass of the jet.
528 // By default the mass is returned in the units as it was stored in the jet
529 // structure. However, the user can select a different mass unit scale by
530 // specification of the scale parameter.
531 // The convention is that scale=1 corresponds to GeV/c**2, so specification
532 // of scale=0.001 will provide the invariant mass in MeV/c**2.
533 // The error can be obtained by invoking GetResultError() after
534 // invokation of GetInvmass().
536 Double_t inv=Dot(*this);
537 Double_t dinv=GetResultError();
541 Double_t m=sqrt(inv);
542 if (m) dm=dinv/(2.*m);
557 ///////////////////////////////////////////////////////////////////////////
558 Float_t AliJet::GetCharge() const
560 // Return the total charge of the jet
563 ///////////////////////////////////////////////////////////////////////////
564 AliTrack* AliJet::GetTrack(Int_t i) const
566 // Return the i-th track of this jet
568 if (!fTracks) return 0;
572 cout << " *AliJet*::GetTrack* Invalid argument i : " << i
573 << " Ntrk = " << fNtrk << endl;
578 return (AliTrack*)fTracks->At(i-1);
581 ///////////////////////////////////////////////////////////////////////////
582 AliTrack* AliJet::GetIdTrack(Int_t id) const
584 // Return the track with user identifier "id" of this jet
585 if (!fTracks) return 0;
588 for (Int_t i=0; i<fNtrk; i++)
590 tx=(AliTrack*)fTracks->At(i);
591 if (id == tx->GetId()) return tx;
593 return 0; // No matching id found
595 ///////////////////////////////////////////////////////////////////////////
596 TObjArray* AliJet::GetTracks(Int_t idmode,Int_t chmode,Int_t pcode)
598 // Provide references to user selected tracks based on the idmode, chmode
599 // and pcode selections as specified by the user.
601 // The following selection combinations are available :
602 // ----------------------------------------------------
603 // idmode = -1 ==> Select tracks with negative user identifier "id"
604 // 0 ==> No selection on user identifier
605 // 1 ==> Select tracks with positive user identifier "id"
607 // chmode = -1 ==> Select tracks with negative charge
608 // 0 ==> Select neutral tracks
609 // 1 ==> Select tracks with positive charge
610 // 2 ==> No selection on charge
611 // 3 ==> Select all charged tracks
613 // pcode = 0 ==> No selection on particle code
614 // X ==> Select tracks with particle code +X or -X
615 // This allows selection of both particles and anti-particles
616 // in case of PDG particle codes.
617 // Selection of either particles or anti-particles can be
618 // obtained in combination with the "chmode" selector.
622 // idmode=-1 chmode=0 pcode=0 : Selection of all neutral tracks with negative id.
623 // idmode=0 chmode=2 pcode=211 : Selection of all charged pions (PDG convention).
624 // idmode=0 chmode=1 pcode=321 : Selection of all positive kaons (PDG convention).
626 // The default values are idmode=0 chmode=2 pcode=0 (i.e. no selections applied).
630 // 1) In case the user has labeled simulated tracks with negative id and
631 // reconstructed tracks with positive id, this memberfunction provides
632 // easy access to either all simulated or reconstructed tracks.
633 // 2) Subsequent invokations of this memberfunction with e.g. chmode=-1 and chmode=1
634 // provides a convenient way to investigate particle pairs with opposite charge
635 // (e.g. for invariant mass analysis).
636 // 3) The selected track pointers are returned via a multi-purpose array,
637 // which will be overwritten by subsequent selections.
638 // In case the selected track list is to be used amongst other selections,
639 // the user is advised to store the selected track pointers in a local
640 // TObjArray or TRefArray.
648 fSelected=new TObjArray();
651 if (!fTracks) return fSelected;
657 for (Int_t i=0; i<fNtrk; i++)
659 tx=(AliTrack*)fTracks->At(i);
662 code=tx->GetParticleCode();
663 if (pcode && abs(pcode)!=abs(code)) continue;
666 if (idmode==-1 && id>=0) continue;
667 if (idmode==1 && id<=0) continue;
670 if (chmode==-1 && q>=0) continue;
671 if (chmode==0 && fabs(q)>1e-10) continue;
672 if (chmode==1 && q<=0) continue;
673 if (chmode==3 && fabs(q)<1e-10) continue;
680 ///////////////////////////////////////////////////////////////////////////
681 TObjArray* AliJet::GetTracks(TString name)
683 // Provide references to all tracks with the specified name.
687 // 1) In case the user has labeled reconstructed tracks with the name of
688 // the applied reconstruction algorithm, this memberfunction provides
689 // easy access to all tracks reconstructed by a certain method.
690 // 2) The selected track pointers are returned via a multi-purpose array,
691 // which will be overwritten by subsequent selections.
692 // In case the selected track list is to be used amongst other selections,
693 // the user is advised to store the selected track pointers in a local
694 // TObjArray or TRefArray.
702 fSelected=new TObjArray();
705 if (!fTracks) return fSelected;
709 for (Int_t i=0; i<fNtrk; i++)
711 tx=(AliTrack*)fTracks->At(i);
715 if (s == name) fSelected->Add(tx);
720 ///////////////////////////////////////////////////////////////////////////
721 void AliJet::RemoveTracks(TString name)
723 // Remove all tracks with the specified name.
724 // If name="*" all tracks will be removed.
728 // In case the user has labeled reconstructed tracks with the name of
729 // the applied reconstruction algorithm, this memberfunction provides
730 // easy removal of all tracks reconstructed by a certain method.
732 if (!fTracks) return;
737 for (Int_t i=0; i<fNtrk; i++)
739 tx=(AliTrack*)fTracks->At(i);
743 if (s==name || name=="*")
745 obj=fTracks->Remove(tx);
746 if (obj && fTracks->IsOwner()) delete tx;
750 fNtrk=fTracks->GetEntries();
752 ///////////////////////////////////////////////////////////////////////////
753 void AliJet::RemoveTracks(Int_t idmode,Int_t chmode,Int_t pcode)
755 // Remove user selected tracks based on the idmode, chmode and pcode
756 // selections as specified by the user.
757 // For defintions of these selections see the corresponding GetTracks()
760 if (!fTracks) return;
762 TObjArray* arr=GetTracks(idmode,chmode,pcode);
765 Int_t ntk=arr->GetEntries();
770 for (Int_t i=0; i<ntk; i++)
772 tx=(AliTrack*)arr->At(i);
775 obj=fTracks->Remove(tx);
776 if (obj && fTracks->IsOwner()) delete tx;
779 fNtrk=fTracks->GetEntries();
782 ///////////////////////////////////////////////////////////////////////////
783 void AliJet::ShowTracks(Int_t mode)
785 // Provide an overview of the available tracks.
786 // The argument mode determines the amount of information as follows :
787 // mode = 0 ==> Only printout of the number of tracks
788 // 1 ==> Provide a listing with 1 line of info for each track
790 // The default is mode=1.
792 Int_t ntk=GetNtracks();
797 cout << " There are " << ntk << " tracks available." << endl;
801 cout << " The following " << ntk << " tracks are available :" << endl;
802 for (Int_t i=1; i<=ntk; i++)
804 AliTrack* tx=GetTrack(i);
807 const char* name=tx->GetName();
808 const char* title=tx->GetTitle();
809 cout << " Track : " << i;
810 cout << " Id : " << tx->GetId();
811 cout << " Q : " << tx->GetCharge() << " m : " << tx->GetMass() << " p : " << tx->GetMomentum();
812 if (strlen(name)) cout << " Name : " << name;
813 if (strlen(title)) cout << " Title : " << title;
821 cout << " No tracks are present." << endl;
824 ///////////////////////////////////////////////////////////////////////////
825 Double_t AliJet::GetPt(Float_t scale)
827 // Provide the transverse momentum value w.r.t. z-axis.
828 // By default the value is returned in the units as it was stored in the jet
829 // structure. However, the user can select a different momentum unit scale by
830 // specification of the scale parameter.
831 // The convention is that scale=1 corresponds to GeV/c, so specification
832 // of scale=0.001 will provide the transverse momentum in MeV/c.
833 // The error on the value can be obtained by GetResultError()
834 // after invokation of GetPt().
837 Double_t norm=v.GetNorm();
838 fDresult=v.GetResultError();
842 fDresult*=fEscale/scale;
847 ///////////////////////////////////////////////////////////////////////////
848 Double_t AliJet::GetPl(Float_t scale)
850 // Provide the longitudinal 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 longitudinal momentum in MeV/c.
856 // Note : the returned value can also be negative.
857 // The error on the value can be obtained by GetResultError()
858 // after invokation of GetPl().
863 Double_t pl=v.GetNorm();
864 fDresult=v.GetResultError();
867 v.GetVector(a,"sph");
868 if (cos(a[1])<0) pl=-pl;
872 fDresult*=fEscale/scale;
877 ///////////////////////////////////////////////////////////////////////////
878 Double_t AliJet::GetEt(Float_t scale)
880 // Provide transverse energy value w.r.t. z-axis.
881 // By default the value is returned in the units as it was stored in the jet
882 // structure. However, the user can select a different energy unit scale by
883 // specification of the scale parameter.
884 // The convention is that scale=1 corresponds to GeV, so specification
885 // of scale=0.001 will provide the transverse energy in MeV.
886 // The error on the value can be obtained by GetResultError()
887 // after invokation of GetEt().
889 Double_t et=GetScaTrans();
893 fDresult*=fEscale/scale;
898 ///////////////////////////////////////////////////////////////////////////
899 Double_t AliJet::GetEl(Float_t scale)
901 // Provide longitudinal energy value w.r.t. z-axis.
902 // By default the value is returned in the units as it was stored in the jet
903 // structure. However, the user can select a different energy unit scale by
904 // specification of the scale parameter.
905 // The convention is that scale=1 corresponds to GeV, so specification
906 // of scale=0.001 will provide the longitudinal energy in MeV.
907 // Note : the returned value can also be negative.
908 // The error on the value can be obtained by GetResultError()
909 // after invokation of GetEl().
911 Double_t el=GetScaLong();
915 fDresult*=fEscale/scale;
920 ///////////////////////////////////////////////////////////////////////////
921 Double_t AliJet::GetMt(Float_t scale)
923 // Provide transverse mass value w.r.t. z-axis.
924 // By default the value is returned in the units as it was stored in the jet
925 // structure. However, the user can select a different energy unit scale by
926 // specification of the scale parameter.
927 // The convention is that scale=1 corresponds to GeV, so specification
928 // of scale=0.001 will provide the transverse mass in MeV.
929 // The error on the value can be obtained by GetResultError()
930 // after invokation of GetMt().
932 Double_t dpt=GetResultError();
933 Double_t m=GetInvmass();
934 Double_t dm=GetResultError();
936 Double_t mt=sqrt(pt*pt+m*m);
938 if (mt) dmt2=(pow((pt*dpt),2)+pow((m*dm),2))/(mt*mt);
944 fDresult*=fEscale/scale;
948 ///////////////////////////////////////////////////////////////////////////
949 Double_t AliJet::GetRapidity()
951 // Provide rapidity value w.r.t. z-axis.
952 // The error on the value can be obtained by GetResultError()
953 // after invokation of GetRapidity().
954 // Note : Also GetPseudoRapidity() is available since this class is
955 // derived from Ali4Vector.
956 Double_t e=GetEnergy();
957 Double_t de=GetResultError();
959 Double_t dpl=GetResultError();
963 Double_t y=9999,dy2=0;
964 if (sum && dif) y=0.5*log(sum/dif);
966 if (sum*dif) dy2=(1./(sum*dif))*(pow((pl*de),2)+pow((e*dpl),2));
971 ///////////////////////////////////////////////////////////////////////////
972 void AliJet::SetTrackCopy(Int_t j)
974 // (De)activate the creation of private copies of the added tracks.
975 // j=0 ==> No private copies are made; pointers of original tracks are stored.
976 // j=1 ==> Private copies of the tracks are made and these pointers are stored.
978 // Note : Once the storage contains pointer(s) to AliTrack(s) one cannot
979 // change the TrackCopy mode anymore.
980 // To change the TrackCopy mode for an existing AliJet containing
981 // tracks one first has to invoke Reset().
990 cout << "*AliJet::SetTrackCopy* Invalid argument : " << j << endl;
995 cout << "*AliJet::SetTrackCopy* Storage already contained tracks."
996 << " ==> TrackCopy mode not changed." << endl;
999 ///////////////////////////////////////////////////////////////////////////
1000 Int_t AliJet::GetTrackCopy() const
1002 // Provide value of the TrackCopy mode.
1003 // 0 ==> No private copies are made; pointers of original tracks are stored.
1004 // 1 ==> Private copies of the tracks are made and these pointers are stored.
1007 ///////////////////////////////////////////////////////////////////////////
1008 void AliJet::SetId(Int_t id)
1010 // Set a user defined identifier for this jet.
1013 ///////////////////////////////////////////////////////////////////////////
1014 Int_t AliJet::GetId() const
1016 // Provide the user defined identifier of this jet.
1019 ///////////////////////////////////////////////////////////////////////////
1020 void AliJet::SetReferencePoint(AliPosition& p)
1022 // Store the position of the jet reference-point.
1023 // The reference-point of a jet provides a means to define a generic
1024 // space-time location for the jet as a whole.
1025 // This doesn't have to be necessarily the location where all the constituent
1026 // tracks originate (e.g. a bundle of parallel tracks doesn't have such
1027 // a location). As such the meaning of this reference-point is different from
1028 // a normal vertex position and allows to provide complimentary information.
1029 // This reference point is the preferable point to start e.g. extrapolations
1030 // and investigate coincidences in space and/or time.
1031 if (fRef) delete fRef;
1032 fRef=new AliPositionObj(p);
1034 ///////////////////////////////////////////////////////////////////////////
1035 AliPosition* AliJet::GetReferencePoint()
1037 // Provide the position of the jet reference-point.
1038 // The reference-point of a jet provides a means to define a generic
1039 // space-time location for the jet as a whole.
1040 // This doesn't have to be necessarily the location where all the constituent
1041 // tracks originate (e.g. a bundle of parallel tracks doesn't have such
1042 // a location). As such the meaning of this reference-point is different from
1043 // a normal vertex position and allows to provide complimentary information.
1044 // This reference point is the preferable point to start e.g. extrapolations
1045 // and investigate coincidences in space and/or time.
1048 ///////////////////////////////////////////////////////////////////////////
1049 TObjArray* AliJet::SortTracks(Int_t mode,TObjArray* tracks)
1051 // Order the references to an array of tracks by looping over the input array "tracks"
1052 // and checking the value of a certain observable.
1053 // The ordered array is returned as a TObjArray.
1054 // In case tracks=0 (default), the registered tracks of the current jet are used.
1055 // Note that the original track array is not modified.
1056 // Via the "mode" argument the user can specify the observable to be checked upon
1057 // and specify whether sorting should be performed in decreasing order (mode<0)
1058 // or in increasing order (mode>0).
1060 // The convention for the observable selection is the following :
1061 // mode : 1 ==> Number of signals associated to the track
1062 // 2 ==> Track energy
1063 // 3 ==> Track momentum
1064 // 4 ==> Mass of the track
1065 // 5 ==> Transverse momentum of the track
1066 // 6 ==> Longitudinal momentum of the track
1067 // 7 ==> Transverse energy of the track
1068 // 8 ==> Longitudinal energy of the track
1069 // 9 ==> Transverse mass of the track
1070 // 10 ==> Track rapidity
1071 // 11 ==> Pseudo-rapidity of the track
1072 // 12 ==> Charge of the track
1073 // 13 ==> Probability of the track hypothesis
1075 // The default is mode=-1.
1077 // Note : This sorting routine uses a common area in memory, which is used
1078 // by various other sorting facilities as well.
1079 // This means that the resulting sorted TObjArray may be overwritten
1080 // when another sorting is invoked.
1081 // To retain the sorted list of pointers, the user is advised to copy
1082 // the pointers contained in the returned TObjArray into a private
1083 // TObjArray instance.
1091 if (!tracks) tracks=fTracks;
1093 if (!mode || abs(mode)>13 || !tracks) return fSelected;
1095 Int_t ntracks=tracks->GetEntries();
1102 fSelected=new TObjArray(ntracks);
1105 Double_t val1,val2; // Values of the observable to be tested upon
1108 for (Int_t i=0; i<ntracks; i++) // Loop over all tracks of the array
1110 AliTrack* tx=(AliTrack*)tracks->At(i);
1114 if (nord == 0) // store the first track at the first ordered position
1117 fSelected->AddAt(tx,nord-1);
1121 for (Int_t j=0; j<=nord; j++) // put track in the right ordered position
1123 if (j == nord) // track has smallest (mode<0) or largest (mode>0) observable value seen so far
1126 fSelected->AddAt(tx,j); // add track at the end
1127 break; // go for next track
1136 val1=tx->GetNsignals();
1137 val2=((AliTrack*)fSelected->At(j))->GetNsignals();
1140 val1=tx->GetEnergy(1);
1141 val2=((AliTrack*)fSelected->At(j))->GetEnergy(1);
1144 val1=tx->GetMomentum(1);
1145 val2=((AliTrack*)fSelected->At(j))->GetMomentum(1);
1148 val1=tx->GetMass(1);
1149 val2=((AliTrack*)fSelected->At(j))->GetMass(1);
1153 val2=((AliTrack*)fSelected->At(j))->GetPt(1);
1157 val2=((AliTrack*)fSelected->At(j))->GetPl(1);
1161 val2=((AliTrack*)fSelected->At(j))->GetEt(1);
1165 val2=((AliTrack*)fSelected->At(j))->GetEl(1);
1169 val2=((AliTrack*)fSelected->At(j))->GetMt(1);
1172 val1=tx->GetRapidity();
1173 val2=((AliTrack*)fSelected->At(j))->GetRapidity();
1176 val1=tx->GetPseudoRapidity();
1177 val2=((AliTrack*)fSelected->At(j))->GetPseudoRapidity();
1180 val1=tx->GetCharge();
1181 val2=((AliTrack*)fSelected->At(j))->GetCharge();
1185 val2=((AliTrack*)fSelected->At(j))->GetProb();
1189 if (mode<0 && val1 <= val2) continue;
1190 if (mode>0 && val1 >= val2) continue;
1193 for (Int_t k=nord-1; k>j; k--) // create empty position
1195 fSelected->AddAt(fSelected->At(k-1),k);
1197 fSelected->AddAt(tx,j); // put track at empty position
1198 break; // go for next track
1203 ///////////////////////////////////////////////////////////////////////////
1204 Double_t AliJet::GetDistance(AliPosition* p,Float_t scale)
1206 // Provide distance of the current jet to the position p.
1207 // The error on the result can be obtained as usual by invoking
1208 // GetResultError() afterwards.
1210 // By default the distance will be provided in the metric unit scale of
1211 // the AliPosition p.
1212 // However, the user can select a different metric unit scale by
1213 // specification of the scale parameter.
1214 // The convention is that scale=1 corresponds to meter, so specification
1215 // of scale=0.01 will provide the distance in cm.
1216 // As such it is possible to obtain a correctly computed distance even in case
1217 // the jet parameters have a different unit scale.
1218 // However, it is recommended to work always with one single unit scale.
1220 // Note : In case of incomplete information, a distance value of -1 is
1226 if (!p) return dist;
1228 // Obtain a defined position on this jet
1229 AliPosition* rx=fRef;
1231 if (!rx) return dist;
1233 Ali3Vector pj=Get3Momentum();
1235 if (pj.GetNorm() <= 0.) return dist;
1238 tj.Set3Momentum(pj);
1239 tj.SetReferencePoint(*rx);
1240 dist=tj.GetDistance(p,scale);
1241 fDresult=tj.GetResultError();
1244 ///////////////////////////////////////////////////////////////////////////
1245 Double_t AliJet::GetDistance(AliTrack* t,Float_t scale)
1247 // Provide distance of the current jet to the track t.
1248 // The error on the result can be obtained as usual by invoking
1249 // GetResultError() afterwards.
1251 // By default the distance will be provided in the metric unit scale of
1253 // However, the user can specify a required metric unit scale by specification
1254 // of the scale parameter.
1255 // The convention is that scale=1 corresponds to meter, so specification
1256 // of scale=0.01 will provide the distance in cm.
1257 // As such it is possible to obtain a correctly computed distance even in case
1258 // the jet and track parameters have a different unit scale.
1259 // However, it is recommended to work always with one single unit scale.
1261 // Note : In case of incomplete information, a distance value of -1 is
1267 if (!t) return dist;
1269 // Obtain a defined position on this jet
1270 AliPosition* rx=fRef;
1272 if (!rx) return dist;
1274 Ali3Vector pj=Get3Momentum();
1276 if (pj.GetNorm() <= 0.) return dist;
1279 tj.Set3Momentum(pj);
1280 tj.SetReferencePoint(*rx);
1281 dist=tj.GetDistance(t,scale);
1282 fDresult=tj.GetResultError();
1285 ///////////////////////////////////////////////////////////////////////////
1286 Double_t AliJet::GetDistance(AliJet* j,Float_t scale)
1288 // Provide distance of the current jet to the jet j.
1289 // The error on the result can be obtained as usual by invoking
1290 // GetResultError() afterwards.
1292 // By default the distance will be provided in the metric unit scale of
1294 // This implies that the results of j1.GetDistance(j2) and j2.GetDistance(j1)
1295 // may be numerically different in case j1 and j2 have different metric units.
1296 // However, the user can specify a required metric unit scale by specification
1297 // of the scale parameter.
1298 // The convention is that scale=1 corresponds to meter, so specification
1299 // of scale=0.01 will provide the distance in cm.
1300 // As such it is possible to obtain a correctly computed distance even in case
1301 // the jet parameters have a different unit scale.
1302 // However, it is recommended to work always with one single unit scale.
1304 // Note : In case of incomplete information, a distance value of -1 is
1310 if (!j) return dist;
1312 // Obtain a defined position on jet j
1313 AliPosition* rx=j->GetReferencePoint();
1315 if (!rx) return dist;
1317 Ali3Vector pj=j->Get3Momentum();
1319 if (pj.GetNorm() <= 0.) return dist;
1322 tj.Set3Momentum(pj);
1323 tj.SetReferencePoint(*rx);
1324 dist=GetDistance(tj,scale);
1327 ///////////////////////////////////////////////////////////////////////////
1328 Int_t AliJet::GetNsignals() const
1330 // Provide the number of signals associated to the jet tracks.
1331 // Note : Multiple occurrences of the same signal are only counted once.
1333 if (fNtrk<1) return 0;
1340 for (Int_t i=1; i<=fNtrk; i++)
1343 for (Int_t j=1; j<=tx->GetNsignals(); j++)
1345 AliSignal* sx=tx->GetSignal(j);
1348 for (Int_t k=0; k<arr.GetEntries(); k++)
1350 if (sx==(AliSignal*)arr.At(k))
1356 if (!exists) arr.Add(sx);
1362 ///////////////////////////////////////////////////////////////////////////
1363 void AliJet::SetEscale(Float_t scale)
1365 // Indicate the energy/momentum scale as used by the user.
1366 // The convention is that scale=1 indicates values in units
1367 // of GeV, GeV/c or GeV/c**2.
1368 // So, in case one decides to store values in units of MeV, MeV/c or MeV/c**2
1369 // the scale indicator should be set to scale=0.001.
1371 // By default scale=1 is set in the constructor.
1379 cout << " *AliJet::SetEscale* Invalid scale value : " << scale << endl;
1382 ///////////////////////////////////////////////////////////////////////////
1383 Float_t AliJet::GetEscale() const
1385 // Provide the energy/momentum scale as used by the user.
1386 // The convention is that scale=1 indicates values in units
1387 // of GeV, GeV/c or GeV/c**2.
1388 // So, a value of scale=0.001 indicates that energy/momentum values are
1389 // stored in units of MeV, MeV/c or MeV/c**2.
1392 ///////////////////////////////////////////////////////////////////////////
1393 TObject* AliJet::Clone(const char* name) const
1395 // Make a deep copy of the current object and provide the pointer to the copy.
1396 // This memberfunction enables automatic creation of new objects of the
1397 // correct type depending on the object type, a feature which may be very useful
1398 // for containers when adding objects in case the container owns the objects.
1399 // This feature allows e.g. AliVertex to store either AliJet objects or
1400 // objects derived from AliJet via the AddJet memberfunction, provided
1401 // these derived classes also have a proper Clone memberfunction.
1403 AliJet* jet=new AliJet(*this);
1406 if (strlen(name)) jet->SetName(name);
1410 ///////////////////////////////////////////////////////////////////////////