KaliVeda  1.13/01
Heavy-Ion Analysis Toolkit
KVParticle.h
Go to the documentation of this file.
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 "KVBase.h"
26 #include "TRef.h"
27 #include "TMath.h"
28 #include "KVList.h"
29 #include "KVUniqueNameList.h"
30 #include "TObjString.h"
31 #include "KVNameValueList.h"
32 #include "KVFrameTransform.h"
33 
34 class KVKinematicalFrame;
35 
105 
334 
398 class KVParticle: public TLorentzVector {
399 
400  friend class KVKinematicalFrame;
401 
402 private:
403  void print_frames(TString fmt = "") const;
404  KVKinematicalFrame* get_frame(const Char_t*) const;
405  KVKinematicalFrame* get_parent_frame(const Char_t*, KVKinematicalFrame* F = nullptr) const;
406 
409 
410  class FrameList : public KVList {
412  public:
414  void Add(TObject*);
416  void Clear(Option_t* = "");
417  void AddAll(const TCollection*);
418  };
419 
420  mutable FrameList fBoosted{this};
423 
425  void SetVect(const TVector3& vect3)
426  {
428  }
429  void SetVectM(const TVector3& spatial, Double_t mass)
430  {
431  TLorentzVector::SetVectM(spatial, mass);
432  }
433  void SetVectMag(const TVector3& spatial, Double_t magnitude)
434  {
435  TLorentzVector::SetVectMag(spatial, magnitude);
436  }
438  {
440  }
442  {
443  TLorentzVector::SetPtEtaPhiE(pt, eta, phi, e);
444  }
446  {
447  TLorentzVector::SetPtEtaPhiM(pt, eta, phi, m);
448  }
449  void SetPx(Double_t a)
450  {
452  }
453  void SetPy(Double_t a)
454  {
456  }
457  void SetPz(Double_t a)
458  {
460  }
462  {
463  TLorentzVector::SetPxPyPzE(px, py, pz, e);
464  }
465  void SetRho(Double_t rho)
466  {
468  }
469  void SetT(Double_t a)
470  {
472  }
473  void SetX(Double_t a)
474  {
476  }
477  void SetY(Double_t a)
478  {
480  }
481  void SetZ(Double_t a)
482  {
484  }
486  {
488  }
490  {
491  TLorentzVector::SetXYZT(x, y, z, t);
492  }
493 
494 protected:
495 
498 
499 private:
503 
505  void SetFrameCopyOnly() const
506  {
508  }
509  void ResetFrameCopyOnly() const
510  {
512  }
514  {
515  return fParentFrame;
516  }
518  {
519  fParentFrame = p;
520  }
522  {
525 
526  auto top_node = this;
527  while (top_node->GetParentFrame()) {
528  top_node = top_node->GetParentFrame();
529  }
530  return top_node;
531  }
532  const KVParticle* GetOriginal() const
533  {
534  return fOriginal ? fOriginal : this;
535  }
537  {
538  return fOriginal ? fOriginal : this;
539  }
541  {
542  fOriginal = p;
543  }
544 
545 protected:
546  void SetGroups(KVUniqueNameList* un);
547  void AddGroups(KVUniqueNameList* un);
548 
549 public:
550 
551  enum {
552  kIsOK = BIT(14), //acceptation/rejection flag
553  kIsOKSet = BIT(15), //flag to indicate flag is set
554  kIsDetected = BIT(16) //flag set when particle is slowed by some KVMaterial
555  };
556 
557  static Double_t C();
558 
559  KVParticle();
562  KVParticle(const KVParticle&);
563  virtual ~ KVParticle();
564  void init();
565  virtual void Copy(TObject&) const;
566  virtual void Clear(Option_t* opt = "");
567 
568  virtual void SetMass(Double_t m)
569  {
570  SetXYZM(Px(), Py(), Pz(), m);
571  };
573  {
574  return M();
575  };
576  void SetMomentum(const TVector3& v)
577  {
578  SetXYZM(v(0), v(1), v(2), M());
579  };
580  void SetMomentum(const TVector3* v)
581  {
582  SetXYZM((*v)(0), (*v)(1), (*v)(2), M());
583  };
584  void SetMomentum(Double_t px, Double_t py, Double_t pz, Option_t* opt =
585  "cart");
586  void SetMomentum(Double_t T, TVector3 dir);
587  void SetRandomMomentum(Double_t T, Double_t thmin, Double_t thmax,
588  Double_t phmin, Double_t phmax,
589  Option_t* opt = "isotropic");
590  virtual void Print(Option_t* t = "") const;
591  void Set4Mom(const TLorentzVector& p)
592  {
593  SetVect(p.Vect());
594  SetT(p.E());
595  }
596  void SetE(Double_t a)
597  {
598  SetKE(a);
599  };
600  void SetKE(Double_t ecin);
602  {
603  SetKE(e);
604  };
605  void SetVelocity(const TVector3&);
607  {
608  return Vect();
609  };
611  {
612  TVector3 perp = GetMomentum();
613  perp.SetZ(0);
614  return perp;
615  }
616  Double_t GetKE() const
617  {
618  Double_t e = E();
619  Double_t m = M();
621  return e - m;
622  };
624  {
625  return GetKE();
626  };
628  {
629  return GetKE() * TMath::Power(TMath::Sin(Theta()), 2.0);
630  }
632  {
633  return GetTransverseEnergy();
634  }
635 
637  {
638  return GetKE() - GetTransverseEnergy();
639  }
641  {
642  return GetLongitudinalEnergy();
643  }
644 
646  {
647  Double_t etran = Mt() - GetMass();
648  return etran;
649  }
651  {
652  return GetRTransverseEnergy();
653  }
654  Double_t GetE() const
655  {
656  return GetKE();
657  };
659  {
661  if (GetMomentum().Mag() == 0)
662  return 0;
663  Double_t h = TMath::H() * TMath::C() / TMath::Qe() * 1e9; //in MeV.fm
664  return h / GetMomentum().Mag();
665  };
667  {
669  Double_t h = TMath::H() * TMath::C() / TMath::Qe() * 1e9; //in MeV.fm
670  return h / TMath::Sqrt(TMath::TwoPi() * temp * GetMass());
671  };
672  TVector3 GetVelocity() const;
673  TVector3 GetV() const
674  {
675  return GetVelocity();
676  };
678  {
679  return GetV().z();
680  };
681  Double_t GetVperp() const;
683  {
684  return TMath::RadToDeg() * Theta();
685  };
687  {
688  return TMath::Cos(Theta());
689  }
690  Double_t GetPhi() const
691  {
692  Double_t phi = TMath::RadToDeg() * Phi();
693  return (phi < 0 ? 360. + phi : phi);
694  };
695  void SetTheta(Double_t theta)
696  {
698  }
699  void SetPhi(Double_t phi)
700  {
702  }
703 
704  Bool_t IsOK() const;
705  void SetIsOK(Bool_t flag = kTRUE);
706  void ResetIsOK()
707  {
709  }
710 
712  {
713  return &fBoosted;
714  }
715  void ls(Option_t* option = "") const;
716 
717  void SetE0(TVector3* e = 0)
718  {
719  if (!fE0)
720  fE0 = new TVector3;
721  if (!e) {
722  *fE0 = GetMomentum();
723  }
724  else {
725  *fE0 = *e;
726  }
727  }
729  {
730  return fE0;
731  }
733  {
735  }
737  {
738  return GetOriginal()->TestBit(kIsDetected);
739  }
740  KVParticle& operator=(const KVParticle& rhs);
741 
742  void ResetEnergy();
743 
744  const Char_t* GetName() const;
745  void SetName(const Char_t* nom);
746 
747  void AddGroup(const Char_t* groupname, const Char_t* from = "") const;
748 
749  Bool_t BelongsToGroup(const Char_t* groupname) const;
750  void RemoveGroup(const Char_t* groupname);
751  void RemoveAllGroups();
752  void ListGroups() const;
754  {
756  return GetGroups()->GetEntries();
757  }
759  {
761  return (KVUniqueNameList*)&GetOriginal()->fGroups;
762  }
763 
764 
766  void ChangeFrame(const KVFrameTransform&, const KVString& = "");
767  void ChangeDefaultFrame(const Char_t*, const Char_t* defname = "");
768  void SetFrame(const Char_t* frame, const KVFrameTransform&);
769  void SetFrame(const Char_t* newframe, const Char_t* oldframe, const KVFrameTransform&);
770  void UpdateAllFrames();
771 
773  {
775  return GetTopmostParentFrame();
776  }
778  {
780  return GetTopmostParentFrame() == this;
781  }
782 
783  Bool_t HasFrame(const Char_t* frame) const
784  {
791 
792  return (GetTopmostParentFrame()->get_frame(frame) != nullptr);
793  }
795  {
798  }
799  KVParticle const* GetFrame(const Char_t* frame, Bool_t warn_and_return_null_if_unknown = kTRUE) const;
800 
801  const Char_t* GetFrameName(void) const
802  {
804  return fFrameName;
805  }
806  void SetFrameName(const Char_t* framename)
807  {
810 
811  fFrameName = framename;
812  if (fFrameName != "") SetParameter("frameName", framename);
813  else GetParameters()->RemoveParameter("frameName");
814  }
815 
817  {
818  return (KVNameValueList*)&fParameters;
819  }
820  template<typename ValType> void SetParameter(const Char_t* name, ValType value) const
821  {
822  GetParameters()->SetValue(name, value);
823  }
824 
825  ClassDef(KVParticle, 8) //General base class for all massive particles
826 };
827 
828 #endif
int Int_t
#define e(i)
char Char_t
const Bool_t kFALSE
bool Bool_t
double Double_t
const Bool_t kTRUE
const char Option_t
#define ClassDef(name, id)
#define BIT(n)
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:27
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)
FrameList(KVParticle *p)
Definition: KVParticle.h:413
void AddAll(const TCollection *)
When all frames in a list are added to this one, this particle becomes the parent of all frames in th...
void Clear(Option_t *="")
When the frame list is cleared, this particle is no longer the parent of any frames in the list.
TObject * Remove(TObject *)
When a kinematical frame is removed, this particle is no longer the parent frame.
void Add(TObject *)
When a kinematical frame is added, this particle becomes the parent frame.
KVParticle * parent
Definition: KVParticle.h:411
Base class for relativistic kinematics of massive particles.
Definition: KVParticle.h:398
Double_t GetWaveLength() const
Definition: KVParticle.h:658
Int_t GetNumberOfDefinedGroups() const
Definition: KVParticle.h:753
void AddGroups(KVUniqueNameList *un)
list of groups added to the current one
Definition: KVParticle.cpp:469
KVParticle * GetOriginal()
Definition: KVParticle.h:536
Double_t GetThermalWaveLength(Double_t temp) const
Definition: KVParticle.h:666
void SetIsOK(Bool_t flag=kTRUE)
Definition: KVParticle.cpp:339
void SetTheta(Double_t theta)
Definition: KVParticle.h:695
TVector3 * GetPInitial() const
Definition: KVParticle.h:728
Double_t GetREtran() const
Definition: KVParticle.h:650
virtual void SetMass(Double_t m)
Definition: KVParticle.h:568
void RemoveGroup(const Char_t *groupname)
Remove group from list of groups.
Definition: KVParticle.cpp:530
void SetPy(Double_t a)
Definition: KVParticle.h:453
virtual ~ KVParticle()
void AddGroup(const Char_t *groupname, const Char_t *from="") const
Definition: KVParticle.cpp:436
void SetPhi(Double_t phi)
Definition: KVParticle.h:699
void UpdateAllFrames()
Definition: KVParticle.cpp:913
const KVParticle * GetOriginal() const
Definition: KVParticle.h:532
TVector3 GetMomentum() const
Definition: KVParticle.h:606
void ResetIsOK()
Definition: KVParticle.h:706
void RemoveAllGroups()
Remove all groups.
Definition: KVParticle.cpp:546
KVUniqueNameList fGroups
list of momenta of the particle in different Lorentz-boosted frames
Definition: KVParticle.h:421
KVNameValueList * GetParameters() const
Definition: KVParticle.h:816
void SetXYZT(Double_t x, Double_t y, Double_t z, Double_t t)
Definition: KVParticle.h:489
Double_t GetTheta() const
Definition: KVParticle.h:682
KVParticle & operator=(const KVParticle &rhs)
KVParticle assignment operator.
Definition: KVParticle.cpp:371
void SetVectM(const TVector3 &spatial, Double_t mass)
Definition: KVParticle.h:429
const Char_t * GetFrameName(void) const
Definition: KVParticle.h:801
void ResetEnergy()
Definition: KVParticle.cpp:393
void SetMomentum(const TVector3 &v)
Definition: KVParticle.h:576
void SetVelocity(const TVector3 &)
Set velocity of particle (in cm/ns units)
const Char_t * GetName() const
return the field fName
Definition: KVParticle.cpp:423
KVParticle * GetParentFrame() const
Definition: KVParticle.h:513
void SetRho(Double_t rho)
Definition: KVParticle.h:465
Double_t GetEnergy() const
Definition: KVParticle.h:623
void ResetFrameCopyOnly() const
Definition: KVParticle.h:509
Double_t GetTransverseEnergy() const
Definition: KVParticle.h:627
void SetVectMag(const TVector3 &spatial, Double_t magnitude)
Definition: KVParticle.h:433
KVUniqueNameList * GetGroups() const
Definition: KVParticle.h:758
Double_t GetLongitudinalEnergy() const
Definition: KVParticle.h:636
static Double_t C()
Definition: KVParticle.cpp:117
Bool_t IsDetected() const
Definition: KVParticle.h:736
void SetKE(Double_t ecin)
Definition: KVParticle.cpp:230
virtual void Clear(Option_t *opt="")
Reset particle properties i.e. before creating/reading a new event.
Definition: KVParticle.cpp:295
void SetPx(Double_t a)
Definition: KVParticle.h:449
void ls(Option_t *option="") const
Definition: KVParticle.cpp:353
void init()
default initialisation
Definition: KVParticle.cpp:50
Double_t GetPhi() const
Definition: KVParticle.h:690
KVNameValueList fParameters
a general-purpose list of parameters associated with this particle
Definition: KVParticle.h:497
Double_t GetCosTheta() const
Definition: KVParticle.h:686
void SetPxPyPzE(Double_t px, Double_t py, Double_t pz, Double_t e)
Definition: KVParticle.h:461
Bool_t HasFrame(const Char_t *frame) const
Definition: KVParticle.h:783
Double_t GetE() const
Definition: KVParticle.h:654
Double_t GetRTransverseEnergy() const
Definition: KVParticle.h:645
void SetMomentum(const TVector3 *v)
Definition: KVParticle.h:580
Double_t GetElong() const
Definition: KVParticle.h:640
virtual void Copy(TObject &) const
Definition: KVParticle.cpp:269
void Set4Mom(const TLorentzVector &p)
Definition: KVParticle.h:591
KVKinematicalFrame * get_parent_frame(const Char_t *, KVKinematicalFrame *F=nullptr) const
Definition: KVParticle.cpp:965
KVParticle * fOriginal
parent kinematical frame
Definition: KVParticle.h:501
virtual void Print(Option_t *t="") const
print out characteristics of particle
Definition: KVParticle.cpp:212
void SetE(Double_t a)
Definition: KVParticle.h:596
KVParticle InFrame(const KVFrameTransform &)
Definition: KVParticle.cpp:582
void SetE0(TVector3 *e=0)
Definition: KVParticle.h:717
void SetY(Double_t a)
Definition: KVParticle.h:477
void SetIsDetected()
Definition: KVParticle.h:732
TString fFrameName
non-persistent frame name field, sets when calling SetFrame method
Definition: KVParticle.h:408
void SetZ(Double_t a)
Definition: KVParticle.h:481
Bool_t IsOK() const
Definition: KVParticle.cpp:318
void SetParameter(const Char_t *name, ValType value) const
Definition: KVParticle.h:820
void ListGroups() const
List all stored groups.
Definition: KVParticle.cpp:558
Bool_t fFrameCopyOnly
if != nullptr, this object is a representation of the particle in a kinematical frame
Definition: KVParticle.h:502
void SetX(Double_t a)
Definition: KVParticle.h:473
void SetFrameName(const Char_t *framename)
Definition: KVParticle.h:806
static Double_t kSpeedOfLight
speed of light in cm/ns
Definition: KVParticle.h:422
void SetParentFrame(KVParticle *p)
Definition: KVParticle.h:517
Double_t GetKE() const
Definition: KVParticle.h:616
Double_t GetVpar() const
Definition: KVParticle.h:677
void SetXYZM(Double_t x, Double_t y, Double_t z, Double_t m)
Definition: KVParticle.h:485
void SetName(const Char_t *nom)
Set Name of the particle.
Definition: KVParticle.cpp:412
Bool_t IsDefaultKinematics() const
Definition: KVParticle.h:777
void SetFrame(const Char_t *frame, const KVFrameTransform &)
Definition: KVParticle.cpp:741
Int_t GetNumberOfDefinedFrames() const
Definition: KVParticle.h:794
KVParticle const * GetFrame(const Char_t *frame, Bool_t warn_and_return_null_if_unknown=kTRUE) const
Definition: KVParticle.cpp:855
void ChangeFrame(const KVFrameTransform &, const KVString &="")
Definition: KVParticle.cpp:601
FrameList fBoosted
Definition: KVParticle.h:420
void SetOriginal(KVParticle *p)
Definition: KVParticle.h:540
void SetGroups(KVUniqueNameList *un)
Define for the particle a new list of groups.
Definition: KVParticle.cpp:457
Bool_t BelongsToGroup(const Char_t *groupname) const
Definition: KVParticle.cpp:509
TVector3 * fE0
the momentum of the particle before it is slowed/stopped by an absorber
Definition: KVParticle.h:496
void ChangeDefaultFrame(const Char_t *, const Char_t *defname="")
Definition: KVParticle.cpp:622
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:445
const KVParticle * GetCurrentDefaultKinematics() const
Definition: KVParticle.h:772
KVParticle * fParentFrame
Definition: KVParticle.h:500
TVector3 GetV() const
Definition: KVParticle.h:673
TString fName
non-persistent name field - Is useful
Definition: KVParticle.h:407
KVKinematicalFrame * get_frame(const Char_t *) const
Definition: KVParticle.cpp:937
void SetEnergy(Double_t e)
Definition: KVParticle.h:601
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
const KVParticle * GetTopmostParentFrame() const
Definition: KVParticle.h:521
void SetPz(Double_t a)
Definition: KVParticle.h:457
KVList * GetListOfFrames() const
Definition: KVParticle.h:711
void SetPerp(Double_t p)
Definition: KVParticle.h:437
void SetVect(const TVector3 &vect3)
TLorentzVector setters should not be used.
Definition: KVParticle.h:425
void SetFrameCopyOnly() const
Definition: KVParticle.h:505
TVector3 GetTransverseMomentum() const
Definition: KVParticle.h:610
void SetT(Double_t a)
Definition: KVParticle.h:469
Int_t _GetNumberOfDefinedFrames() const
used to inhibit full copying of particles in different kinematical frames
Definition: KVParticle.cpp:487
void SetPtEtaPhiE(Double_t pt, Double_t eta, Double_t phi, Double_t e)
Definition: KVParticle.h:441
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:72
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)
Double_t z() const
Double_t Mag() const
void SetZ(Double_t)
TPaveText * pt
Double_t y[n]
Double_t x[n]
TH1 * h
const long double m
Definition: KVUnits.h:70
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
auto * a