1
0
Fork 0
mirror of https://gitlab.rlp.net/mobitar/ReCo.jl.git synced 2024-12-21 00:51:21 +00:00
ReCo.jl/src/shape.jl

68 lines
2 KiB
Julia
Raw Normal View History

2021-12-03 13:41:01 +00:00
using StaticArrays: SVector, SMatrix
2021-12-04 15:19:14 +00:00
using LinearAlgebra: eigvals, Hermitian
2021-12-03 13:16:35 +00:00
2021-12-10 02:16:18 +00:00
function project_to_unit_circle(x::Float64, half_box_len::Float64)
φ = (x + half_box_len) * π / half_box_len
2021-12-08 20:53:15 +00:00
si, co = sincos(φ)
return SVector(co, si)
2021-12-03 13:16:35 +00:00
end
2021-12-10 02:16:18 +00:00
function project_back_from_unit_circle(θ::T, half_box_len::Float64) where {T<:Real}
x = θ * half_box_len / π - half_box_len
return restrict_coordinate(x, half_box_len)
2021-12-03 13:16:35 +00:00
end
2021-12-10 02:16:18 +00:00
function center_of_mass(particles::Vector{Particle}, half_box_len::Float64)
2021-12-07 04:52:10 +00:00
x_proj_sum = SVector(0.0, 0.0)
y_proj_sum = SVector(0.0, 0.0)
2021-12-03 13:16:35 +00:00
2021-12-08 20:53:15 +00:00
for p in particles
2021-12-10 02:16:18 +00:00
x_proj_sum += project_to_unit_circle(p.c[1], half_box_len)
y_proj_sum += project_to_unit_circle(p.c[2], half_box_len)
2021-12-03 13:16:35 +00:00
end
# Prevent for example atan(1e-16, 1e-15) != 0 with rounding
digits = 5
2021-12-10 02:16:18 +00:00
# No need for 1/n_particles with atan
2021-12-03 13:16:35 +00:00
# If proj is (0, 0) then COM is 0 or L or -L. Here, 0 is choosen with θ = π
2021-12-04 15:19:14 +00:00
if round(x_proj_sum[1]; digits=digits) == round(x_proj_sum[2]; digits=digits) == 0
2021-12-03 13:16:35 +00:00
x_θ = π
else
x_θ = atan(x_proj_sum[2], x_proj_sum[1])
end
2021-12-04 15:19:14 +00:00
if round(y_proj_sum[1]; digits=digits) == round(y_proj_sum[2]; digits=digits) == 0
2021-12-03 13:16:35 +00:00
y_θ = π
else
y_θ = atan(y_proj_sum[2], y_proj_sum[1])
end
2021-12-10 02:16:18 +00:00
COM_x = project_back_from_unit_circle(x_θ, half_box_len)
COM_y = project_back_from_unit_circle(y_θ, half_box_len)
2021-12-08 20:53:15 +00:00
return SVector(COM_x, COM_y)
2021-12-03 13:41:01 +00:00
end
2021-12-10 02:16:18 +00:00
function gyration_tensor(particles::Vector{Particle}, half_box_len::Float64)
COM = center_of_mass(particles, half_box_len)
2021-12-03 13:41:01 +00:00
S11 = 0.0
S12 = 0.0
S22 = 0.0
2021-12-08 20:53:15 +00:00
for p in particles
2021-12-15 03:45:15 +00:00
shifted_c = restrict_coordinates(p.c - COM, half_box_len)
2021-12-03 13:41:01 +00:00
2021-12-08 20:53:15 +00:00
S11 += shifted_c[1]^2
S12 += shifted_c[1] * shifted_c[2]
S22 += shifted_c[2]^2
2021-12-03 13:41:01 +00:00
end
2021-12-04 15:19:14 +00:00
return Hermitian(SMatrix{2,2}(S11, S12, S12, S22))
end
2021-12-10 02:16:18 +00:00
function gyration_tensor_eigvals_ratio(particles::Vector{Particle}, half_box_len::Float64)
ev = eigvals(gyration_tensor(particles, half_box_len)) # Eigenvalues are sorted
2021-12-04 15:19:14 +00:00
return ev[1] / ev[2]
2021-12-03 13:16:35 +00:00
end