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