KaliVeda
Toolkit for HIC analysis
KVWilckeReactionParameters.cpp
1 //Created by KVClassFactory on Tue Feb 19 17:16:59 2013
2 //Author: John Frankland,,,
3 
4 #include "KVWilckeReactionParameters.h"
5 
6 #include <TF1Derivative.h>
7 
9 
11 Double_t KVWilckeReactionParameters::mu_Wilcke = 931.5; //atomic mass unit
12 /* Numerical constants from W.D. Myers, Phys. Lett. B 30, 451 (1969) */
13 Double_t KVWilckeReactionParameters::a1_Myers = 15.677; //volume energy coefficient
14 Double_t KVWilckeReactionParameters::a2_Myers = 22.0; //surface energy coefficient
15 Double_t KVWilckeReactionParameters::J_Myers = 35.0; //symmetry energy coefficient
16 Double_t KVWilckeReactionParameters::K_Myers = 300.0; //compressibility coefficient
17 Double_t KVWilckeReactionParameters::L_Myers = 99.0; //density-symmetry coefficient
18 Double_t KVWilckeReactionParameters::M_Myers = 4.5; //symmetry anharmonicity coefficient
19 Double_t KVWilckeReactionParameters::Q_Myers = 25.0; //effective surface stiffness
20 Double_t KVWilckeReactionParameters::r0_Myers = 1.16; //nuclear radius constant
21 Double_t KVWilckeReactionParameters::c1_Myers = 0.745; //Coulomb energy coefficient
22 
23 
25 
27 {
28  ZC = ZP + ZT;
29  NP = AP - ZP;
30  NT = AT - ZT;
31  NC = NP + NT;
32  AC = AP + AT;
33  MU = AP * AT / (1.*AC);
35  R0 = r0_Wilcke(AP, AT);
36  RP = SharpRadius(AP);
37  RT = SharpRadius(AT);
40  CBAR = CP * CT / (CP + CT);
43  RCTOTAL = RCP + RCT;
44  V0 = BSS_V0(ZP, AP, ZT, AT);
45  BSS_K = V0 - e2_Wilcke * ZP * ZT / RCTOTAL;
46  BSS_N = e2_Wilcke * ZP * ZT / RCTOTAL / BSS_K;
47  BSS_K /= pow(RCTOTAL, BSS_N);
48 
50  ASYMMFISSIONTKE = FISSIONTKE * 4.*ZP * ZT / (1.*ZC * ZC);
52  // Factor to convert the dimensionless proximity potential function of
53  // J. Blocki et al., Ann. Phys. (N. Y.) 105, 427 (1977)
54  // into a nuclear potential V_N
55  PROXFACTOR = 4 * TMath::Pi() * GAMMA * CBAR;
57  // Maximum critical angular momentum for fusion
58  Double_t S = CT + CP + 0.3;
59  Double_t PHI = -0.96;
60  LCRIT = PROXFACTOR * PHI + e2_Wilcke * ZP * ZT / S / S;
62  LCRIT = 7. / 5.*pow(LCRIT, 0.5);
63 
64  if (!fBSS) fBSS = new TF1("BSS", this, &KVWilckeReactionParameters::VC, 0, 25, 0, "KVWilckeReactionParameters", "VC");
65  VC_RINT = fBSS->Eval(RINT);
66  if (!fCMThetaQuart) fCMThetaQuart = new TF1("ThetaQuart", this, &KVWilckeReactionParameters::QuarterPointAngle, 0, 100, 0, "KVWilckeReactionParameters", "QuarterPointAngle");
67  if (!fLmax) fLmax = new TF1("LMAX", this, &KVWilckeReactionParameters::Lmax, 0, 100, 0, "KVWilckeReactionParameters", "Lmax");
68  if (!fSigmaR) fSigmaR = new TF1("SIGMA-R", this, &KVWilckeReactionParameters::SigmaR, 0, 100, 0, "KVWilckeReactionParameters", "SigmaR");
69  if (!fProx) fProx = new TF1("PROX", this, &KVWilckeReactionParameters::ProxPot, 0, 25, 0, "KVWilckeReactionParameters", "ProxPot");
70  if (!fPotential) fPotential = new TF1("HIPOT", this, &KVWilckeReactionParameters::Potential, 0, 25, 0, "KVWilckeReactionParameters", "Potential");
71  if (!fCentPot) {
72  fCentPot = new TF1("CENTPOT", this, &KVWilckeReactionParameters::CentrifugalPotential, 0, 25, 1, "KVWilckeReactionParameters", "CentrifugalPotential");
73  fCentPot->SetNpx(1000);
74  fCentPot->SetParName(0, "l [hbar]");
75  }
77 
78  if (!fSigmaFus) fSigmaFus = new TF1("SIGMA-FUS", this, &KVWilckeReactionParameters::SigmaFus, 0, 1.e+02, 0, "KVWilckeReactionParameters", "SigmaFus");
79 }
80 
81 
82 
85 
87  : fBSS(0), fProx(0), fPotential(0), fCMThetaQuart(0), fLmax(0), fSigmaR(0), fSigmaFus(0), fCentPot(0)
88 {
89  // Default constructor
90 }
91 
92 
93 
95 
97  : fBSS(0), fProx(0), fPotential(0), fCMThetaQuart(0), fLmax(0), fSigmaR(0), fSigmaFus(0), fCentPot(0)
98 {
99  SetEntranceChannel(proj, targ);
100 }
101 
102 
103 
106 
108 {
109  // Destructor
110 }
111 
112 
113 
116 
118 {
119  // (Re)set entrance channel to calculate
120  ZP = proj.GetZ();
121  AP = proj.GetA();
122  ZT = targ.GetZ();
123  AT = targ.GetA();
124  init();
125 }
126 
127 
128 
129 
132 
134 {
135  // epsilon_bar, Eq.(7) in W.D. Myers, Phys. Lett. B 30, 451 (1969)
136  Double_t eps = (-2.*a2_Myers * pow(A, -THIRD) +
137  L_Myers * pow(delta_bar_Myers(Z, A), 2) +
138  c1_Myers * Z * Z * pow(A, -4. / 3.)) / K_Myers;
139  return eps;
140 }
141 
142 
143 
146 
148 {
149  // delta_bar, Eq.(8) in W.D. Myers, Phys. Lett. B 30, 451 (1969)
150  Double_t I = (A - 2.*Z) / (1.*A);
151  Double_t delta = (I + (3 * c1_Myers / 8. / Q_Myers) * Z * Z * pow(A, -5. / 3.)) /
152  (1. + (9.*J_Myers / 4. / Q_Myers) * pow(A, -THIRD));
153  return delta;
154 }
155 
156 
157 
163 
165 {
166  // S. Cohen, F. Plasil and W.J. Swiatecki, Ann. Phys. (N. Y.) 82, 557 (1974)
167  // Rotating Liquid Drop model critical angular momentum
168  // Only y_I is returned (normally valid for fissility parameter > 0.81)
169  // as no parameterisation for y_II is given in the paper ???
170 
171  Double_t Esurf = (17.9439 / 0.9517) * NLDSurfaceTensionCoefficient(z, a) * pow(a, 2. / 3.);
172  Double_t Ecoul = 0.7053 * z * z / pow(a, THIRD);
173  Double_t x = Ecoul / 2. / Esurf; // fissility parameter
174  Double_t y_I;
175  if (x < 0.75) y_I = 0.2829 - 0.3475 * x - 0.0016 * x * x + 0.0501 * x * x * x;
176  else y_I = (7. / 5.) * pow((1. - x), 2) - 4.5660 * pow(1. - x, 3.) + 6.7443 * pow(1 - x, 4.);
177  return y_I;
178 }
179 
180 
181 
183 
185 {
186  printf("-------------------------------------------\n");
187  printf("PARAMETERS INDEPENDENT OF BOMBARDING ENERGY\n");
188  printf("-------------------------------------------\n\n");
189  KVNucleus fus(ZC, AC);
190  printf("ATOMIC NUMBERS: ZP=%3d ZT=%3d ZC=%3d (%s)\n", ZP, ZT, ZC, fus.GetSymbol());
191  printf("NEUTRON NUMBERS: NP=%3d NT=%3d NC=%3d\n\n", NP, NT, NC);
192  printf("AP**1/3=%7.3f AT**1/3=%7.3f\n", pow(AP, THIRD), pow(AT, THIRD));
193  printf("REDUCED MASS NUMBER=%6.2f AP+AT=AC=%3d\n\n", MU, AC);
194  printf("INTERACTION RADIUS RINT=%5.2f fm R0=%5.2f fm\n\n", RINT, R0);
195  printf("MATTER HALF-DENSITY RADII [fm]:\n");
196  printf("CP=%5.2f CT=%5.2f CT+CP=%5.2f CBAR=%5.2f\n\n", CP, CT, (CP + CT), CBAR);
197  printf("EQUIVALENT SHARP SURFACE RADII [fm]:\n");
198  printf("RP=%5.2f RT=%5.2f\n\n", RP, RT);
199  printf("COULOMB RADII [fm]:\n");
200  printf("RCP=%5.2f RCT=%5.2f RC=RCP+RCT=%5.2f\n\n", RCP, RCT, RCTOTAL);
201  printf("BSS-COULOMB POTENTIAL [MeV]:\n");
202  printf("VC(r)=%5.3f*ZP*ZT/r for r>RC\n", e2_Wilcke);
203  printf("VC(r)=V0-K*r**n for r<RC\n");
204  printf("V0=%7.2f MeV K=%7.5f n=%5.3f\n", V0, BSS_K, BSS_N);
205  printf("VC(RINT)=%6.1f MeV\n\n", VC_RINT);
206  printf("FISSION-TKE=%5.0f MeV\n", FISSIONTKE);
207  printf("ASYMM. FISSION-TKE=%5.0f MeV\n\n", ASYMMFISSIONTKE);
208  printf("LIQUID DROP PARAMETERS:\n");
209  printf("GAMMA=%6.3f MeV/fm**2 PROX-FACTOR=%6.2f MeV\n\n", GAMMA, PROXFACTOR);
210  printf("FUSION RELATED PARAMETERS:\n");
211  printf("R-BARRIER=%5.2f fm V(RB)=%6.1f MeV\n", RBARRIER, VRB);
212  printf("L-CRITICAL=%4.0f HBAR\n", LCRIT);
213 
214  printf("-------------------------------------------------------------------------\n");
215  printf("EL/u ELAB ECM ECM/VC k ETA LMAX SGMAR SGFUS QP-CM QP-LP EP-QP\n");
216  printf("-------------------------------------------------------------------------\n");
217  Double_t e = 1;
218  Int_t nlines = 0;
219  while (e < 51) {
220  printf(" %4.1f %4.0f %4.0f %5.2f %4.1f %5.1f %4.0f %4.0f %4.0f %5.1f %5.1f %4.0f\n",
221  e, AP * e, ECM(e), ECM(e) / VC_RINT, k(e), Eta(e), Lmax(&e, 0),
222  SigmaR(&e, 0), SigmaFus(&e, 0), TMath::RadToDeg()*QuarterPointAngle(&e, 0),
224  nlines++;
225  if (e < 4) e += 1;
226  else if (e < 12) e += 0.5;
227  else if (e < 20) e += 1;
228  else e += 5;
229  if (!(nlines % 5)) printf("\n");
230  }
231  printf("************************************************************************\n");
232 }
233 
234 
235 
239 
241 {
242  // Find position (radial distance between centres) at which total potential has a
243  // minimum for the given angular momentum. Returns -1.0 if no pocket.
244 
247 
248  // look for extrema
249  Double_t rmax = 20.;
250  Double_t rmin = -1.0;
251  for (int i = 2; i; --i) {
252  rmin = df.GetX(0., 1., rmax);
253  if (rmin == rmax) return -1.0; // no extrema
254  // Test min/max ?
255  if (df2.Eval(rmin) > 0) return rmin; // found minimum
256  // That was a maximum. Try again.
257  rmax = rmin - 0.1;
258  }
259  return -1.0;
260 }
261 
262 
263 
283 
285 {
286  // Retuns maximum angular momentum for which a pocket exists in the interaction potential
287  // which corresponds to a lower potential energy than the maximum
288 
289 // Double_t lmin, lmax;
290 // lmin = 0;
291 // lmax = 200;
292 // Double_t l = lmin;
293 // while (lmax > lmin) {
294 // if (PotentialPocketRadius(l) > 0.) {
295 // if (lmax - lmin == 2.0) return l;
296 // lmin = l;
297 // } else {
298 // if (l < 1) return 0;
299 // lmax = l;
300 // }
301 // if (lmax - lmin == 1.0) return lmin;
302 // l = TMath::Nint((lmin + lmax) / 2.);
303 // }
304 // return 0;
305  Double_t l_last_pock = -1;
306  Double_t l = 0;
307  while (1) {
309  if (Rmax < 0) break;
311  if (Rmin < 0) break;
312  if (GetCentrifugalPotential(l)->Eval(Rmin) >= GetCentrifugalPotential(l)->Eval(Rmax)) break;
313  else l_last_pock = l;
314  ++l;
315  }
316  return l_last_pock;
317 }
318 
319 
320 
324 
326 {
327  // Find position (radial distance between centres) at which total potential has a
328  // maximum for the given angular momentum. Returns -1.0 if no maximum.
329 
332 
333  // look for extrema
334  Double_t rmax = 20.;
335  Double_t rmin = 5.0;
336  for (int i = 2; i; --i) {
337  rmin = df.GetX(0., rmin, rmax);
338  if (rmin == rmax) return -1.0; // no extrema
339  // Test min/max ?
340  if (df2.Eval(rmin) < 0) return rmin; // found minimum
341  // That was a minimum. Try again.
342  rmin += 0.1;
343  }
344  return -1.0;
345 }
346 
347 
348 
351 
353 {
354  // Bass reaction cross-section [mb] for incident energy [MeV/nucleon]
355  Double_t R1 = 1.12 * pow(AP, 1. / 3.) - 0.94 * pow(AP, -1. / 3.);
356  Double_t R2 = 1.12 * pow(AT, 1. / 3.) - 0.94 * pow(AT, -1. / 3.);
357  Double_t Rint = R1 + R2 + 3.2;
358  Double_t VRINT = 1.44 * ZP * ZT / Rint - R1 * R2 / (R1 + R2);
359  return 10 * TMath::Pi() * pow(Rint, 2) * (1. - VRINT / ECM(e_sur_a));
360 }
361 
362 
int Int_t
#define e(i)
double Double_t
Description of properties and kinematics of atomic nuclei.
Definition: KVNucleus.h:126
static Double_t hbar
hbar*c in MeV.fm
Definition: KVNucleus.h:173
const Char_t * GetSymbol(Option_t *opt="") const
Definition: KVNucleus.cpp:81
Int_t GetA() const
Definition: KVNucleus.cpp:802
Int_t GetZ() const
Return the number of proton / atomic number.
Definition: KVNucleus.cpp:773
Reaction parameters for heavy-ion collisions from systematics of Wilcke et al.
Double_t LCRIT
The maximum critical angular momentum for fusion.
Double_t QuarterPointAngle(Double_t *x, Double_t *) const
static Double_t SWaveFusionBarrierRadius(Int_t zp, Int_t ap, Int_t zt, Int_t at)
static Double_t L_Myers
density-symmetry coefficient
static Double_t InteractionRadius(Int_t aproj, Int_t atarg)
Double_t SigmaFus(Double_t *e_sur_a, Double_t *) const
Double_t Potential(Double_t *r, Double_t *)
Double_t Eta(Double_t e_sur_a) const
static Double_t e2_Wilcke
e**2 = 1.438 is value used by Wilcke et al.
Double_t ECM(Double_t e_sur_a) const
TF1 * fLmax
Grazing angular momentum.
static Double_t RLDCriticalAngularMomentum(Int_t z, Int_t a)
virtual ~KVWilckeReactionParameters()
Destructor.
static Double_t MatterHalfDensityRadius(Int_t A)
TF1 * fSigmaR
Reaction cross section.
Int_t NT
Neutron number of the projectile, target.
Double_t GetBassReactionCrossSection(Double_t e_sur_a)
Bass reaction cross-section [mb] for incident energy [MeV/nucleon].
TF1 * fPotential
total (nuclear+coulomb) potential for heavy-ions
Double_t ProjectileLabQP(Double_t e) const
static Double_t epsilon_bar_Myers(Int_t Z, Int_t A)
epsilon_bar, Eq.(7) in W.D. Myers, Phys. Lett. B 30, 451 (1969)
static Double_t M_Myers
symmetry anharmonicity coefficient
static Double_t NLDSurfaceTensionCoefficient(Int_t Z, Int_t A)
static Double_t SharpRadius(Int_t A)
Double_t FISSIONTKE
TKE for symmetric fission of combined system.
TF1 * fSigmaFus
Fusion cross section.
Double_t Lmax(Double_t *x, Double_t *) const
Double_t V0
BSS potential at r=0.
Double_t RBARRIER
Fusion barrier radius RB for s-waves.
Int_t AT
Mass number of the projectile, target.
TF1 * fCMThetaQuart
CM quarter point angle.
static Double_t r0_Myers
nuclear radius constant
Int_t ZT
Atomic number of the projectile, target.
Double_t GAMMA
Nuclear liquid drop surface-tension coefficient.
Double_t k(Double_t e_sur_a) const
Double_t VC_RINT
BSS Coulomb potential at Rint.
static Double_t J_Myers
symmetry energy coefficient
static Double_t r0_Wilcke(Int_t aproj, Int_t atarg)
Double_t SigmaR(Double_t *x, Double_t *) const
static Double_t a2_Myers
surface energy coefficient
static Double_t mu_Wilcke
mu = 931.5 is value used by Wilcke et al.
static Double_t Q_Myers
effective surface stiffness
Double_t ProjectileLabEQP(Double_t e) const
TF1 * GetCentrifugalPotential(Double_t e_sur_a, Double_t b) const
Double_t ASYMMFISSIONTKE
TKE of completely relaxed events in strongly damped collisions.
Double_t CentrifugalPotential(Double_t *x, Double_t *l)
static Double_t ChargeRadius_Myers(Int_t Z, Int_t A)
static Double_t BSS_V0(Int_t zp, Int_t ap, Int_t zt, Int_t at)
static Double_t a1_Myers
volume energy coefficient
void SetEntranceChannel(const KVNucleus &proj, const KVNucleus &targ)
(Re)set entrance channel to calculate
static Double_t delta_bar_Myers(Int_t Z, Int_t A)
delta_bar, Eq.(8) in W.D. Myers, Phys. Lett. B 30, 451 (1969)
static Double_t TKESymFiss(Int_t Z, Int_t A)
Double_t PROXFACTOR
Proximity potential factor.
Double_t ProxPot(Double_t *r, Double_t *)
Double_t VC(Double_t *r, Double_t *)
static Double_t K_Myers
compressibility coefficient
Double_t CT
Matter half-density radii.
static Double_t c1_Myers
Coulomb energy coefficient.
TF1 * fProx
Nuclear proximity potential for heavy-ions.
TF1 * fBSS
BSS Coulomb potential for heavy-ions.
Double_t VRB
The total conservative potential at r=RB for s-waves.
KVWilckeReactionParameters()
Default constructor.
Numerical derivative of a TF1.
Definition: TF1Derivative.h:16
virtual void SetNpx(Int_t npx=100)
virtual Double_t GetX(Double_t y, Double_t xmin=0, Double_t xmax=0, Double_t epsilon=1.E-10, Int_t maxiter=100, Bool_t logx=false) const
virtual void SetParName(Int_t ipar, const char *name)
virtual Double_t Eval(Double_t x, Double_t y=0, Double_t z=0, Double_t t=0) const
RVec< PromoteTypes< T0, T1 > > pow(const T0 &x, const RVec< T1 > &v)
Double_t x[n]
#define I(x, y, z)
RooArgSet S(Args_t &&... args)
Eval
constexpr Double_t Pi()
constexpr Double_t RadToDeg()
#define R1(v, w, x, y, z, i)
#define R2(v, w, x, y, z, i)
TLine l
TArc a
ClassImp(TPyArg)