{ "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "## Julia Implementation of Host SEIR + Vector SEI\n", "\n", "*Author*: Carl A. B. Pearson @pearsonca\n", "\n", "*Date*: 2018-10-02\n", "\n", "This version considers multiple hosts, which introduces some indexing complications. This approach adopts the indexing scheme of:\n", "\n", " - all of host $i$ compartments ($S_H$, $E_H$, etc), for $i \\in 1\\ldots N$\n", " - vector compartments last\n", "\n", "With this approach, we can re-use the solution for one host, with slight modifications:" ] }, { "cell_type": "code", "execution_count": 1, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "F (generic function with 1 method)" ] }, "execution_count": 1, "metadata": {}, "output_type": "execute_result" } ], "source": [ "H_comps = 4\n", "V_comps = 3\n", "\n", "# in sub functions, du / u are the particular relevant slices only\n", "function F1H(du,u,p,t,I_V,N_H)\n", " S_H, E_H, I_H, R_H = u\n", " \n", " # host dynamics\n", " host_infection = (p.β*S_H*I_V)/N_H\n", " host_mortality = p.μ_H .* u # include S_H, so easier to remove mortality\n", " host_births = sum(host_mortality)\n", " host_progression = p.σ_H*E_H\n", " recovery = p.λ*I_H\n", " \n", " du[1] = -host_infection + host_births\n", " du[2] = host_infection - host_progression\n", " du[3] = host_progression - recovery\n", " du[4] = recovery\n", " du[1:end] -= host_mortality \n", "end\n", "\n", "# in sub functions, du / u are the particular relevant slices only\n", "function FV(du,u,p,t,sum_β_I_H)\n", " S_V, E_V, I_V = u\n", " vec_infection = sum_β_I_H*S_V/p.N_H\n", " vec_mortality = p.μ_V .* u # include S_V, so easier to remove mortality\n", " vec_births = sum(vec_mortality)\n", " vec_progression = p.σ_V*E_V\n", " \n", " du[1] = -vec_infection + vec_births\n", " du[2] = vec_infection - vec_progression\n", " du[3] = vec_progression\n", " du[1:end] -= vec_mortality\n", "end\n", "\n", "function F(du,u,p,t)\n", " uvec = @view(u[(p.nHosts*H_comps+1):end]) # grab the vector compartments\n", " S_V, E_V, I_V = uvec\n", " sum_β_I_H = 0.0\n", " for host in 0:(p.nHosts-1)\n", " slice = (1:H_comps).+(H_comps*host)\n", " F1H(@view(du[slice]), @view(u[slice]), p.host[host+1], t, I_V, p.vec.N_H)\n", " # must use @view here, so that these arrays can be modified in F1H\n", " sum_β_I_H += p.host[host+1].β * u[slice[3]] # this host's I compartment\n", " end\n", " FV(@view(du[(p.nHosts*4+1):end]), uvec, p.vec,t,sum_β_I_H)\n", "end" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "First, state initial conditions. This code generates them randomly for convenience, though they could be assigned based on data, desired parameter space, or algorithmically as part of a fitting process:" ] }, { "cell_type": "code", "execution_count": 2, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "23×1 Array{Float64,2}:\n", " 100.0\n", " 0.0\n", " 1.0\n", " 0.0\n", " 100.0\n", " 0.0\n", " 0.0\n", " 0.0\n", " 100.0\n", " 0.0\n", " 0.0\n", " 0.0\n", " 100.0\n", " 0.0\n", " 0.0\n", " 0.0\n", " 100.0\n", " 0.0\n", " 0.0\n", " 0.0\n", " 10000.0\n", " 0.0\n", " 0.0" ] }, "execution_count": 2, "metadata": {}, "output_type": "execute_result" } ], "source": [ "nH = 5\n", "srand(0)\n", "\n", "S_Hs = ones(nH) .* 100.0\n", "E_Hs = zeros(nH)\n", "I_Hs = shuffle(vcat(zeros(nH-1),[1.0]))\n", "R_Hs = zeros(nH)\n", "host0 = reshape(hcat(S_Hs,E_Hs,I_Hs,R_Hs)', nH*4, 1)\n", "vec0 = [10000.0, 0.0, 0.0]\n", "u0 = vcat(host0, vec0)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Now, generate dynamic parameters. Again: this code generates them randomly for convenience, though they could be assigned based on data, desired parameter space, or algorithmically as part of a fitting process:" ] }, { "cell_type": "code", "execution_count": 3, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "(nHosts = 5, vec = (μ_V = 0.03333333333333333, σ_V = 0.14285714285714285, N_H = 501.0), host = NamedTuples._NT_μ__H_σ__H_λ_β{Float64,Float64,Float64,Float64}[(μ_H = 0.0117686, σ_H = 0.790008, λ = 0.0642631, β = 0.0209472), (μ_H = 0.00801628, σ_H = 0.175085, λ = 0.0817059, β = 0.0251379), (μ_H = 0.00888301, σ_H = 0.166683, λ = 0.0840894, β = 0.00203749), (μ_H = 0.351205, σ_H = 0.662263, λ = 0.0461889, β = 0.0287702), (μ_H = 0.00568503, σ_H = 0.168919, λ = 0.127011, β = 0.0859512)])" ] }, "execution_count": 3, "metadata": {}, "output_type": "execute_result" } ], "source": [ "srand(1)\n", "\n", "μs = 1 ./ (rand(nH) .* 360)\n", "σs = 1 ./ (rand(nH) .* 6)\n", "λs = 1 ./ (rand(nH) .* 28)\n", "βs = rand(nH) ./ 10.0\n", "\n", "using NamedTuples\n", "# nb: in >= Julia v0.7, can eliminate this import\n", "# and the @NT syntax\n", "p = @NT(\n", " nHosts = nH,\n", " vec = @NT(μ_V=1/30, σ_V=1/7, N_H = sum(host0)),\n", " host = [@NT(μ_H=μs[i], σ_H=σs[i], λ=λs[i], β=βs[i]) for i in 1:nH]\n", " # just building up a random collection of params for demonstration\n", ")" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Now these values can be used with the ODE solver:" ] }, { "cell_type": "code", "execution_count": 4, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "retcode: Success\n", "Interpolation: 1st order linear\n", "t: 3651-element Array{Float64,1}:\n", " 0.0\n", " 0.1\n", " 0.2\n", " 0.3\n", " 0.4\n", " 0.5\n", " 0.6\n", " 0.7\n", " 0.8\n", " 0.9\n", " 1.0\n", " 1.1\n", " 1.2\n", " ⋮ \n", " 363.9\n", " 364.0\n", " 364.1\n", " 364.2\n", " 364.3\n", " 364.4\n", " 364.5\n", " 364.6\n", " 364.7\n", " 364.8\n", " 364.9\n", " 365.0\n", "u: 3651-element Array{Array{Float64,2},1}:\n", " [100.0; 0.0; … ; 0.0; 0.0] \n", " [100.001; 4.05106e-8; … ; 0.041287; 0.000295821]\n", " [100.002; 3.15487e-7; … ; 0.08154; 0.00117209] \n", " [100.004; 1.03672e-6; … ; 0.120779; 0.00261227] \n", " [100.005; 2.39309e-6; … ; 0.159026; 0.00460017] \n", " [100.006; 4.55256e-6; … ; 0.196298; 0.00711996] \n", " [100.007; 7.66377e-6; … ; 0.232616; 0.0101561] \n", " [100.008; 1.18579e-5; … ; 0.267998; 0.0136934] \n", " [100.009; 1.72497e-5; … ; 0.302464; 0.017717] \n", " [100.011; 2.39395e-5; … ; 0.336032; 0.0222123] \n", " [100.012; 3.20139e-5; … ; 0.36872; 0.0271651] \n", " [100.013; 4.15473e-5; … ; 0.400545; 0.0325615] \n", " [100.014; 5.26028e-5; … ; 0.431526; 0.0383877] \n", " ⋮ \n", " [45.2552; 0.790354; … ; 78.2356; 335.113] \n", " [45.2574; 0.790394; … ; 78.2385; 335.114] \n", " [45.2595; 0.790434; … ; 78.2413; 335.114] \n", " [45.2617; 0.790474; … ; 78.2442; 335.115] \n", " [45.2639; 0.790514; … ; 78.2471; 335.116] \n", " [45.2661; 0.790553; … ; 78.25; 335.117] \n", " [45.2682; 0.790593; … ; 78.2528; 335.117] \n", " [45.2704; 0.790633; … ; 78.2557; 335.118] \n", " [45.2725; 0.790673; … ; 78.2586; 335.119] \n", " [45.2747; 0.790713; … ; 78.2615; 335.12] \n", " [45.2768; 0.790753; … ; 78.2643; 335.121] \n", " [45.279; 0.790792; … ; 78.2672; 335.122] " ] }, "execution_count": 4, "metadata": {}, "output_type": "execute_result" } ], "source": [ "using DifferentialEquations\n", "using IterableTables, DataFrames\n", "\n", "tspan = (0.0, 365.0)\n", "prob = ODEProblem(F, u0, tspan, p)\n", "sol = solve(prob,Tsit5(),reltol=1e-8,abstol=1e-8,saveat=linspace(0,365,365*10+1))" ] }, { "cell_type": "code", "execution_count": 5, "metadata": {}, "outputs": [ { "data": { "text/html": [ "
variablevaluetimestampindexnamefacetcompartment
1value1100.00.01S_H1H1S
2value1100.00117612465690.11S_H1H1S
3value1100.002350620499740.21S_H1H1S
4value1100.003523250351390.31S_H1H1S
5value1100.004693784142530.41S_H1H1S
6value1100.005861998757510.51S_H1H1S
7value1100.007027677894070.61S_H1H1S
8value1100.008190611906170.71S_H1H1S
9value1100.009350597680180.81S_H1H1S
10value1100.010507438471120.91S_H1H1S
11value1100.011660943787361.01S_H1H1S
12value1100.012810929250071.11S_H1H1S
13value1100.013957216439921.21S_H1H1S
14value1100.015099632805471.31S_H1H1S
15value1100.016238011509741.41S_H1H1S
16value1100.017372191289471.51S_H1H1S
17value1100.018502016372151.61S_H1H1S
18value1100.019627336345351.71S_H1H1S
19value1100.020748005986771.81S_H1H1S
20value1100.021863885203461.91S_H1H1S
21value1100.022974838928352.01S_H1H1S
22value1100.02408073694862.11S_H1H1S
23value1100.025181453807082.21S_H1H1S
24value1100.026276868728022.31S_H1H1S
25value1100.027366865528522.41S_H1H1S
26value1100.02845133242412.51S_H1H1S
27value1100.029530161968172.61S_H1H1S
28value1100.030603250975522.71S_H1H1S
29value1100.031670500446422.81S_H1H1S
30value1100.032731815377872.91S_H1H1S
" ], "text/plain": [ "83973×7 DataFrames.DataFrame\n", "│ Row │ variable │ value │ timestamp │ index │ name │ facet │\n", "├───────┼──────────┼─────────┼───────────┼───────┼────────┼───────┤\n", "│ 1 │ value1 │ 100.0 │ 0.0 │ 1 │ \"S_H1\" │ \"H1\" │\n", "│ 2 │ value1 │ 100.001 │ 0.1 │ 1 │ \"S_H1\" │ \"H1\" │\n", "│ 3 │ value1 │ 100.002 │ 0.2 │ 1 │ \"S_H1\" │ \"H1\" │\n", "│ 4 │ value1 │ 100.004 │ 0.3 │ 1 │ \"S_H1\" │ \"H1\" │\n", "│ 5 │ value1 │ 100.005 │ 0.4 │ 1 │ \"S_H1\" │ \"H1\" │\n", "│ 6 │ value1 │ 100.006 │ 0.5 │ 1 │ \"S_H1\" │ \"H1\" │\n", "│ 7 │ value1 │ 100.007 │ 0.6 │ 1 │ \"S_H1\" │ \"H1\" │\n", "│ 8 │ value1 │ 100.008 │ 0.7 │ 1 │ \"S_H1\" │ \"H1\" │\n", "│ 9 │ value1 │ 100.009 │ 0.8 │ 1 │ \"S_H1\" │ \"H1\" │\n", "│ 10 │ value1 │ 100.011 │ 0.9 │ 1 │ \"S_H1\" │ \"H1\" │\n", "│ 11 │ value1 │ 100.012 │ 1.0 │ 1 │ \"S_H1\" │ \"H1\" │\n", "⋮\n", "│ 83962 │ value23 │ 335.113 │ 363.9 │ 23 │ \"I_V\" │ \"V\" │\n", "│ 83963 │ value23 │ 335.114 │ 364.0 │ 23 │ \"I_V\" │ \"V\" │\n", "│ 83964 │ value23 │ 335.114 │ 364.1 │ 23 │ \"I_V\" │ \"V\" │\n", "│ 83965 │ value23 │ 335.115 │ 364.2 │ 23 │ \"I_V\" │ \"V\" │\n", "│ 83966 │ value23 │ 335.116 │ 364.3 │ 23 │ \"I_V\" │ \"V\" │\n", "│ 83967 │ value23 │ 335.117 │ 364.4 │ 23 │ \"I_V\" │ \"V\" │\n", "│ 83968 │ value23 │ 335.117 │ 364.5 │ 23 │ \"I_V\" │ \"V\" │\n", "│ 83969 │ value23 │ 335.118 │ 364.6 │ 23 │ \"I_V\" │ \"V\" │\n", "│ 83970 │ value23 │ 335.119 │ 364.7 │ 23 │ \"I_V\" │ \"V\" │\n", "│ 83971 │ value23 │ 335.12 │ 364.8 │ 23 │ \"I_V\" │ \"V\" │\n", "│ 83972 │ value23 │ 335.121 │ 364.9 │ 23 │ \"I_V\" │ \"V\" │\n", "│ 83973 │ value23 │ 335.122 │ 365.0 │ 23 │ \"I_V\" │ \"V\" │\n", "\n", "│ Row │ compartment │\n", "├───────┼─────────────┤\n", "│ 1 │ \"S\" │\n", "│ 2 │ \"S\" │\n", "│ 3 │ \"S\" │\n", "│ 4 │ \"S\" │\n", "│ 5 │ \"S\" │\n", "│ 6 │ \"S\" │\n", "│ 7 │ \"S\" │\n", "│ 8 │ \"S\" │\n", "│ 9 │ \"S\" │\n", "│ 10 │ \"S\" │\n", "│ 11 │ \"S\" │\n", "⋮\n", "│ 83962 │ \"I\" │\n", "│ 83963 │ \"I\" │\n", "│ 83964 │ \"I\" │\n", "│ 83965 │ \"I\" │\n", "│ 83966 │ \"I\" │\n", "│ 83967 │ \"I\" │\n", "│ 83968 │ \"I\" │\n", "│ 83969 │ \"I\" │\n", "│ 83970 │ \"I\" │\n", "│ 83971 │ \"I\" │\n", "│ 83972 │ \"I\" │\n", "│ 83973 │ \"I\" │" ] }, "execution_count": 5, "metadata": {}, "output_type": "execute_result" } ], "source": [ "# rename!(df, Dict(:timestamp => :t,\n", "# :value1 => :S_H, :value2 => :E_H, :value3 => :I_H, :value4 => :R_H,\n", "# :value5 => :S_V, :value6 => :E_V, :value7 => :I_V\n", "# ))\n", "# mlt[:host] = contains.(string.(mlt[:variable]),\"H\"); # tag which entries are host vs vector\n", "# df\n", "df = DataFrame(sol)\n", "mlt = melt(df,:timestamp) # convert results into long format for plotting\n", "mlt[:index] = parse.(Int,replace.(string.(mlt[:variable]),r\"[^\\d]+\"=>\"\"))\n", "mlt[:name] = [\n", " \"S_H1\",\"E_H1\",\"I_H1\",\"R_H1\",\n", " \"S_H2\",\"E_H2\",\"I_H2\",\"R_H2\",\n", " \"S_H3\",\"E_H3\",\"I_H3\",\"R_H3\",\n", " \"S_H4\",\"E_H4\",\"I_H4\",\"R_H4\",\n", " \"S_H5\",\"E_H5\",\"I_H5\",\"R_H5\",\n", " \"S_V\",\"E_V\",\"I_V\"\n", "][mlt[:index]]\n", "mlt[:facet] = replace.(string.(mlt[:name]),r\"\\w+_\"=>\"\")\n", "mlt[:compartment] = replace.(string.(mlt[:name]),r\"_\\w+\"=>\"\")\n", "mlt" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Now that we have a solution, we want to view what is happening in host vs mosquito population:" ] }, { "cell_type": "code", "execution_count": 6, "metadata": { "scrolled": true }, "outputs": [ { "name": "stderr", "output_type": "stream", "text": [ "WARNING: Method definition promote_rule(Type{T}, Type{Any}) in module Missings at /home/simon/.julia/v0.6/Missings/src/Missings.jl:52 overwritten in module Nulls at /home/simon/.julia/v0.6/Nulls/src/Nulls.jl:29.\n", "WARNING: Method definition ==(Base.Nullable{S}, Base.Nullable{T}) in module Base at nullable.jl:238 overwritten in module NullableArrays at /home/simon/.julia/v0.6/NullableArrays/src/operators.jl:99.\n" ] }, { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAbAAAAFoCAMAAAAB/V5aAAAAG1BMVEUAAAAAv8QaGhpNTU18rgDHfP/r6+v4dm3///+XVUguAAAOtUlEQVR4nO2djXqjOAxF03Zmlvd/4m0CBAP+k9G1r4jON9vdnarE8YlAMY76mBxTPEYPwJHhwozhwozhwoyRF/bz+2f558i/8rE/MgRNVtivqFXZ6XtMk8QUgqaUYTth/wL+1vHvE+ghakUkLKTu9Xj013YUSyFowMIOnLwxTfU9hEGKjs0a01TfQViGy0/v5Yxpql1YMSR2VYM80KcKAxRYt6gkUXJiDM2wOaSQZkzpYy3DQjRnIKuMyYYBYT8/z+WOSJmoPANpZ0w2DAhb3oudUZ+BlDEmG1aErQmGvXTbLT/AjnYUhf0EX0MQL9noeZEpfSxk2M/65wBmBuwvN6KpyrCe98NOyphsWBCWAjcDB2NMNlxYNGSfZEw2rAnrVnYZKxhRcmIwZtiTLcuY0icb8j1TPsolWIVtlzIOG3UhaFtTTZX4Wpo6/z1+BprucX69UR1LZQiBsJ+p6/uwPS9l+ZCvI8lvxIJ0h8shbNnq9vqf/pfz9O2yde5bjhq32D7KeapormGD9yVGFj+OWaLyQNmUXB7OQFnPsJF0t5MgdkrrMJa6s2qHc2JF0cGwVXt2lpop1bFU7pA9/3yXMyJxWX8MyezYubCkUqVC8kBorAibUysxoXVHKbnxpSk1djXce7YLP3TWAxsfSk4MAxmWuMLLrjKfufgb0mkG0gVZ/7HYENZp11Scpy0mGwaEJYr6Xu99+jyQZgiaqpWOIe/DBKsLTCFoara5Tf3XEi8t7XUHq2hPefF3wGr9rtBgSh8DGfbT/X7YsS5ksmFAWBLQDIS2/hzoPZamEDRUwpbk2ht6h6T1ubAE0Cv3a/1pdlETftAHHVoJlJwYHBn2um+SOOtVHeV49jzl4WdmWIjaDKyyIA+UFHmZ8qgwVAkDVYnLvdvCk2e6QBnJMN33YcFt9j9UU31zYb9nhcK2lRR//vx3A/RV1CESFlZGz6mXPcVRNRwerKI9/U+JNw9B01x0OGPwFrLGcGHGcGHGcGHG8K7ayiFovKu2cgiavj1/PyAEDboNOnqXNAU9RK30ybDkR0OYcuMeGaZcdBytMU31HYRluPD0AmlMU+3CciGLM6apvqEw3Uv1fYoRlJwYY/d0nKsQ0AN9aIaFaM1A0RmTDQPCOnTVLihjsmFAWJeu2lc/uc8UgkbUVTtEdwYyyphsWBAWdtUGFlqmK0asoj08XbX9Fw1UQdRVO1F8MNmwICwFYgaixphsuLBjSERZ+SiNO5AxPeHQ8Hw+bEZUe6wz3/ZQWZmiI6HkxCDLsPrfDBHmiP5YxOk5EXUkjYO7KJS7RB1njOYaRtDzd0RX7b/Z3wwBu/rcQ9iYrtrJm9LRc5ELC9h11Q4Bz8CmbAuJXzrwY6kMobmGxXdNoTmtVTVXg3jmeXq6Gt+keeC+xCXJ5pBUcnUaS03IK7nGN2keuFX7ve0jeS7sN5aKkO+pR4oRlvUBT2X/pmJPUhZhk7dBfyrzFrIhbEtTR76+DNwpQ8mJwZ1hX6+ev5Et3gPGUhmChlnY19bzN+/MhQUM66r9dej5m3EGGYukK35HWLtqR5fiU1PV/EBSJRYybEhX7fB91/ke52k6Kx6oQspNTokDumrXLEIVTlhnkAPGKtrD11X7uKrBVFEYyLDeXbWJb3bVhaDhKutja4ZMNlzYcR9FlweChqChWZoKK41cr13kGFpByYlBkWHvTTWxDsj7o0T7JHuGVaA2A7OsdK/q5FFEfa4/SxiqSpx3+M3TfaXX7vxAVxubyyg/fQxDWsi+d2MGqiqP0qr1ozLsurDjttm3KNA8yrpHN1EeFYZLXbXreekZVcTBwSra4121lUPQeFdtY3gLWWO4MGO4MGO4MGM0N7h0xuBdtZVD0HhXbeUQNJ26ao/aHtOHHqJWxmfY0R/sgfqEoKH7VR4ndUw26IVlwM9AJuu6j6U+BA2xsC2kpM6FJRh9cT+VOqOHs4CSE8NEhqUo1Z+VDyQrauvO0jBMCyuECCf/HqfEDl217xVygYSLR01QQIeu2ncKuYCeMP+VivUhSR6P15fH/OX172n9v+d/Ln/9WAN33w8OU3qcTl21TVPla7ayfN10zN+bJW1/wu8LM6z758OMh6R4rF8eoZVpTbjlu5vLx/b983Fy9OuqfY+QFElh7/88ZNjUKCwF0yQxhaRInhJTwloz7MLYPzIkSVB0vK2s9cVcVzwegbPHlmGyomM3XidKq8JuP/SE6VXNFCLAhY0NoWkhG4djkshCCJo0j2iDbjeEQNiYNuhWQxiEhW3QR1djrMxTRXMN832JlSGf3gbdWsj3+wuSiqLD1xLrQr6XTuhYvKzXC/m9fnkbdHMhaHwtUYEr83+4oVz+gdZHYnpVM4UIEc+/C1MOEaIvzLe5yUKSxBvpqJ8SR7VBNxsiBHFK3Nqgj764syKddNn8S38gaIMewvSqZgoRon9K9G1ushA0lXvrz99gmiSmEDRe1iuHoPkQYbLujsnfCWhN2OhqLElZgfLh9gdEyYlhNMNqU6LHWPpiSFjaENG5F06VsHFVYl0aubA9Hd+HnS4VWl3mP1OY9MLcQDDJij3LlbSHMBcdvtIhCBECuR8WFeaogLmB6R1JYbQKEy4ZO1q0Z5h0md8RcihbF66cEl3ZANpPiWunCKcrykXHgI6k1kPQZIV5G3R5CJpShml31RbeN7FBD1ErImEhba/HhlaRTOljLcNC1Gag5K63DVFTzP5QFR3R+RAdpfFEXauCPsMyYF/46rP9IafEDEzXDaYQNPfYhDMYlJwYnmHKIWhcmHIImsqt2r7SURsiBPCRWW+DLgoRcuo4WvEDBYIEG31xZ6XF1Db/ysLCNughTK9qppAk+91h69+Kb0LWCPPVekGIEEyG+f2w+hAhgGtYAqZJYgoR4sJGh6DxpSkFUHJieIaphdB0JI3DMUlkIQQ9f72rtiSEQJh31ZaEMAjzrtpl5qmiuYb5vsTKELytybtqa4ZQtEH3rdrVIV3OiF7WK4a0yNJfrU9BMkl0ITLkHxnypSkFRDPePv2NP7HA9KpmCkkS3/Dqp8TRIXK0b2D6rilRiAz9a5h31RaGCEGcEr2rdgm5p3a8q7ZyCBrvqq0cgsa7aiuHoPGyXjkEjQtTDkHjS1MKoOTE8AxTDkHjwpRD0NB21W5owwsbiyQEDUnP38j0U7Ribnt1INEUlnt+oCbc47g+9W10a4M+qoLrAVbRHpJT4n1C0HhXbWN401hjuDBjuDBjuDBjUDW4vEMIGu+qrRyCRqGr9qi3qzz0ELVyvUlzupcr0wv/MzMspDD2mpbH5aMYDEGDLjoK0pim+g7CMgieXkYa01S7sDAk4Yxpqm8o7FotpVFQXrmDdfnBk6DkxOi6p+OcZtVHyd43rDrKFdcRyg+JofcmnIOx+qkGjAUSgqZ7V+19kpWPUvNqdmEh6l21Q2WlozwvPc0PNCQEjairdsiVp7cZyx/llVtMNiwIC7tqq5VVVfWinY07WEV7RnXVXpIsHfK+dDGlj5UMQ9wPm40lQ7ZKg8mGBWEpLj+93IJ+WBky2fhoYS9j8ZBdIc9kw5ow5Wt1qvSwU20soOTEGPpxo19j5788vk9mSp9sCE1H0jgqMxBZwT+ta3DYqAsh6EiK7qp9NHZeh6KxURFCIAzeVXtvLLJuSGOjIoRBWNhVO0RNWGgsts5LY6MUQnMNi29z0+NdK0JvMiJZJup76tBWm6Gr9pJjV+5OMoS8mmqPb9LcYav2fFZM3fbisFERQiIsheYMPI0lb1Ny2KgJ+e7RuJ5C2K+x9G1lEhu1IWhIWhfZ3qGPkhODIsN+y8P09mCm9LGWYSGKM/A6HSaNMdmwIAzfVXu9fCWUMdkwIAzfVXsrN+LGmGwYEDaFXbVDlGZgVx3eYfs9GlFXbf366rAYZbRYxCraM7Sr9nlTb+QzLkzpYyDDgF21n7rOISdjTDYMCEty9el9pbb1HpOMycbnCltPhtEQ6eclqELQDFmaOt/4OrTDNNZMAiUnRv8M231+6G1qawT7/OuXsr9zn1am9LGWYSFNT2/30bzV1D5k8za3jahprysYS6KriJzyQ2Lo11V7k7V2U84dZfW2n6WDO92Jvk+GXXwftn4qOOh6LThK2Db7NPsK3bFbqZg3CAhhx89vZ56kwqta3i5dhfLAMUC6ar/sjCrZBoBVtMe7aiuHoPGu2sbwFrLGcGHGcGHGcGHGiAsLGlumtmo7Y4gKC7ppe1dtYQiaVIa5sMYQNCJh4bv7piXugFGrEghohYVcfD1Wr6FTpQ9xhhWLDtUZyLmrPMrVjJdSHhUGhr31O0AzePMMq4DpNMQUgobk82G2QcmJ4RmmHILGhSmHoEkJe3fT9l/4JgtBk86wn3fbIl/pEISgyQrz39AnD0GTFLZ0004uTTkbfVTN5IT54m9DCJpshnVbmrpRCBov65VD0Lgw5RA0vjSlAEpODM8w5RA0Lkw5BE26SgyWpbxKrA4Z1vM3XJby92F1IR1aNGcXf12YLKSDrSl7DcvvmnI21hnrkWKZU6LvS2wJGda3vuuuqRuFoPGyXjkEjQtTDkHjS1MKoOTE8AxTDkHjwpRD0Lgw5RA0hbVE3+YmDUFTeOPsu6akIWiyS1PHbW6jqzFWOrl6kduXOPnSlDwETXq13m+vNIWgKeyt97VEaQgaL+uVQ9C4MOUQNL6WqABKTgzPMOUQNNzCKtsku7ADvarEqJTkUbQbY99HGOp92HkKG1tZN45F03Y3mrtqX3q2ET01F3dYT/MGmIsOX+kQhKBpFuaMwdugG8N7/hrDhRnDhRnDhRmjVZhCHaLSYz23T6j3UbrQKEyh0s/1RRKNg+UofRgnbNIQNilNtdJYOmBeGJP2HhgXFu4VujKS2wtjKTp+tIoOjQKoC17WG8OFGcOFGcOFGcOFGcOFGcOqMKvjvozVJ2513Jcx+sQfD6MDv4zV52113Jex+sStjvsyVp+41XFfxuoTtzruy1h94l50ODZwYcZwYcZwYcZwYcZwYcZwYcZwYcZwYcb4HwtcFz0kasMDAAAAAElFTkSuQmCC" }, "metadata": {}, "output_type": "display_data" }, { "data": { "text/plain": [ "RCall.RObject{RCall.VecSxp}\n" ] }, "execution_count": 6, "metadata": {}, "output_type": "execute_result" } ], "source": [ "using RCall\n", "# current version RCall supports better transfers, which would simplify this mess\n", "# but requires Julia v >= 0.7\n", "vals = mlt[:value]\n", "tstamps = mlt[:timestamp]\n", "fcts = mlt[:facet]\n", "comps = mlt[:compartment]\n", "@rput vals tstamps fcts comps\n", "R\"\n", "library(ggplot2)\n", "suppressPackageStartupMessages(library(data.table))\n", "dt <- data.table(t=tstamps, y=vals, species=fcts, compartment=comps)\n", "ggplot(dt) + aes(x=t, y=y, color=compartment) + facet_grid(species ~ ., scale = 'free_y') +\n", " theme_minimal() +\n", " geom_line()\n", "\"" ] } ], "metadata": { "kernelspec": { "display_name": "Julia 0.6.3", "language": "julia", "name": "julia-0.6" }, "language_info": { "file_extension": ".jl", "mimetype": "application/julia", "name": "julia", "version": "0.6.3" } }, "nbformat": 4, "nbformat_minor": 2 }