Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Patch for new SArray to MArray issue in IMUDeltaFactor #714

Merged
merged 2 commits into from
Nov 9, 2023
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Jump to
Jump to file
Failed to load files.
Diff view
Diff view
14 changes: 14 additions & 0 deletions src/factors/Inertial/IMUDeltaFactor.jl
Original file line number Diff line number Diff line change
Expand Up @@ -369,6 +369,20 @@ function (cf::CalcFactor{<:IMUDeltaFactor})(
return cf(Δmeas, p, q, b)
end

# ArrayPartition{Float64, Tuple{MMatrix{3, 3, Float64, 9}, MVector{3, Float64}, MVector{3, Float64}}}
function (cf::CalcFactor{<:IMUDeltaFactor})(
Δmeas, # point on IMUDeltaGroup
_p::ArrayPartition{Float64, Tuple{MMatrix{3, 3, Float64, 9}, MVector{3, Float64}, MVector{3, Float64}}},
_q::ArrayPartition{Float64, Tuple{MMatrix{3, 3, Float64, 9}, MVector{3, Float64}, MVector{3, Float64}}},
b::AbstractVector = zeros(SVector{6,Float64})
)
p_t = Dates.value(cf.cache.timestams[1])*1e-9
q_t = Dates.value(cf.cache.timestams[2])*1e-9
p = ArrayPartition(SMatrix{3,3,Float64,9}(_p.x[1]), SVector{3}(_p.x[2]), SVector{3}(_p.x[3]), p_t)
q = ArrayPartition(SMatrix{3,3,Float64,9}(_q.x[1]), SVector{3}(_q.x[2]), SVector{3}(_q.x[3]), q_t)
return cf(Δmeas, p, q, SVector{6}(b))
end

function (cf::CalcFactor{<:IMUDeltaFactor})(
Δmeas,
p_SE3,
Expand Down
38 changes: 38 additions & 0 deletions test/inertial/testIMUDeltaFactor.jl
Original file line number Diff line number Diff line change
Expand Up @@ -7,6 +7,7 @@ using DistributedFactorGraphs
using RoME
using RoME: IMUDeltaGroup
using Dates
using StaticArrays
# using ManifoldDiff

##
Expand Down Expand Up @@ -130,6 +131,11 @@ Y = hat(M, SA[0.9, 0.8, 0.7, 0.6, 0.5, 0.4, 0.3, 0.2, 0.1, 1] * 0.1)
# Z1 = ManifoldDiff.differential_exp_argument_lie_approx(M, p, X, Y)
Z2 = RoME.Jr(M, X) * vee(M, Y)

#test right jacobian with [Ad(g)] = Jl*Jr⁻¹ - Chirikjian p29
jr = RoME.Jr(M, X)
jl = RoME.Jr(M, -X)
@test isapprox(jl*inv(jr), Adₚ)

θ=asin(0.1)*10 # for precicely 0.1
X = hat(M, SA[1,0,0, 0,0,0, 0,0,θ, 1] * 0.1)
p = exp(M, ϵ, X)
Expand Down Expand Up @@ -193,6 +199,38 @@ q = ArrayPartition(SMatrix{3,3}(ΔR), SA[0.,0,-1], SA[1.,0,-0.5], 1.0)
a_b = SA[0.,0,0]
ω_b = SA[0.,0,0]

@defVariable RotVelPos Manifolds.ProductGroup(ProductManifold(SpecialOrthogonal(3), TranslationGroup(3), TranslationGroup(3))) ArrayPartition(SA[1 0 0; 0 1 0; 0 0 1.0], SA[0; 0; 0.0], SA[0;0;0.0])
Base.convert(::Type{<:Tuple}, ::IIF.InstanceType{typeof(getManifold(RotVelPos))}) = (:Circular,:Circular,:Circular,:Euclid,:Euclid,:Euclid,:Euclid,:Euclid,:Euclid)

fg = initfg()
fg.solverParams.graphinit = false

foreach(enumerate(Nanosecond.(timestamps[[1,end]] * 10^9))) do (i, nanosecondtime)
addVariable!(fg, Symbol("x",i-1), RotVelPos; nanosecondtime)
end

addFactor!(
fg,
[:x0],
ManifoldPrior(
getManifold(RotVelPos),
ArrayPartition(SA[1.0 0.0 0.0; 0.0 1.0 0.0; 0.0 0.0 1.0], SA[10.0, 0.0, 0.0], SA[0.0, 0.0, 0.0]),
MvNormal(diagm(ones(9)*1e-3))
)
)

addFactor!(fg, [:x0, :x1], fac)

@time IIF.solveGraphParametric!(fg; is_sparse=false, damping_term_min=1e-12, expect_zero_residual=true);
# @time IIF.solveGraphParametric!(fg; stopping_criterion, debug, is_sparse=false, damping_term_min=1e-12, expect_zero_residual=true);

getVariableSolverData(fg, :x0, :parametric).val[1] ≈ ArrayPartition(SA[1.0 0.0 0.0; 0.0 1.0 0.0; 0.0 0.0 1.0], SA[10.0, 0.0, 0.0], SA[0.0, 0.0, 0.0])
x1 = getVariableSolverData(fg, :x1, :parametric).val[1]
@test isapprox(SpecialOrthogonal(3), x1.x[1], ΔR, atol=1e-5)
@test isapprox(x1.x[2], [10, 0, -1], atol=1e-3)
@test isapprox(x1.x[3], [10, 0, -0.5], atol=1e-3)


dt = 0.01
N = 11
dT = (N-1)*dt
Expand Down