KaliVeda
Toolkit for HIC analysis
Loading...
Searching...
No Matches
KVWilckeReactionParameters.h
1
3
4#ifndef __KVWILCKEREACTIONPARAMETERS_H
5#define __KVWILCKEREACTIONPARAMETERS_H
6
7#define THIRD (1./3.)
8#include "KVNucleus.h"
9#include "TF1.h"
10#include "TH1.h"
11#include "TPad.h"
12#include <TF2.h>
13
22 /* The following variables are exactly copied from pp. 395-400 of Wilcke et al. */
29
49
50 /* Numerical constants from W.D. Myers, Phys. Lett. B 30, 451 (1969) */
60
69
70 void init();
71
72public:
74 KVWilckeReactionParameters(const KVNucleus& proj, const KVNucleus& targ);
76
77 void SetEntranceChannel(const KVNucleus& proj, const KVNucleus& targ);
78
80 {
82 return MatterHalfDensityRadius(atarg) + MatterHalfDensityRadius(aproj) + 4.49 -
83 (MatterHalfDensityRadius(atarg) + MatterHalfDensityRadius(aproj)) / 6.35;
84 }
85 static Double_t r0_Wilcke(Int_t aproj, Int_t atarg)
86 {
88 return InteractionRadius(aproj, atarg) / (pow(aproj, THIRD) + pow(atarg, THIRD));
89 }
91 {
93 return (1.28 * pow(A, THIRD) - 0.76 + 0.8 * pow(A, -THIRD));
94 }
96 {
99 return (SharpRadius(A) * (1. - pow(SharpRadius(A), -2.)));
100 }
102 static Double_t delta_bar_Myers(Int_t Z, Int_t A);
104 {
106 return r0_Myers * pow(2.*Z / (1. - 3.*epsilon_bar_Myers(Z, A)) / (1. - delta_bar_Myers(Z, A)), THIRD);
107 }
108 static Double_t BSS_V0(Int_t zp, Int_t ap, Int_t zt, Int_t at)
109 {
115 Double_t v0 = pow((zp + zt), 2.) / pow(pow(ChargeRadius_Myers(zt, at), 3.) + pow(ChargeRadius_Myers(zp, ap), 3.), THIRD);
116 v0 -= (pow(zt, 2.) / ChargeRadius_Myers(zt, at) + pow(zp, 2.) / ChargeRadius_Myers(zp, ap));
117 v0 *= 3.*e2_Wilcke / 5.;
118 return v0;
119 }
121 {
123 if (*r >= RCTOTAL) {
124 return (ZP * ZT * e2_Wilcke / (*r));
125 }
126 if (*r == 0) return V0;
127 return (V0 - BSS_K * pow(*r, BSS_N));
128 }
130 {
134 static Double_t zeta1 = 1.2511;
135 static Double_t zeta0 = 2.54;
136 static Double_t Kprox = 0.0852;
137
138 Double_t Phi;
139 Double_t zeta = *r - CP - CT;
140 Double_t dzeta = zeta - zeta0;
141 if (zeta <= zeta1) {
142 Phi = -0.5 * pow(dzeta, 2) - Kprox * pow(dzeta, 3);
143 }
144 else {
145 Phi = -3.437 * exp(-zeta / 0.75);
146 }
147 return PROXFACTOR * Phi;
148 }
150 {
152 return ProxPot(r, 0) + VC(r, 0);
153 }
155 {
159
160 Double_t R = TMath::Max(0.1, x[0]);
161 Double_t Vcent = l[0] * (l[0] + 1.) * pow(KVNucleus::hbar, 2) / (2.*MU * mu_Wilcke * R * R);
162 return ProxPot(&x[0], 0) + VC(&x[0], 0) + Vcent;
163 }
165 {
168 return b * k(e_sur_a);
169 }
171 {
174 return l / k(e_sur_a);
175 }
177 {
180 return 10.*(TMath::Pi() / pow(k(e_sur_a), 2)) * pow(lmax + 0.5, 2);
181 }
183 {
187 return sqrt(sigma / (10.*(TMath::Pi() / pow(k(e_sur_a), 2)))) - 0.5;
188 }
190 {
193 return 10.*TMath::Pi() * pow(bmax, 2);
194 }
196 {
199 return sqrt(sigma / (10.*TMath::Pi()));
200 }
202 {
207 return fCentPot;
208 }
209
211 {
214 return fCentPot;
215 }
217 {
218 return fBSS;
219 }
221 {
222 return fProx;
223 }
225 {
226 return fPotential;
227 }
229 {
232 return (0.1071 * Z * Z / pow(A, THIRD) + 22.3);
233 }
235 {
237 return ASYMMFISSIONTKE;
238 }
240 {
243 Double_t I = (A - 2.*Z) / (1.*A);
244 return 0.9517 * (1. - 1.7826 * I * I);
245 }
247 {
249 Double_t zpzt = zp * zt;
250 Double_t D;
251 if (zpzt < 500) {
252 D = 0.3117 * pow(zpzt, 0.2122);
253 }
254 else {
255 D = 1.096 + 1.391e-04 * zpzt;
256 }
257 return InteractionRadius(ap, at) - D;
258 }
260 Double_t Eta(Double_t e_sur_a) const
261 {
263 return 0.15746 * ZP * ZT / pow(e_sur_a, 0.5);
264 }
265 Double_t k(Double_t e_sur_a) const
266 {
268 return 0.2187 * AT * AP * pow(e_sur_a, 0.5) / (1.*AC);
269 }
271 {
273 Double_t krint = k(*x) * RINT;
274 Double_t eta = Eta(*x);
275 if (krint < 2.*eta) return TMath::Pi();
276 return 2.*asin(eta / (krint - eta));
277 }
279 {
281 Double_t QP = QuarterPointAngle(&e, 0);
282 return atan(sin(QP) / (cos(QP) + AP / (1.*AT)));
283 }
285 {
287 return AP * e * (AP * AP + 2.*AP * AT * cos(QuarterPointAngle(&e, 0)) + AT * AT) / (1.*AC * AC);
288 }
289
291 {
292 return Eta(*x) / tan(QuarterPointAngle(x, 0) / 2.);
293 }
295 {
296 return fCMThetaQuart;
297 }
298 TF1* GetLmax() const
299 {
300 return fLmax;
301 }
303 {
306 }
307 Double_t ECM(Double_t e_sur_a) const
308 {
310 return MU * e_sur_a;
311 }
313 {
314 return fSigmaR;
315 }
317 {
319 if (*e_sur_a <= 0) return 0.;
320 Double_t e = VRB / ECM(*e_sur_a);
323 Double_t S = TMath::Min(S1, S2);
324 return S;
325 }
327 {
328 return fSigmaFus;
329 }
330 void DrawAllPotentials(Double_t l = 0) const
331 {
332 TF1* totpot;
333 if (l > 0) {
334 totpot = GetCentrifugalPotential(l);
335 totpot->SetTitle(Form("HIPOT l=%3.0f", l));
336 }
337 else totpot = GetTotalPotential();
338 totpot->SetNpx(1000);
341 totpot->SetLineColor(kBlack);
342 totpot->Draw();
344 GetBSSCoulombPotential()->Draw("same");
345 Double_t max = TMath::Max(GetBSSCoulombPotential()->GetMaximum() * 2., totpot->GetMinimum() * 10);
349 totpot->GetHistogram()->GetYaxis()->SetRangeUser(min, max);
350 ((TPad*)gPad)->BuildLegend();
351 }
352
353 void Print() const;
357
359 {
361 return MU;
362 }
363 Int_t GetAP() const
364 {
365 return AP;
366 }
367 Int_t GetAT() const
368 {
369 return AT;
370 }
371 Int_t GetAC() const
372 {
373 return AC;
374 }
376 {
377 return CP;
378 }
380 {
381 return CT;
382 }
384 {
385 return RINT;
386 }
388
389 ClassDef(KVWilckeReactionParameters, 1) //Reaction parameters for heavy-ion collisions (Wilcke et al)
390};
391
392#endif
int Int_t
ROOT::R::TRInterface & r
#define S1(x)
#define e(i)
double Double_t
#define ClassDef(name, id)
kRed
kBlack
kBlue
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t b
char * Form(const char *fmt,...)
#define gPad
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
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
Double_t GetCrossSectionFromMaxAngularMomentum(Double_t e_sur_a, Double_t lmax) 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)
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.
Double_t GetCrossSectionFromMaxImpactParameter(Double_t bmax) const
static Double_t J_Myers
symmetry energy coefficient
static Double_t r0_Wilcke(Int_t aproj, Int_t atarg)
void DrawAllPotentials(Double_t l=0) const
Double_t GetImpactParameterFromAngularMomentum(Double_t e_sur_a, Double_t l) const
Double_t SigmaR(Double_t *x, Double_t *) const
Double_t GetMaxAngularMomentumFromCrossSection(Double_t e_sur_a, Double_t sigma) 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
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
TF1 * GetCentrifugalPotential(Double_t e_sur_a, Double_t b) const
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 GetMaxImpactParameterFromCrossSection(Double_t sigma) const
Double_t PROXFACTOR
Proximity potential factor.
Double_t GetAngularMomentumFromImpactParameter(Double_t e_sur_a, Double_t b) const
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 * GetCentrifugalPotential(Double_t l) const
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.
virtual void SetLineColor(Color_t lcolor)
virtual void SetRangeUser(Double_t ufirst, Double_t ulast)
virtual Double_t GetMinimum(Double_t xmin=0, Double_t xmax=0, Double_t epsilon=1.E-10, Int_t maxiter=100, Bool_t logx=false) const
virtual TH1 * GetHistogram() const
void SetTitle(const char *title="") override
virtual void SetNpx(Int_t npx=100)
void Draw(Option_t *option="") override
virtual void SetParameter(const TString &name, Double_t value)
TAxis * GetYaxis()
RVec< PromoteType< T > > tan(const RVec< T > &v)
RVec< PromoteType< T > > cos(const RVec< T > &v)
RVec< PromoteType< T > > asin(const RVec< T > &v)
RVec< PromoteTypes< T0, T1 > > pow(const RVec< T0 > &v, const T1 &y)
RVec< PromoteType< T > > exp(const RVec< T > &v)
RVec< PromoteType< T > > atan(const RVec< T > &v)
RVec< PromoteType< T > > sin(const RVec< T > &v)
const Double_t sigma
Double_t x[n]
#define I(x, y, z)
Expr< UnaryOp< Sqrt< T >, Expr< A, T, D, D2, R >, T >, T, D, D2, R > sqrt(const Expr< A, T, D, D2, R > &rhs)
RooArgSet S(Args_t &&... args)
double min(double x, double y)
double max(double x, double y)
Double_t Min(Double_t a, Double_t b)
constexpr Double_t Pi()
constexpr Double_t R()
Double_t Max(Double_t a, Double_t b)
v0
TLine l