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