Skip to content

Commit cb421ce

Browse files
authored
Merge pull request #97 from JuliaComputing/jointrender
update Universal and Spherical rendering
2 parents 4832876 + 44dc7ca commit cb421ce

File tree

5 files changed

+51
-33
lines changed

5 files changed

+51
-33
lines changed

docs/src/examples/spherical_pendulum.md

Lines changed: 2 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -17,9 +17,9 @@ D = Differential(t)
1717
world = Multibody.world
1818
1919
systems = @named begin
20-
joint = Spherical(state=true, isroot=true, phi = 1, radius=0.2, color=[1,1,0,1])
20+
joint = Spherical(state=true, isroot=true, phi = 1, phi_d = 3, radius=0.1, color=[1,1,0,1])
2121
bar = FixedTranslation(r = [0, -1, 0])
22-
body = Body(; m = 1, isroot = false)
22+
body = Body(; m = 1, isroot = false, r_cm=[0.1, 0, 0])
2323
end
2424
2525
connections = [connect(world.frame_b, joint.frame_a)

docs/src/rendering.md

Lines changed: 3 additions & 3 deletions
Original file line numberDiff line numberDiff line change
@@ -2,16 +2,16 @@
22

33
Multibody.jl has an automatic 3D-rendering feature that draws a mechanism in 3D. This can be used to create animations of the mechanism's motion from a solution trajectory, as well as to create interactive applications where the evolution of time can be controlled by the user.
44

5-
The functionality requires the user to load any of the Makie frontend packages, e.g.,
5+
The functionality requires the user to install and load one of the [Makie backend packages](https://docs.makie.org/), e.g.,
66
```julia
7-
using GLMakie # Preferred when running locally
7+
using GLMakie # Preferred
88
```
99
or
1010
```julia
1111
using CairoMakie
1212
```
1313
!!! note "Backend choice"
14-
GLMakie and WGLMakie produce much nicer-looking animations and is also significantly faster than CairoMakie. CairoMakie may be used to produce the graphics in some web environments if constraints imposed by the web environment do not allow any of the GL alternatives. CairoMakie struggles with the Z-order of drawn objects, sometimes making bodies that should have been visible hidden behind bodies that are further back in the scene.
14+
GLMakie and WGLMakie produce much nicer-looking animations and are also significantly faster than CairoMakie. CairoMakie may be used to produce the graphics in some web environments if constraints imposed by the web environment do not allow any of the GL alternatives. CairoMakie struggles with the Z-order of drawn objects, sometimes making bodies that should have been visible hidden behind bodies that are further back in the scene.
1515

1616
After that, the [`render`](@ref) function is the main entry point to create 3D renderings. This function has the following methods:
1717

ext/Render.jl

Lines changed: 10 additions & 5 deletions
Original file line numberDiff line numberDiff line change
@@ -309,8 +309,7 @@ function render!(scene, ::typeof(World), sys, sol, t)
309309
true
310310
end
311311

312-
function render!(scene, ::Union{typeof(Revolute), typeof(RevolutePlanarLoopConstraint)}, sys, sol, t)
313-
# TODO: change to cylinder
312+
function render!(scene, T::Union{typeof(Revolute), typeof(RevolutePlanarLoopConstraint)}, sys, sol, t)
314313
r_0 = get_fun(sol, collect(sys.frame_a.r_0))
315314
n = get_fun(sol, collect(sys.n))
316315
color = get_color(sys, sol, :red)
@@ -321,20 +320,26 @@ function render!(scene, ::Union{typeof(Revolute), typeof(RevolutePlanarLoopConst
321320
catch
322321
0.05f0
323322
end |> Float32
323+
length = try
324+
sol(sol.t[1], idxs=sys.length)
325+
catch
326+
radius
327+
end |> Float32
324328
thing = @lift begin
325-
# radius = sol($t, idxs=sys.radius)
326329
O = r_0($t)
327330
n_a = n($t)
328331
R_w_a = rotfun($t)
329332
n_w = R_w_a*n_a # Rotate to the world frame
330-
p1 = Point3f(O + radius*n_w)
331-
p2 = Point3f(O - radius*n_w)
333+
p1 = Point3f(O + length*n_w)
334+
p2 = Point3f(O - length*n_w)
332335
Makie.GeometryBasics.Cylinder(p1, p2, radius)
333336
end
334337
mesh!(scene, thing; color, specular = Vec3f(1.5), shininess=20f0, diffuse=Vec3f(1))
335338
true
336339
end
337340

341+
render!(scene, ::typeof(Universal), sys, sol, t) = false # To recurse down to rendering the revolute joints
342+
338343
function render!(scene, ::typeof(Spherical), sys, sol, t)
339344
vars = get_fun(sol, collect(sys.frame_a.r_0))
340345
color = get_color(sys, sol, :yellow)

src/joints.jl

Lines changed: 31 additions & 13 deletions
Original file line numberDiff line numberDiff line change
@@ -19,10 +19,11 @@ If `axisflange`, flange connectors for ModelicaStandardLibrary.Mechanics.Rotatio
1919
2020
# Rendering options
2121
- `radius = 0.05`: Radius of the joint in animations
22+
- `length = radius`: Length of the joint in animations
2223
- `color`: Color of the joint in animations, a vector of length 4 with values between [0, 1] providing RGBA values
2324
"""
2425
@component function Revolute(; name, phi0 = 0, w0 = 0, n = Float64[0, 0, 1], axisflange = false,
25-
isroot = true, iscut = false, radius = 0.05, color = [0.5019608f0,0.0f0,0.5019608f0,1.0f0], state_priority = 3.0)
26+
isroot = true, iscut = false, radius = 0.05, length = radius, color = [0.5019608f0,0.0f0,0.5019608f0,1.0f0], state_priority = 3.0)
2627
if !(eltype(n) <: Num)
2728
norm(n) 1 || error("Axis of rotation must be a unit vector")
2829
end
@@ -31,6 +32,7 @@ If `axisflange`, flange connectors for ModelicaStandardLibrary.Mechanics.Rotatio
3132
@parameters n[1:3]=n [description = "axis of rotation"]
3233
pars = @parameters begin
3334
radius = radius, [description = "radius of the joint in animations"]
35+
length = length, [description = "length of the joint in animations"]
3436
color[1:4] = color, [description = "color of the joint in animations (RGBA)"]
3537
end
3638
@variables tau(t)=0 [
@@ -263,16 +265,32 @@ Joint with 3 constraints that define that the origin of `frame_a` and the origin
263265
add_params(sys, pars; name)
264266
end
265267

268+
269+
"""
270+
Universal(; name, n_a, n_b, phi_a = 0, phi_b = 0, w_a = 0, w_b = 0, a_a = 0, a_b = 0, state_priority=10)
271+
272+
Joint where `frame_a` rotates around axis `n_a` which is fixed in `frame_a` and `frame_b` rotates around axis `n_b` which is fixed in `frame_b`. The two frames coincide when `revolute_a.phi=0` and `revolute_b.phi=0`. This joint has the following potential states;
273+
274+
- The relative angle `phi_a = revolute_a.phi` [rad] around axis `n_a`
275+
- the relative angle `phi_b = revolute_b.phi` [rad] around axis `n_b`
276+
- the relative angular velocity `w_a = D(phi_a)`
277+
- the relative angular velocity `w_b = D(phi_b)`
278+
"""
266279
@component function Universal(; name, n_a = [1, 0, 0], n_b = [0, 1, 0], phi_a = 0,
267280
phi_b = 0,
281+
state_priority = 10,
268282
w_a = 0,
269283
w_b = 0,
270284
a_a = 0,
271-
a_b = 0)
285+
a_b = 0,
286+
radius = 0.05f0,
287+
length = radius,
288+
color = [1,0,0,1]
289+
)
272290
@named begin
273291
ptf = PartialTwoFrames()
274-
revolute_a = Revolute(n = n_a, isroot = false)
275-
revolute_b = Revolute(n = n_b, isroot = false)
292+
revolute_a = Revolute(n = n_a, isroot = false, radius, length, color)
293+
revolute_b = Revolute(n = n_b, isroot = false, radius, length, color)
276294
end
277295
@unpack frame_a, frame_b = ptf
278296
@parameters begin
@@ -288,32 +306,32 @@ end
288306
@variables begin
289307
(phi_a(t) = phi_a),
290308
[
291-
state_priority = 10,
309+
state_priority = state_priority,
292310
description = "Relative rotation angle from frame_a to intermediate frame",
293311
]
294312
(phi_b(t) = phi_b),
295313
[
296-
state_priority = 10,
314+
state_priority = state_priority,
297315
description = "Relative rotation angle from intermediate frame to frame_b",
298316
]
299317
(w_a(t) = w_a),
300318
[
301-
state_priority = 10,
319+
state_priority = state_priority,
302320
description = "First derivative of angle phi_a (relative angular velocity a)",
303321
]
304322
(w_b(t) = w_b),
305323
[
306-
state_priority = 10,
324+
state_priority = state_priority,
307325
description = "First derivative of angle phi_b (relative angular velocity b)",
308326
]
309327
(a_a(t) = a_a),
310328
[
311-
state_priority = 10,
329+
state_priority = state_priority,
312330
description = "Second derivative of angle phi_a (relative angular acceleration a)",
313331
]
314332
(a_b(t) = a_b),
315333
[
316-
state_priority = 10,
334+
state_priority = state_priority,
317335
description = "Second derivative of angle phi_b (relative angular acceleration b)",
318336
]
319337
end
@@ -855,12 +873,12 @@ LinearAlgebra.normalize(a::Vector{Num}) = a / norm(a)
855873
"""
856874
Planar(; n = [0,0,1], n_x = [1,0,0], cylinderlength = 0.1, cylinderdiameter = 0.05, cylindercolor = [1, 0, 1, 1], boxwidth = 0.3*cylinderdiameter, boxheight = boxwidth, boxcolor = [0, 0, 1, 1])
857875
858-
Joint where frame_b can move in a plane and can rotate around an
876+
Joint where `frame_b` can move in a plane and can rotate around an
859877
axis orthogonal to the plane. The plane is defined by
860878
vector `n` which is perpendicular to the plane and by vector `n_x`,
861879
which points in the direction of the x-axis of the plane.
862-
frame_a and frame_b coincide when s_x=prismatic_x.s=0,
863-
s_y=prismatic_y.s=0 and phi=revolute.phi=0.
880+
`frame_a` and `frame_b` coincide when `s_x=prismatic_x.s=0,
881+
s_y=prismatic_y.s=0` and `phi=revolute.phi=0`.
864882
"""
865883
@mtkmodel Planar begin
866884
@structural_parameters begin

test/runtests.jl

Lines changed: 5 additions & 10 deletions
Original file line numberDiff line numberDiff line change
@@ -545,17 +545,13 @@ ssys = structural_simplify(IRSystem(model))#, alias_eliminate = true)
545545
# ssys = structural_simplify(model, allow_parameters = false)
546546
prob = ODEProblem(ssys,
547547
[collect(body.body.w_a .=> 0);
548-
collect(body.body.v_0 .=> 0);
549-
collect(D.(body.body.phi)) .=> 1;
550-
# collect((body.body.r_0)) .=> collect((body.r_0));
551-
collect(D.(D.(body.body.phi))) .=> 1], (0, 10))
548+
collect(body.body.v_0 .=> 0);], (0, 10))
552549

553550
# @test_skip begin # The modelica example uses angles_fixed = true, which causes the body component to run special code for variable initialization. This is not yet supported by MTK
554551
# Without proper initialization, the example fails most of the time. Random perturbation of u0 can make it work sometimes.
555552
sol = solve(prob, Rodas4())
556553
@test SciMLBase.successful_retcode(sol)
557554

558-
@info "Initialization broken, initial value for body.r_0 not respected, add tests when MTK has a working initialization"
559555
doplot() && plot(sol, idxs = [body.r_0...]) |> display
560556
# end
561557

@@ -650,29 +646,28 @@ end
650646
# ==============================================================================
651647
## universal pendulum
652648
# ==============================================================================
653-
649+
using LinearAlgebra
654650
@testset "Universal pendulum" begin
655651
t = Multibody.t
656652
D = Differential(t)
657653
world = Multibody.world
658654
@named begin
659-
joint = Universal()
655+
joint = Universal(length=0.1)
660656
bar = FixedTranslation(r = [0, -1, 0])
661-
body = Body(; m = 1, isroot = false)
657+
body = Body(; m = 1, isroot = false, r_cm=[0.1, 0, 0])
662658
end
663659
connections = [connect(world.frame_b, joint.frame_a)
664660
connect(joint.frame_b, bar.frame_a)
665661
connect(bar.frame_b, body.frame_a)]
666662
@named model = ODESystem(connections, t, systems = [world, joint, bar, body])
663+
model = complete(model)
667664
ssys = structural_simplify(IRSystem(model))
668665

669666
prob = ODEProblem(ssys,
670667
[
671668
joint.phi_b => sqrt(2);
672669
joint.revolute_a.phi => sqrt(2);
673-
D(joint.revolute_a.phi) => 0;
674670
joint.revolute_b.phi => 0;
675-
D(joint.revolute_b.phi) => 0
676671
], (0, 10))
677672
sol2 = solve(prob, Rodas4())
678673
@test SciMLBase.successful_retcode(sol2)

0 commit comments

Comments
 (0)