GRSISort "v4.0.0.5"
An extension of the ROOT analysis Framework
Loading...
Searching...
No Matches
TKinematics.h
Go to the documentation of this file.
1#ifndef TKINEMATICS_H
2#define TKINEMATICS_H
3
4/** \addtogroup Fitting Fitting & Analysis
5 * @{
6 */
7
8#include "TNamed.h"
9#include "TSpline.h"
10#include "TGraph.h"
11#include "TLorentzVector.h"
12
13#include "TNucleus.h"
14
15#ifndef PI
16#define PI (TMath::Pi())
17#endif
18
19/////////////////////////////////////////////////////////////////
20///
21/// \class TKinematics
22///
23/// This class calculates 2 body kinematics from a beam, target,
24/// recoil, ejectile, and beam energy.
25///
26/////////////////////////////////////////////////////////////////
27
28class TKinematics : public TNamed {
29public:
30 TKinematics(double eBeam, const char* beam, const char* targ, const char* ejec = nullptr, const char* reco = nullptr, const char* name = "");
31 TKinematics(TNucleus* projectile, TNucleus* target, double eBeam, const char* name = "");
32 TKinematics(TNucleus* projectile, TNucleus* target, TNucleus* recoil, TNucleus* ejectile, double eBeam, const char* name = "");
33 TKinematics(TNucleus* projectile, TNucleus* target, TNucleus* recoil, TNucleus* ejectile, double eBeam, double ex3, const char* name = "");
34 TKinematics(const char* beam, const char* targ, const char* ejec, const char* reco, double eBeam, double ex3 = 0.0, const char* name = "");
35
36 void Initial();
37 void FinalCm();
38 void Final(double angle, int part);
39 // void SetAngles(double angle, int part);
40 void SetAngles(double angle, int part, bool upper = false);
41 // get
42 TSpline3* Evslab(double thmin, double thmax, double size, int part = 2);
43 TSpline3* Evscm(double thmin, double thmax, double size, int part = 2);
44 TSpline3* labvscm(double thmin, double thmax, double size, int part = 2) const;
45 TSpline3* cmvslab(double thmin, double thmax, double size, int part = 2) const;
46 TSpline3* Steffen_labvscminverse(double thmin, double thmax, double size, int part = 2) const;
47
48 TGraph* Evslab_graph(double thmin, double thmax, double size, int part = 2);
49
50 double GetCmEnergy(double eBeam) const;
51 double GetCmEnergy() const;
52 double GetBeamEnergy(double LabAngle, double LabEnergy) const;
53 double GetMaxAngle(double vcm) const;
54 double GetMaxAngle(int part) const;
55 double NormalkinEnergy() const;
56 bool CheckMaxAngle(double angle, int part) const;
57 double GetExcEnergy(TLorentzVector recoil);
58 double GetExcEnergy(TVector3 position, double KinE);
59 double GetExcEnergy(double theta, double KinE);
60 double ELab(double angle_lab, int part);
61
62 double GetQValue() const { return fQValue; }
63 double GetElab(int i) const { return fE[i]; }
64 double GetM(int i) const { return fM[i]; }
65 double GetTlab(int i) const { return fT[i]; }
66 double GetEcm(int i) const { return fEcm[i]; }
67 double GetTcm(int i) const { return fTcm[i]; }
68
69 double GetThetalab(int i) const
70 {
71 if(fTheta[i] < 1e-5) {
72 return 0;
73 }
74 return fTheta[i];
75 }
76 double GetThetacm(int i) const
77 {
78 if(fThetacm[i] < 1e-5) {
79 return 0;
80 }
81 return fThetacm[i];
82 }
83
84 double GetBetacm() const { return fBeta_cm; }
85 double GetGammacm() const { return fGamma_cm; }
86 double GetBetacm(int i) const { return fBetacm[i]; }
87 double GetVcm(int i) const { return fVcm[i]; }
88 double GetV(int i) const { return fV[i]; }
89
90 TSpline3* Ruthvscm(double thmin, double thmax, double size) const;
91 TSpline3* Ruthvslab(double thmin, double thmax, double size, int part) const;
92 double Angle_lab2cm(double vcm, double angle_lab) const;
93 double Angle_lab2cminverse(double vcm, double angle_lab, bool upper = true) const;
94 double Steffen_cm2labinverse(double theta_cm, int part = 2) const; // NEW FUNCTIN+
95 double Steffen_lab2cminverse(double theta_lab); // assumes part = 2;
96 double Angle_cm2lab(double vcm, double angle_cm) const;
97 double Sigma_cm2lab(double angle_cm, double sigma_cm) const;
98 double Sigma_lab2cm(double angle_cm, double sigma_lab) const;
99 void Transform2cm(double& angle, double& sigma) const;
100 void Transform2cm(double& angle, double& errangle, double& sigma, double& errsigma) const;
101 void AngleErr_lab2cm(double angle, double& err) const;
102 void SigmaErr_lab2cm(double angle, double err, double& sigma, double& errsigma) const;
103 double Rutherford(double angle_cm) const;
104
105private:
106 double fTCm_i{0.};
107 double fTCm_f{0.};
108
109 // I don't think it makes sense to have these as arrays, we never use these as arrays (apart from setting the masses)
110 // and it would make it easier to understand to see fProjectile or fRecoil instead of fParticle[0] or fParticle[2]
111 TNucleus* fParticle[4]{nullptr}; // NOLINT(*-avoid-c-arrays)
112 double fM[4]{0.}; // NOLINT(*-avoid-c-arrays)
113 double fEBeam{0.};
114 double fQValue{0.};
115 double fT[4]{0.}; // NOLINT(*-avoid-c-arrays)
116 double fE[4]{0.}; // NOLINT(*-avoid-c-arrays)
117 double fP[4]{0.}; // NOLINT(*-avoid-c-arrays)
118 double fV[4]{0.}; // NOLINT(*-avoid-c-arrays)
119 double fTheta[4]{0.}; // NOLINT(*-avoid-c-arrays)
120
121 double fTcm[4]{0.}; // NOLINT(*-avoid-c-arrays)
122 double fEcm[4]{0.}; // NOLINT(*-avoid-c-arrays)
123 double fPcm[4]{0.}; // NOLINT(*-avoid-c-arrays)
124 double fVcm[4]{0.}; // NOLINT(*-avoid-c-arrays)
125 double fBetacm[4]{0.}; // NOLINT(*-avoid-c-arrays)
126 double fThetacm[4]{0.}; // NOLINT(*-avoid-c-arrays)
127
128 double fBeta_cm{0.};
129 double fGamma_cm{0.};
130
131 TSpline3* fCm2LabSpline{nullptr};
132
133 double Pcm_em(double, double) const;
134 double P_tm(double, double) const;
135 double E_tm(double, double) const;
136 double T_em(double, double) const;
137 double betacm_tm(double, double) const;
138 double V_pe(double, double) const;
139 double E_final(int) const;
140 double T_final(int) const;
141
142 /// \cond CLASSIMP
143 ClassDefOverride(TKinematics, 1) // NOLINT(readability-else-after-return)
144 /// \endcond
145};
146/*! @} */
147#endif
double fTCm_i
double fQValue
double fTheta[4]
double Angle_cm2lab(double vcm, double angle_cm) const
double Sigma_lab2cm(double angle_cm, double sigma_lab) const
double GetVcm(int i) const
Definition TKinematics.h:87
double GetCmEnergy() const
double GetQValue() const
Definition TKinematics.h:62
double GetThetacm(int i) const
Definition TKinematics.h:76
double GetM(int i) const
Definition TKinematics.h:64
double GetExcEnergy(TLorentzVector recoil)
double NormalkinEnergy() const
double GetGammacm() const
Definition TKinematics.h:85
double fVcm[4]
void AngleErr_lab2cm(double angle, double &err) const
double GetMaxAngle(double vcm) const
double GetElab(int i) const
Definition TKinematics.h:63
double fTcm[4]
double fM[4]
double fThetacm[4]
TSpline3 * fCm2LabSpline
double GetEcm(int i) const
Definition TKinematics.h:66
double fV[4]
TSpline3 * Ruthvslab(double thmin, double thmax, double size, int part) const
double Angle_lab2cminverse(double vcm, double angle_lab, bool upper=true) const
TSpline3 * labvscm(double thmin, double thmax, double size, int part=2) const
TSpline3 * cmvslab(double thmin, double thmax, double size, int part=2) const
TSpline3 * Evscm(double thmin, double thmax, double size, int part=2)
bool CheckMaxAngle(double angle, int part) const
TKinematics(double eBeam, const char *beam, const char *targ, const char *ejec=nullptr, const char *reco=nullptr, const char *name="")
double fT[4]
double Sigma_cm2lab(double angle_cm, double sigma_cm) const
double Pcm_em(double, double) const
double GetBetacm(int i) const
Definition TKinematics.h:86
void Transform2cm(double &angle, double &sigma) const
double P_tm(double, double) const
double betacm_tm(double, double) const
double fBetacm[4]
double V_pe(double, double) const
TNucleus * fParticle[4]
double fPcm[4]
double Rutherford(double angle_cm) const
double GetBetacm() const
Definition TKinematics.h:84
void SigmaErr_lab2cm(double angle, double err, double &sigma, double &errsigma) const
double fGamma_cm
double GetBeamEnergy(double LabAngle, double LabEnergy) const
double Steffen_cm2labinverse(double theta_cm, int part=2) const
double GetTcm(int i) const
Definition TKinematics.h:67
TSpline3 * Ruthvscm(double thmin, double thmax, double size) const
double fTCm_f
TGraph * Evslab_graph(double thmin, double thmax, double size, int part=2)
double fBeta_cm
double fP[4]
TSpline3 * Steffen_labvscminverse(double thmin, double thmax, double size, int part=2) const
void SetAngles(double angle, int part, bool upper=false)
TSpline3 * Evslab(double thmin, double thmax, double size, int part=2)
void Final(double angle, int part)
double E_final(int) const
double GetTlab(int i) const
Definition TKinematics.h:65
double fE[4]
double E_tm(double, double) const
double Steffen_lab2cminverse(double theta_lab)
double fEcm[4]
double T_final(int) const
double ELab(double angle_lab, int part)
double GetV(int i) const
Definition TKinematics.h:88
double GetThetalab(int i) const
Definition TKinematics.h:69
double T_em(double, double) const
double Angle_lab2cm(double vcm, double angle_lab) const
double fEBeam