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 ///////////////////////////////////////////////////////////////////////////
99 #include "Riostream.h"
101 ClassImp(AliJet) // Class implementation to enable ROOT I/O
103 AliJet::AliJet() : TNamed(),Ali4Vector()
105 // Default constructor
106 // All variables initialised to 0
107 // Initial maximum number of tracks is set to the default value
112 ///////////////////////////////////////////////////////////////////////////
115 // Initialisation of pointers etc...
123 ///////////////////////////////////////////////////////////////////////////
124 AliJet::AliJet(Int_t n) : TNamed(),Ali4Vector()
126 // Create a jet to hold initially a maximum of n tracks
127 // All variables initialised to 0
137 cout << " *AliJet* Initial max. number of tracks entered : " << n << endl;
138 cout << " This is invalid. Default initial maximum will be used." << endl;
143 ///////////////////////////////////////////////////////////////////////////
146 // Default destructor
163 ///////////////////////////////////////////////////////////////////////////
164 void AliJet::SetOwner(Bool_t own)
166 // Set ownership of all added objects.
167 // The default parameter is own=kTRUE.
169 // Invokation of this memberfunction also sets all the copy modes
170 // (e.g. TrackCopy & co.) according to the value of own.
172 // This function (with own=kTRUE) is particularly useful when reading data
173 // from a tree/file, since Reset() will then actually remove all the
174 // added objects from memory irrespective of the copy mode settings
175 // during the tree/file creation process. In this way it provides a nice way
176 // of preventing possible memory leaks in the reading/analysis process.
178 // In addition this memberfunction can also be used as a shortcut to set all
179 // copy modes in one go during a tree/file creation process.
180 // However, in this case the user has to take care to only set/change the
181 // ownership (and copy mode) for empty objects (e.g. newly created objects
182 // or after invokation of the Reset() memberfunction) otherwise it will
183 // very likely result in inconsistent destructor behaviour.
187 if (fTracks) fTracks->SetOwner(own);
190 ///////////////////////////////////////////////////////////////////////////
191 AliJet::AliJet(const AliJet& j) : TNamed(j),Ali4Vector(j)
199 fTrackCopy=j.fTrackCopy;
201 if (j.fRef) fRef=new AliPositionObj(*(j.fRef));
208 fTracks=new TObjArray(fNtmax);
209 if (fTrackCopy) fTracks->SetOwner();
212 for (Int_t i=1; i<=fNtrk; i++)
214 AliTrack* tx=j.GetTrack(i);
217 fTracks->Add(tx->Clone());
225 ///////////////////////////////////////////////////////////////////////////
226 void AliJet::SetNtinit(Int_t n)
228 // Set the initial maximum number of tracks for this jet
243 ///////////////////////////////////////////////////////////////////////////
246 // Reset all variables to 0
247 // The max. number of tracks is set to the initial value again
248 // Note : The scale for the energy/momentum units will not be changed.
252 Double_t a[4]={0,0,0,0};
254 if (fNtinit > 0) SetNtinit(fNtinit);
256 ///////////////////////////////////////////////////////////////////////////
257 void AliJet::AddTrack(AliTrack& t)
259 // Add a track to the jet.
260 // In case the maximum number of tracks has been reached
261 // space will be extended to hold an additional amount of tracks as
262 // was initially reserved.
263 // See SetTrackCopy() to tailor the functionality of the stored structures.
267 // In case a private copy is made, this is performed via the Clone() memberfunction.
268 // All AliTrack and derived classes have the default TObject::Clone() memberfunction.
269 // However, derived classes generally contain an internal data structure which may
270 // include pointers to other objects. Therefore it is recommended to provide
271 // for all derived classes a specific copy constructor and override the default Clone()
272 // memberfunction using this copy constructor.
273 // An example for this may be seen from AliTrack.
275 // In case NO private copy is made, a check will be performed if this
276 // specific track is already present in the jet.
277 // If this is the case, no action is performed to prevent multiple
278 // additions of the same track.
283 ///////////////////////////////////////////////////////////////////////////
284 void AliJet::AddTrack(AliTrack& t,Int_t copy)
286 // Internal memberfunction to actually add a track to the jet.
287 // In case the maximum number of tracks has been reached
288 // space will be extended to hold an additional amount of tracks as
289 // was initially reserved.
291 // If copy=0 NO copy of the track will be made, irrespective of the setting
292 // of the TrackCopy flag.
293 // This allows a proper treatment of automatically generated connecting
294 // tracks between vertices.
296 // In case NO copy of the track is made, a check will be performed if this
297 // specific track is already present in the jet.
298 // If this is the case, no action is performed to prevent multiple
299 // additions of the same track.
302 // In case a private copy is made, this is performed via the Clone() memberfunction.
306 fTracks=new TObjArray(fNtmax);
307 if (fTrackCopy) fTracks->SetOwner();
309 else if (!fTrackCopy || !copy) // Check if this track is already present
311 for (Int_t i=0; i<fNtrk; i++)
313 AliTrack* tx=(AliTrack*)fTracks->At(i);
314 if (tx == &t) return;
318 if (fNtrk == fNtmax) // Check if maximum track number is reached
321 fTracks->Expand(fNtmax);
324 // Add the track to this jet
326 if (fTrackCopy && copy)
328 fTracks->Add(t.Clone());
337 Ali4Vector p4=(Ali4Vector)t;
338 Float_t tscale=t.GetEscale();
339 if ((tscale/fEscale > 1.1) || (fEscale/tscale > 1.1)) p4=p4*(tscale/fEscale);
343 ///////////////////////////////////////////////////////////////////////////
344 void AliJet::Data(TString f,TString u)
346 // Provide jet information within the coordinate frame f
348 // The string argument "u" allows to choose between different angular units
349 // in case e.g. a spherical frame is selected.
350 // u = "rad" : angles provided in radians
351 // "deg" : angles provided in degrees
353 // The defaults are f="car" and u="rad".
355 const char* name=GetName();
356 const char* title=GetTitle();
358 cout << " *AliJet::Data*";
359 if (strlen(name)) cout << " Name : " << GetName();
360 if (strlen(title)) cout << " Title : " << GetTitle();
362 cout << " Id : " << fUserId << " Invmass : " << GetInvmass() << " Charge : " << fQ
363 << " Momentum : " << GetMomentum() << " Energy scale : " << fEscale << " GeV" << endl;
367 Ali4Vector::Data(f,u);
369 ///////////////////////////////////////////////////////////////////////////
370 void AliJet::List(TString f,TString u)
372 // Provide jet and primary track information within the coordinate frame f
374 // The string argument "u" allows to choose between different angular units
375 // in case e.g. a spherical frame is selected.
376 // u = "rad" : angles provided in radians
377 // "deg" : angles provided in degrees
379 // The defaults are f="car" and u="rad".
381 Data(f,u); // Information of the current jet
382 if (fRef) { cout << " Ref-point :"; fRef->Data(f,u); }
384 // The tracks of this jet
386 for (Int_t it=1; it<=fNtrk; it++)
391 cout << " ---Track no. " << it << endl;
397 cout << " *AliJet::List* Error : No track present." << endl;
401 ///////////////////////////////////////////////////////////////////////////
402 void AliJet::ListAll(TString f,TString u)
404 // Provide jet and prim.+sec. track information within the coordinate frame f
406 // The string argument "u" allows to choose between different angular units
407 // in case e.g. a spherical frame is selected.
408 // u = "rad" : angles provided in radians
409 // "deg" : angles provided in degrees
411 // The defaults are f="car" and u="rad".
413 Data(f,u); // Information of the current jet
414 if (fRef) { cout << " Ref-point :"; fRef->Data(f,u); }
416 // The tracks of this jet
418 for (Int_t it=1; it<=fNtrk; it++)
423 cout << " ---Track no. " << it << endl;
429 cout << " *AliJet::List* Error : No track present." << endl;
433 ///////////////////////////////////////////////////////////////////////////
434 Int_t AliJet::GetNtracks(Int_t idmode,Int_t chmode,Int_t pcode)
436 // Provide the number of user selected tracks in this jet based on the
437 // idmode, chmode and pcode selections as specified by the user.
438 // For specification of the selection parameters see GetTracks().
439 // The default parameters correspond to no selection, which implies
440 // that invokation of GetNtracks() just returns the total number of
441 // tracks registered in this jet.
443 // Note : In case certain selections are specified, this function
444 // invokes GetTracks(idmode,chmode,pcode) to determine the
445 // number of tracks corresponding to the selections.
446 // When the jet contains a large number of tracks, invokation
447 // of GetTracks(idmode,chmode,pcode) and subsequently invoking
448 // GetEntries() for the resulting TObjArray* might be slightly
452 if (idmode==0 && chmode==2 && pcode==0)
458 TObjArray* arr=GetTracks(idmode,chmode,pcode);
459 if (arr) n=arr->GetEntries();
463 ///////////////////////////////////////////////////////////////////////////
464 Int_t AliJet::GetNtracks(TString name)
466 // Provide the number of tracks with the specified name.
470 // This facility invokes the corresponding GetTracks memberfunction
471 // and as such may result in overwriting existing track selection
472 // arrays. Please refer to the docs of GetTracks for further details.
476 TObjArray* arr=GetTracks(name);
477 if (arr) n=arr->GetEntries();
480 ///////////////////////////////////////////////////////////////////////////
481 Double_t AliJet::GetEnergy(Float_t scale)
483 // Return the total energy of the jet.
484 // By default the energy is returned in the units as it was stored in the jet
485 // structure. However, the user can select a different energy unit scale by
486 // specification of the scale parameter.
487 // The convention is that scale=1 corresponds to GeV, so specification
488 // of scale=0.001 will provide the energy in MeV.
489 // The error can be obtained by invoking GetResultError() after
490 // invokation of GetEnergy().
491 Double_t E=GetScalar();
497 fDresult*=fEscale/scale;
506 ///////////////////////////////////////////////////////////////////////////
507 Double_t AliJet::GetMomentum(Float_t scale)
509 // Return the value of the total jet 3-momentum
510 // By default the momentum is returned in the units as it was stored in the jet
511 // structure. However, the user can select a different momentum unit scale by
512 // specification of the scale parameter.
513 // The convention is that scale=1 corresponds to GeV/c, so specification
514 // of scale=0.001 will provide the momentum in MeV/c.
515 // The error can be obtained by invoking GetResultError() after
516 // invokation of GetMomentum().
518 Double_t norm=fV.GetNorm();
519 fDresult=fV.GetResultError();
523 fDresult*=fEscale/scale;
527 ///////////////////////////////////////////////////////////////////////////
528 Ali3Vector AliJet::Get3Momentum(Float_t scale) const
530 // Return the the total jet 3-momentum
531 // By default the components of the 3-momentum are returned in the units
532 // as they were stored in the jet structure.
533 // However, the user can select a different momentum unit scale for the
534 // components by specification of the scale parameter.
535 // The convention is that scale=1 corresponds to GeV/c, so specification
536 // of scale=0.001 will provide the 3-momentum in MeV/c.
538 Ali3Vector p=Get3Vector();
539 if (scale>0) p*=fEscale/scale;
542 ///////////////////////////////////////////////////////////////////////////
543 Double_t AliJet::GetInvmass(Float_t scale)
545 // Return the invariant mass of the jet.
546 // By default the mass is returned in the units as it was stored in the jet
547 // structure. However, the user can select a different mass unit scale by
548 // specification of the scale parameter.
549 // The convention is that scale=1 corresponds to GeV/c**2, so specification
550 // of scale=0.001 will provide the invariant mass in MeV/c**2.
551 // The error can be obtained by invoking GetResultError() after
552 // invokation of GetInvmass().
554 Double_t inv=Dot(*this);
555 Double_t dinv=GetResultError();
559 Double_t m=sqrt(inv);
560 if (m) dm=dinv/(2.*m);
575 ///////////////////////////////////////////////////////////////////////////
576 Float_t AliJet::GetCharge() const
578 // Return the total charge of the jet
581 ///////////////////////////////////////////////////////////////////////////
582 AliTrack* AliJet::GetTrack(Int_t i) const
584 // Return the i-th track of this jet
586 if (!fTracks) return 0;
590 cout << " *AliJet*::GetTrack* Invalid argument i : " << i
591 << " Ntrk = " << fNtrk << endl;
596 return (AliTrack*)fTracks->At(i-1);
599 ///////////////////////////////////////////////////////////////////////////
600 AliTrack* AliJet::GetIdTrack(Int_t id) const
602 // Return the track with user identifier "id" of this jet
603 if (!fTracks) return 0;
606 for (Int_t i=0; i<fNtrk; i++)
608 tx=(AliTrack*)fTracks->At(i);
609 if (id == tx->GetId()) return tx;
611 return 0; // No matching id found
613 ///////////////////////////////////////////////////////////////////////////
614 TObjArray* AliJet::GetTracks(Int_t idmode,Int_t chmode,Int_t pcode)
616 // Provide references to user selected tracks based on the idmode, chmode
617 // and pcode selections as specified by the user.
619 // The following selection combinations are available :
620 // ----------------------------------------------------
621 // idmode = -1 ==> Select tracks with negative user identifier "id"
622 // 0 ==> No selection on user identifier
623 // 1 ==> Select tracks with positive user identifier "id"
625 // chmode = -1 ==> Select tracks with negative charge
626 // 0 ==> Select neutral tracks
627 // 1 ==> Select tracks with positive charge
628 // 2 ==> No selection on charge
629 // 3 ==> Select all charged tracks
631 // pcode = 0 ==> No selection on particle code
632 // X ==> Select tracks with particle code +X or -X
633 // This allows selection of both particles and anti-particles
634 // in case of PDG particle codes.
635 // Selection of either particles or anti-particles can be
636 // obtained in combination with the "chmode" selector.
640 // idmode=-1 chmode=0 pcode=0 : Selection of all neutral tracks with negative id.
641 // idmode=0 chmode=2 pcode=211 : Selection of all charged pions (PDG convention).
642 // idmode=0 chmode=1 pcode=321 : Selection of all positive kaons (PDG convention).
644 // The default values are idmode=0 chmode=2 pcode=0 (i.e. no selections applied).
648 // 1) In case the user has labeled simulated tracks with negative id and
649 // reconstructed tracks with positive id, this memberfunction provides
650 // easy access to either all simulated or reconstructed tracks.
651 // 2) Subsequent invokations of this memberfunction with e.g. chmode=-1 and chmode=1
652 // provides a convenient way to investigate particle pairs with opposite charge
653 // (e.g. for invariant mass analysis).
654 // 3) The selected track pointers are returned via a multi-purpose array,
655 // which will be overwritten by subsequent selections.
656 // In case the selected track list is to be used amongst other selections,
657 // the user is advised to store the selected track pointers in a local
658 // TObjArray or TRefArray.
666 fSelected=new TObjArray();
669 if (!fTracks) return fSelected;
675 for (Int_t i=0; i<fNtrk; i++)
677 tx=(AliTrack*)fTracks->At(i);
680 code=tx->GetParticleCode();
681 if (pcode && abs(pcode)!=abs(code)) continue;
684 if (idmode==-1 && id>=0) continue;
685 if (idmode==1 && id<=0) continue;
688 if (chmode==-1 && q>=0) continue;
689 if (chmode==0 && fabs(q)>1e-10) continue;
690 if (chmode==1 && q<=0) continue;
691 if (chmode==3 && fabs(q)<1e-10) continue;
698 ///////////////////////////////////////////////////////////////////////////
699 TObjArray* AliJet::GetTracks(TString name)
701 // Provide references to all tracks with the specified name.
705 // 1) In case the user has labeled reconstructed tracks with the name of
706 // the applied reconstruction algorithm, this memberfunction provides
707 // easy access to all tracks reconstructed by a certain method.
708 // 2) The selected track pointers are returned via a multi-purpose array,
709 // which will be overwritten by subsequent selections.
710 // In case the selected track list is to be used amongst other selections,
711 // the user is advised to store the selected track pointers in a local
712 // TObjArray or TRefArray.
720 fSelected=new TObjArray();
723 if (!fTracks) return fSelected;
727 for (Int_t i=0; i<fNtrk; i++)
729 tx=(AliTrack*)fTracks->At(i);
733 if (s == name) fSelected->Add(tx);
738 ///////////////////////////////////////////////////////////////////////////
739 void AliJet::RemoveTracks(TString name)
741 // Remove all tracks with the specified name.
742 // If name="*" all tracks will be removed.
746 // In case the user has labeled reconstructed tracks with the name of
747 // the applied reconstruction algorithm, this memberfunction provides
748 // easy removal of all tracks reconstructed by a certain method.
750 if (!fTracks) return;
755 for (Int_t i=0; i<fNtrk; i++)
757 tx=(AliTrack*)fTracks->At(i);
761 if (s==name || name=="*")
763 obj=fTracks->Remove(tx);
764 if (obj && fTracks->IsOwner()) delete tx;
768 fNtrk=fTracks->GetEntries();
770 ///////////////////////////////////////////////////////////////////////////
771 void AliJet::RemoveTracks(Int_t idmode,Int_t chmode,Int_t pcode)
773 // Remove user selected tracks based on the idmode, chmode and pcode
774 // selections as specified by the user.
775 // For defintions of these selections see the corresponding GetTracks()
778 if (!fTracks) return;
780 TObjArray* arr=GetTracks(idmode,chmode,pcode);
783 Int_t ntk=arr->GetEntries();
788 for (Int_t i=0; i<ntk; i++)
790 tx=(AliTrack*)arr->At(i);
793 obj=fTracks->Remove(tx);
794 if (obj && fTracks->IsOwner()) delete tx;
797 fNtrk=fTracks->GetEntries();
800 ///////////////////////////////////////////////////////////////////////////
801 void AliJet::ShowTracks(Int_t mode)
803 // Provide an overview of the available tracks.
804 // The argument mode determines the amount of information as follows :
805 // mode = 0 ==> Only printout of the number of tracks
806 // 1 ==> Provide a listing with 1 line of info for each track
808 // The default is mode=1.
810 Int_t ntk=GetNtracks();
815 cout << " There are " << ntk << " tracks available." << endl;
819 cout << " The following " << ntk << " tracks are available :" << endl;
820 for (Int_t i=1; i<=ntk; i++)
822 AliTrack* tx=GetTrack(i);
825 const char* name=tx->GetName();
826 const char* title=tx->GetTitle();
827 cout << " Track : " << i;
828 cout << " Id : " << tx->GetId();
829 cout << " Q : " << tx->GetCharge() << " m : " << tx->GetMass() << " p : " << tx->GetMomentum();
830 if (strlen(name)) cout << " Name : " << name;
831 if (strlen(title)) cout << " Title : " << title;
839 cout << " No tracks are present." << endl;
842 ///////////////////////////////////////////////////////////////////////////
843 Double_t AliJet::GetPt(Float_t scale)
845 // Provide the transverse momentum value w.r.t. z-axis.
846 // By default the value is returned in the units as it was stored in the jet
847 // structure. However, the user can select a different momentum unit scale by
848 // specification of the scale parameter.
849 // The convention is that scale=1 corresponds to GeV/c, so specification
850 // of scale=0.001 will provide the transverse momentum in MeV/c.
851 // The error on the value can be obtained by GetResultError()
852 // after invokation of GetPt().
855 Double_t norm=v.GetNorm();
856 fDresult=v.GetResultError();
860 fDresult*=fEscale/scale;
865 ///////////////////////////////////////////////////////////////////////////
866 Double_t AliJet::GetPl(Float_t scale)
868 // Provide the longitudinal momentum value w.r.t. z-axis.
869 // By default the value is returned in the units as it was stored in the jet
870 // structure. However, the user can select a different momentum unit scale by
871 // specification of the scale parameter.
872 // The convention is that scale=1 corresponds to GeV/c, so specification
873 // of scale=0.001 will provide the longitudinal momentum in MeV/c.
874 // Note : the returned value can also be negative.
875 // The error on the value can be obtained by GetResultError()
876 // after invokation of GetPl().
881 Double_t pl=v.GetNorm();
882 fDresult=v.GetResultError();
885 v.GetVector(a,"sph");
886 if (cos(a[1])<0) pl=-pl;
890 fDresult*=fEscale/scale;
895 ///////////////////////////////////////////////////////////////////////////
896 Double_t AliJet::GetEt(Float_t scale)
898 // Provide transverse energy value w.r.t. z-axis.
899 // By default the value is returned in the units as it was stored in the jet
900 // structure. However, the user can select a different energy unit scale by
901 // specification of the scale parameter.
902 // The convention is that scale=1 corresponds to GeV, so specification
903 // of scale=0.001 will provide the transverse energy in MeV.
904 // The error on the value can be obtained by GetResultError()
905 // after invokation of GetEt().
907 Double_t et=GetScaTrans();
911 fDresult*=fEscale/scale;
916 ///////////////////////////////////////////////////////////////////////////
917 Double_t AliJet::GetEl(Float_t scale)
919 // Provide longitudinal energy value w.r.t. z-axis.
920 // By default the value is returned in the units as it was stored in the jet
921 // structure. However, the user can select a different energy unit scale by
922 // specification of the scale parameter.
923 // The convention is that scale=1 corresponds to GeV, so specification
924 // of scale=0.001 will provide the longitudinal energy in MeV.
925 // Note : the returned value can also be negative.
926 // The error on the value can be obtained by GetResultError()
927 // after invokation of GetEl().
929 Double_t el=GetScaLong();
933 fDresult*=fEscale/scale;
938 ///////////////////////////////////////////////////////////////////////////
939 Double_t AliJet::GetMt(Float_t scale)
941 // Provide transverse mass value w.r.t. z-axis.
942 // By default the value is returned in the units as it was stored in the jet
943 // structure. However, the user can select a different energy unit scale by
944 // specification of the scale parameter.
945 // The convention is that scale=1 corresponds to GeV, so specification
946 // of scale=0.001 will provide the transverse mass in MeV.
947 // The error on the value can be obtained by GetResultError()
948 // after invokation of GetMt().
950 Double_t dpt=GetResultError();
951 Double_t m=GetInvmass();
952 Double_t dm=GetResultError();
954 Double_t mt=sqrt(pt*pt+m*m);
956 if (mt) dmt2=(pow((pt*dpt),2)+pow((m*dm),2))/(mt*mt);
962 fDresult*=fEscale/scale;
966 ///////////////////////////////////////////////////////////////////////////
967 Double_t AliJet::GetRapidity()
969 // Provide rapidity value w.r.t. z-axis.
970 // The error on the value can be obtained by GetResultError()
971 // after invokation of GetRapidity().
972 // Note : Also GetPseudoRapidity() is available since this class is
973 // derived from Ali4Vector.
974 Double_t e=GetEnergy();
975 Double_t de=GetResultError();
977 Double_t dpl=GetResultError();
981 Double_t y=9999,dy2=0;
982 if (sum && dif) y=0.5*log(sum/dif);
984 if (sum*dif) dy2=(1./(sum*dif))*(pow((pl*de),2)+pow((e*dpl),2));
989 ///////////////////////////////////////////////////////////////////////////
990 void AliJet::SetTrackCopy(Int_t j)
992 // (De)activate the creation of private copies of the added tracks.
993 // j=0 ==> No private copies are made; pointers of original tracks are stored.
994 // j=1 ==> Private copies of the tracks are made and these pointers are stored.
996 // Note : Once the storage contains pointer(s) to AliTrack(s) one cannot
997 // change the TrackCopy mode anymore.
998 // To change the TrackCopy mode for an existing AliJet containing
999 // tracks one first has to invoke Reset().
1008 cout << "*AliJet::SetTrackCopy* Invalid argument : " << j << endl;
1013 cout << "*AliJet::SetTrackCopy* Storage already contained tracks."
1014 << " ==> TrackCopy mode not changed." << endl;
1017 ///////////////////////////////////////////////////////////////////////////
1018 Int_t AliJet::GetTrackCopy() const
1020 // Provide value of the TrackCopy mode.
1021 // 0 ==> No private copies are made; pointers of original tracks are stored.
1022 // 1 ==> Private copies of the tracks are made and these pointers are stored.
1025 ///////////////////////////////////////////////////////////////////////////
1026 void AliJet::SetId(Int_t id)
1028 // Set a user defined identifier for this jet.
1031 ///////////////////////////////////////////////////////////////////////////
1032 Int_t AliJet::GetId() const
1034 // Provide the user defined identifier of this jet.
1037 ///////////////////////////////////////////////////////////////////////////
1038 void AliJet::SetReferencePoint(AliPosition& p)
1040 // Store the position of the jet reference-point.
1041 // The reference-point of a jet provides a means to define a generic
1042 // space-time location for the jet as a whole.
1043 // This doesn't have to be necessarily the location where all the constituent
1044 // tracks originate (e.g. a bundle of parallel tracks doesn't have such
1045 // a location). As such the meaning of this reference-point is different from
1046 // a normal vertex position and allows to provide complimentary information.
1047 // This reference point is the preferable point to start e.g. extrapolations
1048 // and investigate coincidences in space and/or time.
1049 if (fRef) delete fRef;
1050 fRef=new AliPositionObj(p);
1052 ///////////////////////////////////////////////////////////////////////////
1053 AliPosition* AliJet::GetReferencePoint()
1055 // Provide the position of the jet reference-point.
1056 // The reference-point of a jet provides a means to define a generic
1057 // space-time location for the jet as a whole.
1058 // This doesn't have to be necessarily the location where all the constituent
1059 // tracks originate (e.g. a bundle of parallel tracks doesn't have such
1060 // a location). As such the meaning of this reference-point is different from
1061 // a normal vertex position and allows to provide complimentary information.
1062 // This reference point is the preferable point to start e.g. extrapolations
1063 // and investigate coincidences in space and/or time.
1066 ///////////////////////////////////////////////////////////////////////////
1067 TObjArray* AliJet::SortTracks(Int_t mode,TObjArray* tracks)
1069 // Order the references to an array of tracks by looping over the input array "tracks"
1070 // and checking the value of a certain observable.
1071 // The ordered array is returned as a TObjArray.
1072 // In case tracks=0 (default), the registered tracks of the current jet are used.
1073 // Note that the original track array is not modified.
1074 // Via the "mode" argument the user can specify the observable to be checked upon
1075 // and specify whether sorting should be performed in decreasing order (mode<0)
1076 // or in increasing order (mode>0).
1078 // The convention for the observable selection is the following :
1079 // mode : 1 ==> Number of signals associated to the track
1080 // 2 ==> Track energy
1081 // 3 ==> Track momentum
1082 // 4 ==> Mass of the track
1083 // 5 ==> Transverse momentum of the track
1084 // 6 ==> Longitudinal momentum of the track
1085 // 7 ==> Transverse energy of the track
1086 // 8 ==> Longitudinal energy of the track
1087 // 9 ==> Transverse mass of the track
1088 // 10 ==> Track rapidity
1089 // 11 ==> Pseudo-rapidity of the track
1090 // 12 ==> Charge of the track
1091 // 13 ==> Probability of the track hypothesis
1093 // The default is mode=-1.
1095 // Note : This sorting routine uses a common area in memory, which is used
1096 // by various other sorting facilities as well.
1097 // This means that the resulting sorted TObjArray may be overwritten
1098 // when another sorting is invoked.
1099 // To retain the sorted list of pointers, the user is advised to copy
1100 // the pointers contained in the returned TObjArray into a private
1101 // TObjArray instance.
1109 if (!tracks) tracks=fTracks;
1111 if (!mode || abs(mode)>13 || !tracks) return fSelected;
1113 Int_t ntracks=tracks->GetEntries();
1120 fSelected=new TObjArray(ntracks);
1123 Double_t val1,val2; // Values of the observable to be tested upon
1126 for (Int_t i=0; i<ntracks; i++) // Loop over all tracks of the array
1128 AliTrack* tx=(AliTrack*)tracks->At(i);
1132 if (nord == 0) // store the first track at the first ordered position
1135 fSelected->AddAt(tx,nord-1);
1139 for (Int_t j=0; j<=nord; j++) // put track in the right ordered position
1141 if (j == nord) // track has smallest (mode<0) or largest (mode>0) observable value seen so far
1144 fSelected->AddAt(tx,j); // add track at the end
1145 break; // go for next track
1154 val1=tx->GetNsignals();
1155 val2=((AliTrack*)fSelected->At(j))->GetNsignals();
1158 val1=tx->GetEnergy(1);
1159 val2=((AliTrack*)fSelected->At(j))->GetEnergy(1);
1162 val1=tx->GetMomentum(1);
1163 val2=((AliTrack*)fSelected->At(j))->GetMomentum(1);
1166 val1=tx->GetMass(1);
1167 val2=((AliTrack*)fSelected->At(j))->GetMass(1);
1171 val2=((AliTrack*)fSelected->At(j))->GetPt(1);
1175 val2=((AliTrack*)fSelected->At(j))->GetPl(1);
1179 val2=((AliTrack*)fSelected->At(j))->GetEt(1);
1183 val2=((AliTrack*)fSelected->At(j))->GetEl(1);
1187 val2=((AliTrack*)fSelected->At(j))->GetMt(1);
1190 val1=tx->GetRapidity();
1191 val2=((AliTrack*)fSelected->At(j))->GetRapidity();
1194 val1=tx->GetPseudoRapidity();
1195 val2=((AliTrack*)fSelected->At(j))->GetPseudoRapidity();
1198 val1=tx->GetCharge();
1199 val2=((AliTrack*)fSelected->At(j))->GetCharge();
1203 val2=((AliTrack*)fSelected->At(j))->GetProb();
1207 if (mode<0 && val1 <= val2) continue;
1208 if (mode>0 && val1 >= val2) continue;
1211 for (Int_t k=nord-1; k>j; k--) // create empty position
1213 fSelected->AddAt(fSelected->At(k-1),k);
1215 fSelected->AddAt(tx,j); // put track at empty position
1216 break; // go for next track
1221 ///////////////////////////////////////////////////////////////////////////
1222 Double_t AliJet::GetDistance(AliPosition* p,Float_t scale)
1224 // Provide distance of the current jet to the position p.
1225 // The error on the result can be obtained as usual by invoking
1226 // GetResultError() afterwards.
1228 // By default the distance will be provided in the metric unit scale of
1229 // the AliPosition p.
1230 // However, the user can select a different metric unit scale by
1231 // specification of the scale parameter.
1232 // The convention is that scale=1 corresponds to meter, so specification
1233 // of scale=0.01 will provide the distance in cm.
1234 // As such it is possible to obtain a correctly computed distance even in case
1235 // the jet parameters have a different unit scale.
1236 // However, it is recommended to work always with one single unit scale.
1238 // Note : In case of incomplete information, a distance value of -1 is
1244 if (!p) return dist;
1246 // Obtain a defined position on this jet
1247 AliPosition* rx=fRef;
1249 if (!rx) return dist;
1251 Ali3Vector pj=Get3Momentum();
1253 if (pj.GetNorm() <= 0.) return dist;
1256 tj.Set3Momentum(pj);
1257 tj.SetReferencePoint(*rx);
1258 dist=tj.GetDistance(p,scale);
1259 fDresult=tj.GetResultError();
1262 ///////////////////////////////////////////////////////////////////////////
1263 Double_t AliJet::GetDistance(AliTrack* t,Float_t scale)
1265 // Provide distance of the current jet to the track t.
1266 // The error on the result can be obtained as usual by invoking
1267 // GetResultError() afterwards.
1269 // By default the distance will be provided in the metric unit scale of
1271 // However, the user can specify a required metric unit scale by specification
1272 // of the scale parameter.
1273 // The convention is that scale=1 corresponds to meter, so specification
1274 // of scale=0.01 will provide the distance in cm.
1275 // As such it is possible to obtain a correctly computed distance even in case
1276 // the jet and track parameters have a different unit scale.
1277 // However, it is recommended to work always with one single unit scale.
1279 // Note : In case of incomplete information, a distance value of -1 is
1285 if (!t) return dist;
1287 // Obtain a defined position on this jet
1288 AliPosition* rx=fRef;
1290 if (!rx) return dist;
1292 Ali3Vector pj=Get3Momentum();
1294 if (pj.GetNorm() <= 0.) return dist;
1297 tj.Set3Momentum(pj);
1298 tj.SetReferencePoint(*rx);
1299 dist=tj.GetDistance(t,scale);
1300 fDresult=tj.GetResultError();
1303 ///////////////////////////////////////////////////////////////////////////
1304 Double_t AliJet::GetDistance(AliJet* j,Float_t scale)
1306 // Provide distance of the current jet to the jet j.
1307 // The error on the result can be obtained as usual by invoking
1308 // GetResultError() afterwards.
1310 // By default the distance will be provided in the metric unit scale of
1312 // This implies that the results of j1.GetDistance(j2) and j2.GetDistance(j1)
1313 // may be numerically different in case j1 and j2 have different metric units.
1314 // However, the user can specify a required metric unit scale by specification
1315 // of the scale parameter.
1316 // The convention is that scale=1 corresponds to meter, so specification
1317 // of scale=0.01 will provide the distance in cm.
1318 // As such it is possible to obtain a correctly computed distance even in case
1319 // the jet parameters have a different unit scale.
1320 // However, it is recommended to work always with one single unit scale.
1322 // Note : In case of incomplete information, a distance value of -1 is
1328 if (!j) return dist;
1330 // Obtain a defined position on jet j
1331 AliPosition* rx=j->GetReferencePoint();
1333 if (!rx) return dist;
1335 Ali3Vector pj=j->Get3Momentum();
1337 if (pj.GetNorm() <= 0.) return dist;
1340 tj.Set3Momentum(pj);
1341 tj.SetReferencePoint(*rx);
1342 dist=GetDistance(tj,scale);
1345 ///////////////////////////////////////////////////////////////////////////
1346 Int_t AliJet::GetNsignals() const
1348 // Provide the number of signals associated to the jet tracks.
1349 // Note : Multiple occurrences of the same signal are only counted once.
1351 if (fNtrk<1) return 0;
1358 for (Int_t i=1; i<=fNtrk; i++)
1361 for (Int_t j=1; j<=tx->GetNsignals(); j++)
1363 AliSignal* sx=tx->GetSignal(j);
1366 for (Int_t k=0; k<arr.GetEntries(); k++)
1368 if (sx==(AliSignal*)arr.At(k))
1374 if (!exists) arr.Add(sx);
1380 ///////////////////////////////////////////////////////////////////////////
1381 void AliJet::SetEscale(Float_t scale)
1383 // Indicate the energy/momentum scale as used by the user.
1384 // The convention is that scale=1 indicates values in units
1385 // of GeV, GeV/c or GeV/c**2.
1386 // So, in case one decides to store values in units of MeV, MeV/c or MeV/c**2
1387 // the scale indicator should be set to scale=0.001.
1389 // By default scale=1 is set in the constructor.
1397 cout << " *AliJet::SetEscale* Invalid scale value : " << scale << endl;
1400 ///////////////////////////////////////////////////////////////////////////
1401 Float_t AliJet::GetEscale() const
1403 // Provide the energy/momentum scale as used by the user.
1404 // The convention is that scale=1 indicates values in units
1405 // of GeV, GeV/c or GeV/c**2.
1406 // So, a value of scale=0.001 indicates that energy/momentum values are
1407 // stored in units of MeV, MeV/c or MeV/c**2.
1410 ///////////////////////////////////////////////////////////////////////////
1411 TObject* AliJet::Clone(const char* name) const
1413 // Make a deep copy of the current object and provide the pointer to the copy.
1414 // This memberfunction enables automatic creation of new objects of the
1415 // correct type depending on the object type, a feature which may be very useful
1416 // for containers when adding objects in case the container owns the objects.
1417 // This feature allows e.g. AliVertex to store either AliJet objects or
1418 // objects derived from AliJet via the AddJet memberfunction, provided
1419 // these derived classes also have a proper Clone memberfunction.
1421 AliJet* jet=new AliJet(*this);
1424 if (strlen(name)) jet->SetName(name);
1428 ///////////////////////////////////////////////////////////////////////////