KaliVeda  1.12/06
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;
36 
106 
141 
171 
201 
328 
398 class KVParticle: public TLorentzVector {
399 
400 private:
401  void print_frames(TString fmt = "") const;
402  KVKinematicalFrame* get_frame(const Char_t*) const;
404 
410 
412  void SetVect(const TVector3& vect3)
413  {
415  }
416  void SetVectM(const TVector3& spatial, Double_t mass)
417  {
418  TLorentzVector::SetVectM(spatial, mass);
419  }
420  void SetVectMag(const TVector3& spatial, Double_t magnitude)
421  {
422  TLorentzVector::SetVectMag(spatial, magnitude);
423  }
425  {
427  }
429  {
430  TLorentzVector::SetPtEtaPhiE(pt, eta, phi, e);
431  }
433  {
434  TLorentzVector::SetPtEtaPhiM(pt, eta, phi, m);
435  }
436  void SetPx(Double_t a)
437  {
439  }
440  void SetPy(Double_t a)
441  {
443  }
444  void SetPz(Double_t a)
445  {
447  }
449  {
450  TLorentzVector::SetPxPyPzE(px, py, pz, e);
451  }
452  void SetRho(Double_t rho)
453  {
455  }
456  void SetT(Double_t a)
457  {
459  }
460  void SetX(Double_t a)
461  {
463  }
464  void SetY(Double_t a)
465  {
467  }
468  void SetZ(Double_t a)
469  {
471  }
473  {
475  }
477  {
478  TLorentzVector::SetXYZT(x, y, z, t);
479  }
480 
481 protected:
482 
485 
486  virtual void AddGroup_Withcondition(const Char_t*, KVParticleCondition*);
487  virtual void AddGroup_Sanscondition(const Char_t* groupname, const Char_t* from = "");
488  void CreateGroups();
489  void SetGroups(KVUniqueNameList* un);
490  void AddGroups(KVUniqueNameList* un);
491 
492 public:
493 
494  Bool_t HasFrame(const Char_t* frame)
495  {
499 
500  return (get_frame(frame) != nullptr);
501  }
504  KVUniqueNameList* GetGroups() const;
505 
506  enum {
507  kIsOK = BIT(14), //acceptation/rejection flag
508  kIsOKSet = BIT(15), //flag to indicate flag is set
509  kIsDetected = BIT(16) //flag set when particle is slowed by some KVMaterial
510  };
511 
512  static Double_t C();
513 
514  KVParticle();
517  KVParticle(const KVParticle&);
518  virtual ~ KVParticle();
519  void init();
520 #if ROOT_VERSION_CODE >= ROOT_VERSION(3,4,0)
521  virtual void Copy(TObject&) const;
522 #else
523  virtual void Copy(TObject&);
524 #endif
525  virtual void Clear(Option_t* opt = "");
526 
527  virtual void SetMass(Double_t m)
528  {
529  SetXYZM(Px(), Py(), Pz(), m);
530  };
532  {
533  return M();
534  };
535  void SetMomentum(const TVector3& v)
536  {
537  SetXYZM(v(0), v(1), v(2), M());
538  };
539  void SetMomentum(const TVector3* v)
540  {
541  SetXYZM((*v)(0), (*v)(1), (*v)(2), M());
542  };
543  void SetMomentum(Double_t px, Double_t py, Double_t pz, Option_t* opt =
544  "cart");
545  void SetMomentum(Double_t T, TVector3 dir);
546  void SetRandomMomentum(Double_t T, Double_t thmin, Double_t thmax,
547  Double_t phmin, Double_t phmax,
548  Option_t* opt = "isotropic");
549  virtual void Print(Option_t* t = "") const;
550  void Set4Mom(const TLorentzVector& p)
551  {
552  SetVect(p.Vect());
553  SetT(p.E());
554  }
555  void SetE(Double_t a)
556  {
557  SetKE(a);
558  };
559  void SetKE(Double_t ecin);
561  {
562  SetKE(e);
563  };
564  void SetVelocity(const TVector3&);
566  {
567  return Vect();
568  };
570  {
571  TVector3 perp = GetMomentum();
572  perp.SetZ(0);
573  return perp;
574  }
575  Double_t GetKE() const
576  {
577  Double_t e = E();
578  Double_t m = M();
580  return e - m;
581  };
583  {
584  return GetKE();
585  };
587  {
588  return GetKE() * TMath::Power(TMath::Sin(Theta()), 2.0);
589  }
591  {
592  return GetTransverseEnergy();
593  }
594 
596  {
597  return GetKE() - GetTransverseEnergy();
598  }
600  {
601  return GetLongitudinalEnergy();
602  }
603 
605  {
606  Double_t etran = Mt() - GetMass();
607  return etran;
608  }
610  {
611  return GetRTransverseEnergy();
612  }
613  Double_t GetE() const
614  {
615  return GetKE();
616  };
618  {
620  if (GetMomentum().Mag() == 0)
621  return 0;
622  Double_t h = TMath::H() * TMath::C() / TMath::Qe() * 1e9; //in MeV.fm
623  return h / GetMomentum().Mag();
624  };
626  {
628  Double_t h = TMath::H() * TMath::C() / TMath::Qe() * 1e9; //in MeV.fm
629  return h / TMath::Sqrt(TMath::TwoPi() * temp * GetMass());
630  };
631  TVector3 GetVelocity() const;
632  TVector3 GetV() const
633  {
634  return GetVelocity();
635  };
637  {
638  return GetV().z();
639  };
640  Double_t GetVperp() const;
642  {
643  return TMath::RadToDeg() * Theta();
644  };
646  {
647  return TMath::Cos(Theta());
648  }
649  Double_t GetPhi() const
650  {
651  Double_t phi = TMath::RadToDeg() * Phi();
652  return (phi < 0 ? 360. + phi : phi);
653  };
654  void SetTheta(Double_t theta)
655  {
657  }
658  void SetPhi(Double_t phi)
659  {
661  }
662 
663  virtual Bool_t IsOK();
664  void SetIsOK(Bool_t flag = kTRUE);
665  void ResetIsOK()
666  {
668  }
669 
671  {
672  return (KVList*)&fBoosted;
673  }
674  void ls(Option_t* option = "") const;
675 
676  void SetE0(TVector3* e = 0)
677  {
678  if (!fE0)
679  fE0 = new TVector3;
680  if (!e) {
681  *fE0 = GetMomentum();
682  }
683  else {
684  *fE0 = *e;
685  }
686  };
688  {
689  return fE0;
690  };
692  {
694  };
696  {
697  return TestBit(kIsDetected);
698  }
699  KVParticle& operator=(const KVParticle& rhs);
700 
701  void ResetEnergy();
702 
703  const Char_t* GetName() const;
704  void SetName(const Char_t* nom);
705 
706  void AddGroup(const Char_t* groupname, const Char_t* from = "");
707  void AddGroup(const Char_t* groupname, KVParticleCondition*);
708 
709  Bool_t BelongsToGroup(const Char_t* groupname) const;
710  void RemoveGroup(const Char_t* groupname);
711  void RemoveAllGroups();
712  void ListGroups(void) const;
713 
715  void ChangeFrame(const KVFrameTransform&, const KVString& = "");
716  void ChangeDefaultFrame(const Char_t*, const Char_t* defname = "");
717  void SetFrame(const Char_t* frame, const KVFrameTransform&);
718  void SetFrame(const Char_t* newframe, const Char_t* oldframe, const KVFrameTransform&);
719 
720  KVParticle const* GetFrame(const Char_t* frame, Bool_t warn_and_return_null_if_unknown = kTRUE) const;
721 
722  const Char_t* GetFrameName(void) const
723  {
724  return fFrameName;
725  }
726  void SetFrameName(const Char_t* framename)
727  {
730  fFrameName = framename;
731  if (fFrameName != "") SetParameter("frameName", framename);
732  else GetParameters()->RemoveParameter("frameName");
733  }
734 
736  {
737  return (KVNameValueList*)&fParameters;
738  }
739  template<typename ValType> void SetParameter(const Char_t* name, ValType value) const
740  {
741  GetParameters()->SetValue(name, value);
742  }
743 
744  void UpdateAllFrames();
745 
746  ClassDef(KVParticle, 8) //General base class for all massive particles
747 };
748 
749 #endif
int Int_t
#define e(i)
char Char_t
bool Bool_t
double Double_t
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)
Handles particle selection criteria for data analysis classes ,.
Base class for relativistic kinematics of massive particles.
Definition: KVParticle.h:398
Double_t GetWaveLength() const
Definition: KVParticle.h:617
void AddGroups(KVUniqueNameList *un)
list of groups added to the current one
Definition: KVParticle.cpp:535
Double_t GetThermalWaveLength(Double_t temp) const
Definition: KVParticle.h:625
void SetIsOK(Bool_t flag=kTRUE)
Definition: KVParticle.cpp:333
void SetTheta(Double_t theta)
Definition: KVParticle.h:654
TVector3 * GetPInitial() const
Definition: KVParticle.h:687
Double_t GetREtran() const
Definition: KVParticle.h:609
virtual void SetMass(Double_t m)
Definition: KVParticle.h:527
void RemoveGroup(const Char_t *groupname)
Definition: KVParticle.cpp:619
void SetPy(Double_t a)
Definition: KVParticle.h:440
void AddGroup(const Char_t *groupname, const Char_t *from="")
Definition: KVParticle.cpp:485
virtual ~ KVParticle()
void SetPhi(Double_t phi)
Definition: KVParticle.h:658
KVUniqueNameList * GetGroups() const
return the KVUniqueNameList pointeur where list of groups are stored
Definition: KVParticle.cpp:584
void UpdateAllFrames()
TVector3 GetMomentum() const
Definition: KVParticle.h:565
void ResetIsOK()
Definition: KVParticle.h:665
void RemoveAllGroups()
Definition: KVParticle.cpp:646
KVUniqueNameList fGroups
list of TObjString for manage different group name
Definition: KVParticle.h:408
KVNameValueList * GetParameters() const
Definition: KVParticle.h:735
void SetXYZT(Double_t x, Double_t y, Double_t z, Double_t t)
Definition: KVParticle.h:476
virtual Bool_t IsOK()
Definition: KVParticle.cpp:312
Double_t GetTheta() const
Definition: KVParticle.h:641
KVParticle & operator=(const KVParticle &rhs)
KVParticle assignment operator.
Definition: KVParticle.cpp:365
void SetVectM(const TVector3 &spatial, Double_t mass)
Definition: KVParticle.h:416
const Char_t * GetFrameName(void) const
Definition: KVParticle.h:722
void ResetEnergy()
Definition: KVParticle.cpp:385
void SetMomentum(const TVector3 &v)
Definition: KVParticle.h:535
virtual void AddGroup_Withcondition(const Char_t *, KVParticleCondition *)
Definition: KVParticle.cpp:428
void SetVelocity(const TVector3 &)
Set velocity of particle (in cm/ns units)
Bool_t HasFrame(const Char_t *frame)
Definition: KVParticle.h:494
const Char_t * GetName() const
return the field fName
Definition: KVParticle.cpp:416
void SetRho(Double_t rho)
Definition: KVParticle.h:452
Bool_t IsDetected()
Definition: KVParticle.h:695
Double_t GetEnergy() const
Definition: KVParticle.h:582
Double_t GetTransverseEnergy() const
Definition: KVParticle.h:586
void SetVectMag(const TVector3 &spatial, Double_t magnitude)
Definition: KVParticle.h:420
Double_t GetLongitudinalEnergy() const
Definition: KVParticle.h:595
static Double_t C()
Definition: KVParticle.cpp:117
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:289
KVList fBoosted
list of momenta of the particle in different Lorentz-boosted frames
Definition: KVParticle.h:407
void SetPx(Double_t a)
Definition: KVParticle.h:436
void ls(Option_t *option="") const
Definition: KVParticle.cpp:347
void init()
default initialisation
Definition: KVParticle.cpp:50
Double_t GetPhi() const
Definition: KVParticle.h:649
KVNameValueList fParameters
a general-purpose list of parameters associated with this particle
Definition: KVParticle.h:484
Double_t GetCosTheta() const
Definition: KVParticle.h:645
void SetPxPyPzE(Double_t px, Double_t py, Double_t pz, Double_t e)
Definition: KVParticle.h:448
Double_t GetE() const
Definition: KVParticle.h:613
Double_t GetRTransverseEnergy() const
Definition: KVParticle.h:604
void SetMomentum(const TVector3 *v)
Definition: KVParticle.h:539
Double_t GetElong() const
Definition: KVParticle.h:599
virtual void Copy(TObject &) const
Definition: KVParticle.cpp:266
KVList * GetListOfFrames(void)
Definition: KVParticle.h:670
void Set4Mom(const TLorentzVector &p)
Definition: KVParticle.h:550
KVKinematicalFrame * get_parent_frame(const Char_t *, KVKinematicalFrame *F=nullptr) const
virtual void Print(Option_t *t="") const
print out characteristics of particle
Definition: KVParticle.cpp:212
void SetE(Double_t a)
Definition: KVParticle.h:555
KVParticle InFrame(const KVFrameTransform &)
Definition: KVParticle.cpp:688
void SetE0(TVector3 *e=0)
Definition: KVParticle.h:676
void SetY(Double_t a)
Definition: KVParticle.h:464
void SetIsDetected()
Definition: KVParticle.h:691
TString fFrameName
non-persistent frame name field, sets when calling SetFrame method
Definition: KVParticle.h:406
void SetZ(Double_t a)
Definition: KVParticle.h:468
void SetParameter(const Char_t *name, ValType value) const
Definition: KVParticle.h:739
Int_t GetNumberOfDefinedFrames(void)
Definition: KVParticle.cpp:553
void SetX(Double_t a)
Definition: KVParticle.h:460
void SetFrameName(const Char_t *framename)
Definition: KVParticle.h:726
static Double_t kSpeedOfLight
speed of light in cm/ns
Definition: KVParticle.h:409
Double_t GetKE() const
Definition: KVParticle.h:575
Double_t GetVpar() const
Definition: KVParticle.h:636
void SetXYZM(Double_t x, Double_t y, Double_t z, Double_t m)
Definition: KVParticle.h:472
void SetName(const Char_t *nom)
Set Name of the particle.
Definition: KVParticle.cpp:404
void SetFrame(const Char_t *frame, const KVFrameTransform &)
Definition: KVParticle.cpp:840
KVParticle const * GetFrame(const Char_t *frame, Bool_t warn_and_return_null_if_unknown=kTRUE) const
Definition: KVParticle.cpp:952
void ChangeFrame(const KVFrameTransform &, const KVString &="")
Definition: KVParticle.cpp:707
void SetGroups(KVUniqueNameList *un)
Definition: KVParticle.cpp:522
Bool_t BelongsToGroup(const Char_t *groupname) const
Definition: KVParticle.cpp:597
TVector3 * fE0
the momentum of the particle before it is slowed/stopped by an absorber
Definition: KVParticle.h:483
void ChangeDefaultFrame(const Char_t *, const Char_t *defname="")
Definition: KVParticle.cpp:726
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
Int_t GetNumberOfDefinedGroups(void)
return the number of defined groups for the particle
Definition: KVParticle.cpp:573
void SetPtEtaPhiM(Double_t pt, Double_t eta, Double_t phi, Double_t m)
Definition: KVParticle.h:432
TVector3 GetV() const
Definition: KVParticle.h:632
TString fName
non-persistent name field - Is useful
Definition: KVParticle.h:405
KVKinematicalFrame * get_frame(const Char_t *) const
void SetEnergy(Double_t e)
Definition: KVParticle.h:560
Double_t GetMass() const
Definition: KVParticle.h:531
virtual void AddGroup_Sanscondition(const Char_t *groupname, const Char_t *from="")
Definition: KVParticle.cpp:448
void ListGroups(void) const
List all stored groups.
Definition: KVParticle.cpp:666
TVector3 GetVelocity() const
returns velocity vector in cm/ns units
Double_t GetVperp() const
Double_t GetEtran() const
Definition: KVParticle.h:590
void SetPz(Double_t a)
Definition: KVParticle.h:444
void SetPerp(Double_t p)
Definition: KVParticle.h:424
void SetVect(const TVector3 &vect3)
TLorentzVector setters should not be used.
Definition: KVParticle.h:412
void CreateGroups()
TVector3 GetTransverseMomentum() const
Definition: KVParticle.h:569
void SetT(Double_t a)
Definition: KVParticle.h:456
void SetPtEtaPhiE(Double_t pt, Double_t eta, Double_t phi, Double_t e)
Definition: KVParticle.h:428
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.
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
#define F(x, y, z)
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