NuclearRMF/nucleons.jl

209 lines
7.5 KiB
Julia
Raw Permalink Blame History

This file contains ambiguous Unicode characters

This file contains Unicode characters that might be confused with other characters. If you think that this is intentional, you can safely ignore this warning. Use the Escape button to reveal them.

using LinearAlgebra, DifferentialEquations, Interpolations
include("bisection.jl")
include("common.jl")
include("system.jl")
const M_n = 939.0 # MeV/c2
const M_p = 939.0 # MeV/c2
"The spherical Dirac equation that returns du=[dg, df] in-place where
u=[g, f] are the reduced radial components evaluated at r,
κ is the generalized angular momentum,
p is true for proton and false for neutron,
E in the energy in MeV,
f1(r) = M-Φ0(r)-W0(r)-(p-0.5)*2B0(r)-p*A0(r) is a function of r in MeV (see optimized_dirac_potentials()),
f2(r) = -M+Φ0(r)-W0(r)-(p-0.5)*2B0(r)-p*A0(r) is a function of r in MeV (see optimized_dirac_potentials()),
r is the radius in fm.
Reference: P. Giuliani, K. Godbey, E. Bonilla, F. Viens, and J. Piekarewicz, Frontiers in Physics 10, (2023)"
function dirac!(du::Vector{Float64}, u::Vector{Float64}, (κ, E, f1, f2), r::Float64) # TODO: Static typing
(g, f) = u
@inbounds du[1] = -(κ/(r + r_reg)) * g + (E + f1(r)) * f / ħc
@inbounds du[2] = (κ/(r + r_reg)) * f - (E + f2(r)) * g / ħc
end
"Get the potentials f1 and f2 that goes into the Dirac equation, defined as
f1(r) = M-Φ0(r)-W0(r)-(p-0.5)*2B0(r)-p*A0(r),
f2(r) = -M+Φ0(r)-W0(r)-(p-0.5)*2B0(r)-p*A0(r)."
function optimized_dirac_potentials(p::Bool, s::system)
M = p ? M_p : M_n
f1s = zero_array(s)
f2s = zero_array(s)
@. f1s = M - s.Φ0 - s.W0 - (p - 0.5) * 2 * s.B0 - p * s.A0
@. f2s = -M + s.Φ0 - s.W0 - (p - 0.5) * 2 * s.B0 - p * s.A0
f1 = linear_interpolation(rs(s), f1s)
f2 = linear_interpolation(rs(s), f2s)
return (f1, f2)
end
"Approximate boundary condition for u(r)=[g(r), f(r)] at r -> ∞ where
κ is the generalized angular momentum,
p is true for proton and false for neutron,
E is the energy in MeV,
r is the radius in fm."
function asymp_BC(κ::Int, p::Bool, E::Float64, r::Float64)
M = p ? M_p : M_n
g = 1.0
f = ħc / (E + M) * (-(M^2 - E^2) + κ/r) * g
return [g, f]
end
"Initial boundary condition for u(r)=[g(r), f(r)] at r=0"
init_BC() = [1.0, 1.0] # TODO: Why not [0.0, 1.0]?
"Solve the Dirac equation and return the wave function u(r)=[g(r), f(r)] where
divs is the number of mesh divisions so solution would be returned as a 2×(1+divs) matrix,
the other parameters are the same from dirac!(...)."
function solveNucleonWf(κ, p::Bool, E, s::system; normalize=true, algo=Vern9())
(f1, f2) = optimized_dirac_potentials(p, s)
# partitioning
mid_idx = s.divs ÷ 2
r_mid = rs(s)[mid_idx]
left_r = rs(s)[1:mid_idx]
right_r = rs(s)[mid_idx:end]
# left partition
prob = ODEProblem(dirac!, init_BC(), (0, r_mid))
sol = solve(prob, algo, p=(κ, E, f1, f2), saveat=left_r)
wf_left = hcat(sol.u...)
# right partition
prob = ODEProblem(dirac!, asymp_BC(κ, p, E, s.r_max), (s.r_max, r_mid))
sol = solve(prob, algo, p=(κ, E, f1, f2), saveat=right_r)
wf_right = reverse(hcat(sol.u...); dims=2)
# join two segments
u1 = wf_left[:, end]
u2 = wf_right[:, 1]
if norm(u2) < 1e-10
@warn "Right partition too small to rescale, setting to zero"
wf_right .= 0.0
else
proj = only(u1' * u2) / norm(u2)^2
wf_right .*= proj
end
wf = hcat(wf_left[:, 1:(end - 1)], wf_right)
if normalize
wf ./= norm(wf) * sqrt(Δr(s)) # integration by Reimann sum
end
return wf
end
"Returns a function that solves the Dirac equation in two partitions and returns
the determinant of [g_left(r) g_right(r); f_left(r) f_right(r)],
where is r is in fm."
function determinantFunc(κ, p::Bool, s::system, r::Float64=s.r_max/2, algo=Vern9())
(f1, f2) = optimized_dirac_potentials(p, s)
prob_left = ODEProblem(dirac!, init_BC(), (0, r))
function func(E)
prob_right = ODEProblem(dirac!, asymp_BC(κ, p, E, s.r_max), (s.r_max, r))
u_left = solve(prob_left, algo, p=(κ, E, f1, f2), saveat=[r])
u_right = solve(prob_right, algo, p=(κ, E, f1, f2), saveat=[r])
return u_left[1, 1] * u_right[2, 1] - u_right[1, 1] * u_left[2, 1]
end
return func
end
"Find all bound energies between E_min (=850.0) and E_max (=938.0) where
tol_digits determines the precision for root finding and the threshold for identifying duplicates,
the other parameters are the same from dirac!(...)."
function findEs(κ, p::Bool, s::system, E_min=850.0, E_max=938.0; tol_digits=8)
func = determinantFunc(κ, p, s)
Es = find_all_zeros(func, E_min, E_max; partitions=20, tol=1/10^tol_digits)
return unique(E -> round(E; digits=tol_digits), Es)
end
"Find all orbitals and return two lists containing κ values and corresponding energies for a single species where
the other parameters are defined above"
function findAllOrbitals(p::Bool, s::system, E_min=850.0, E_max=938.0)
κs = Int[]
Es = Float64[]
# start from κ=-1 and go both up and down
for direction in [-1, 1]
for κ in direction * (1:100) # cutoff is 100
new_Es = findEs(κ, p, s, E_min, E_max)
if isempty(new_Es); break; end
append!(Es, new_Es)
append!(κs, fill(κ, length(new_Es)))
end
end
return (κs, Es)
end
"For a given list of κ values with corresponding energies, attempt to fill Z_or_N lowest lying orbitals and return the spectrum"
function fillNucleons(Z_or_N::Int, κs, Es)::spectrum
sort_i = sortperm(Es)
occ = zeros(Int, length(κs))
for i in sort_i
if Z_or_N 0; break; end;
max_occ = 2 * j_κ(κs[i]) + 1
occ[i] = min(max_occ, Z_or_N)
Z_or_N -= occ[i]
end
@assert Z_or_N == 0 "All orbitals could not be filled"
return spectrum(κs, Es, occ)
end
"Total angular momentum j for a given κ value"
j_κ(κ::Int) = abs(κ) - 1/2
"Orbital angular momentum l for a given κ value"
l_κ(κ::Int) = abs(κ) - (κ < 0) # since true = 1 and false = 0
"Calculate scalar and vector densities of a nucleon species on [0,r_max] divided into (divs+1) points and returns them as vectors (ρ_s, ρ_v) where
the arrays κs, Es, occs tabulate the energies and occupation numbers corresponding to each κ,
the other parameters are defined above"
function calculateNucleonDensity(p::Bool, s::system)::Tuple{Vector{Float64}, Vector{Float64}}
spectrum = p ? s.p_spectrum : s.n_spectrum
(κs, Es, occs) = (spectrum.κ, spectrum.E, spectrum.occ)
ρr2 = zeros(2, s.divs + 1) # ρ×r² values
for (κ, E, occ) in zip(κs, Es, occs)
wf = solveNucleonWf(κ, p, E, s; normalize=true)
wf2 = wf .* wf
ρr2 += (occ / (4 * pi)) * wf2 # 2j+1 factor is accounted in the occupancy number
end
r2s = rs(s).^2
ρ = ρr2 ./ transpose(r2s)
ρ[:, 1] .= ρ[:, 2] # dirty fix for NaN at r=0
ρ_s = ρ[1, :] - ρ[2, :] # g^2 - f^2
ρ_v = ρ[1, :] + ρ[2, :] # g^2 + f^2
return (ρ_s, ρ_v)
end
"For a nucleon species, solve the Dirac equation and save the spectrum & densities in-place where
the other parameters are defined above"
function solveNucleonDensity!(p::Bool, s::system, E_min=850.0, E_max=938.0)
κs, Es = findAllOrbitals(p, s, E_min, E_max)
spec = fillNucleons(Z_or_N(s, p), κs, Es)
if p
s.p_spectrum = spec
else
s.n_spectrum = spec
end
(ρ_s, ρ_v) = calculateNucleonDensity(p, s)
if p
s.ρ_sp = ρ_s
s.ρ_vp = ρ_v
else
s.ρ_sn = ρ_s
s.ρ_vn = ρ_v
end
end
"Total energy of filled nucleons in the system"
nucleon_E(s::system) = sum(s.p_spectrum.occ .* (s.p_spectrum.E .- M_p)) + sum(s.n_spectrum.occ .* (s.n_spectrum.E .- M_n))