GRSISort "v4.0.0.5"
An extension of the ROOT analysis Framework
Loading...
Searching...
No Matches
TReaction.cxx
Go to the documentation of this file.
1// g++ -c -fPIC TReaction.cxx -I./ `root-config --cflags`
2
3#include "TReaction.h"
4
5#include <iostream>
6
7#include "TStyle.h"
8#include "TMath.h"
9
10TReaction::TReaction(const char* beam, const char* targ, const char* ejec, const char* reco, double eBeam, double ex3, bool inverse)
11 : fTBeam(eBeam), fInverse(inverse), fExc(ex3)
12{
13 Clear();
14 // I THINK INVERSE KINEMATICS NECESSITATES THE BEAM<->TARGET ENTIRELY ?
15 fNuc[0] = new TNucleus(beam);
16 fNuc[1] = new TNucleus(targ);
17 fNuc[2] = new TNucleus(ejec);
18 fNuc[3] = new TNucleus(reco);
19
20 for(int i = 0; i < 4; i++) {
21 fM[i] = fNuc[i]->GetMass();
22 }
23
24 fQVal = (fM[0] + fM[1]) - (fM[2] + fM[3]) - ex3; // effective Q value (includes excitation)
25
26 if(fInverse) {
27 SetName(Form("%s(%s,%s)%s", beam, targ, ejec, reco));
28 } else {
29 SetName(Form("%s(%s,%s)%s", targ, beam, ejec, reco));
30 }
31
33}
34
36{
37
38 // An initializing function that sets the energies and momenta of the beam and target in the lab and CM frame,
39 // as well as a few basic calculations.
40
41 fTLab[0] = fTBeam; // target kinetic energy is ebeam
42 fELab[0] = fTLab[0] + fM[0]; // total E of beam in lab
43 fPLab[0] = sqrt(pow(fTLab[0], 2) + 2 * fTLab[0] * fM[0]); // momentum of beam in lab
44 fVLab[0] = fPLab[0] / fELab[0]; // velocity of beam in lab
45 fGLab[0] = 1 / sqrt(1 - pow(fVLab[0], 2)); // gamma factor of beam
46
47 fTLab[1] = 0; // target kinetic energy is always 0
48 fELab[1] = fM[1]; // total E of target in lab
49 fPLab[1] = 0; // momentum of target in lab
50 fVLab[1] = 0; // velocity of target in lab
51 fGLab[1] = 1; // gamma factor of target
52
53 fS = pow(fM[0], 2) + pow(fM[1], 2) + 2 * fELab[0] * fM[1]; // fELab[1] = fM[1]
54 fInvariantMass = sqrt(fS); // √s
55
56 // CM motion
58 fCmTi = fCmE - fM[0] - fM[1];
59 fCmTf = fCmTi + fQVal;
60 fCmV = fPLab[0] / (fELab[0] + fM[1]);
61 fCmP = fCmV * fCmE;
62 fCmG = 1 / sqrt(1 - pow(fCmV, 2));
63
64 // take care of particles in CM frame using the reaction excitation energy
66
67 // options to make graphs draw nicely
68 gStyle->SetTitleYOffset(1.5);
69 gStyle->SetTitleXOffset(1.2);
70 gStyle->SetDrawOption("AC");
71}
72
73double TReaction::GetTBeam(bool inverse) const
74{
75 if(fInverse || inverse) {
76 return (fGLab[0] - 1) * fM[1];
77 }
78 return fTLab[0];
79}
80
81void TReaction::SetCmFrame(double exc)
82{
83
84 // particles in CM frame
85 fPCm[0] = sqrt((fS - pow(fM[0] - fM[1], 2)) * (fS - pow(fM[0] + fM[1], 2))) / (2 * sqrt(fS));
86 fPCm[1] = fPCm[0];
87 fPCm[2] = sqrt((fS - pow(exc + fM[2] - fM[3], 2)) * (fS - pow(exc + fM[2] + fM[3], 2))) / (2 * sqrt(fS));
88 fPCm[3] = fPCm[2];
89
90 for(int i = 0; i < 4; i++) {
91 fECm[i] = sqrt(pow(fM[i], 2) + pow(fPCm[i], 2));
92 fTCm[i] = fECm[i] - fM[i];
93 fVCm[i] = fPCm[i] / fECm[i];
94 fGCm[i] = 1 / sqrt(1 - pow(fVCm[i], 2));
95
96 // max lab frame theta
97 if(i < 2) {
98 fThetaMax[i] = 0;
99 } else {
100 double val = fPCm[i] / (fM[i] * fCmV * fCmG);
101 if(val < 1) {
102 fThetaMax[i] = asin(val);
103 } else if(val < 1.001) { // catches elastic channels with small numerical rounding errors
104 fThetaMax[i] = PI / 2;
105 } else {
106 fThetaMax[i] = PI;
107 }
108 }
109 }
110}
111
112// particle properties in LAB frame
113// These guys do the math to get lab frame values using CM angles
114double TReaction::GetELabFromThetaCm(double theta_cm, int part) const
115{
116 if(part == 0 || part == 1) {
117 return fTLab[part];
118 }
119
120 return fCmG * (fECm[part] - fCmV * fPCm[part] * cos(theta_cm));
121}
122
123double TReaction::GetTLabFromThetaCm(double theta_cm, int part) const
124{
125 if(part == 0 || part == 1) {
126 return fTLab[part];
127 }
128
129 double ELab = GetELabFromThetaCm(theta_cm, part);
130 return ELab - fM[part]; // T = E - M
131}
132
133double TReaction::GetVLabFromThetaCm(double theta_cm, int part) const
134{
135 if(part == 0 || part == 1) {
136 return fTLab[part];
137 }
138
139 double PLab = GetPLabFromThetaCm(theta_cm, part);
140 double ELab = GetELabFromThetaCm(theta_cm, part);
141 return PLab / ELab; // V = P / E
142}
143
144double TReaction::GetPLabFromThetaCm(double theta_cm, int part) const
145{
146 if(part == 0 || part == 1) {
147 return fTLab[part];
148 }
149
150 double Pz = fCmG * (fPCm[part] * cos(theta_cm) - fCmV * fECm[part]);
151 double Pperp = fPCm[part] * sin(theta_cm);
152 return sqrt(pow(Pperp, 2) + pow(Pz, 2));
153}
154
155double TReaction::GetGLabFromThetaCm(double theta_cm, int part) const
156{
157 if(part == 0 || part == 1) {
158 return fTLab[part];
159 }
160
161 double VLab = GetVLabFromThetaCm(theta_cm, part);
162 return 1 / sqrt(1 - pow(VLab, 2));
163}
164
165double TReaction::GetExcEnergy(double ekin, double theta_lab, int) const
166{
167 if(ekin == 0.00 && theta_lab == 0.00) {
168 return fExc;
169 }
170
171 double val1 = pow(fM[0] + fM[1], 2) + pow(fM[2], 2) + 2 * fM[1] * fTLab[0];
172 double val2 = 2 * fCmG * sqrt(pow(fM[0] + fM[1], 2) + 2 * fM[1] * fTLab[0]);
173 double val3 = fM[2] + ekin - fCmV * sqrt(pow(ekin, 2) + 2 * fM[2] * ekin) * TMath::Cos(theta_lab);
174
175 return sqrt(val1 - val2 * val3) - fM[3];
176}
177
178void TReaction::AnalysisAngDist(double ekin, double theta_lab, int part, double& exc, double& theta_cm, double& omega_lab2cm)
179{
180 exc = GetExcEnergy(ekin, theta_lab, part);
181
182 // adjust cm frame to include the excited state
183 SetCmFrame(exc);
184
185 theta_cm = ConvertThetaLabToCm(theta_lab, part);
186 omega_lab2cm = ConvertOmegaLabToCm(theta_lab, part);
187
188 // reset the cm frame to normal
190}
191
192double TReaction::AnalysisBeta(double ekin, int part) const
193{
194 return sqrt(pow(ekin, 2) + 2 * fM[part] * ekin) / (ekin + fM[part]);
195}
196
197// THIS IS ACTUALLY MOTT SCATTERING (RELATIVISTIC RUTHERFORD)
198// taken from http://www7b.biglobe.ne.jp/~kcy05t/rathef.html (eqn. 61)
199// alpha obtained from http://www.pas.rochester.edu/~cline/Gosia/Gosia_Manual_20120510.pdf
200double TReaction::GetRutherfordCm(double theta_cm, int, bool Units_mb) const
201{
202 static const double alpha = 1.29596; // 0.359994;
203 double scale = 1;
204 if(!Units_mb) {
205 scale = 0.1; // fm^2 = 10^-30 m^2, mb = 10^-31 m^2
206 }
207 // motion is described in lab frame (so TLab instead of TCm) because
208 // Rutherford scattering approximates lab frame == cm frame
209 double a = pow(fNuc[0]->GetZ() * fNuc[1]->GetZ() / fTLab[0], 2);
210 return scale * alpha * a / pow(sin(theta_cm / 2), 4); // ?
211}
212
213double TReaction::GetRutherfordLab(double theta_lab, int part, bool Units_mb) const
214{
215 double theta_cm = ConvertThetaLabToCm(theta_lab, part);
216 double jacobian = ConvertOmegaCmToLab(theta_cm, part);
217
218 return GetRutherfordCm(theta_cm, part, Units_mb) * jacobian;
219}
220
221// Conversion from LAB frame to CM frame
222double TReaction::ConvertThetaLabToCm(double theta_lab, int part) const
223{
224 theta_lab = std::max(theta_lab, fThetaMax[part]);
225
226 // Uses the particle velocity in the CM frame, which makes it more complex
227 double gtan2 = pow(tan(theta_lab) * fCmG, 2);
228 double x = fCmV / fVCm[part];
229 double expr = sqrt(1 + gtan2 * (1 - pow(x, 2)));
230 double theta_cm = 0.;
231
232 // deals with double valued thetas in lab frame
233 if(tan(theta_lab) >= 0) {
234 theta_cm = acos((-x * gtan2 + expr) / (1 + gtan2));
235 } else {
236 theta_cm = acos((-x * gtan2 - expr) / (1 + gtan2));
237 }
238
239 if(fInverse) {
240 theta_cm = PI - theta_cm;
241 }
242
243 if(part == 3) {
244 theta_cm = -theta_cm;
245 }
246
247 return theta_cm;
248}
249
250// dOmegaLab/dOmegaCm[ThetaLab,ThetaCm]
251double TReaction::ConvertOmegaLabToCm(double theta_lab, int part) const
252{
253 // the way to test this function is to use the known 4*cos(theta_lab) for elastics
254 double theta_cm = ConvertThetaLabToCm(theta_lab, part);
255 return 1 / ConvertOmegaCmToLab(theta_cm, part);
256}
257
258void TReaction::ConvertLabToCm(double theta_lab, double omega_lab, double& theta_cm, double& omega_cm, int part) const
259{
260 theta_cm = ConvertThetaLabToCm(theta_lab, part);
261 omega_cm = ConvertOmegaLabToCm(omega_lab, part);
262}
263
264// Conversion from CM frame to LAB frame
265double TReaction::ConvertThetaCmToLab(double theta_cm, int part) const
266{
267 if(fInverse) {
268 theta_cm = PI - theta_cm;
269 }
270
271 double theta_lab = TMath::ATan2(sin(theta_cm), fCmG * (cos(theta_cm) + fCmV / fVCm[part]));
272
273 if(theta_lab > fThetaMax[part]) {
274 return fThetaMax[part];
275 }
276 return theta_lab;
277}
278
279double TReaction::ConvertOmegaCmToLab(double theta_cm, int part) const
280{
281 // the way to test this function is to use the known 4*cos(theta_lab) for elastics
282 double x = fCmV / fVCm[part];
283 if(fInverse) {
284 theta_cm = PI - theta_cm;
285 }
286 double val1 = pow(pow(sin(theta_cm), 2) + pow(fCmG * (x + cos(theta_cm)), 2), 1.5);
287 double val2 = (fCmG * (1 + x * cos(theta_cm)));
288
289 return val1 / val2;
290}
291
292void TReaction::ConvertCmToLab(double theta_cm, double omega_cm, double& theta_lab, double& omega_lab, int part) const
293{
294 theta_lab = ConvertThetaCmToLab(theta_cm, part);
295 omega_lab = ConvertOmegaCmToLab(omega_cm, part);
296}
297
298////////////////////////////////////////////////////////////////////////////////////
299// // // // // // // GRAPHS // // // // // // // // // // // // // // //
300////////////////////////////////////////////////////////////////////////////////////
301
302// Kinetic energy (lab frame) versus theta (either frame)
303TGraph* TReaction::KinVsTheta(double thmin, double thmax, int part, bool Frame_Lab, bool Units_keV) const
304{
305 auto* g = new TGraph();
306 const char* frame = Form("%s", Frame_Lab ? "Lab" : "Cm");
307
308 g->SetName(Form("KinVsTheta%s_%s", frame, GetName()));
309 g->SetTitle(
310 Form("Kinematics for %s; Theta_{%s} [deg]; Kinetic energy [%s]", GetName(), frame, Units_keV ? "keV" : "MeV"));
311
312 double theta = 0.;
313 double T = 0.;
314
315 for(int i = 0; i <= 180; i++) {
316 theta = static_cast<double>(i); // always in CM frame since function is continuous
317
318 T = GetTLabFromThetaCm(theta * D2R, part);
319 if(Units_keV) {
320 T *= 1e3;
321 }
322
323 if(Frame_Lab) { // this is now converted to specified frame (from Frame_Lab)
324 theta = ConvertThetaCmToLab(theta * D2R, part) * R2D;
325 // if(theta==g->GetX[g->GetN()-1]) // if angle is the same
326 // continue;
327 }
328
329 if(theta < thmin || theta > thmax) {
330 continue; // set angular range
331 }
332 g->SetPoint(i, theta, T);
333 }
334
335 return g;
336}
337
338// Frame_Lab -> ThetaCm[ThetaLab] and Frame_Cm -> ThetaLab[ThetaCm]
339TGraph* TReaction::ThetaVsTheta(double thmin, double thmax, int part, bool Frame_Lab) const
340{
341 auto* g = new TGraph();
342 const char* frame = Form("%s", Frame_Lab ? "Lab" : "Cm");
343 const char* other = Form("%s", !Frame_Lab ? "Lab" : "Cm");
344
345 g->SetName(Form("ThetaVsTheta%s_%s", frame, GetName()));
346 g->SetTitle(Form("Angle conversion for %s; Theta_{%s} [deg]; Theta_{%s} [deg]", GetName(), frame, other));
347
348 double theta_cm = 0.;
349 double theta_lab = 0.;
350
351 for(int i = 0; i <= 180; i++) {
352 theta_cm = static_cast<double>(i); // always in CM frame
353 theta_lab = ConvertThetaCmToLab(theta_cm * D2R, part) * R2D;
354
355 if((Frame_Lab && (theta_lab < thmin || theta_lab > thmax))) {
356 continue;
357 }
358 if(!Frame_Lab && (theta_cm < thmin || theta_cm > thmax)) {
359 continue;
360 }
361 // set angular range
362
363 if(Frame_Lab) { // this is now converted to specified frame (from Frame_Lab)
364 g->SetPoint(i, theta_lab, theta_cm);
365 } else {
366 g->SetPoint(i, theta_cm, theta_lab);
367 }
368 }
369
370 return g;
371}
372
373// Frame_Lab -> dOmegaCm/dOmegaLab[ThetaLab] and Frame_Cm -> dOmegaLab/dOmegaCm[ThetaCm]
374TGraph* TReaction::OmegaVsTheta(double thmin, double thmax, int part, bool Frame_Lab) const
375{
376 auto* g = new TGraph();
377 const char* frame = Form("%s", Frame_Lab ? "Lab" : "Cm");
378 const char* other = Form("%s", !Frame_Lab ? "Lab" : "Cm");
379
380 g->SetName(Form("%s_OmegaVsTheta%s", GetName(), frame));
381 g->SetTitle(Form("Solid angle conversion for %s; Theta_{%s} [deg]; dOmega_{%s} / dOmega_{%s}", GetName(), frame,
382 other, frame));
383
384 double theta = 0.;
385 double Om = 0.;
386 for(int i = 0; i <= 180; i++) {
387
388 theta = static_cast<double>(i); // always in CM frame
389 Om = 1 / ConvertOmegaCmToLab(theta * D2R, part);
390
391 if(Frame_Lab) { // this is now converted to specified frame (from Frame_Lab)
392 theta = ConvertThetaCmToLab(theta * D2R, part) * R2D;
393 Om = 1 / Om;
394 }
395
396 if(theta < thmin || theta > thmax || Om > 1e3) { //|| Om<=0 || isnan(Om) || isinf(Om))
397 continue; // set angular range and remove singularities
398 }
399 g->SetPoint(g->GetN(), theta, Om);
400 }
401
402 return g;
403}
404
405// Frame_Lab -> dSigma/dOmegaLab[ThetaLab] and Frame_Cm -> dSigma/dOmegaCm[ThetaCm]
406TGraph* TReaction::RutherfordVsTheta(double thmin, double thmax, int part, bool Frame_Lab, bool Units_mb) const
407{
408 auto* g = new TGraph();
409 const char* frame = Form("%s", Frame_Lab ? "Lab" : "Cm");
410
411 g->SetName(Form("%s_RutherfordVsTheta%s", GetName(), frame));
412 g->SetTitle(Form("Rutherford cross section for %s; Theta_{%s} [deg]; dSigma / dOmega_{%s} [%s/sr]", GetName(), frame,
413 frame, Units_mb ? "mb" : "fm^2"));
414
415 double theta = 0.;
416 double R = 0.;
417
418 for(int i = 1; i <= 180; i++) {
419 theta = static_cast<double>(i); // always in CM frame
420
421 R = GetRutherfordCm(theta * D2R, part, Units_mb); //*ConvertOmega?
422
423 if(Frame_Lab) { // this is now converted to specified frame (from Frame_Lab)
424 R *= ConvertOmegaCmToLab(theta * D2R, part);
425 theta = ConvertThetaCmToLab(theta * D2R, part) * R2D;
426 }
427
428 if(theta < thmin || theta > thmax) {
429 continue; // set angular range
430 }
431 g->SetPoint(g->GetN(), theta, R);
432 }
433
434 return g;
435}
436
437void TReaction::Print(Option_t* opt) const
438{
439 std::string pstring;
440 pstring.assign(opt);
441
442 std::cout << std::endl
443 << std::endl
444 << " * * * * * * * * * * * * * * * * * * * * * * * * *" << std::endl;
445 std::cout << std::endl
446 << std::endl
447 << "\tTReaction '" << GetName() << "' :" << std::endl
448 << std::endl;
449
450 std::cout << " -> Beam '" << fNuc[0]->GetName() << "' kinetic energy = " << fTBeam << " [MeV]" << std::endl;
451 std::cout << " -> Reaction Q value (total) = " << fQVal << " [MeV]" << std::endl;
452 std::cout << " -> Reaction kinematics type = '" << (fInverse ? "INVERSE" : "NORMAL") << "'" << std::endl;
453 if(fInverse) {
454 std::cout << std::endl
455 << " Inverse beam '" << fNuc[1]->GetName() << "' [lab frame] :- " << std::endl;
456 std::cout << "\t Kinetic energy = " << (fGLab[0] - 1) * fM[1] << " [MeV]" << std::endl;
457 std::cout << "\t Velocity = " << fVLab[0] << " [/c] " << std::endl;
458 }
459
460 std::cout << std::endl
461 << " Center of mass motion :- " << std::endl;
462 std::cout << "\t CmE = " << fCmE << " [MeV]" << std::endl;
463 std::cout << "\t CmTi = " << fCmTi << " [MeV]" << std::endl;
464 std::cout << "\t CmTf = " << fCmTf << " [MeV]" << std::endl;
465 std::cout << "\t CmV = " << fCmV << " [/c]" << std::endl;
466 std::cout << "\t CmP = " << fCmP << " [MeV/c]" << std::endl;
467 std::cout << "\t CmG = " << fCmG << std::endl;
468 std::cout << std::endl;
469
470 if(pstring.find("all") != std::string::npos) {
471 for(int i = 0; i < 4; i++) {
472 std::cout << std::endl
473 << " Particle " << i << " : '" << fNuc[i]->GetName() << "' : \t A = " << fNuc[i]->GetA() << ", Z = " << fNuc[i]->GetZ() << ", Mass = " << fM[i] << " [MeV]" << std::endl;
474
475 if(i < 2) {
476 std::cout << "\t ECm = " << fECm[i] << " [MeV]\t\t ELab = " << fELab[i] << " [MeV]" << std::endl;
477 std::cout << "\t TCm = " << fTCm[i] << " [MeV]\t\t TLab = " << fTLab[i] << " [MeV]" << std::endl;
478 std::cout << "\t VCm = " << fVCm[i] << " [/c]\t\t VLab = " << fVLab[i] << " [/c]" << std::endl;
479 std::cout << "\t PCm = " << fPCm[i] << " [MeV/c]\t PLab = " << fPLab[i] << " [MeV/c]" << std::endl;
480 std::cout << "\t GCm = " << fGCm[i] << " \t\t GLab = " << fGLab[i] << std::endl;
481 } else {
482 std::cout << "\t ECm = " << fECm[i] << " [MeV]\t\t ELab = N/A" << std::endl;
483 std::cout << "\t TCm = " << fTCm[i] << " [MeV]\t\t TLab = N/A" << std::endl;
484 std::cout << "\t VCm = " << fVCm[i] << " [/c]\t\t VLab = N/A" << std::endl;
485 std::cout << "\t PCm = " << fPCm[i] << " [MeV/c]\t PLab = N/A" << std::endl;
486 std::cout << "\t GCm = " << fGCm[i] << " \t\t GLab = N/A" << std::endl;
487 std::cout << "\t\t ThetaLab_max = " << fThetaMax[i] * R2D << " [deg]" << std::endl;
488 }
489 }
490 }
491 std::cout << std::endl
492 << std::endl
493 << " * * * * * * * * * * * * * * * * * * * * * * * * *" << std::endl
494 << std::endl;
495}
496
497void TReaction::Clear(Option_t*)
498{
499 fQVal = 0.;
500 fS = 0.;
501 fInvariantMass = 0.;
502 fTBeam = 0.;
503 fInverse = false;
504
505 fCmTi = 0.;
506 fCmTf = 0.;
507 fCmE = 0.;
508 fCmV = 0.;
509 fCmP = 0.;
510 fCmG = 0.;
511
512 for(int i = 0; i < 4; i++) {
513 fNuc[i] = nullptr;
514 fM[i] = 0.;
515
516 fTCm[i] = 0.;
517 fECm[i] = 0.;
518 fVCm[i] = 0.;
519 fPCm[i] = 0.;
520 fGCm[i] = 0.;
521
522 fThetaMax[i] = 0.;
523
524 if(i < 2) {
525 fTLab[i] = 0.;
526 fELab[i] = 0.;
527 fVLab[i] = 0.;
528 fPLab[i] = 0.;
529 fGLab[i] = 0.;
530 }
531 }
532}
int GetZ() const
Gets the Z (# of protons) of the nucleus.
Definition TNucleus.h:57
int GetA() const
Gets the A (Z + N) of the nucleus.
Definition TNucleus.h:59
double GetMass() const
Gets the mass of the nucleus (in MeV)
Definition TNucleus.h:61
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
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 fECm[4]
Definition TReaction.h:175
void Clear(Option_t *opt="") override
double ConvertThetaCmToLab(double theta_cm, int part=2) const
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 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 fVCm[4]
Definition TReaction.h:177
double GetPLabFromThetaCm(double theta_cm=0.0, int part=0) const
void Print(Option_t *opt="") const override
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
void AnalysisAngDist(double ekin, double theta_lab, int part, double &exc, double &theta_cm, double &omega_lab2cm)
double fInvariantMass
Definition TReaction.h:165
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
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
bool fInverse
Definition TReaction.h:158
double fTCm[4]
Definition TReaction.h:174
double fCmP
Definition TReaction.h:170
double fExc
Definition TReaction.h:159
double fCmE
Definition TReaction.h:168
double GetRutherfordLab(double theta_lab, int part=2, bool Units_mb=true) const
void SetCmFrame(double exc)
Definition TReaction.cxx:81
double GetExcEnergy(double ekin=0.00, double theta_lab=0.00, int part=2) const
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
#define D2R
Definition TReaction.h:24
#define R2D
Definition TReaction.h:20
#define PI
Definition TReaction.h:16