KaliVeda  1.13/01
Heavy-Ion Analysis Toolkit
List of all members | Classes | Public Types | Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | Private Attributes | Static Private Attributes | Friends | List of all members
KVParticle Class Reference

Base class for relativistic kinematics of massive particles.

Implements all kinematical manipulations necessary for obtaining angles, kinetic energies, transverse energies etc. etc., defining and changing relativistic reference frames.

Unless otherwise stated,

This class derives from TLorentzVector - a particle is basically a Lorentz 4-vector with added attributes - and therefore any methods which are not documented here will be found in that class

Note
most 'Set'ter methods of TLorentzVector should not be used directly and have been made 'private' in KVParticle in order to prevent their misuse

This also means that all particles and kinematics are relativistic in KaliVeda.

Methods defined/redefined in this class:

GetMass() return mass of particle in MeV/c**2 (same as TLorentzVector::M())
GetMomentum() returns momentum 3-vector (TVector3) of particle (units MeV/c) (same as TLorentzVector::Vect())
GetKE()/GetEnergy() returns kinetic energy in MeV (same as TLorentzVector::E()-TLorentzVector::M())
GetVelocity()/GetV() returns velocity 3-vector (TVector3) of particle (units cm/ns)
GetVpar() returns velocity component in beam (z) direction (cm/ns) (same as GetV().Z())
GetVperp() same as GetV().Perp(), but sign is same as y-component of velocity, GetV().Y()
GetTheta() same as TLorentzVector::Theta() but in degrees, not radians
GetPhi() same as TLorentzVector::Phi() but in degrees, not radians, and always positive, between 0 and TMath::TwoPi()
static Double_t energy[]
Double_t GetREtran() const
Definition: KVParticle.h:650
TVector3 GetMomentum() const
Definition: KVParticle.h:606
Double_t GetTheta() const
Definition: KVParticle.h:682
Double_t GetEnergy() const
Definition: KVParticle.h:623
Double_t GetTransverseEnergy() const
Definition: KVParticle.h:627
Double_t GetPhi() const
Definition: KVParticle.h:690
Double_t GetRTransverseEnergy() const
Definition: KVParticle.h:645
Double_t GetKE() const
Definition: KVParticle.h:616
Double_t GetVpar() const
Definition: KVParticle.h:677
TVector3 GetV() const
Definition: KVParticle.h:673
Double_t GetMass() const
Definition: KVParticle.h:572
TVector3 GetVelocity() const
returns velocity vector in cm/ns units
Double_t GetVperp() const
Double_t GetEtran() const
Definition: KVParticle.h:631
Double_t Y() const
TVector3 Vect() const
Double_t M() const
Double_t Theta() const
Double_t Phi() const
Double_t Perp() const
Double_t E() const
Double_t Z() const
const long double MeV
energies
Definition: KVUnits.h:88
const long double cm
Definition: KVUnits.h:66

Particle properties can be defined either using one of the constructors :

KVParticle(Double_t m, TVector3 & p); // rest mass and momentum 3-vector
KVParticle(Double_t m, Double_t px, Double_t py, Double_t pz); // rest mass and Cartesian components of momentum 3-vector
double Double_t
const long double m
Definition: KVUnits.h:70

or with the usual 'Set' methods:

SetMass(Double_t m); // changes rest mass, leaves momentum unchanged
SetEnergy(Double_t e); // set kinetic energy (in MeV)
SetMomentum(const TVector3 & v); // changes momentum, leaves rest mass unchanged
SetMomentum(const TVector3 * v);
virtual void SetMass(Double_t m)
Definition: KVParticle.h:568
void SetMomentum(const TVector3 &v)
Definition: KVParticle.h:576
void SetEnergy(Double_t e)
Definition: KVParticle.h:601

or, for the following method:

SetMomentum(Double_t px, Double_t py, Double_t pz, Option_t * opt ="cart");
const char Option_t

Set momentum 3-vector, leaving rest mass unchanged, with :

Other methods:

SetMomentum(Double_t T, TVector3 dir); // set kinetic energy to T [MeV] and direction given by 'dir' unit vector
SetRandomMomentum(Double_t T, Double_t thmin, Double_t thmax, Double_t phmin, Double_t phmax, Option_t * opt = "isotropic");
/// a handy tool for giving random momenta to particles.
void SetRandomMomentum(Double_t T, Double_t thmin, Double_t thmax, Double_t phmin, Double_t phmax, Option_t *opt="isotropic")
Definition: KVParticle.cpp:136
Double_t T() const

GetEtran() and GetTransverseEnergy() return the non-relativistic transverse energy defined as the kinetic energy multiplied by the squared sinus of the polar angle: i.e. \(E_{tran} = E\sin^2(\theta)\) where \(\theta\) is the polar angle of the particle measured with respect to the beam (z-) axis.

GetREtran() and GetRTransverseEnergy() return the relativistic transverse energy which is frame invariant.

Kinematical reference frames

1. Accessing particle kinematics in different frames

Particle kinematics in different frames can be accessed with method InFrame(). This method does not modify the particle's kinematics or default reference frame (see 2 below). If you want to access lots of information from this frame, it is probably more efficient to define and store it in the particle's list of reference frames (see 3 below).

Example: inspect kinematics of particle (accessed through pointer KVParticle* p)

p->InFrame(TVector3(0,0,5)).GetVpar();
p->InFrame(rot).GetPhi();
TRotation & RotateZ(Double_t)
constexpr Double_t PiOver2()
p->InFrame(KVFrameTransform(TVector3(0,0,0.1),kTRUE)).GetKE();
p->InFrame({{0,0,0.1},kTRUE}).GetKE(); // C++11 equivalent of above
const Bool_t kTRUE
Utility class for kinematical transformations of KVParticle class.

In this case we need to pass two arguments to the KVFrameTransform constructor: before C++11, this required an explicit call to the constructor (top line); since C++11, uniform initialization allows to just give the expected arguments, using {} to enclose the arguments for any aggregate (class) types.

2. Modifying particle kinematics

Particle kinematics can be modified using method ChangeFrame():

Example: change kinematics of particle (accessed through pointer KVParticle* p)

p->ChangeFrame(TVector3(0,0,5));
p->ChangeFrame(rot);
p->ChangeFrame(KVFrameTransform(TVector3(0,0,0.1),kTRUE));
p->ChangeFrame({{0,0,0.1},kTRUE}); // C++11 equivalent of above

3. Using several reference frames

Rather than changing the reference frame of the particle, you can define and use several different reference frames while keeping the original kinematics as the default. Each frame can be used independently, and new frames can be defined based on any of the existing frames:

Example: (for a particle accessed through pointer KVParticle* p{.cpp}):

p->SetFrame("moving_frame", TVector3(0,0,5));
p->SetFrame("rotated_moving_frame", "moving_frame", rot);

Note that the same frame can be defined directly from the original particle by using a combined boost-then-rotation transform:

p->SetFrame("rotated_moving_frame", KVFrameTransform(TVector3(0,0,5),rot));
p->SetFrame("rotated_moving_frame", {{0,0,5},rot}); // C++11 equivalent of above
p->SetFrame("rotated_frame", rot);
p->GetFrame("moving_frame")->GetVpar();
p->GetFrame("rotated_frame")->GetPhi();
p->GetFrame("rotated_moving_frame")->GetTransverseEnergy();

Note that the frame "rotated_moving_frame" is directly accessible even if it is defined in two steps as a rotation of the "moving_frame".

The total number of defined frames (including the default frame) is given by GetNumberOfDefinedFrames(). Calling Print() will show all reference frames defined for the particle:

p->Print()
KVParticle mass=0 Theta=0 Phi=0 KE=0 Vpar=0
moving_frame: Theta=0 Phi=0 KE=0 Vpar=0
rotated_moving_frame: Theta=0 Phi=0 KE=0 Vpar=0
rotated_frame: Theta=0 Phi=0 KE=0 Vpar=0
Base class for relativistic kinematics of massive particles.
Definition: KVParticle.h:398

Indentation indicates the relationships between frames: "rotated_moving_frame" is a child frame of "moving_frame". The first line is the default kinematics, which by default has no name (a name can be set using SetFrameName()).

The current default kinematics can always be accessed using GetCurrentDefaultKinematics(): this method returns the default kinematical frame for the particle when accessed from any kinematical frame:

p->SetFrameName("lab"); // set name for default kinematics
auto show_default_frame = [](const KVParticle* p){ std::cout << p->GetCurrentDefaultKinematics()->GetFrameName() << std::endl; };
show_default_frame( p->GetFrame("moving_frame") )
"lab"

4. Changing the default kinematics

Let us consider a particle for which the different reference frames in the previous paragraph have been defined. For an example, imagine that the default kinematics are that of particle with rest mass 939 MeV moving with a momentum of 250 MeV/c at an angle of \(45^o\) to the +ve z-direction:

KVParticle p(939,0,0,250);
p.SetTheta(45);
p.SetFrame("moving_frame", TVector3(0,0,5));
p.SetFrame("rotated_moving_frame", "moving_frame", rot);
p.SetFrame("rotated_frame", rot);

Now if we want to change the default kinematical frame for this particle:

p.ChangeDefaultFrame("rotated_moving_frame", "lab"); // second argument is name for previous default frame, if not defined before
p.Print();
KVParticle mass=939 Theta=85.1751 Phi=270 KE=16.6117 Vpar=0.468125
moving_frame: Theta=85.1751 Phi=0 KE=16.6117 Vpar=0.468125
lab: Theta=45 Phi=0 KE=32.7103 Vpar=5.45392
rotated_frame: Theta=45 Phi=270 KE=32.7103 Vpar=5.45392
KVNameValueList::ParticleParameters : Parameters associated with a particle in an event (0x7f5a1ff8b1b8)
<frameName=rotated_moving_frame>
RooCmdArg Parameters(const RooArgSet &params)
auto * a

The relationships between frames are preserved, i.e. if we present the frames as graphs:

with "lab" as default frame:

lab
|
+--moving_frame
| |
| +--rotated_moving_frame
|
+--rotated_frame

with "rotated_moving_frame" as default frame:

rotated_moving_frame
|
+--moving_frame
|
+--lab
|
+--rotated_frame

5. Updating stored kinematical frames

If you call SetFrame() several times with the same frame name [note that frame names are case insensitive], the existing reference frame will be updated to use the new transformation, which will be applied to the kinematics of the particle in the 'parent' frame used to define the frame. Any frames which were defined based on the frame will be updated too.

Example: for the previous particle & frame definitions, after resetting the default kinematics:

p.ChangeDefaultFrame("lab");
p.SetFrame("rotated_moving_frame", rot);
p.Print();
KVParticle mass=939 Theta=45 Phi=0 KE=32.7103 Vpar=5.45392
rotated_frame: Theta=45 Phi=270 KE=32.7103 Vpar=5.45392
moving_frame: Theta=85.1751 Phi=0 KE=16.6117 Vpar=0.468125
rotated_moving_frame: Theta=45 Phi=315 KE=32.7103 Vpar=5.45392
KVNameValueList::ParticleParameters : Parameters associated with a particle in an event (0x7faae41241b8)
<frameName=lab>
virtual void Print(Option_t *option="") const
constexpr Double_t PiOver4()
p.SetFrame("moving_frame", KVFrameTransform(TVector3(0,0,0.1),kTRUE));
p.SetFrame("moving_frame", {{0,0,0.1},kTRUE});
p.Print();
KVParticle mass=939 Theta=45 Phi=0 KE=32.7103 Vpar=5.45392
rotated_frame: Theta=45 Phi=270 KE=32.7103 Vpar=5.45392
moving_frame: Theta=65.6491 Phi=0 KE=19.8389 Vpar=2.50151
rotated_moving_frame: Theta=65.6491 Phi=315 KE=19.8389 Vpar=2.50151
KVNameValueList::ParticleParameters : Parameters associated with a particle in an event (0x7faae41241b8)
<frameName=lab>

Note that in this case, the "rotated_moving_frame" is updated automatically to take account of the new velocity of "moving_frame".

However, if you change the kinematics of the particle in its original (default) frame, you have to update the kinematics in all defined frames by hand, it does not occur automatically:

p.SetTheta(30);
p.Print();
KVParticle mass=939 Theta=30 Phi=0 KE=32.7103 Vpar=6.67966
rotated_frame: Theta=45 Phi=270 KE=32.7103 Vpar=5.45392
moving_frame: Theta=65.6491 Phi=0 KE=19.8389 Vpar=2.50151
rotated_moving_frame: Theta=65.6491 Phi=315 KE=19.8389 Vpar=2.50151
KVNameValueList::ParticleParameters : Parameters associated with a particle in an event (0x7faae41241b8)
<frameName=lab>
p.UpdateAllFrames();
p.Print();
KVParticle mass=939 Theta=30 Phi=0 KE=32.7103 Vpar=6.67966
rotated_frame: Theta=30 Phi=270 KE=32.7103 Vpar=6.67966
moving_frame: Theta=46.1843 Phi=0 KE=15.8459 Vpar=3.76564
rotated_moving_frame: Theta=46.1843 Phi=315 KE=15.8459 Vpar=3.76564
KVNameValueList::ParticleParameters : Parameters associated with a particle in an event (0x7faae41241b8)
<frameName=lab>

Definition of groups

AddGroup() and BelongsToGroup() methods allow to sort particles in an event into subsets based on various criteria such as "emitted by QP", "backwards emitted", or "include in calorimetry". The group attribution can be done 'by hand' or using a KVParticleCondition object to define a selection. The number of such groups for a given particle is unlimited.

Group names are case insensitive. Example:

part.AddGroup("forward");
part.BelongsToGroup("ForWaRD");// -> return kTRUE
void AddGroup(const Char_t *groupname, const Char_t *from="") const
Definition: KVParticle.cpp:436
Bool_t BelongsToGroup(const Char_t *groupname) const
Definition: KVParticle.cpp:509

Groups are valid for all kinematical frames, and can be defined using a pointer to the nucleus in any of its defined kinematical frames:

part.SetFrame("CM", [some transformation]);
part.AddGroup("toto");
part.GetFrame("CM")->BelongsToGroup("toto");// returns kTRUE
part.GetFrame("CM")->AddGroup("titi");
part.BelongsToGroup("titi");// returns kRUE
void SetFrame(const Char_t *frame, const KVFrameTransform &)
Definition: KVParticle.cpp:741
KVParticle const * GetFrame(const Char_t *frame, Bool_t warn_and_return_null_if_unknown=kTRUE) const
Definition: KVParticle.cpp:855

Definition at line 398 of file KVParticle.h.

Classes

class  FrameList
 

Public Types

enum  { kIsOK = BIT(14) , kIsOKSet = BIT(15) , kIsDetected = BIT(16) }
 
- Public Types inherited from TLorentzVector
typedef Double_t Scalar
 
- Public Types inherited from TObject
enum  EDeprecatedStatusBits
 
enum  EStatusBits
 

Public Member Functions

 KVParticle ()
 
 KVParticle (const KVParticle &)
 copy ctor More...
 
 KVParticle (Double_t m, Double_t px, Double_t py, Double_t pz)
 create particle with given mass and momentum vector More...
 
 KVParticle (Double_t m, TVector3 &p)
 create particle with given mass and momentum vector More...
 
virtual ~ KVParticle ()
 
void AddGroup (const Char_t *groupname, const Char_t *from="") const
 
Bool_t BelongsToGroup (const Char_t *groupname) const
 
void ChangeDefaultFrame (const Char_t *, const Char_t *defname="")
 
void ChangeFrame (const KVFrameTransform &, const KVString &="")
 
virtual void Clear (Option_t *opt="")
 Reset particle properties i.e. before creating/reading a new event. More...
 
virtual void Copy (TObject &) const
 
Double_t GetCosTheta () const
 
const KVParticleGetCurrentDefaultKinematics () const
 
Double_t GetE () const
 
Double_t GetElong () const
 
Double_t GetEnergy () const
 
Double_t GetEtran () const
 
KVParticle const * GetFrame (const Char_t *frame, Bool_t warn_and_return_null_if_unknown=kTRUE) const
 
const Char_tGetFrameName (void) const
 
KVUniqueNameListGetGroups () const
 
Double_t GetKE () const
 
KVListGetListOfFrames () const
 
Double_t GetLongitudinalEnergy () const
 
Double_t GetMass () const
 
TVector3 GetMomentum () const
 
const Char_tGetName () const
 return the field fName More...
 
Int_t GetNumberOfDefinedFrames () const
 
Int_t GetNumberOfDefinedGroups () const
 
KVNameValueListGetParameters () const
 
Double_t GetPhi () const
 
TVector3GetPInitial () const
 
Double_t GetREtran () const
 
Double_t GetRTransverseEnergy () const
 
Double_t GetThermalWaveLength (Double_t temp) const
 
Double_t GetTheta () const
 
Double_t GetTransverseEnergy () const
 
TVector3 GetTransverseMomentum () const
 
TVector3 GetV () const
 
TVector3 GetVelocity () const
 returns velocity vector in cm/ns units More...
 
Double_t GetVpar () const
 
Double_t GetVperp () const
 
Double_t GetWaveLength () const
 
Bool_t HasFrame (const Char_t *frame) const
 
KVParticle InFrame (const KVFrameTransform &)
 
void init ()
 default initialisation More...
 
Bool_t IsDefaultKinematics () const
 
Bool_t IsDetected () const
 
Bool_t IsOK () const
 
void ListGroups () const
 List all stored groups. More...
 
void ls (Option_t *option="") const
 
KVParticleoperator= (const KVParticle &rhs)
 KVParticle assignment operator. More...
 
virtual void Print (Option_t *t="") const
 print out characteristics of particle More...
 
void RemoveAllGroups ()
 Remove all groups. More...
 
void RemoveGroup (const Char_t *groupname)
 Remove group from list of groups. More...
 
void ResetEnergy ()
 
void ResetIsOK ()
 
void Set4Mom (const TLorentzVector &p)
 
void SetE (Double_t a)
 
void SetE0 (TVector3 *e=0)
 
void SetEnergy (Double_t e)
 
void SetFrame (const Char_t *frame, const KVFrameTransform &)
 
void SetFrame (const Char_t *newframe, const Char_t *oldframe, const KVFrameTransform &)
 
void SetFrameName (const Char_t *framename)
 
void SetIsDetected ()
 
void SetIsOK (Bool_t flag=kTRUE)
 
void SetKE (Double_t ecin)
 
virtual void SetMass (Double_t m)
 
void SetMomentum (const TVector3 &v)
 
void SetMomentum (const TVector3 *v)
 
void SetMomentum (Double_t px, Double_t py, Double_t pz, Option_t *opt="cart")
 
void SetMomentum (Double_t T, TVector3 dir)
 
void SetName (const Char_t *nom)
 Set Name of the particle. More...
 
template<typename ValType >
void SetParameter (const Char_t *name, ValType value) const
 
void SetPhi (Double_t phi)
 
void SetRandomMomentum (Double_t T, Double_t thmin, Double_t thmax, Double_t phmin, Double_t phmax, Option_t *opt="isotropic")
 
void SetTheta (Double_t theta)
 
void SetVelocity (const TVector3 &)
 Set velocity of particle (in cm/ns units) More...
 
void UpdateAllFrames ()
 
- Public Member Functions inherited from TLorentzVector
 TLorentzVector ()
 
 TLorentzVector (const Double_t *carray)
 
 TLorentzVector (const Float_t *carray)
 
 TLorentzVector (const TLorentzVector &lorentzvector)
 
 TLorentzVector (const TVector3 &vector3, Double_t t)
 
 TLorentzVector (Double_t x, Double_t y, Double_t z, Double_t t)
 
virtual ~TLorentzVector ()
 
Double_t Angle (const TVector3 &v) const
 
Double_t Beta () const
 
void Boost (const TVector3 &)
 
void Boost (Double_t, Double_t, Double_t)
 
TVector3 BoostVector () const
 
Double_t CosTheta () const
 
Double_t DeltaPhi (const TLorentzVector &) const
 
Double_t DeltaR (const TLorentzVector &, Bool_t useRapidity=kFALSE) const
 
Double_t Dot (const TLorentzVector &) const
 
Double_t DrEtaPhi (const TLorentzVector &) const
 
Double_t DrRapidityPhi (const TLorentzVector &) const
 
Double_t E () const
 
Double_t Energy () const
 
Double_t Et () const
 
Double_t Et (const TVector3 &) const
 
Double_t Et2 () const
 
Double_t Et2 (const TVector3 &) const
 
Double_t Eta () const
 
TVector2 EtaPhiVector ()
 
Double_t Gamma () const
 
void GetXYZT (Double_t *carray) const
 
void GetXYZT (Float_t *carray) const
 
Double_t M () const
 
Double_t M2 () const
 
Double_t Mag () const
 
Double_t Mag2 () const
 
Double_t Minus () const
 
Double_t Mt () const
 
Double_t Mt2 () const
 
Bool_t operator!= (const TLorentzVector &) const
 
Double_toperator() (int i)
 
Double_t operator() (int i) const
 
Double_t operator* (const TLorentzVector &) const
 
TLorentzVector operator* (Double_t a) const
 
TLorentzVectoroperator*= (const TLorentzRotation &)
 
TLorentzVectoroperator*= (const TRotation &)
 
TLorentzVectoroperator*= (Double_t a)
 
TLorentzVector operator+ (const TLorentzVector &) const
 
TLorentzVectoroperator+= (const TLorentzVector &)
 
TLorentzVector operator- () const
 
TLorentzVector operator- (const TLorentzVector &) const
 
TLorentzVectoroperator-= (const TLorentzVector &)
 
TLorentzVectoroperator= (const TLorentzVector &)
 
Bool_t operator== (const TLorentzVector &) const
 
Double_toperator[] (int i)
 
Double_t operator[] (int i) const
 
Double_t P () const
 
Double_t Perp () const
 
Double_t Perp (const TVector3 &v) const
 
Double_t Perp2 () const
 
Double_t Perp2 (const TVector3 &v) const
 
Double_t Phi () const
 
Double_t Plus () const
 
Double_t PseudoRapidity () const
 
Double_t Pt () const
 
Double_t Pt (const TVector3 &v) const
 
Double_t Px () const
 
Double_t Py () const
 
Double_t Pz () const
 
Double_t Rapidity () const
 
Double_t Rho () const
 
void Rotate (Double_t, const TVector3 &)
 
void RotateUz (const TVector3 &newUzVector)
 
void RotateX (Double_t angle)
 
void RotateY (Double_t angle)
 
void RotateZ (Double_t angle)
 
void SetE (Double_t a)
 
void SetPerp (Double_t)
 
void SetPhi (Double_t phi)
 
void SetPtEtaPhiE (Double_t pt, Double_t eta, Double_t phi, Double_t e)
 
void SetPtEtaPhiM (Double_t pt, Double_t eta, Double_t phi, Double_t m)
 
void SetPx (Double_t a)
 
void SetPxPyPzE (Double_t px, Double_t py, Double_t pz, Double_t e)
 
void SetPy (Double_t a)
 
void SetPz (Double_t a)
 
void SetRho (Double_t rho)
 
void SetT (Double_t a)
 
void SetTheta (Double_t theta)
 
void SetVect (const TVector3 &vect3)
 
void SetVectM (const TVector3 &spatial, Double_t mass)
 
void SetVectMag (const TVector3 &spatial, Double_t magnitude)
 
void SetX (Double_t a)
 
void SetXYZM (Double_t x, Double_t y, Double_t z, Double_t m)
 
void SetXYZT (Double_t x, Double_t y, Double_t z, Double_t t)
 
void SetY (Double_t a)
 
void SetZ (Double_t a)
 
Double_t T () const
 
Double_t Theta () const
 
TLorentzVectorTransform (const TLorentzRotation &)
 
TLorentzVectorTransform (const TRotation &)
 
TVector3 Vect () const
 
Double_t X () const
 
Double_t Y () const
 
Double_t Z () const
 
- Public Member Functions inherited from TObject
 TObject ()
 
 TObject (const TObject &object)
 
virtual ~TObject ()
 
void AbstractMethod (const char *method) const
 
virtual void AppendPad (Option_t *option="")
 
virtual void Browse (TBrowser *b)
 
ULong_t CheckedHash ()
 
virtual const char * ClassName () const
 
virtual TObjectClone (const char *newname="") const
 
virtual Int_t Compare (const TObject *obj) const
 
virtual void Delete (Option_t *option="")
 
virtual Int_t DistancetoPrimitive (Int_t px, Int_t py)
 
virtual void Draw (Option_t *option="")
 
virtual void DrawClass () const
 
virtual TObjectDrawClone (Option_t *option="") const
 
virtual void Dump () const
 
virtual void Error (const char *method, const char *msgfmt,...) const
 
virtual void Execute (const char *method, const char *params, Int_t *error=0)
 
virtual void Execute (TMethod *method, TObjArray *params, Int_t *error=0)
 
virtual void ExecuteEvent (Int_t event, Int_t px, Int_t py)
 
virtual void Fatal (const char *method, const char *msgfmt,...) const
 
virtual TObjectFindObject (const char *name) const
 
virtual TObjectFindObject (const TObject *obj) const
 
virtual Option_tGetDrawOption () const
 
virtual const char * GetIconName () const
 
virtual char * GetObjectInfo (Int_t px, Int_t py) const
 
virtual Option_tGetOption () const
 
virtual const char * GetTitle () const
 
virtual UInt_t GetUniqueID () const
 
virtual Bool_t HandleTimer (TTimer *timer)
 
virtual ULong_t Hash () const
 
Bool_t HasInconsistentHash () const
 
virtual void Info (const char *method, const char *msgfmt,...) const
 
virtual Bool_t InheritsFrom (const char *classname) const
 
virtual Bool_t InheritsFrom (const TClass *cl) const
 
virtual void Inspect () const
 
void InvertBit (UInt_t f)
 
virtual Bool_t IsEqual (const TObject *obj) const
 
virtual Bool_t IsFolder () const
 
R__ALWAYS_INLINE Bool_t IsOnHeap () const
 
virtual Bool_t IsSortable () const
 
R__ALWAYS_INLINE Bool_t IsZombie () const
 
void MayNotUse (const char *method) const
 
virtual Bool_t Notify ()
 
void Obsolete (const char *method, const char *asOfVers, const char *removedFromVers) const
 
void operator delete (void *ptr)
 
void operator delete[] (void *ptr)
 
voidoperator new (size_t sz)
 
voidoperator new (size_t sz, void *vp)
 
voidoperator new[] (size_t sz)
 
voidoperator new[] (size_t sz, void *vp)
 
TObjectoperator= (const TObject &rhs)
 
virtual void Paint (Option_t *option="")
 
virtual void Pop ()
 
virtual Int_t Read (const char *name)
 
virtual void RecursiveRemove (TObject *obj)
 
void ResetBit (UInt_t f)
 
virtual void SaveAs (const char *filename="", Option_t *option="") const
 
virtual void SavePrimitive (std::ostream &out, Option_t *option="")
 
void SetBit (UInt_t f)
 
void SetBit (UInt_t f, Bool_t set)
 
virtual void SetDrawOption (Option_t *option="")
 
virtual void SetUniqueID (UInt_t uid)
 
virtual void SysError (const char *method, const char *msgfmt,...) const
 
R__ALWAYS_INLINE Bool_t TestBit (UInt_t f) const
 
Int_t TestBits (UInt_t f) const
 
virtual void UseCurrentStyle ()
 
virtual void Warning (const char *method, const char *msgfmt,...) const
 
virtual Int_t Write (const char *name=0, Int_t option=0, Int_t bufsize=0)
 
virtual Int_t Write (const char *name=0, Int_t option=0, Int_t bufsize=0) const
 

Static Public Member Functions

static Double_t C ()
 
- Static Public Member Functions inherited from TObject
static Longptr_t GetDtorOnly ()
 
static Bool_t GetObjectStat ()
 
static void SetDtorOnly (void *obj)
 
static void SetObjectStat (Bool_t stat)
 

Protected Member Functions

void AddGroups (KVUniqueNameList *un)
 list of groups added to the current one More...
 
void SetGroups (KVUniqueNameList *un)
 Define for the particle a new list of groups. More...
 
- Protected Member Functions inherited from TObject
virtual void DoError (int level, const char *location, const char *fmt, va_list va) const
 
void MakeZombie ()
 

Protected Attributes

TVector3fE0
 the momentum of the particle before it is slowed/stopped by an absorber More...
 
KVNameValueList fParameters
 a general-purpose list of parameters associated with this particle More...
 
- Protected Attributes inherited from TObject
 kOnlyPrepStep
 

Private Member Functions

Int_t _GetNumberOfDefinedFrames () const
 used to inhibit full copying of particles in different kinematical frames More...
 
KVKinematicalFrameget_frame (const Char_t *) const
 
KVKinematicalFrameget_parent_frame (const Char_t *, KVKinematicalFrame *F=nullptr) const
 
KVParticleGetOriginal ()
 
const KVParticleGetOriginal () const
 
KVParticleGetParentFrame () const
 
const KVParticleGetTopmostParentFrame () const
 
void print_frames (TString fmt="") const
 recursive print out of all defined kinematical frames More...
 
void ResetFrameCopyOnly () const
 
void SetFrameCopyOnly () const
 
void SetOriginal (KVParticle *p)
 
void SetParentFrame (KVParticle *p)
 
void SetPerp (Double_t p)
 
void SetPtEtaPhiE (Double_t pt, Double_t eta, Double_t phi, Double_t e)
 
void SetPtEtaPhiM (Double_t pt, Double_t eta, Double_t phi, Double_t m)
 
void SetPx (Double_t a)
 
void SetPxPyPzE (Double_t px, Double_t py, Double_t pz, Double_t e)
 
void SetPy (Double_t a)
 
void SetPz (Double_t a)
 
void SetRho (Double_t rho)
 
void SetT (Double_t a)
 
void SetVect (const TVector3 &vect3)
 TLorentzVector setters should not be used. More...
 
void SetVectM (const TVector3 &spatial, Double_t mass)
 
void SetVectMag (const TVector3 &spatial, Double_t magnitude)
 
void SetX (Double_t a)
 
void SetXYZM (Double_t x, Double_t y, Double_t z, Double_t m)
 
void SetXYZT (Double_t x, Double_t y, Double_t z, Double_t t)
 
void SetY (Double_t a)
 
void SetZ (Double_t a)
 

Private Attributes

FrameList fBoosted {this}
 
Bool_t fFrameCopyOnly {kFALSE}
 if != nullptr, this object is a representation of the particle in a kinematical frame More...
 
TString fFrameName
 non-persistent frame name field, sets when calling SetFrame method More...
 
KVUniqueNameList fGroups
 list of momenta of the particle in different Lorentz-boosted frames More...
 
TString fName
 non-persistent name field - Is useful More...
 
KVParticlefOriginal {nullptr}
 parent kinematical frame More...
 
KVParticlefParentFrame {nullptr}
 

Static Private Attributes

static Double_t kSpeedOfLight = TMath::C() * 1.e-07
 speed of light in cm/ns More...
 

Friends

class KVKinematicalFrame
 

Additional Inherited Members

- Public Attributes inherited from TLorentzVector
 kNUM_COORDINATES
 
 kSIZE
 
 kT
 
 kX
 
 kY
 
 kZ
 
- Public Attributes inherited from TObject
 kBitMask
 
 kCanDelete
 
 kCannotPick
 
 kHasUUID
 
 kInconsistent
 
 kInvalidObject
 
 kIsOnHeap
 
 kIsReferenced
 
 kMustCleanup
 
 kNoContextMenu
 
 kNotDeleted
 
 kObjInCanvas
 
 kOverwrite
 
 kSingleKey
 
 kWriteDelete
 
 kZombie
 

#include <KVParticle.h>

Inheritance diagram for KVParticle:
Inheritance graph
[legend]

Member Enumeration Documentation

◆ anonymous enum

anonymous enum
Enumerator
kIsOK 
kIsOKSet 
kIsDetected 

Definition at line 551 of file KVParticle.h.

Constructor & Destructor Documentation

◆ KVParticle() [1/4]

KVParticle::KVParticle ( )

Definition at line 40 of file KVParticle.cpp.

◆ KVParticle() [2/4]

KVParticle::KVParticle ( Double_t  m,
TVector3 p 
)

create particle with given mass and momentum vector

Definition at line 79 of file KVParticle.cpp.

◆ KVParticle() [3/4]

KVParticle::KVParticle ( Double_t  m,
Double_t  px,
Double_t  py,
Double_t  pz 
)

create particle with given mass and momentum vector

Definition at line 92 of file KVParticle.cpp.

◆ KVParticle() [4/4]

KVParticle::KVParticle ( const KVParticle obj)

copy ctor

Definition at line 63 of file KVParticle.cpp.

◆ ~ KVParticle()

virtual KVParticle::~ KVParticle ( )
virtual

Member Function Documentation

◆ _GetNumberOfDefinedFrames()

Int_t KVParticle::_GetNumberOfDefinedFrames ( ) const
private

used to inhibit full copying of particles in different kinematical frames

Returns
the number of kinematical frames defined by transformations of this frame
Note
this is always at least equal to one (we count the present frame of the particle)

Definition at line 487 of file KVParticle.cpp.

◆ AddGroup()

void KVParticle::AddGroup ( const Char_t groupname,
const Char_t from = "" 
) const

Associate this particle with the given named group. Optional argument "from" allows to put a condition on the already stored group list, is set to "" by default

Definition at line 436 of file KVParticle.cpp.

◆ AddGroups()

void KVParticle::AddGroups ( KVUniqueNameList un)
protected

list of groups added to the current one

Definition at line 469 of file KVParticle.cpp.

◆ BelongsToGroup()

Bool_t KVParticle::BelongsToGroup ( const Char_t groupname) const

Check if particle belong to a given group return kTRUE if groupname="". return kFALSE if no group has be defined

Definition at line 509 of file KVParticle.cpp.

◆ C()

Double_t KVParticle::C ( )
static

Static function. Returns speed of light in cm/ns units.

Definition at line 117 of file KVParticle.cpp.

◆ ChangeDefaultFrame()

void KVParticle::ChangeDefaultFrame ( const Char_t newdef,
const Char_t defname = "" 
)

Make existing reference frame 'newdef' the new default frame for particle kinematics. The current default frame will then be accessible from the list of frames using its name (if set with SetFrameName).

You can change/set the name of the previous default frame with 'defname'. If no name was set and none given, it will be renamed "default" by default.

Definition at line 622 of file KVParticle.cpp.

◆ ChangeFrame()

void KVParticle::ChangeFrame ( const KVFrameTransform t,
const KVString name = "" 
)

Permanently modify kinematics of particle according to the given transformation. You can optionally set the name of this new default kinematics. NB the current kinematics will be lost. If you want to keep it after changing the default kinematics, define the new frame with SetFrame and then use ChangeDefaultFrame.

Definition at line 601 of file KVParticle.cpp.

◆ Clear()

void KVParticle::Clear ( Option_t opt = "")
virtual

Reset particle properties i.e. before creating/reading a new event.

Reimplemented from TObject.

Reimplemented in KVINDRAReconNuc, KVFAZIAReconNuc, KVReconstructedNucleus, and KVNucleus.

Definition at line 295 of file KVParticle.cpp.

◆ Copy()

void KVParticle::Copy ( TObject obj) const
virtual

Copy this to obj Particle kinematics are copied using operator=(const KVParticle&) List of particle's groups is copied The particle's name is copied The list of parameters associated with the particle is copied

Reimplemented from TObject.

Reimplemented in KVSimNucleus, KVNucleus, KVReconstructedNucleus, KVINDRAReconNuc, and KVFAZIAReconNuc.

Definition at line 269 of file KVParticle.cpp.

◆ get_frame()

KVKinematicalFrame * KVParticle::get_frame ( const Char_t frame) const
private
Parameters
[in]framename of a kinematical frame
Returns
pointer to the kinematical frame if previously defined, nullptr otherwise
Note
PRIVATE method for internal use only. This method allows to modify the returned frame, i.e. in order to define new frames in SetFrame()

Definition at line 937 of file KVParticle.cpp.

◆ get_parent_frame()

KVKinematicalFrame * KVParticle::get_parent_frame ( const Char_t f,
KVKinematicalFrame F = nullptr 
) const
private

PRIVATE method for internal use only Returns pointer to parent frame of 'f'

Definition at line 965 of file KVParticle.cpp.

◆ GetCosTheta()

Double_t KVParticle::GetCosTheta ( ) const
inline

Definition at line 686 of file KVParticle.h.

◆ GetCurrentDefaultKinematics()

const KVParticle* KVParticle::GetCurrentDefaultKinematics ( ) const
inline
Returns
pointer to current default kinematics for particle

Definition at line 772 of file KVParticle.h.

◆ GetE()

Double_t KVParticle::GetE ( ) const
inline

Definition at line 654 of file KVParticle.h.

◆ GetElong()

Double_t KVParticle::GetElong ( ) const
inline

Definition at line 640 of file KVParticle.h.

◆ GetEnergy()

Double_t KVParticle::GetEnergy ( ) const
inline

Definition at line 623 of file KVParticle.h.

◆ GetEtran()

Double_t KVParticle::GetEtran ( ) const
inline

Definition at line 631 of file KVParticle.h.

◆ GetFrame()

KVParticle const * KVParticle::GetFrame ( const Char_t frame,
Bool_t  warn_and_return_null_if_unknown = kTRUE 
) const

Return the momentum of the particle in the Lorentz-boosted frame corresponding to the name "frame" given as argument (see SetFrame() for definition of different frames). If the default frame name has been set (see KVEvent::SetFrameName) and 'frame' is the name of this default frame (KVParticle::fFrameName), we return the address of the particle itself.

This frame may have been defined by a direct transformation of the original kinematics of the particle (using SetFrame(newframe,...)) or by a transformation of the kinematics in another user-defined frame (using SetFrame(newframe,oldframe,...)).

Note that frames are not "dynamic": if any changes are made to the original particle's kinematics after definition of a frame, if you want these changes to affect also the other frames you need to update them by hand by calling KVParticle::UpdateAllFrames().

Frame names are case insensitive: "CM" or "cm" or "Cm" are all good...

By default, if no frame with the given name is found, we return nullptr and print a warning. If warn_and_return_null_if_unknown=kFALSE, we return the address of the particle itself, i.e. the original/default kinematics. [Note that this is an inversion of the previous default behaviour]

Note that the properties of the particle returned by this method can not be modified: this is deliberate, as any modifications e.g. to kinematics will have no effect in any other frames.

The returned pointer corresponds to a "pseudoparticle" in the desired frame, therefore you can use any KVParticle method in order to access the kinematics of the particle in the boosted frame, e.g.

(...supposing a valid pointer KVParticle* my_part...)
my_part->GetFrame("cm_frame")->GetVpar();// //el velocity in "cm_frame"
my_part->GetFrame("QP_frame")->GetTheta();// polar angle in "QP_frame"
etc. etc.

Definition at line 855 of file KVParticle.cpp.

◆ GetFrameName()

const Char_t* KVParticle::GetFrameName ( void  ) const
inline
Returns
name of the current kinematical frame

Definition at line 801 of file KVParticle.h.

◆ GetGroups()

KVUniqueNameList* KVParticle::GetGroups ( ) const
inline
Returns
the list of groups

Definition at line 758 of file KVParticle.h.

◆ GetKE()

Double_t KVParticle::GetKE ( ) const
inline

return (E() - M());

Definition at line 616 of file KVParticle.h.

◆ GetListOfFrames()

KVList* KVParticle::GetListOfFrames ( ) const
inline

Definition at line 711 of file KVParticle.h.

◆ GetLongitudinalEnergy()

Double_t KVParticle::GetLongitudinalEnergy ( ) const
inline

Definition at line 636 of file KVParticle.h.

◆ GetMass()

Double_t KVParticle::GetMass ( ) const
inline

Definition at line 572 of file KVParticle.h.

◆ GetMomentum()

TVector3 KVParticle::GetMomentum ( ) const
inline

Definition at line 606 of file KVParticle.h.

◆ GetName()

const Char_t * KVParticle::GetName ( ) const
virtual

return the field fName

Reimplemented from TObject.

Definition at line 423 of file KVParticle.cpp.

◆ GetNumberOfDefinedFrames()

Int_t KVParticle::GetNumberOfDefinedFrames ( ) const
inline
Returns
the total number of kinematic frames defined for this particle

Definition at line 794 of file KVParticle.h.

◆ GetNumberOfDefinedGroups()

Int_t KVParticle::GetNumberOfDefinedGroups ( ) const
inline
Returns
the number of defined groups for the particle

Definition at line 753 of file KVParticle.h.

◆ GetOriginal() [1/2]

KVParticle* KVParticle::GetOriginal ( )
inlineprivate

Definition at line 536 of file KVParticle.h.

◆ GetOriginal() [2/2]

const KVParticle* KVParticle::GetOriginal ( ) const
inlineprivate

Definition at line 532 of file KVParticle.h.

◆ GetParameters()

KVNameValueList* KVParticle::GetParameters ( ) const
inline

Definition at line 816 of file KVParticle.h.

◆ GetParentFrame()

KVParticle* KVParticle::GetParentFrame ( ) const
inlineprivate

Definition at line 513 of file KVParticle.h.

◆ GetPhi()

Double_t KVParticle::GetPhi ( ) const
inline
Examples
ExampleCorrelationAnalysis.cpp.

Definition at line 690 of file KVParticle.h.

◆ GetPInitial()

TVector3* KVParticle::GetPInitial ( ) const
inline

Definition at line 728 of file KVParticle.h.

◆ GetREtran()

Double_t KVParticle::GetREtran ( ) const
inline

Definition at line 650 of file KVParticle.h.

◆ GetRTransverseEnergy()

Double_t KVParticle::GetRTransverseEnergy ( ) const
inline

Definition at line 645 of file KVParticle.h.

◆ GetThermalWaveLength()

Double_t KVParticle::GetThermalWaveLength ( Double_t  temp) const
inline

should be in fm

Definition at line 666 of file KVParticle.h.

◆ GetTheta()

Double_t KVParticle::GetTheta ( ) const
inline

Definition at line 682 of file KVParticle.h.

◆ GetTopmostParentFrame()

const KVParticle* KVParticle::GetTopmostParentFrame ( ) const
inlineprivate
Returns
pointer to current default kinematics for particle, i.e. the one at the top of the graph of kinematical frames which has no parent frame

Definition at line 521 of file KVParticle.h.

◆ GetTransverseEnergy()

Double_t KVParticle::GetTransverseEnergy ( ) const
inline

Definition at line 627 of file KVParticle.h.

◆ GetTransverseMomentum()

TVector3 KVParticle::GetTransverseMomentum ( ) const
inline

Definition at line 610 of file KVParticle.h.

◆ GetV()

TVector3 KVParticle::GetV ( ) const
inline

Definition at line 673 of file KVParticle.h.

◆ GetVelocity()

TVector3 KVParticle::GetVelocity ( ) const

returns velocity vector in cm/ns units

Definition at line 1013 of file KVParticle.cpp.

◆ GetVpar()

Double_t KVParticle::GetVpar ( ) const
inline

◆ GetVperp()

Double_t KVParticle::GetVperp ( ) const

returns transverse velocity in cm/ns units sign is +ve if py>=0, -ve if py<0

Examples
ExampleE789ReconAnalysis.cpp, and ExampleReconAnalysis.cpp.

Definition at line 1033 of file KVParticle.cpp.

◆ GetWaveLength()

Double_t KVParticle::GetWaveLength ( ) const
inline

should be in fm

Definition at line 658 of file KVParticle.h.

◆ HasFrame()

Bool_t KVParticle::HasFrame ( const Char_t frame) const
inline
Parameters
[in]framename of a previously defined kinematical frame
Returns
kTRUE if a kinematical frame called frame has been defined for this particle
Note
Frame names are case-insensitive
We search the entire hierarchy of kinematical frames, starting from the current default kinematics

Definition at line 783 of file KVParticle.h.

◆ InFrame()

KVParticle KVParticle::InFrame ( const KVFrameTransform t)

Use this method to obtain 'on-the-fly' some information on particle kinematics in a different reference frame. The default kinematics of the particle are unchanged.

Definition at line 582 of file KVParticle.cpp.

◆ init()

void KVParticle::init ( void  )

default initialisation

Definition at line 50 of file KVParticle.cpp.

◆ IsDefaultKinematics()

Bool_t KVParticle::IsDefaultKinematics ( ) const
inline
Returns
kTRUE if this particle is in its default kinematic frame

Definition at line 777 of file KVParticle.h.

◆ IsDetected()

Bool_t KVParticle::IsDetected ( ) const
inline

Definition at line 736 of file KVParticle.h.

◆ IsOK()

Bool_t KVParticle::IsOK ( ) const

Determine whether this particle is considered "good" or not for analysis, depending on a previous call to SetIsOK(Bool_t flag). If SetIsOK has not been called, IsOK returns kTRUE by default.

Definition at line 318 of file KVParticle.cpp.

◆ ListGroups()

void KVParticle::ListGroups ( void  ) const

List all stored groups.

Definition at line 558 of file KVParticle.cpp.

◆ ls()

void KVParticle::ls ( Option_t option = "") const
virtual

Reimplemented from TObject.

Definition at line 353 of file KVParticle.cpp.

◆ operator=()

KVParticle & KVParticle::operator= ( const KVParticle rhs)

KVParticle assignment operator.

Definition at line 371 of file KVParticle.cpp.

◆ Print()

void KVParticle::Print ( Option_t t = "") const
virtual

print out characteristics of particle

Reimplemented from TLorentzVector.

Reimplemented in KVSimNucleus, KVNucleus, KVReconstructedNucleus, KVINDRAReconNuc, and KVFAZIAReconNuc.

Definition at line 212 of file KVParticle.cpp.

◆ print_frames()

void KVParticle::print_frames ( TString  fmt = "") const
private

recursive print out of all defined kinematical frames

Definition at line 189 of file KVParticle.cpp.

◆ RemoveAllGroups()

void KVParticle::RemoveAllGroups ( )

Remove all groups.

Definition at line 546 of file KVParticle.cpp.

◆ RemoveGroup()

void KVParticle::RemoveGroup ( const Char_t groupname)

Remove group from list of groups.

Definition at line 530 of file KVParticle.cpp.

◆ ResetEnergy()

void KVParticle::ResetEnergy ( )

Used for simulated particles after passage through some absorber/detector. The passage of a particle through the different absorbers modifies its kinetic energies, indeed the particle may be stopped in the detector. Calling this method will reset the particle's momentum to its initial value i.e. before it entered the first absorber. Particles which have not encountered any absorbers/detectors are left as they are.

Definition at line 393 of file KVParticle.cpp.

◆ ResetFrameCopyOnly()

void KVParticle::ResetFrameCopyOnly ( ) const
inlineprivate

Definition at line 509 of file KVParticle.h.

◆ ResetIsOK()

void KVParticle::ResetIsOK ( )
inline

Definition at line 706 of file KVParticle.h.

◆ Set4Mom()

void KVParticle::Set4Mom ( const TLorentzVector p)
inline

Definition at line 591 of file KVParticle.h.

◆ SetE()

void KVParticle::SetE ( Double_t  a)
inline

Definition at line 596 of file KVParticle.h.

◆ SetE0()

void KVParticle::SetE0 ( TVector3 e = 0)
inline

Definition at line 717 of file KVParticle.h.

◆ SetEnergy()

void KVParticle::SetEnergy ( Double_t  e)
inline

Definition at line 601 of file KVParticle.h.

◆ SetFrame() [1/2]

void KVParticle::SetFrame ( const Char_t frame,
const KVFrameTransform ft 
)

Define a Lorentz-boosted and/or rotated frame in which to calculate this particle's momentum and energy.

The new frame will have the name given in the string "frame", which can then be used to access the kinematics of the particle in different frames using GetFrame() (frame names are case-insensitive). Calling this method with the name of an existing frame will update the kinematics of the particle in that frame using the given transform (note that kinematics in all frames defined as 'subframes' of this frame will also be updated).

USING BOOSTS The boost velocity vector is that of the boosted frame with respect to the original frame of the particles in the event. The velocity vector can be given either in cm/ns units (default) or in units of 'c' (beta=kTRUE).

E.g. to define a frame moving at 0.1c in the +ve z-direction with respect to the original event frame:

(...supposing a valid pointer KVParticle* my_part...) TVector3 vframe(0,0,0.1); my_part->SetFrame("my_frame", KVFrameTransform(vframe, kTRUE));

or with velocity in cm/ns units (default): TVector3 vframe(0,0,3); my_part->SetFrame("my_frame", KVFrameTransform(vframe)); OR my_part->SetFrame("my_frame", vframe);

USING ROTATIONS According to the conventions adopted for the TRotation and TLorentzRotation classes, we actually use the inverse of the TLorentzRotation to make the transformation, to get the coordinates corresponding to a rotated coordinate system, not the coordinates of a rotated vector in the same coordinate system => you do not need to invert the transformation matrix

E.g. if you want to define a new frame whose coordinate axes are rotated with respect to the original axes, you can set up a TRotation like so:

TRotation rot; TVector3 newX, newY, newZ; // the new coordinate axes rot.RotateAxes(newX, newY, newZ);

If you are using one of the two global variables which calculate the event tensor (KVTensP and KVTensPCM) you can obtain the transformation to the tensor frame using:

TRotation rot; KVTensP* tens_gv;// pointer to tensor global variable tens_gv->GetTensor()->GetRotation(rot);// see KVTenseur3::GetRotation

Then the new frame can be defined by my_part->SetFrame("my_frame", KVFrameTransform(rot)); OR my_part->SetFrame("my_frame", rot);

USING COMBINED BOOST AND ROTATION You can define a frame using both a boost and a rotation like so: my_part->SetFrame("my_frame", KVFrameTransform(vframe,rot,kTRUE));

ACCESSING KINEMATICS IN NEW FRAMES In order to access the kinematics of the particle in the new frame:

my_part->GetFrame("my_frame")->GetTransverseEnergy();// transverse energy in "my_frame"

Definition at line 741 of file KVParticle.cpp.

◆ SetFrame() [2/2]

void KVParticle::SetFrame ( const Char_t newframe,
const Char_t oldframe,
const KVFrameTransform ft 
)

Define new kinematical frame by transformation from existing frame See SetFrame(const Char_t*,const KVFrameTransform&) for details on defining kinematically-transformed frames.

Definition at line 993 of file KVParticle.cpp.

◆ SetFrameCopyOnly()

void KVParticle::SetFrameCopyOnly ( ) const
inlineprivate

Definition at line 505 of file KVParticle.h.

◆ SetFrameName()

void KVParticle::SetFrameName ( const Char_t framename)
inline

Set the (non-persistent) name of the reference frame for this particle's kinematics. Also sets a (persistent) parameter "frameName"

Definition at line 806 of file KVParticle.h.

◆ SetGroups()

void KVParticle::SetGroups ( KVUniqueNameList un)
protected

Define for the particle a new list of groups.

Definition at line 457 of file KVParticle.cpp.

◆ SetIsDetected()

void KVParticle::SetIsDetected ( )
inline

Definition at line 732 of file KVParticle.h.

◆ SetIsOK()

void KVParticle::SetIsOK ( Bool_t  flag = kTRUE)

Set acceptation/rejection status for the particle.

In order to 'forget' this status (accept all particles) use ResetIsOK()

Definition at line 339 of file KVParticle.cpp.

◆ SetKE()

void KVParticle::SetKE ( Double_t  ecin)

Change particle KE, keeping momentum direction constant If momentum is zero (i.e. no direction defined) the particle will be given a velocity in the positive z-direction

Definition at line 230 of file KVParticle.cpp.

◆ SetMass()

virtual void KVParticle::SetMass ( Double_t  m)
inlinevirtual

Definition at line 568 of file KVParticle.h.

◆ SetMomentum() [1/4]

void KVParticle::SetMomentum ( const TVector3 v)
inline

Definition at line 576 of file KVParticle.h.

◆ SetMomentum() [2/4]

void KVParticle::SetMomentum ( const TVector3 v)
inline

Definition at line 580 of file KVParticle.h.

◆ SetMomentum() [3/4]

void KVParticle::SetMomentum ( Double_t  px,
Double_t  py,
Double_t  pz,
Option_t opt = "cart" 
)

Set Momentum components (in MeV/c) if option is "cart" or "cartesian" we give cartesian components (x,y,z) if option is "spher" or "spherical" we give components (rho,theta,phi) in spherical system with theta, phi in DEGREES

Definition at line 1049 of file KVParticle.cpp.

◆ SetMomentum() [4/4]

void KVParticle::SetMomentum ( Double_t  T,
TVector3  dir 
)

set momentum with kinetic energy t and unit direction vector d (d is normalised first in case it is not a unit vector)

Definition at line 169 of file KVParticle.cpp.

◆ SetName()

void KVParticle::SetName ( const Char_t nom)

Set Name of the particle.

Definition at line 412 of file KVParticle.cpp.

◆ SetOriginal()

void KVParticle::SetOriginal ( KVParticle p)
inlineprivate

Definition at line 540 of file KVParticle.h.

◆ SetParameter()

template<typename ValType >
void KVParticle::SetParameter ( const Char_t name,
ValType  value 
) const
inline

Definition at line 820 of file KVParticle.h.

◆ SetParentFrame()

void KVParticle::SetParentFrame ( KVParticle p)
inlineprivate

Definition at line 517 of file KVParticle.h.

◆ SetPerp()

void KVParticle::SetPerp ( Double_t  p)
inlineprivate

Definition at line 437 of file KVParticle.h.

◆ SetPhi()

void KVParticle::SetPhi ( Double_t  phi)
inline

Definition at line 699 of file KVParticle.h.

◆ SetPtEtaPhiE()

void KVParticle::SetPtEtaPhiE ( Double_t  pt,
Double_t  eta,
Double_t  phi,
Double_t  e 
)
inlineprivate

Definition at line 441 of file KVParticle.h.

◆ SetPtEtaPhiM()

void KVParticle::SetPtEtaPhiM ( Double_t  pt,
Double_t  eta,
Double_t  phi,
Double_t  m 
)
inlineprivate

Definition at line 445 of file KVParticle.h.

◆ SetPx()

void KVParticle::SetPx ( Double_t  a)
inlineprivate

Definition at line 449 of file KVParticle.h.

◆ SetPxPyPzE()

void KVParticle::SetPxPyPzE ( Double_t  px,
Double_t  py,
Double_t  pz,
Double_t  e 
)
inlineprivate

Definition at line 461 of file KVParticle.h.

◆ SetPy()

void KVParticle::SetPy ( Double_t  a)
inlineprivate

Definition at line 453 of file KVParticle.h.

◆ SetPz()

void KVParticle::SetPz ( Double_t  a)
inlineprivate

Definition at line 457 of file KVParticle.h.

◆ SetRandomMomentum()

void KVParticle::SetRandomMomentum ( Double_t  T,
Double_t  thmin,
Double_t  thmax,
Double_t  phmin,
Double_t  phmax,
Option_t opt = "isotropic" 
)

Give randomly directed momentum to particle with kinetic energy T Direction will be between (thmin,thmax) [degrees] limits in polar angle, and (phmin,phmax) [degrees] limits in azimuthal angle.

If opt = "" or "isotropic" (default) : direction is isotropically distributed over the solid angle If opt = "random" : direction is randomly distributed over solid angle

Based on KVPosition::GetRandomDirection().

Definition at line 136 of file KVParticle.cpp.

◆ SetRho()

void KVParticle::SetRho ( Double_t  rho)
inlineprivate

Definition at line 465 of file KVParticle.h.

◆ SetT()

void KVParticle::SetT ( Double_t  a)
inlineprivate

Definition at line 469 of file KVParticle.h.

◆ SetTheta()

void KVParticle::SetTheta ( Double_t  theta)
inline

Definition at line 695 of file KVParticle.h.

◆ SetVect()

void KVParticle::SetVect ( const TVector3 vect3)
inlineprivate

TLorentzVector setters should not be used.

Definition at line 425 of file KVParticle.h.

◆ SetVectM()

void KVParticle::SetVectM ( const TVector3 spatial,
Double_t  mass 
)
inlineprivate

Definition at line 429 of file KVParticle.h.

◆ SetVectMag()

void KVParticle::SetVectMag ( const TVector3 spatial,
Double_t  magnitude 
)
inlineprivate

Definition at line 433 of file KVParticle.h.

◆ SetVelocity()

void KVParticle::SetVelocity ( const TVector3 vel)

Set velocity of particle (in cm/ns units)

Definition at line 1080 of file KVParticle.cpp.

◆ SetX()

void KVParticle::SetX ( Double_t  a)
inlineprivate

Definition at line 473 of file KVParticle.h.

◆ SetXYZM()

void KVParticle::SetXYZM ( Double_t  x,
Double_t  y,
Double_t  z,
Double_t  m 
)
inlineprivate

Definition at line 485 of file KVParticle.h.

◆ SetXYZT()

void KVParticle::SetXYZT ( Double_t  x,
Double_t  y,
Double_t  z,
Double_t  t 
)
inlineprivate

Definition at line 489 of file KVParticle.h.

◆ SetY()

void KVParticle::SetY ( Double_t  a)
inlineprivate

Definition at line 477 of file KVParticle.h.

◆ SetZ()

void KVParticle::SetZ ( Double_t  a)
inlineprivate

Definition at line 481 of file KVParticle.h.

◆ UpdateAllFrames()

void KVParticle::UpdateAllFrames ( )

Call this method to update particle kinematics in all defined frames if you change the kinematics of the particle in its original/default frame.

Definition at line 913 of file KVParticle.cpp.

Friends And Related Function Documentation

◆ KVKinematicalFrame

friend class KVKinematicalFrame
friend

Definition at line 400 of file KVParticle.h.

Member Data Documentation

◆ fBoosted

FrameList KVParticle::fBoosted {this}
mutableprivate

Definition at line 420 of file KVParticle.h.

◆ fE0

TVector3* KVParticle::fE0
protected

the momentum of the particle before it is slowed/stopped by an absorber

Definition at line 496 of file KVParticle.h.

◆ fFrameCopyOnly

Bool_t KVParticle::fFrameCopyOnly {kFALSE}
mutableprivate

if != nullptr, this object is a representation of the particle in a kinematical frame

Definition at line 502 of file KVParticle.h.

◆ fFrameName

TString KVParticle::fFrameName
private

non-persistent frame name field, sets when calling SetFrame method

Definition at line 408 of file KVParticle.h.

◆ fGroups

KVUniqueNameList KVParticle::fGroups
mutableprivate

list of momenta of the particle in different Lorentz-boosted frames

list of TObjString for manage different group name

Definition at line 421 of file KVParticle.h.

◆ fName

TString KVParticle::fName
private

non-persistent name field - Is useful

Definition at line 407 of file KVParticle.h.

◆ fOriginal

KVParticle* KVParticle::fOriginal {nullptr}
private

parent kinematical frame

Definition at line 501 of file KVParticle.h.

◆ fParameters

KVNameValueList KVParticle::fParameters
protected

a general-purpose list of parameters associated with this particle

Definition at line 497 of file KVParticle.h.

◆ fParentFrame

KVParticle* KVParticle::fParentFrame {nullptr}
private

Definition at line 500 of file KVParticle.h.

◆ kSpeedOfLight

Double_t KVParticle::kSpeedOfLight = TMath::C() * 1.e-07
staticprivate

speed of light in cm/ns

Definition at line 422 of file KVParticle.h.


The documentation for this class was generated from the following files: