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

Replaced goal ratio with a_b_ratio and added abs to eigvals ratio

This commit is contained in:
Mo8it 2022-01-26 00:20:53 +01:00
parent 1003ccf851
commit 9c26795790
4 changed files with 17 additions and 23 deletions

View file

@ -7,7 +7,7 @@ struct EnvHelperSharedProps{H<:AbstractHook}
n_steps_before_actions_update::Int64
goal_gyration_tensor_eigvals_ratio::Float64
elliptical_a_b_ratio::Float64
n_particles::Int64
@ -22,7 +22,7 @@ struct EnvHelperSharedProps{H<:AbstractHook}
agent::Agent,
hook::H,
n_steps_before_actions_update::Int64,
goal_gyration_tensor_eigvals_ratio::Float64,
elliptical_a_b_ratio::Float64,
n_particles::Int64,
) where {H<:AbstractHook}
return new{H}(
@ -30,7 +30,7 @@ struct EnvHelperSharedProps{H<:AbstractHook}
agent,
hook,
n_steps_before_actions_update,
goal_gyration_tensor_eigvals_ratio,
elliptical_a_b_ratio,
n_particles,
fill(0, n_particles),
fill(0, n_particles),

View file

@ -51,8 +51,7 @@ mutable struct LocalCOMEnvHelper <: EnvHelper
max_elliptical_distance::Float64
function LocalCOMEnvHelper(shared::EnvHelperSharedProps, half_box_len::Float64, skin_r)
max_elliptical_distance =
sqrt(2) * half_box_len / shared.goal_gyration_tensor_eigvals_ratio
max_elliptical_distance = sqrt(2) * half_box_len / shared.elliptical_a_b_ratio
max_distance_to_local_center_of_mass = skin_r
@ -193,7 +192,7 @@ function update_reward!(env::LocalCOMEnv, env_helper::LocalCOMEnvHelper, particl
env_helper.center_of_mass,
env_helper.gyration_tensor_eigvec_to_smaller_eigval,
env_helper.gyration_tensor_eigvec_to_bigger_eigval,
env_helper.shared.goal_gyration_tensor_eigvals_ratio,
env_helper.shared.elliptical_a_b_ratio,
env_helper.half_box_len,
)

View file

@ -53,7 +53,7 @@ end
function run_rl(;
EnvType::Type{E},
parent_dir_appendix::String,
goal_gyration_tensor_eigvals_ratio::Float64,
elliptical_a_b_ratio::Float64,
n_episodes::Int64=200,
episode_duration::Float64=50.0,
update_actions_at::Float64=0.1,
@ -64,7 +64,7 @@ function run_rl(;
packing_ratio::Float64=0.22,
show_progress::Bool=true,
) where {E<:Env}
@assert 0.0 <= goal_gyration_tensor_eigvals_ratio <= 1.0
@assert 0.0 <= elliptical_a_b_ratio <= 1.0
@assert n_episodes > 0
@assert episode_duration > 0
@assert update_actions_at in 0.001:0.001:episode_duration
@ -93,12 +93,7 @@ function run_rl(;
hook = TotalRewardPerEpisode()
env_helper_shared = EnvHelperSharedProps(
env,
agent,
hook,
n_steps_before_actions_update,
goal_gyration_tensor_eigvals_ratio,
n_particles,
env, agent, hook, n_steps_before_actions_update, elliptical_a_b_ratio, n_particles
)
env_helper_args = (half_box_len=sim_consts.half_box_len, skin_r=sim_consts.skin_r)

View file

@ -4,7 +4,7 @@ export center_of_mass,
gyration_tensor_eigvals_ratio, gyration_tensor_eigvecs, elliptical_distance
using StaticArrays: SVector, SMatrix
using LinearAlgebra: eigvals, eigvecs, Hermitian, dot
using LinearAlgebra: LinearAlgebra as LA
using ..ReCo: ReCo, Particle
@ -88,7 +88,7 @@ function gyration_tensor(
S22 += shifted_c[2]^2
end
return Hermitian(SMatrix{2,2}(S11, S12, S12, S22))
return LA.Hermitian(SMatrix{2,2}(S11, S12, S12, S22))
end
function gyration_tensor(particles::AbstractVector{Particle}, half_box_len::Real)
@ -101,15 +101,15 @@ function gyration_tensor_eigvals_ratio(
particles::AbstractVector{Particle}, half_box_len::Real
)
g_tensor = gyration_tensor(particles, half_box_len)
ev = eigvals(g_tensor) # Eigenvalues are sorted
return ev[1] / ev[2]
ev = LA.eigvals(g_tensor) # Eigenvalues are sorted
return abs(ev[1] / ev[2])
end
function gyration_tensor_eigvecs(
particles::AbstractVector{Particle}, half_box_len::R, COM::SVector{2,R}
) where {R<:Real}
g_tensor = gyration_tensor(particles, half_box_len, COM)
eig_vecs = eigvecs(g_tensor)
eig_vecs = LA.eigvecs(g_tensor)
v1 = eig_vecs[:, 1]
v2 = eig_vecs[:, 2]
@ -122,15 +122,15 @@ function elliptical_distance(
COM::SVector{2,R},
gyration_tensor_eigvec_to_smaller_eigval::SVector{2,R},
gyration_tensor_eigvec_to_bigger_eigval::SVector{2,R},
goal_gyration_tensor_eigvals_ratio::R,
elliptical_a_b_ratio::R,
half_box_len::R,
) where {R<:Real}
v = ReCo.minimum_image(v - COM, half_box_len)
x = dot(v, gyration_tensor_eigvec_to_bigger_eigval)
y = dot(v, gyration_tensor_eigvec_to_smaller_eigval)
x = LA.dot(v, gyration_tensor_eigvec_to_bigger_eigval)
y = LA.dot(v, gyration_tensor_eigvec_to_smaller_eigval)
return sqrt(x^2 + (y / goal_gyration_tensor_eigvals_ratio)^2)
return sqrt(x^2 + (y / elliptical_a_b_ratio)^2)
end
end # module