KaliVeda
Toolkit for HIC analysis
KVParticle.h
1 /***************************************************************************
2  kvparticle.h - description
3  -------------------
4  begin : Sun May 19 2002
5  copyright : (C) 2002 by J.D. Frankland
6  email : frankland@ganil.fr
7 
8 $Id: KVParticle.h,v 1.41 2008/05/21 13:19:56 ebonnet Exp $
9  ***************************************************************************/
10 
11 /***************************************************************************
12  * *
13  * This program is free software; you can redistribute it and/or modify *
14  * it under the terms of the GNU General Public License as published by *
15  * the Free Software Foundation; either version 2 of the License, or *
16  * (at your option) any later version. *
17  * *
18  ***************************************************************************/
19 
20 #ifndef KVPARTICLE_H
21 #define KVPARTICLE_H
22 
23 #include "TVector3.h"
24 #include "TLorentzVector.h"
25 #include "TMath.h"
26 #include "KVList.h"
27 #include "KVUniqueNameList.h"
28 #include "KVNameValueList.h"
29 #include "KVFrameTransform.h"
30 
31 class KVKinematicalFrame;
32 
103 
332 
396 class KVParticle: public TLorentzVector {
397 
398  friend class KVKinematicalFrame;
399 
400 private:
401  void print_frames(TString fmt = "") const;
402  KVKinematicalFrame* get_frame(const Char_t*) const;
403  KVKinematicalFrame* get_parent_frame(const Char_t*, KVKinematicalFrame* F = nullptr) const;
404 
407 
408 public:
409  class FrameList : public KVList {
411  public:
413  void Add(TObject*) override;
414  TObject* Remove(TObject*) override;
415  void Clear(Option_t* = "") override;
416  void AddAll(const TCollection*) override;
417  void Copy(TObject& _new_list) const override;
419  };
420 private:
421  mutable FrameList fBoosted{this};
424 
426  void SetVect(const TVector3& vect3)
427  {
429  }
430  void SetVectM(const TVector3& spatial, Double_t mass)
431  {
432  TLorentzVector::SetVectM(spatial, mass);
433  }
434  void SetVectMag(const TVector3& spatial, Double_t magnitude)
435  {
436  TLorentzVector::SetVectMag(spatial, magnitude);
437  }
439  {
441  }
443  {
444  TLorentzVector::SetPtEtaPhiE(pt, eta, phi, e);
445  }
447  {
448  TLorentzVector::SetPtEtaPhiM(pt, eta, phi, m);
449  }
450  void SetPx(Double_t a)
451  {
453  }
454  void SetPy(Double_t a)
455  {
457  }
458  void SetPz(Double_t a)
459  {
461  }
463  {
464  TLorentzVector::SetPxPyPzE(px, py, pz, e);
465  }
466  void SetRho(Double_t rho)
467  {
469  }
470  void SetT(Double_t a)
471  {
473  }
474  void SetX(Double_t a)
475  {
477  }
478  void SetY(Double_t a)
479  {
481  }
482  void SetZ(Double_t a)
483  {
485  }
487  {
489  }
491  {
492  TLorentzVector::SetXYZT(x, y, z, t);
493  }
494 
495 protected:
496 
499 
500 private:
504 
506  void SetFrameCopyOnly() const
507  {
509  }
510  void ResetFrameCopyOnly() const
511  {
513  }
515  {
516  return fParentFrame;
517  }
519  {
520  fParentFrame = p;
521  }
523  {
526 
527  auto top_node = this;
528  while (top_node->GetParentFrame()) {
529  top_node = top_node->GetParentFrame();
530  }
531  return top_node;
532  }
533  const KVParticle* GetOriginal() const
534  {
535  return fOriginal ? fOriginal : this;
536  }
538  {
539  return fOriginal ? fOriginal : this;
540  }
542  {
543  fOriginal = p;
544  }
545  void SetMomentum(const TVector3* v)
546  {
547  SetMomentum(*v);
548  }
549 
550 protected:
551  void SetGroups(KVUniqueNameList* un);
552  void AddGroups(KVUniqueNameList* un);
553 
554 public:
555 
556  enum {
557  kIsOK = BIT(14), //acceptation/rejection flag
558  kIsOKSet = BIT(15), //flag to indicate flag is set
559  kIsDetected = BIT(16) //flag set when particle is slowed by some KVMaterial
560  };
561 
562  static Double_t C();
563 
564  KVParticle();
565  KVParticle(Double_t m, const TVector3& p);
567  KVParticle(const KVParticle&);
568  virtual ~ KVParticle();
569  void init();
570  void Copy(TObject&) const override;
571  void Clear(Option_t* opt = "") override;
572 
573  virtual void SetMass(Double_t m)
574  {
575  SetXYZM(Px(), Py(), Pz(), m);
576  }
578  {
579  return M();
580  }
581  void SetMomentum(const TVector3& v)
582  {
583  SetXYZM(v(0), v(1), v(2), M());
584  }
585  void SetMomentum(Double_t px, Double_t py, Double_t pz, Option_t* opt =
586  "cart");
587  void SetMomentum(Double_t T, const TVector3& dir);
588  void SetRandomMomentum(Double_t T, Double_t thmin, Double_t thmax,
589  Double_t phmin, Double_t phmax,
590  Option_t* opt = "isotropic");
591  void Print(Option_t* t = "") const override;
592  void Set4Mom(const TLorentzVector& p)
593  {
594  SetVect(p.Vect());
595  SetT(p.E());
596  }
597  void SetE(Double_t a)
598  {
599  SetKE(a);
600  }
601  void SetKE(Double_t ecin);
603  {
604  SetKE(e);
605  }
606  void SetVelocity(const TVector3&);
608  {
609  return Vect();
610  }
612  {
613  TVector3 perp = GetMomentum();
614  perp.SetZ(0);
615  return perp;
616  }
617  Double_t GetKE() const
618  {
619  Double_t e = E();
620  Double_t m = M();
622  return e - m;
623  }
625  {
626  return GetKE();
627  }
629  {
630  return GetKE() * TMath::Power(TMath::Sin(Theta()), 2.0);
631  }
633  {
634  return GetTransverseEnergy();
635  }
636 
638  {
639  return GetKE() - GetTransverseEnergy();
640  }
642  {
643  return GetLongitudinalEnergy();
644  }
645 
647  {
648  Double_t etran = Mt() - GetMass();
649  return etran;
650  }
652  {
653  return GetRTransverseEnergy();
654  }
655  Double_t GetE() const
656  {
657  return GetKE();
658  };
660  {
662  if (GetMomentum().Mag() == 0)
663  return 0;
664  Double_t h = TMath::H() * TMath::C() / TMath::Qe() * 1e9; //in MeV.fm
665  return h / GetMomentum().Mag();
666  };
668  {
670  Double_t h = TMath::H() * TMath::C() / TMath::Qe() * 1e9; //in MeV.fm
671  return h / TMath::Sqrt(TMath::TwoPi() * temp * GetMass());
672  };
673  TVector3 GetVelocity() const;
674  TVector3 GetV() const
675  {
676  return GetVelocity();
677  };
679  {
680  return GetV().z();
681  };
682  Double_t GetVperp() const;
684  {
685  return TMath::RadToDeg() * Theta();
686  };
688  {
689  return TMath::Cos(Theta());
690  }
691  Double_t GetPhi() const
692  {
693  Double_t phi = TMath::RadToDeg() * Phi();
694  return (phi < 0 ? 360. + phi : phi);
695  };
696  void SetTheta(Double_t theta)
697  {
699  }
700  void SetPhi(Double_t phi)
701  {
703  }
704 
705  Bool_t IsOK() const;
706  void SetIsOK(Bool_t flag = kTRUE);
707  void ResetIsOK()
708  {
710  }
711 
713  {
714  return &fBoosted;
715  }
716  void ls(Option_t* option = "") const override;
717 
718  void SetE0(TVector3* e = 0)
719  {
720  if (!fE0)
721  fE0 = new TVector3;
722  if (!e) {
723  *fE0 = GetMomentum();
724  }
725  else {
726  *fE0 = *e;
727  }
728  }
730  {
731  return fE0;
732  }
734  {
736  }
738  {
739  return GetOriginal()->TestBit(kIsDetected);
740  }
741  KVParticle& operator=(const KVParticle& rhs);
742 
743  void ResetEnergy();
744 
745  const Char_t* GetName() const override;
746  void SetName(const Char_t* nom);
747 
748  void AddGroup(const Char_t* groupname, const Char_t* from = "") const;
749 
750  Bool_t BelongsToGroup(const Char_t* groupname) const;
751  void RemoveGroup(const Char_t* groupname);
752  void RemoveAllGroups();
753  void ListGroups() const;
755  {
757  return GetGroups()->GetEntries();
758  }
760  {
762  return (KVUniqueNameList*)&GetOriginal()->fGroups;
763  }
764 
765 
767  void ChangeFrame(const KVFrameTransform&, const KVString& = "");
768  void ChangeDefaultFrame(const Char_t*, const Char_t* defname = "");
769  void SetFrame(const Char_t* frame, const KVFrameTransform&);
770  void SetFrame(const Char_t* newframe, const Char_t* oldframe, const KVFrameTransform&);
771  void UpdateAllFrames();
772 
774  {
776  return GetTopmostParentFrame();
777  }
779  {
781  return GetTopmostParentFrame() == this;
782  }
783 
784  Bool_t HasFrame(const Char_t* frame) const
785  {
792 
793  if (!fFrameName.CompareTo(frame, TString::kIgnoreCase)) return true; // case where 'frame'=default frame
794  return (GetTopmostParentFrame()->get_frame(frame) != nullptr);
795  }
797  {
800  }
801  KVParticle const* GetFrame(const Char_t* frame, Bool_t warn_and_return_null_if_unknown = kTRUE) const;
802 
803  const Char_t* GetFrameName(void) const
804  {
806  return fFrameName;
807  }
808  void SetFrameName(const Char_t* framename)
809  {
812 
813  fFrameName = framename;
814  if (fFrameName != "") SetParameter("frameName", framename);
815  else GetParameters()->RemoveParameter("frameName");
816  }
817 
819  {
821  }
822  template<typename ValType> void SetParameter(const Char_t* name, ValType value) const
823  {
825  }
826 
827  ClassDefOverride(KVParticle, 8) //General base class for all massive particles
828 };
829 
830 #endif
int Int_t
#define e(i)
bool Bool_t
char Char_t
constexpr Bool_t kFALSE
double Double_t
constexpr Bool_t kTRUE
const char Option_t
#define BIT(n)
#define ClassDefOverride(name, id)
winID h TVirtualViewer3D TVirtualGLPainter p
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void value
char name[80]
Utility class for kinematical transformations of KVParticle class.
Kinematical representation of a particle in different reference frames.
Extended TList class which owns its objects by default.
Definition: KVList.h:28
Handles lists of named parameters with different types, a list of KVNamedParameter objects.
void SetValue(const Char_t *name, value_type value)
void RemoveParameter(const Char_t *name)
void Clear(Option_t *="") override
When the frame list is cleared, this particle is no longer the parent of any frames in the list.
TObject * Remove(TObject *) override
When a kinematical frame is removed, this particle is no longer the parent frame.
FrameList(KVParticle *p)
Definition: KVParticle.h:412
void Add(TObject *) override
When a kinematical frame is added, this particle becomes the parent frame.
void Copy(TObject &_new_list) const override
KVParticle * parent
Definition: KVParticle.h:410
void AddAll(const TCollection *) override
When all frames in a list are added to this one, this particle becomes the parent of all frames in th...
Base class for relativistic kinematics of massive particles.
Definition: KVParticle.h:396
Double_t GetWaveLength() const
Definition: KVParticle.h:659
Int_t GetNumberOfDefinedGroups() const
Definition: KVParticle.h:754
void AddGroups(KVUniqueNameList *un)
list of groups added to the current one
Definition: KVParticle.cpp:501
KVParticle * GetOriginal()
Definition: KVParticle.h:537
Double_t GetThermalWaveLength(Double_t temp) const
Definition: KVParticle.h:667
void SetIsOK(Bool_t flag=kTRUE)
Definition: KVParticle.cpp:371
void SetTheta(Double_t theta)
Definition: KVParticle.h:696
TVector3 * GetPInitial() const
Definition: KVParticle.h:729
Double_t GetREtran() const
Definition: KVParticle.h:651
virtual void SetMass(Double_t m)
Definition: KVParticle.h:573
void RemoveGroup(const Char_t *groupname)
Remove group from list of groups.
Definition: KVParticle.cpp:562
void SetPy(Double_t a)
Definition: KVParticle.h:454
virtual ~ KVParticle()
void AddGroup(const Char_t *groupname, const Char_t *from="") const
Definition: KVParticle.cpp:468
void SetPhi(Double_t phi)
Definition: KVParticle.h:700
void UpdateAllFrames()
Definition: KVParticle.cpp:955
const KVParticle * GetOriginal() const
Definition: KVParticle.h:533
TVector3 GetMomentum() const
Definition: KVParticle.h:607
void ResetIsOK()
Definition: KVParticle.h:707
void RemoveAllGroups()
Remove all groups.
Definition: KVParticle.cpp:578
KVUniqueNameList fGroups
list of momenta of the particle in different Lorentz-boosted frames
Definition: KVParticle.h:422
KVNameValueList * GetParameters() const
Definition: KVParticle.h:818
void SetXYZT(Double_t x, Double_t y, Double_t z, Double_t t)
Definition: KVParticle.h:490
Double_t GetTheta() const
Definition: KVParticle.h:683
KVParticle & operator=(const KVParticle &rhs)
KVParticle assignment operator.
Definition: KVParticle.cpp:403
void SetVectM(const TVector3 &spatial, Double_t mass)
Definition: KVParticle.h:430
const Char_t * GetFrameName(void) const
Definition: KVParticle.h:803
void ResetEnergy()
Definition: KVParticle.cpp:425
void SetMomentum(const TVector3 &v)
Definition: KVParticle.h:581
void ls(Option_t *option="") const override
Definition: KVParticle.cpp:385
void SetVelocity(const TVector3 &)
Set velocity of particle (in cm/ns units)
KVParticle * GetParentFrame() const
Definition: KVParticle.h:514
void SetRho(Double_t rho)
Definition: KVParticle.h:466
Double_t GetEnergy() const
Definition: KVParticle.h:624
void ResetFrameCopyOnly() const
Definition: KVParticle.h:510
Double_t GetTransverseEnergy() const
Definition: KVParticle.h:628
void SetVectMag(const TVector3 &spatial, Double_t magnitude)
Definition: KVParticle.h:434
KVUniqueNameList * GetGroups() const
Definition: KVParticle.h:759
Double_t GetLongitudinalEnergy() const
Definition: KVParticle.h:637
static Double_t C()
Definition: KVParticle.cpp:117
Bool_t IsDetected() const
Definition: KVParticle.h:737
const Char_t * GetName() const override
return the field fName
Definition: KVParticle.cpp:455
void SetKE(Double_t ecin)
Definition: KVParticle.cpp:246
void SetPx(Double_t a)
Definition: KVParticle.h:450
void init()
default initialisation
Definition: KVParticle.cpp:50
Double_t GetPhi() const
Definition: KVParticle.h:691
KVNameValueList fParameters
a general-purpose list of parameters associated with this particle
Definition: KVParticle.h:498
Double_t GetCosTheta() const
Definition: KVParticle.h:687
void SetPxPyPzE(Double_t px, Double_t py, Double_t pz, Double_t e)
Definition: KVParticle.h:462
Bool_t HasFrame(const Char_t *frame) const
Definition: KVParticle.h:784
void Copy(TObject &) const override
Definition: KVParticle.cpp:286
Double_t GetE() const
Definition: KVParticle.h:655
Double_t GetRTransverseEnergy() const
Definition: KVParticle.h:646
void SetMomentum(const TVector3 *v)
Definition: KVParticle.h:545
Double_t GetElong() const
Definition: KVParticle.h:641
void Set4Mom(const TLorentzVector &p)
Definition: KVParticle.h:592
KVKinematicalFrame * get_parent_frame(const Char_t *, KVKinematicalFrame *F=nullptr) const
KVParticle * fOriginal
parent kinematical frame
Definition: KVParticle.h:502
void SetE(Double_t a)
Definition: KVParticle.h:597
KVParticle InFrame(const KVFrameTransform &)
Definition: KVParticle.cpp:614
void SetE0(TVector3 *e=0)
Definition: KVParticle.h:718
void SetY(Double_t a)
Definition: KVParticle.h:478
void SetIsDetected()
Definition: KVParticle.h:733
TString fFrameName
non-persistent frame name field, sets when calling SetFrame method
Definition: KVParticle.h:406
void Clear(Option_t *opt="") override
Reset particle properties i.e. before creating/reading a new event.
Definition: KVParticle.cpp:327
void SetZ(Double_t a)
Definition: KVParticle.h:482
Bool_t IsOK() const
Definition: KVParticle.cpp:350
void SetParameter(const Char_t *name, ValType value) const
Definition: KVParticle.h:822
void ListGroups() const
List all stored groups.
Definition: KVParticle.cpp:590
Bool_t fFrameCopyOnly
if != nullptr, this object is a representation of the particle in a kinematical frame
Definition: KVParticle.h:503
void SetX(Double_t a)
Definition: KVParticle.h:474
void SetFrameName(const Char_t *framename)
Definition: KVParticle.h:808
static Double_t kSpeedOfLight
speed of light in cm/ns
Definition: KVParticle.h:423
void SetParentFrame(KVParticle *p)
Definition: KVParticle.h:518
Double_t GetKE() const
Definition: KVParticle.h:617
Double_t GetVpar() const
Definition: KVParticle.h:678
void SetXYZM(Double_t x, Double_t y, Double_t z, Double_t m)
Definition: KVParticle.h:486
void SetName(const Char_t *nom)
Set Name of the particle.
Definition: KVParticle.cpp:444
Bool_t IsDefaultKinematics() const
Definition: KVParticle.h:778
void SetFrame(const Char_t *frame, const KVFrameTransform &)
Definition: KVParticle.cpp:775
Int_t GetNumberOfDefinedFrames() const
Definition: KVParticle.h:796
KVParticle const * GetFrame(const Char_t *frame, Bool_t warn_and_return_null_if_unknown=kTRUE) const
Definition: KVParticle.cpp:897
void ChangeFrame(const KVFrameTransform &, const KVString &="")
Definition: KVParticle.cpp:633
FrameList fBoosted
Definition: KVParticle.h:421
void SetOriginal(KVParticle *p)
Definition: KVParticle.h:541
void SetGroups(KVUniqueNameList *un)
Define for the particle a new list of groups.
Definition: KVParticle.cpp:489
Bool_t BelongsToGroup(const Char_t *groupname) const
Definition: KVParticle.cpp:541
TVector3 * fE0
the momentum of the particle before it is slowed/stopped by an absorber
Definition: KVParticle.h:497
void ChangeDefaultFrame(const Char_t *, const Char_t *defname="")
Definition: KVParticle.cpp:654
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
void SetPtEtaPhiM(Double_t pt, Double_t eta, Double_t phi, Double_t m)
Definition: KVParticle.h:446
const KVParticle * GetCurrentDefaultKinematics() const
Definition: KVParticle.h:773
KVParticle * fParentFrame
Definition: KVParticle.h:501
TVector3 GetV() const
Definition: KVParticle.h:674
TString fName
non-persistent name field - Is useful
Definition: KVParticle.h:405
KVKinematicalFrame * get_frame(const Char_t *) const
Definition: KVParticle.cpp:979
void SetEnergy(Double_t e)
Definition: KVParticle.h:602
void Print(Option_t *t="") const override
print out characteristics of particle
Definition: KVParticle.cpp:212
Double_t GetMass() const
Definition: KVParticle.h:577
TVector3 GetVelocity() const
returns velocity vector in cm/ns units
Double_t GetVperp() const
Double_t GetEtran() const
Definition: KVParticle.h:632
const KVParticle * GetTopmostParentFrame() const
Definition: KVParticle.h:522
void SetPz(Double_t a)
Definition: KVParticle.h:458
KVList * GetListOfFrames() const
Definition: KVParticle.h:712
void SetPerp(Double_t p)
Definition: KVParticle.h:438
void SetVect(const TVector3 &vect3)
TLorentzVector setters should not be used.
Definition: KVParticle.h:426
void SetFrameCopyOnly() const
Definition: KVParticle.h:506
TVector3 GetTransverseMomentum() const
Definition: KVParticle.h:611
void SetT(Double_t a)
Definition: KVParticle.h:470
Int_t _GetNumberOfDefinedFrames() const
used to inhibit full copying of particles in different kinematical frames
Definition: KVParticle.cpp:519
void SetPtEtaPhiE(Double_t pt, Double_t eta, Double_t phi, Double_t e)
Definition: KVParticle.h:442
void print_frames(TString fmt="") const
recursive print out of all defined kinematical frames
Definition: KVParticle.cpp:189
Extension of ROOT TString class which allows backwards compatibility with ROOT v3....
Definition: KVString.h:73
Optimised list in which named objects can only be placed once.
virtual Int_t GetEntries() const
void SetPy(Double_t a)
void SetXYZT(Double_t x, Double_t y, Double_t z, Double_t t)
void SetY(Double_t a)
void SetPerp(Double_t)
TVector3 Vect() const
void SetT(Double_t a)
Double_t Px() const
Double_t Mag() const
Double_t M() const
void SetPtEtaPhiE(Double_t pt, Double_t eta, Double_t phi, Double_t e)
Double_t Pz() const
Double_t Py() const
void SetPtEtaPhiM(Double_t pt, Double_t eta, Double_t phi, Double_t m)
void SetPhi(Double_t phi)
void SetPz(Double_t a)
void SetRho(Double_t rho)
void SetVectMag(const TVector3 &spatial, Double_t magnitude)
Double_t Theta() const
void SetPx(Double_t a)
void SetPxPyPzE(Double_t px, Double_t py, Double_t pz, Double_t e)
void SetTheta(Double_t theta)
Double_t Phi() const
Double_t Mt() const
void SetZ(Double_t a)
Double_t E() const
void SetVect(const TVector3 &vect3)
void SetVectM(const TVector3 &spatial, Double_t mass)
void SetXYZM(Double_t x, Double_t y, Double_t z, Double_t m)
Double_t T() const
void SetX(Double_t a)
void SetBit(UInt_t f)
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
void ResetBit(UInt_t f)
int CompareTo(const char *cs, ECaseCompare cmp=kExact) const
Double_t z() const
Double_t Mag() const
void SetZ(Double_t)
TPaveText * pt
Double_t y[n]
Double_t x[n]
TH1 * h
constexpr Double_t C()
constexpr Double_t Qe()
Double_t Power(Double_t x, Double_t y)
constexpr Double_t DegToRad()
constexpr Double_t H()
Double_t Sqrt(Double_t x)
Double_t Cos(Double_t)
Double_t Sin(Double_t)
constexpr Double_t RadToDeg()
constexpr Double_t TwoPi()
v
TMarker m
TArc a