Skip to content

Commit a4cd668

Browse files
Fix #22 extra fields
1 parent 0259154 commit a4cd668

File tree

4 files changed

+101
-12
lines changed

4 files changed

+101
-12
lines changed

README.md

Lines changed: 7 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -161,6 +161,7 @@ Each type above then contains three fields:
161161
- `values::Vector{B}`
162162
- `end_idxs::Vector{Int}``
163163

164+
Note that the ordering of the fields matters.
164165
`B` is the `BottomType`, which has to be the same as the eltype for the array
165166
in the leaf types. `T` is another `AbstractMultiScaleArray`. Thus at each level,
166167
an` AbstractMultiScaleArray` contains some information of its own (`values`), the
@@ -201,8 +202,12 @@ struct Cell{B} <: AbstractMultiScaleArrayLeaf{B}
201202
end
202203
```
203204

204-
Then we'd construct cells with `cell3 = Cell([3.0; 2.0; 5.0], :BCell)`, and can
205-
give it a cell type. This information is part of the call, so
205+
Note that the ordering of the fields matters here: the extra fields must come
206+
after the standard fields (so for a leaf it comes after `values`, for a standard
207+
multiscale array it would come after `nodes,values,end_idxs`).
208+
Then we'd construct cells with
209+
`cell3 = Cell([3.0; 2.0; 5.0], :BCell)`, and can give it a cell type.
210+
This information is part of the call, so
206211

207212
```julia
208213
for (cell, y, z) in level_iter_idx(embryo, 2)

src/shape_construction.jl

Lines changed: 16 additions & 10 deletions
Original file line numberDiff line numberDiff line change
@@ -9,17 +9,23 @@ size(m::AbstractMultiScaleArray) = (length(m),)
99
parameterless_type(T::Type) = Base.typename(T).wrapper
1010
parameterless_type(x) = parameterless_type(typeof(x))
1111

12-
similar(m::AbstractMultiScaleArray) =
13-
construct(typeof(m), deepcopy(m.nodes), deepcopy(m.values))
14-
similar(m::AbstractMultiScaleArrayLeaf) =
15-
construct(typeof(m), similar(m.values))
16-
similar(m::AbstractMultiScaleArrayLeaf, T::Type) =
17-
construct(parameterless_type(m), similar(m.values, T))
18-
similar(m::AbstractMultiScaleArray,T::Type) =
19-
construct(parameterless_type(m), [similar(v, T) for v in m.nodes], similar(m.values, T))
12+
@generated function similar(m::AbstractMultiScaleArrayLeaf,::Type{T}=eltype(m)) where T
13+
assignments = [s == :x ? :(similar(m.x, T)) :
14+
(sq = Meta.quot(s); :(deepcopy(getfield(m, $sq))))
15+
for s in fieldnames(m)[2:end]] # 1 is values
16+
:(construct(parameterless_type(m), similar(m.values,T),$(assignments...)))
17+
end
18+
19+
@generated function similar(m::AbstractMultiScaleArray,::Type{T}=eltype(m)) where {T}
20+
assignments = [s == :x ? :(similar(m.x, T, dims)) :
21+
(sq = Meta.quot(s); :(deepcopy(getfield(m, $sq))))
22+
for s in fieldnames(m)[4:end]] # 1:3 is nodes,values,end_idxs
23+
:(construct(parameterless_type(m), recursive_similar(m.nodes,T),similar(m.values, T),$(assignments...)))
24+
end
25+
26+
recursive_similar(x,T) = [similar(y, T) for y in x]
2027

21-
construct(::Type{T}, values, args...) where {T<:AbstractMultiScaleArrayLeaf} =
22-
T(values, args...)
28+
construct(::Type{T}, args...) where {T<:AbstractMultiScaleArrayLeaf} = T(args...)
2329

2430
function __construct(nodes::Vector{<:AbstractMultiScaleArray})
2531
end_idxs = Vector{Int}(length(nodes))

test/additional_fields_test.jl

Lines changed: 77 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,77 @@
1+
using OrdinaryDiffEq, Base.Test
2+
using MultiScaleArrays
3+
4+
struct CellGenotype{B<:Float64} <: AbstractMultiScaleArrayLeaf{B}
5+
values::Vector{B}
6+
Iam::Symbol
7+
end
8+
9+
struct Host{T<:AbstractMultiScaleArray,B<:Float64} <: AbstractMultiScaleArrayHead{B}
10+
nodes::Vector{T}
11+
values::Vector{B}
12+
end_idxs::Vector{Int}
13+
end
14+
15+
host = construct(Host,[
16+
CellGenotype(# Genotype present at start for symbiont 1
17+
[Float64(1)],
18+
:aSymbiontGenotype, # I am
19+
),
20+
CellGenotype(# Genotype present at start for symbiont 1
21+
[Float64(1)],
22+
:aSymbiontGenotype, # I am
23+
)]
24+
)
25+
26+
evolve = function (t,u::Host,du)
27+
for i in 1:length(u)
28+
du[i] = 2*u[i]
29+
end
30+
end
31+
32+
similar(host)
33+
similar(host.nodes[1])
34+
similar(host.nodes[1],Float64)
35+
similar(host,Float64)
36+
37+
## Solve
38+
tspan = (0.0,50)
39+
prob = ODEProblem(evolve, host, tspan)
40+
sol = solve(prob,Tsit5())
41+
42+
43+
44+
struct Host2{T<:AbstractMultiScaleArray,B<:Float64} <: AbstractMultiScaleArrayHead{B}
45+
nodes::Vector{T}
46+
values::Vector{B}
47+
end_idxs::Vector{Int}
48+
Iam::Symbol
49+
end
50+
51+
host2 = construct(Host2,[
52+
CellGenotype(# Genotype present at start for symbiont 1
53+
[Float64(1)],
54+
:aSymbiontGenotype, # I am
55+
),
56+
CellGenotype(# Genotype present at start for symbiont 1
57+
[Float64(1)],
58+
:aSymbiontGenotype, # I am
59+
)],Float64[],
60+
:my_cell_type
61+
)
62+
63+
evolve = function (t,u::Host2,du)
64+
for i in 1:length(u)
65+
du[i] = 2*u[i]
66+
end
67+
end
68+
69+
similar(host2)
70+
similar(host2.nodes[1])
71+
similar(host2.nodes[1],Float64)
72+
similar(host2,Float64)
73+
74+
## Solve
75+
tspan = (0.0,50)
76+
prob = ODEProblem(evolve, host2, tspan)
77+
sol = solve(prob,Tsit5())

test/runtests.jl

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -5,4 +5,5 @@ using Base.Test
55
@time @testset "Indexing and Creation Tests" begin include("indexing_and_creation_tests.jl") end
66
@time @testset "Values Indexing" begin include("values_indexing.jl") end
77
@time @testset "Get Indices Tests" begin include("get_indices.jl") end
8+
@time @testset "Additional Fields Test" begin include("additional_fields_test.jl") end
89
@time @testset "Dynamic DiffEq Tests" begin include("dynamic_diffeq.jl") end

0 commit comments

Comments
 (0)