36 for (
int i = 0; i < 5; i++) {
43 fThetaLabVsThetaCM[i] =
nullptr;
44 fELabVsThetaCM[i] =
nullptr;
45 fELabVsThetaLab[i] =
nullptr;
47 fKoxReactionXSec =
nullptr;
48 fEqbmChargeState =
nullptr;
49 fEqbmChargeStateShSol =
nullptr;
50 fEqbmChargeStateShGas =
nullptr;
51 fSetOutgoing = kFALSE;
53 SetIntegralPrecision(1e-10);
110 syst.ReplaceAll(
" ",
"");
113 if (syst.Contains(
"@")) {
115 scouple = syst.
Next();
119 if (sener.Contains(
"MeV/A")) {
120 sener.ReplaceAll(
"MeV/A",
"");
123 else if (sener.Contains(
"MeV/u")) {
124 sener.ReplaceAll(
"MeV/u",
"");
133 if (nnuc1.IsZombie()) MakeZombie();
137 if (nnuc2.IsZombie()) MakeZombie();
143 SetTitle(Form(
"%s + %s %.1f MeV/A",
fNuclei[1].GetSymbol(),
fNuclei[2].GetSymbol(), ee));
156 Obsolete(
"KV2Body(KVNucleus*, KVNucleus*, KVNucleus*, Double_t)",
"1.11",
"1.12");
159 Info(
"KV2Body",
"Use constructor KV2Body(const KVNucleus& compound) to simulate binary decay of compound");
161 else if (!proj_out) {
162 Info(
"KV2Body",
"Use constructor KV2Body(const KVNucleus& proj, const KVNucleus& targ) to define entrance channel");
165 Info(
"KV2Body",
"Use constructor KV2Body(const KVNucleus& proj, const KVNucleus& targ, const KVNucleus& outgoing) to define entrance & exit channels");
197 : fNuclei(5), fEDiss(-Exx)
243 fNuclei(5), fEDiss(Ediss)
260 SetTitle(Form(
"%s + %s %.1f MeV/A",
fNuclei[1].GetSymbol(),
fNuclei[2].GetSymbol(),
fNuclei[1].GetAMeV()));
283 fNuclei(5), fEDiss(Ediss)
306 SetTitle(Form(
"%s + %s %.1f MeV/A",
fNuclei[1].GetSymbol(),
fNuclei[2].GetSymbol(),
fNuclei[1].GetAMeV()));
366 for (
int i = 0; i < 5; i++) {
385 Double_t p = TMath::Sqrt(E * E - mass * mass);
403 SetTitle(Form(
"%s + %s %.1f MeV/A",
fNuclei[1].GetSymbol(),
fNuclei[2].GetSymbol(),
fNuclei[1].GetAMeV()));
408 for (
unsigned i = 1; i <= 4; ++i)
if (
fNuclei[i].IsDefined())
fNuclei[i].SetExcitEnergy(0);
466 if (i > 0 && i < 5) {
470 Warning(
"GetNucleus(Int_t i)",
"Index i out of bounds, i=%d", i);
485 Warning(
"GetQReaction",
"Parameters for outgoing nuclei not set");
504 Warning(
"GetQGroundStates",
505 "Parameters for outgoing nuclei not set");
540 if (i < 3 || i > 4) {
541 Warning(
"GetMaxAngleLab(Int_t i)",
542 "Returns maximum scattering angle in lab for nuclei i=3 (quasiproj) and i=4 (quasitarget)");
558 if (i < 3 || i > 4) {
559 Warning(
"GetMinAngleLab(Int_t i)",
560 "Returns minimum scattering angle in lab for nuclei i=3 (quasiproj) and i=4 (quasitarget)");
609 if (i < 1 || i > 2) {
610 Warning(
"GetLabGrazingAngle(Int_t i)",
611 "i should be 1 (proj) or 2 (targ)");
615 Warning(
"GetLabGrazingAngle(Int_t i)",
616 "No target defined for reaction");
619 Double_t RP = 1.28 * TMath::Power(
GetNucleus(1)->GetA(), 1. / 3.)
620 - 0.76 + 0.8 * TMath::Power(
GetNucleus(1)->GetA(), -1. / 3.);
621 Double_t RT = 1.28 * TMath::Power(
GetNucleus(2)->GetA(), 1. / 3.)
622 - 0.76 + 0.8 * TMath::Power(
GetNucleus(2)->GetA(), -1. / 3.);
623 Double_t CP = RP * (1. - TMath::Power(1. / RP, 2.));
624 Double_t CT = RT * (1. - TMath::Power(1. / RT, 2.));
625 Double_t RINT = CP + CT + 4.49 - (CT + CP) / 6.35;
643 Double_t GR = 2. * TMath::ASin(U);
644 Double_t GPART = TMath::Pi() - GR;
646 TMath::Sin(GR) / (TMath::Cos(GR) +
649 Double_t GRAZ = TMath::ATan(ARAZ);
655 return (GRAZ * TMath::RadToDeg());
658 return (GPART * TMath::RadToDeg() / 2.);
693 Error(
"CalculateKinematics",
"Set outgoing (decay) product first");
702 for (
int i = 0; i < 5; i++) {
710 VCM.SetXYZ(0., 0., 0.);
755 if (
VC[3] > 0.)
K[3] =
VCM.Mag() /
VC[3];
758 if (TMath::AreEqualAbs(
K[3], 1., 1.e-16))
759 T3MAX = TMath::PiOver2();
763 T3MAX = TMath::ATan((1. / TMath::Sqrt(
K[3] *
K[3] - 1.)) /
GetCMGamma());
764 TETAMAX[3] = T3MAX * TMath::RadToDeg();
773 if (
VC[4] > 0.)
K[4] =
VCM.Mag() /
VC[4];
776 if (TMath::AreEqualAbs(
K[4], 1., 1.e-16))
777 T4MAX = TMath::PiOver2();
781 T4MAX = TMath::ATan((1. / TMath::Sqrt(
K[4] *
K[4] - 1.)) /
GetCMGamma());
783 T4MAX = TMath::PiOver2();
784 TETAMAX[4] = T4MAX * TMath::RadToDeg();
829 cout <<
" ******" << endl;
838 cout <<
" ==> Q-REACTION = " <<
GetQReaction() <<
" MEV";
841 cout <<
" AVAILABLE ENERGY IN C.M. : ECM = " <<
GetCMEnergy() <<
844 GetA() : 0.)) <<
" MEV/A)" << endl;
845 cout <<
" PROJECTILE VELOCITY IN LAB " <<
GetNucleus(1)->
GetV().Mag()
846 <<
" CM/NS ( " <<
GetNucleus(1)->Beta() <<
" * C )" << endl;
848 Mag() <<
" CM/NS" << endl;
850 for (
int i = 1; i <= 4; i++) {
852 cout <<
" ENERGY - VELOCITY OF NUCLEUS " << i <<
" IN CM : " <<
857 cout <<
" MAXIMUM SCATTERING ANGLE IN LABORATORY" << endl;
865 cout <<
" GRAZING ANGLE IN LABORATORY : PROJECTILE " <<
867 cout <<
" GRAZING ANGLE IN LABORATORY : TARGET " <<
874 if (!strcmp(opt,
"lab") || !strcmp(opt,
"LAB")
875 || !strcmp(opt,
"Lab")) {
878 " LABORATORY ANGULAR DISTRIBUTION" << endl
881 " TETA 3 EN.DIF.3 V3 TETA 4 EN.DIF.4 TETA 3"
884 " (LAB) (LAB) (LAB) (LAB) (LAB) (C.M)"
887 for (
int step = 0; step <= nsteps; step++) {
888 Double_t theta = dtheta * step;
890 Double_t elabP1, elabP2;
891 Double_t tcmP1, tcmP2;
894 Double_t elabT1, elabT2;
895 GetELab(3, theta, 3, elabP1, elabP2);
897 GetVLab(3, theta, 3, vP1, vP2);
899 GetELab(4, theta, 3, elabT1, elabT2);
901 (
" %6.2f %7.2f/%7.2f %5.2f/%5.2f %6.2f/%6.2f %7.2f/%7.2f %6.2f/%6.2f\n",
902 theta, elabP1, elabP2, vP1, vP2,
903 tlT1, tlT2, elabT1, elabT2,
907 if (
GetNucleus(2) && (!strcmp(opt,
"ruth") || !strcmp(opt,
"RUTH")
908 || !strcmp(opt,
"Ruth"))) {
911 " RUTHERFORD LABORATORY ANGULAR DISTRIBUTION" <<
914 " TETA 3 TETA 3 EN.DIF.3 SIG.RUT. SIG.RUT."
917 " (LAB) (CM) (LAB) (LAB) (b/sr) (CM) (b/sr)"
921 for (
int step = 0; step < nsteps; step++) {
922 Double_t theta = dtheta * (step + 1);
923 Double_t elabP1, elabP2;
925 GetELab(3, theta, 3, elabP1, elabP2);
928 (
" %6.2f %6.2f %7.2f %10.4g %10.4g\n",
945 Int_t
KV2Body::GetELab(Int_t OfNucleus, Double_t ThetaLab, Int_t AngleNucleus, Double_t& e1, Double_t& e2)
const
954 if (ThetaLab >
TETAMAX[AngleNucleus])
return 0;
957 Int_t nsol =
GetThetaCM(ThetaLab, AngleNucleus, TCM1, TCM2);
958 TCM1 =
GetThetaCM(OfNucleus, TCM1, AngleNucleus);
959 if (nsol > 1) TCM2 =
GetThetaCM(OfNucleus, TCM2, AngleNucleus);
962 if (nsol == 2) e2 =
GetELab(TCM2, OfNucleus);
980 Int_t nsol =
FindRoots(TLF, 0., 180., ThetaLab, t1, t2);
994 Int_t nsol =
GetELab((Int_t)par[0], x[0], (Int_t)par[0], e1, e2);
1010 Int_t OfNucleus = (Int_t)par[0];
1011 Double_t TCM = x[0] * TMath::DegToRad();
1030 Double_t ThetaCM = x[0] * TMath::DegToRad();
1031 Int_t OfNucleus = (Int_t)par[0];
1032 Double_t ThetaL = TMath::ATan2(TMath::Sin(ThetaCM) + 1e-10, (
K[OfNucleus] + TMath::Cos(ThetaCM)) *
GetCMGamma()) * TMath::RadToDeg();
1034 if (ThetaL < 0.) ThetaL += 180.;
1050 fELabVsThetaCM[OfNucleus] =
new TF1(Form(
"KV2Body:ELabVsThetaCM:%d", OfNucleus),
1070 fELabVsThetaLab[OfNucleus] =
new TF1(Form(
"KV2Body:ELabVsThetaLab:%d", OfNucleus),
1088 Int_t
KV2Body::GetVLab(Int_t OfNucleus, Double_t ThetaLab, Int_t AngleNucleus, Double_t& v1, Double_t& v2)
const
1098 Int_t nsol =
GetELab(OfNucleus, ThetaLab, AngleNucleus, e1, e2);
1099 if (!nsol)
return nsol;
1117 Int_t
KV2Body::GetThetaLab(Int_t OfNucleus, Double_t ThetaLab, Int_t AngleNucleus, Double_t& t1, Double_t& t2)
const
1126 if (ThetaLab >
TETAMAX[AngleNucleus])
return 0;
1127 if (!(TMath::Abs(OfNucleus - AngleNucleus) % 2)) {
1133 Double_t TCM1, TCM2;
1134 Int_t nsol =
GetThetaCM(ThetaLab, AngleNucleus, TCM1, TCM2);
1135 TCM1 =
GetThetaCM(OfNucleus, TCM1, AngleNucleus);
1136 if (nsol > 1) TCM2 =
GetThetaCM(OfNucleus, TCM2, AngleNucleus);
1155 fThetaLabVsThetaCM[OfNucleus] =
new TF1(Form(
"KV2Body:ThetaLabVsThetaCM:%d", OfNucleus),
1180 Double_t par = OfNucleus;
1202 Warning(
"GetXSecRuthCM",
"No target defined for reaction");
1210 if (TCM < 0)
return 0;
1211 Double_t D = 1. / (16. * (TMath::Power(TMath::Sin(TCM * TMath::DegToRad() / 2.), 4.)));
1213 return ((TMath::Power(PB, 2.)) * D / 100.);
1231 if (theta >
TETAMAX[OtherNucleus])
return -1;
1233 Int_t nsol =
GetThetaCM(theta, OtherNucleus, t1, t2);
1234 if (!nsol)
return -1.;
1235 Double_t Pt1 =
GetThetaCM(OfNucleus, t1, OtherNucleus);
1236 Double_t Pt2 =
GetThetaCM(OfNucleus, t2, OtherNucleus);
1239 Pt = TMath::Min(Pt1, Pt2);
1257 Warning(
"XSecRuthCMVsThetaCM",
"No target defined for reaction");
1263 Double_t TCM = x[0] * TMath::DegToRad();
1264 Int_t OfNucleus = (Int_t)par[0];
1265 if (OfNucleus == 2 || OfNucleus == 4) TCM = TMath::Pi() - TCM;
1266 Double_t D = 1. / (16. * (TMath::Power(TMath::Sin(TCM / 2.), 4.)));
1267 return ((TMath::Power(PB, 2.)) * D / 100.);
1287 Double_t par = OfNucleus;
1312 Double_t T3L =
GetThetaLab(T3CM, par[0]) * TMath::DegToRad();
1313 T3CM *= TMath::DegToRad();
1314 Double_t RLC = (TMath::Power(TMath::Sin(T3CM), 3.)) /
1315 ((TMath::Power(TMath::Sin(T3L), 3.)) *
GetCMGamma() *
1316 (1. +
K[(int)par[0]] * TMath::Cos(T3CM)));
1317 if (DSIDTB * RLC < 0) {
1318 Warning(
"GetXSecRuthLab",
"negative value for choosen parameters : %lf %d\n", x[0], (Int_t)par[0]);
1339 return XSecRuthLab(x, par) * TMath::Sin(x[0] * TMath::DegToRad());
1399 if (th2 < th1)
return 0;
1402 if (phi1 == -1 || phi2 == -1) dphi = 2 * TMath::Pi();
1403 else if (phi1 == phi2) dphi = 1;
1406 dphi *= TMath::DegToRad();
1408 Double_t theta_min = 1.;
1409 Double_t theta_max = 179.;
1410 if (th1 < theta_min) theta_min = th1;
1411 if (th2 > theta_max) theta_max = th2;
1413 #if ROOT_VERSION_CODE > ROOT_VERSION(5,99,01)
1416 const Double_t* para = 0;
1434 TString name =
"RuthXSecCM: ";
1440 name += Form(
"%6.1f AMeV ", elab);
1441 if (OfNucleus == 3) name +=
"(projectile)";
1442 else name +=
"(target)";
1443 TF1* fXSecRuthCM = (TF1*)gROOT->GetListOfFunctions()->FindObject(name.Data());
1445 fXSecRuthCM =
new TF1(name.Data(),
1447 theta_cm_min, theta_cm_max, 1,
"KV2Body",
"XSecRuthCMVsThetaCM");
1448 fXSecRuthCM->SetParameter(0, OfNucleus);
1449 fXSecRuthCM->SetNpx(1000);
1451 fXSecRuthCM->SetRange(theta_cm_min, theta_cm_max);
1469 const Double_t r0 = 1.07;
1470 const Double_t e2 = 1.44;
1471 Double_t A1third = pow(
GetNucleus(1)->GetA(), 1. / 3.);
1472 Double_t A2third = pow(
GetNucleus(2)->GetA(), 1. / 3.);
1473 Double_t R12 = r0 * (A1third + A2third);
1476 - 2.9 * A1third * A2third / (A1third + A2third);
1497 const Double_t r0 = 1.05;
1498 const Double_t rc = 1.3;
1499 const Double_t e2 = 1.44;
1500 const Double_t a = 1.9;
1508 Double_t c = 11.3315 * TMath::Landau(TMath::Log10(eproj[0]), 2.47498, 0.554937);
1509 Double_t A1third = pow(proj->
GetA(), 1. / 3.);
1510 Double_t A2third = pow(targ->
GetA(), 1. / 3.);
1511 Double_t BC = proj->
GetZ() * targ->
GetZ() * e2 / (rc * (A1third + A2third));
1512 Double_t EFac = 1 - BC / ECM;
1514 Double_t Xsec = TMath::Pi() * pow(r0, 2) *
1515 pow((A1third + A2third + a * A1third * A2third / (A1third + A2third) - c), 2) *
1533 Double_t A1third = pow(
GetNucleus(1)->GetA(), 1. / 3.);
1534 Double_t A2third = pow(
GetNucleus(2)->GetA(), 1. / 3.);
1536 Double_t Xsec = TMath::Pi() * pow(r0, 2) *
1537 pow(A1third + A2third, 2);
1555 return 10 * TMath::Sqrt(ReacXsec / TMath::Pi());
1574 return TMath::Pi() * (TMath::Power(b2, 2.) - TMath::Power(b1, 2.)) / 100;
1664 Double_t Zp = proj->
GetZ();
1666 Double_t beta = proj->Beta();
1670 Double_t q = Zp * (1. - TMath::Exp(-83.275 * beta / pow(Zp, 0.477)));
1673 Double_t g = 0.929 + 0.269 * TMath::Exp(-0.16 * Zt) + (0.022 - 0.249 * TMath::Exp(-0.322 * Zt)) * vp / pow(Zp, 0.477);
1678 Double_t f = 1. - TMath::Exp(-12.905 + 0.2124 * Zp - 0.00122 * pow(Zp, 2));
1733 Double_t Zp = proj->
GetZ();
1736 Double_t Zt = targ->
GetZ();
1737 Double_t Vp = proj->Beta();
1738 const Double_t V0 = 2.19e+06 / TMath::C();
1739 Double_t x = -0.019 * pow(Zp, -0.52) * Vp / V0;
1740 x = Vp / V0 * pow(Zp, -0.52) * pow(Zt, x) / 1.68;
1741 x = pow(x, 1. + 1.8 / Zp);
1742 Double_t q = 0.07 / x + 6. + 0.3 * pow(x, 0.5) + 10.37 * x + pow(x, 4);
1743 q = (Zp * (12 * x + pow(x, 4.))) / q;
1773 name +=
" SHIWIETZ-SOLID";
1793 Double_t Zp = proj->
GetZ();
1796 Double_t Zt = targ->
GetZ();
1797 Double_t Vp = proj->Beta();
1798 const Double_t V0 = 2.19e+06 / TMath::C();
1799 Double_t x = 0.03 - 0.017 * pow(Zp, -0.52) * Vp / V0;
1800 x = Vp / V0 * pow(Zp, -0.52) * pow(Zt, x);
1801 x = pow(x, 1. + 0.4 / Zp);
1802 Double_t q = 1428. - 1206.*pow(x, 0.5) + 690 * x + pow(x, 6.);
1803 q = (Zp * (376.*x + pow(x, 6.))) / q;
1834 name +=
" SHIWIETZ-GAS";
1858 Info(
"GetXSecRuthLabFunc",
"Maximum angle set to %lf degrees", theta_max);
1861 TString name =
"RuthXSec: ";
1867 name += Form(
"%6.1f AMeV ", elab);
1868 if (OfNucleus == 3) name +=
"(projectile)";
1869 else name +=
"(target)";
1870 TF1*
fXSecRuthLab = (TF1*)gROOT->GetListOfFunctions()->FindObject(name.Data());
1910 TString name =
"RuthXSecInt: ";
1916 name += Form(
"%6.1f AMeV ", elab);
1917 if (OfNucleus == 3) name +=
"(projectile)";
1918 else name +=
"(target)";
1920 TF1*
fXSecRuthLab = (TF1*)gROOT->GetListOfFunctions()->FindObject(name.Data());
1944 Int_t
KV2Body::FindRoots(TF1* fonc, Double_t xmin, Double_t xmax, Double_t val, Double_t& x1, Double_t& x2)
const
1955 x1 = x2 = xmin - 1.;
1956 Double_t max = fonc->GetMaximum(xmin, xmax);
1958 if (val > max)
return nRoots;
1959 Double_t maxX = fonc->GetMaximumX(xmin, xmax);
1961 if (TMath::AreEqualAbs(val, max, 1.e-10)) {
1967 if ((maxX < xmax) && (maxX > xmin) && (fonc->Eval(xmin) < val) && (fonc->Eval(xmax) < val)) nRoots = 2;
1968 Double_t xmax1 = (nRoots == 2 ? maxX : xmax);
1969 x1 = fonc->GetX(val, xmin, xmax1);
1970 if (nRoots == 1)
return nRoots;
1972 x2 = fonc->GetX(val, maxX, xmax);
Relativistic binary kinematics calculator.
Double_t XSecRuthLabInt(Double_t *, Double_t *)
Double_t TETAMIN[5]
defined only for nuclei 3 et 4
Double_t eqbm_charge_state_shiwietz_gas(Double_t *t, Double_t *)
Double_t VC[5]
cm velocities
void SetTarget(const KVNucleus &)
Set target for reaction.
TF1 * GetEqbmChargeStateFunc() const
TF1 * GetXSecRuthLabFunc(Int_t OfNucleus=3, Double_t theta_min=1., Double_t theta_max=179.) const
Double_t K[5]
ratio of c.m. velocity to velocity of nucleus in c.m. v_cm/v_i_cm
Double_t GetMaxAngleLab(Int_t i) const
std::vector< KVNucleus > fNuclei
nuclei involved in calculation
void Print(Option_t *opt="") const
TF1 * GetXSecRuthLabIntegralFunc(Int_t OfNucleus=3, Double_t theta_min=1., Double_t theta_max=179.) const
Double_t WCT
total cm energy
Double_t GetSphereDureReactionXSec(Double_t r0=1.05)
Double_t ELabVsThetaLab(Double_t *, Double_t *)
Function calculating lab energy of nucleus par[0] for any lab angle x[0].
TF1 * fThetaLabVsThetaCM[5]
Double_t GetMinAngleLab(Int_t i) const
Double_t GetThetaLab(Double_t ThetaCM, Int_t OfNucleus) const
void SetProjectile(const KVNucleus &)
Set projectile for reaction.
Double_t fIntPrec
Precision of the TF1::Integral method.
Double_t GetXSecRuthLab(Double_t ThetaLab_Proj, Int_t OfNucleus=3) const
Double_t GetQReaction() const
Calculate Q-value for reaction, including dissipated (excitation) energy.
TF1 * fEqbmChargeStateShSol
function equilibrium charge state of projectile vs. E/A projectile (Shiwietz et al solid)
Double_t EqbmChargeState(Double_t *t, Double_t *)
TF1 * GetThetaLabVsThetaCMFunc(Int_t OfNucleus) const
void SetOutgoing(const KVNucleus &proj_out)
TF1 * fKoxReactionXSec
function Kox reaction cross-section [barns] vs. E/A projectile
TF1 * GetELabVsThetaCMFunc(Int_t OfNucleus) const
Double_t XSecRuthCM(Double_t *, Double_t *)
TF1 * GetShiwietzEqbmChargeStateFuncForGasTargets() const
Double_t XSecRuthLab(Double_t *, Double_t *)
TF1 * GetKoxReactionXSecFunc() const
Double_t BCM
beta of centre of mass
TF1 * fEqbmChargeStateShGas
function equilibrium charge state of projectile vs. E/A projectile (Shiwietz et al gas)
Double_t GetELab(Double_t ThetaCM, Int_t OfNucleus) const
Int_t FindRoots(TF1 *, Double_t, Double_t, Double_t, Double_t &, Double_t &) const
Double_t GetEDiss() const
static Double_t GetVelocity(Double_t mass, Double_t E)
Double_t GetIntegratedXSecRuthLab(Float_t th1, Float_t th2, Float_t phi1=-1, Float_t phi2=-1, Int_t OfNucleus=3)
Double_t GetBmaxFromReactionXSec(Double_t ReacXsec)
TVector3 VCM
velocity of centre of mass
Double_t KoxReactionXSec(Double_t *, Double_t *)
Double_t ELabVsThetaCM(Double_t *, Double_t *)
Function calculating lab energy of nucleus par[0] for any CM angle x[0].
Double_t ThetaLabVsThetaCM(Double_t *, Double_t *)
Double_t GetCMGamma() const
Double_t TETAMAX[5]
defined only for nuclei 3 et 4
Bool_t fSetOutgoing
= kTRUE if SetOutgoing is called before CalculateKinematics
Double_t WC[5]
cm energy of each nucleus
Double_t EC[5]
cm energies
Double_t GetExcitEnergy() const
Double_t GetMinThetaCMFromThetaLab(Int_t OfNucleus, Double_t theta, Int_t OtherNucleus) const
Double_t WLT
total lab energy
Int_t GetVLab(Int_t OfNucleus, Double_t ThetaLab, Int_t AngleNucleus, Double_t &e1, Double_t &e2) const
Double_t BassIntBarrier()
TVector3 GetCMVelocity() const
Return vector velocity of centre of mass of reaction (units: cm/ns)
void CalculateKinematics()
void init()
Default initialisations.
Double_t GetIntegratedXsec(Double_t b1, Double_t b2)
TF1 * GetXSecRuthCMFunc(Int_t OfNucleus=3, Double_t theta_cm_min=1., Double_t theta_cm_max=179.) const
Double_t GetLabGrazingAngle(Int_t i=1) const
Double_t GetQGroundStates() const
Calculate Q-value for reaction, assuming all nuclei in ground state.
KVNucleus * GetNucleus(Int_t i) const
Double_t XSecRuthCMVsThetaCM(Double_t *, Double_t *)
Double_t GetCMEnergy() const
Return available kinetic energy in centre of mass.
TF1 * fEqbmChargeState
function equilibrium charge state of projectile vs. E/A projectile (Leon et al)
Double_t eqbm_charge_state_shiwietz_solid(Double_t *t, Double_t *)
Double_t GetXSecRuthCM(Double_t ThetaLab_Proj, Int_t OfNucleus=3) const
Double_t fEDiss
dissipated energy, 0 means elastic scattering
Int_t GetThetaCM(Double_t ThetaLab, Int_t OfNucleus, Double_t &t1, Double_t &t2) const
TF1 * GetShiwietzEqbmChargeStateFuncForSolidTargets() const
TF1 * GetELabVsThetaLabFunc(Int_t OfNucleus) const
Description of properties and kinematics of atomic nuclei.
const Char_t * GetSymbol(Option_t *opt="") const
Double_t GetExcitEnergy() const
Double_t GetMassExcess(Int_t z=-1, Int_t a=-1) const
Int_t GetZ() const
Return the number of proton / atomic number.
TVector3 GetMomentum() const
Double_t GetEnergy() const
void SetFrame(const Char_t *frame, const KVFrameTransform &)
KVParticle const * GetFrame(const Char_t *frame, Bool_t warn_and_return_null_if_unknown=kTRUE) const
void SetEnergy(Double_t e)
TVector3 GetVelocity() const
returns velocity vector in cm/ns units
Extension of ROOT TString class which allows backwards compatibility with ROOT v3....
void Begin(TString delim) const
KVString Next(Bool_t strip_whitespace=kFALSE) const