GRSISort "v4.0.0.5"
An extension of the ROOT analysis Framework
Loading...
Searching...
No Matches
TReaction.h
Go to the documentation of this file.
1#ifndef TREACTION_H
2#define TREACTION_H
3
4/** \addtogroup Fitting Fitting & Analysis
5 * @{
6 */
7
8#include <cmath>
9
10#include "TNamed.h"
11#include "TGraph.h"
12
13#include "TNucleus.h"
14
15#ifndef PI
16#define PI (3.14159265358979312e+00)
17#endif
18
19#ifndef R2D
20#define R2D (5.72957795130823229e+01)
21#endif
22
23#ifndef D2R
24#define D2R (1.74532925199432955e-02)
25#endif
26
27/// *********** *********** *********** *********** *********** *********** *********** //
28///
29/// \class TReaction
30///
31/// A simple, intuitive two body reaction class.
32/// Lab frame variables are generally calculated by boosting from the CM frame
33/// Reactions have the form A(B,C)D where;
34/// [0] - B is "beam"
35/// [1] - A is "targ"
36/// [2] - C is "ejec"
37/// [3] - D is "reco"
38///
39/// Uses planck units (c=1 and is dimensionless)
40/// M = mc^2
41/// E = full mass energy (inc. kinetic energy)
42/// T = kinetic energy
43/// V = magnitude of velocity (=beta)
44/// P = magnitude of momentum
45/// G = gamma factor (=[1-V^2]^-0.5)
46///
47/// * An excitation energy can be included in the final state heavy recoil nucleus using ex3
48/// in the reaction initialisation or by using SetExcEnergy(exc)
49/// * ALL ENERGIES ARE IN MEV - MASSES TOO !
50/// * currently only works for part==2, I think that part==3 just means theta_cm -> -theta_cm
51///
52/// BUGS
53/// THERE ARE PERSISTANT NUMERICAL ERRORS (~0.1-1%) ...??! WHY?!!?!
54/// *********** *********** *********** *********** *********** *********** *********** //
55
56class TReaction : public TNamed {
57public:
58 TReaction(const char* beam, const char* targ, const char* ejec, const char* reco, double eBeam = 0.0, double ex3 = 0.0, bool inverse = false);
59
60 void InitReaction();
61
62 // returns reaction input parameters
63 const char* GetNameFull() const { return Form("%s @ %.3f MeV/u", this->GetName(), fTLab[0] / fNuc[0]->GetA()); }
64 TNucleus* GetNucleus(int part) const { return fNuc[part]; }
65 double GetM(int part) const { return fM[part]; }
66 double GetExc() const { return fExc; }
67 double GetQVal() const { return fQVal; }
68 bool Inverse() const { return fInverse; }
69 double GetTBeam(bool inverse) const;
70 double GetVBeam() const { return fVLab[0]; }
71
72 // CM frame properties
73 double GetInvariantMass() const { return fInvariantMass; }
74 double GetCmE() const { return fCmE; }
75 double GetCmTi() const { return fCmTi; }
76 double GetCmTf() const { return fCmTf; }
77 double GetCmV() const { return fCmV; }
78 double GetCmP() const { return fCmP; }
79 double GetCmG() const { return fCmG; }
80
81 // particle properties in CM frame
82 double GetECm(int part) const { return fECm[part]; }
83 double GetTCm(int part) const { return fTCm[part]; }
84 double GetVCm(int part) const { return fVCm[part]; }
85 double GetPCm(int part) const { return fPCm[part]; }
86 double GetGCm(int part) const { return fGCm[part]; }
87
88 // particle properties in LAB frame (default is beam)
89 double GetThetaMax(int part) const { return fThetaMax[part]; }
90 // this stuff depends on the lab angle as the cm motion and the particle in the cm are coupled
91 double GetELab(double theta_lab = 0.0, int part = 0) const
92 {
93 return GetELabFromThetaCm(ConvertThetaLabToCm(theta_lab, part), part);
94 }
95 double GetTLab(double theta_lab = 0.0, int part = 0) const
96 {
97 return GetTLabFromThetaCm(ConvertThetaLabToCm(theta_lab, part), part);
98 }
99 double GetVLab(double theta_lab = 0.0, int part = 0) const
100 {
101 return GetVLabFromThetaCm(ConvertThetaLabToCm(theta_lab, part), part);
102 }
103 double GetPLab(double theta_lab = 0.0, int part = 0) const
104 {
105 return GetPLabFromThetaCm(ConvertThetaLabToCm(theta_lab, part), part);
106 }
107 double GetGLab(double theta_lab = 0.0, int part = 0) const
108 {
109 return GetGLabFromThetaCm(ConvertThetaLabToCm(theta_lab, part), part);
110 }
111
112 // this stuff depends on the CM angle as the cm motion and the particle in the cm are coupled
113 double GetELabFromThetaCm(double theta_cm = 0.0, int part = 0) const; // FULL MASS+KINETIC ENERGY
114 double GetTLabFromThetaCm(double theta_cm = 0.0, int part = 0) const; // KINETIC ENERGY
115 double GetVLabFromThetaCm(double theta_cm = 0.0, int part = 0) const;
116 double GetPLabFromThetaCm(double theta_cm = 0.0, int part = 0) const;
117 double GetGLabFromThetaCm(double theta_cm = 0.0, int part = 0) const;
118
119 double GetExcEnergy(double ekin = 0.00, double theta_lab = 0.00, int part = 2) const;
120 void AnalysisAngDist(double ekin, double theta_lab, int part, double& exc, double& theta_cm, double& omega_lab2cm);
121 double AnalysisBeta(double ekin, int part) const;
122
123 double GetRutherfordCm(double theta_cm, int part = 2, bool Units_mb = true) const;
124 double GetRutherfordLab(double theta_lab, int part = 2, bool Units_mb = true) const;
125
126 // Conversion from LAB frame to CM frame
127 double ConvertThetaLabToCm(double theta_lab, int part = 2) const;
128 double ConvertOmegaLabToCm(double theta_lab, int part = 2) const;
129 void ConvertLabToCm(double theta_lab, double omega_lab, double& theta_cm, double& omega_cm, int part = 2) const;
130
131 // Conversion from CM frame to LAB frame
132 double ConvertThetaCmToLab(double theta_cm, int part = 2) const;
133 double ConvertOmegaCmToLab(double theta_cm, int part = 2) const;
134 void ConvertCmToLab(double theta_cm, double omega_cm, double& theta_lab, double& omega_lab, int part = 2) const;
135
136 // Graphs for conversions and kinematic/cross-section curves
137 // Frame_Lab -> TLab[ThetaLab] and Frame_Cm -> TLab[ThetaCm]
138 TGraph* KinVsTheta(double thmin = 0.0, double thmax = 180.0, int part = 2, bool Frame_Lab = true, bool Units_keV = true) const;
139 // Frame_Lab -> ThetaCm[ThetaLab] and Frame_Cm -> ThetaLab[ThetaCm]
140 TGraph* ThetaVsTheta(double thmin = 0.0, double thmax = 180.0, int part = 2, bool Frame_Lab = true) const;
141 // Frame_Lab -> dOmegaCm/dOmegaLab[ThetaLab] and Frame_Cm -> dOmegaLab/dOmegaCm[ThetaCm]
142 TGraph* OmegaVsTheta(double thmin = 0.0, double thmax = 180.0, int part = 2, bool Frame_Lab = true) const;
143 // Frame_Lab -> dSigma/dThetaLab[ThetaLab] and Frame_Cm -> dSigma/dThetaCm[ThetaCm]
144 TGraph* RutherfordVsTheta(double thmin = 1.0, double thmax = 179.0, int part = 2, bool Frame_Lab = true, bool Units_mb = true) const;
145
146 void Print(Option_t* opt = "") const override;
147 void Clear(Option_t* opt = "") override;
148
149 void SetExcEnergy(double exc) { SetCmFrame(exc); }
150
151private:
152 void SetCmFrame(double exc); // enables the reaction to be modified using excitation energy
153
154 // same as in TKinematics, should probably change from arrays to individual members
155 // USER INPUTS
156 TNucleus* fNuc[4]{nullptr}; // NOLINT(*-avoid-c-arrays)
157 double fTBeam{0.};
158 bool fInverse{false};
159 double fExc{0.};
160 double fM[4]{0.}; // NOLINT(*-avoid-c-arrays)
161
162 // CM FRAME MOTION
163 double fQVal{0.}; // effective Q value (includes excitation)
164 double fS{0.}; // 'S' = M^2
165 double fInvariantMass{0.};
166 double fCmTi{0.};
167 double fCmTf{0.};
168 double fCmE{0.};
169 double fCmV{0.};
170 double fCmP{0.};
171 double fCmG{0.};
172
173 // PARTICLES IN CM FRAME
174 double fTCm[4]{0.}; // NOLINT(*-avoid-c-arrays)
175 double fECm[4]{0.}; // NOLINT(*-avoid-c-arrays)
176 double fPCm[4]{0.}; // NOLINT(*-avoid-c-arrays)
177 double fVCm[4]{0.}; // NOLINT(*-avoid-c-arrays)
178 double fGCm[4]{0.}; // NOLINT(*-avoid-c-arrays)
179
180 // PARTICLE IN LAB FRAME
181 // Note that in the lab frame only the initial state (beam/targ) is fixed in the reaction
182 double fTLab[2]{0.}; // NOLINT(*-avoid-c-arrays)
183 double fELab[2]{0.}; // NOLINT(*-avoid-c-arrays)
184 double fPLab[2]{0.}; // NOLINT(*-avoid-c-arrays)
185 double fVLab[2]{0.}; // NOLINT(*-avoid-c-arrays)
186 double fGLab[2]{0.}; // NOLINT(*-avoid-c-arrays)
187 double fThetaMax[4]{0.}; // NOLINT(*-avoid-c-arrays)
188
189 /// \cond CLASSIMP
190 ClassDefOverride(TReaction, 1) // NOLINT(readability-else-after-return)
191 /// \endcond
192};
193/*! @} */
194#endif
double fCmG
Definition TReaction.h:171
TGraph * OmegaVsTheta(double thmin=0.0, double thmax=180.0, int part=2, bool Frame_Lab=true) const
double fTLab[2]
Definition TReaction.h:182
double GetCmV() const
Definition TReaction.h:77
TGraph * KinVsTheta(double thmin=0.0, double thmax=180.0, int part=2, bool Frame_Lab=true, bool Units_keV=true) const
double ConvertThetaLabToCm(double theta_lab, int part=2) const
double GetVLabFromThetaCm(double theta_cm=0.0, int part=0) const
void ConvertCmToLab(double theta_cm, double omega_cm, double &theta_lab, double &omega_lab, int part=2) const
double AnalysisBeta(double ekin, int part) const
double GetELabFromThetaCm(double theta_cm=0.0, int part=0) const
double GetGCm(int part) const
Definition TReaction.h:86
double fECm[4]
Definition TReaction.h:175
void Clear(Option_t *opt="") override
double GetELab(double theta_lab=0.0, int part=0) const
Definition TReaction.h:91
double ConvertThetaCmToLab(double theta_cm, int part=2) const
double GetCmE() const
Definition TReaction.h:74
double GetTLab(double theta_lab=0.0, int part=0) const
Definition TReaction.h:95
double GetTBeam(bool inverse) const
Definition TReaction.cxx:73
double ConvertOmegaLabToCm(double theta_lab, int part=2) const
double fM[4]
Definition TReaction.h:160
double fQVal
Definition TReaction.h:163
double GetGLab(double theta_lab=0.0, int part=0) const
Definition TReaction.h:107
double fS
Definition TReaction.h:164
double fCmV
Definition TReaction.h:169
double fCmTi
Definition TReaction.h:166
double fTBeam
Definition TReaction.h:157
double fGCm[4]
Definition TReaction.h:178
double GetInvariantMass() const
Definition TReaction.h:73
double fVCm[4]
Definition TReaction.h:177
double GetTCm(int part) const
Definition TReaction.h:83
double GetPLabFromThetaCm(double theta_cm=0.0, int part=0) const
double GetVBeam() const
Definition TReaction.h:70
void Print(Option_t *opt="") const override
const char * GetNameFull() const
Definition TReaction.h:63
TReaction(const char *beam, const char *targ, const char *ejec, const char *reco, double eBeam=0.0, double ex3=0.0, bool inverse=false)
Definition TReaction.cxx:10
double fGLab[2]
Definition TReaction.h:186
TGraph * RutherfordVsTheta(double thmin=1.0, double thmax=179.0, int part=2, bool Frame_Lab=true, bool Units_mb=true) const
void ConvertLabToCm(double theta_lab, double omega_lab, double &theta_cm, double &omega_cm, int part=2) const
bool Inverse() const
Definition TReaction.h:68
void AnalysisAngDist(double ekin, double theta_lab, int part, double &exc, double &theta_cm, double &omega_lab2cm)
double GetThetaMax(int part) const
Definition TReaction.h:89
double GetExc() const
Definition TReaction.h:66
double fInvariantMass
Definition TReaction.h:165
TNucleus * GetNucleus(int part) const
Definition TReaction.h:64
TGraph * ThetaVsTheta(double thmin=0.0, double thmax=180.0, int part=2, bool Frame_Lab=true) const
double fPLab[2]
Definition TReaction.h:184
double GetRutherfordCm(double theta_cm, int part=2, bool Units_mb=true) const
double GetCmG() const
Definition TReaction.h:79
double GetCmTf() const
Definition TReaction.h:76
double GetCmTi() const
Definition TReaction.h:75
void InitReaction()
Definition TReaction.cxx:35
double ConvertOmegaCmToLab(double theta_cm, int part=2) const
double GetGLabFromThetaCm(double theta_cm=0.0, int part=0) const
double fELab[2]
Definition TReaction.h:183
double fPCm[4]
Definition TReaction.h:176
double GetCmP() const
Definition TReaction.h:78
double GetM(int part) const
Definition TReaction.h:65
bool fInverse
Definition TReaction.h:158
double fTCm[4]
Definition TReaction.h:174
double GetPCm(int part) const
Definition TReaction.h:85
double fCmP
Definition TReaction.h:170
double fExc
Definition TReaction.h:159
double fCmE
Definition TReaction.h:168
double GetQVal() const
Definition TReaction.h:67
void SetExcEnergy(double exc)
Definition TReaction.h:149
double GetRutherfordLab(double theta_lab, int part=2, bool Units_mb=true) const
double GetPLab(double theta_lab=0.0, int part=0) const
Definition TReaction.h:103
void SetCmFrame(double exc)
Definition TReaction.cxx:81
double GetExcEnergy(double ekin=0.00, double theta_lab=0.00, int part=2) const
double GetVLab(double theta_lab=0.0, int part=0) const
Definition TReaction.h:99
double GetVCm(int part) const
Definition TReaction.h:84
double fCmTf
Definition TReaction.h:167
TNucleus * fNuc[4]
Definition TReaction.h:156
double fVLab[2]
Definition TReaction.h:185
double fThetaMax[4]
Definition TReaction.h:187
double GetTLabFromThetaCm(double theta_cm=0.0, int part=0) const
double GetECm(int part) const
Definition TReaction.h:82