1
0
Fork 0
mirror of https://gitlab.rlp.net/mobitar/ReCo.jl.git synced 2024-11-08 22:21:08 +00:00

Added elliptical distance

This commit is contained in:
Mo8it 2022-01-07 00:28:40 +01:00
parent 380c4cc3d5
commit f446c7fbc3

View file

@ -1,20 +1,23 @@
module Shape module Shape
export center_of_mass, gyration_tensor_eigvals_ratio export center_of_mass,
gyration_tensor_eigvals_ratio, gyration_tensor_eigvecs, elliptical_distance
using StaticArrays: SVector, SMatrix using StaticArrays: SVector, SMatrix
using LinearAlgebra: eigvals, Hermitian using LinearAlgebra: eigvals, eigvecs, Hermitian, dot
using ..ReCo: Particle, restrict_coordinate, restrict_coordinates using ..ReCo: Particle, restrict_coordinate, restrict_coordinates
function project_to_unit_circle(x::Float64, half_box_len::Float64) function project_to_unit_circle(x::Float64, half_box_len::Float64)
φ = (x + half_box_len) * π / half_box_len φ = (x + half_box_len) * π / half_box_len
si, co = sincos(φ) si, co = sincos(φ)
return SVector(co, si) return SVector(co, si)
end end
function project_back_from_unit_circle(θ::T, half_box_len::Float64) where {T<:Real} function project_back_from_unit_circle(θ::T, half_box_len::Float64) where {T<:Real}
x = θ * half_box_len / π - half_box_len x = θ * half_box_len / π - half_box_len
return restrict_coordinate(x, half_box_len) return restrict_coordinate(x, half_box_len)
end end
@ -87,8 +90,31 @@ function gyration_tensor(particles::Vector{Particle}, half_box_len::Float64)
end end
function gyration_tensor_eigvals_ratio(particles::Vector{Particle}, half_box_len::Float64) function gyration_tensor_eigvals_ratio(particles::Vector{Particle}, half_box_len::Float64)
ev = eigvals(gyration_tensor(particles, half_box_len)) # Eigenvalues are sorted g_tensor = gyration_tensor(particles, half_box_len)
ev = eigvals(g_tensor) # Eigenvalues are sorted
return ev[1] / ev[2] return ev[1] / ev[2]
end end
function gyration_tensor_eigvecs(particles::Vector{Particle}, half_box_len::Float64)
g_tensor = gyration_tensor(particles, half_box_len)
eig_vecs = eigvecs(g_tensor)
v1 = eig_vecs[:, 1]
v2 = eig_vecs[:, 2]
return (v1, v2)
end
function elliptical_distance(
particle::Particle,
gyration_tensor_eigvec_to_smaller_eigval::SVector{2,Float64},
gyration_tensor_eigvec_to_bigger_eigval::SVector{2,Float64},
goal_gyration_tensor_eigvals_ratio::Float64,
)
cx = dot(particle.c, gyration_tensor_eigvec_to_bigger_eigval)
cy = dot(particle.c, gyration_tensor_eigvec_to_smaller_eigval)
return cx^2 + (cy / goal_gyration_tensor_eigvals_ratio)^2
end
end # module end # module