Решение системы диффренциальных уравнений при различных начальных условиях¶
В этом примере рассматриваются два подхода к решению дифференциального уравнения (д.у.) при различных начальных условиях (н.у.):
- используя цикл
for
с изменением н.у. - используя
EnsembleProblem
для решения в многопоточном режиме
Описание задачи¶
В качестве д.у. рассмотрим уравнение Лотки-Вольтерры, описывающего модель "хищник-жертва":
$$\begin{cases} \frac{d x}{dt} = \alpha x - \beta x y \\ \frac{d y}{dt} = \delta xy - \gamma y \\ \end{cases}$$
где $x$ и $y$ - популяции жертвы и хищника, соответственно
$t$ - время
$\alpha, \beta, \delta, \gamma$ - постоянные параметры, описывающие взаимодействие между представителями видов.
Для определённости возьмём их равными $\alpha = \gamma = 1, \beta =0.01, \delta = 0.02$.
В поставленной задаче н.у. являются размеры популяции видов. Решая д.у. мы сможем наблюдать изменение популяции во времени.
Рассмотрим, какие пакеты нужны для решения нашей задачи:
OrdinaryDiffEq
- пакет, позволяющий решать обыкновенные д.у., который может работать без подключения пакетаDifferentialEquations
Plots
имеет возможность отобразить решение дифференциального уравнения буквально одной строчкойDistributed
позволит разместить функции на разных ядрах процессора, чтобы повысить эффективность вычисления
using OrdinaryDiffEq, Plots
Обыкновенное д.у. имеет вид $$\frac{du}{dt}=f(u,,p,t)$$
Для решения задачи необходимо определить следующие исходные данные:
f
- саму функцию (структуру д.у.)
u0
- начальные условия
tspan
- промежуток времени, на котором отыскивается решение
p
- параметры системы, которые также определяют её вид
kwargs
- ключевые слова, передаваемые в решение
Есть два способа задать функцию f
:
function lotka(du, u, p, t)
du[1] = p[1] * u[1] - p[2] * u[1] * u[2]
du[2] = p[4] * u[1] * u[2] - p[3] * u[2]
end
α = 1; β = 0.01; γ = 1; δ = 0.02;
p = [α, β, γ, δ]
tspan = (0.0, 6.5)
u0 = [50; 50]
prob = ODEProblem(lotka, u0, tspan, p)
После того, как мы поставили задачу, перейдём к её решению. Функция solve
позволяет решать д.у. с применением различных решателей (solver). В нашем случае имеем ОДУ, которое можно решать с постоянным шагом. В качестве параметров будет передаваться поставленная проблема, решатель Tsit5()
и параметр saveat=0.01
, указывающий функции шаг, с которым необходимо сохранять решение.
sol = solve(prob, Tsit5(), saveat=0.01)
Теперь отобразим наше решение, с автоматически указанными легендой и подписью оси t
plot(sol)
Теперь отобразим фазовую траекторию - изменение фазовых координат (x
и y
) от времени. На рисунке начальные положения в оттенках красного, а в конце решения - синие.
scatter(
sol[1, :], sol[2, :], # x и y
legend=false, marker_z=sol.t, color=:RdYlBu, # marker_z позволяет изменять цвет маркера от t
xlabel="Жертва x", ylabel="Хищник y", # RdYlBu - от Красного до Синего через Жёлтый
title="Фазовая траектория системы")
Теперь, когда мы смогли изобразить фазовую траекторию для одной пары н.у., построим фазовый портрет - совокупность этих тракторий.
Для этого будем использовать цикл for
, постепенно увеличивая начальный размер популяции хищника.
tspan=(0., 15.)
initial_conditions = range(10, stop=400, length=40)
plt = plot()
for y0 in initial_conditions
u0 = [50; y0]
prob = ODEProblem(lotka, u0, tspan, p)
sol = solve(prob, Tsit5(), saveat=0.01)
plot!(plt, sol[1, :], sol[2, :], title="Фазовый портрет системы", xlabel="Жертва", ylabel="Хищник" ,legend=false)
end
display(plt)
Многопоточное построение фазового портрета¶
Теперь рассмотрим другой вариант решения задачи построения фазового портрета. Для этого подключим стандартный пакет Distributed
, позволяющий выполнять задачи распределённо. И будем распределять нашу задачу на 4 потока (процесса).
using Distributed
N_PROCS = 4 # зададим количество процессов
addprocs(N_PROCS)
Подключив DifferentialEquations
- пакет для решения д.у. различных типов, Воспользуемся функцией EnsembleProb
, позволяющей задать изменение начальных условий посредством вызова функции prob_func
.
Например, изменять начальные условия можно случайно:
function prob_func(prob, i, repeat)
@. prob.u0 = randn() * prob.u0 prob
end
либо же использовать внешний контейнер данных:
function prob_func(prob, i, repeat)
remake(prob, u0=[50, initial_conditions[i]])
end
using DifferentialEquations
u0 = [50; 10];
prob = ODEProblem(lotka, u0, tspan, p);
function prob_func(prob, i, repeat)
remake(prob, u0=[50, initial_conditions[i]])
end
ensemble_prob = EnsembleProblem(prob, prob_func=prob_func)
Поставив задачу, получим решение, используя solve
. Но в этот раз нам вернёnся не одно решение д.у., а целый набор симуляции sim
. Заметим, что в solve
добавились аргументы:
EnsembleThreads
() - означающий, что распределение будет по потокам процессора, а не, например, по разным компьютерам кластера
trajectories=N_TRAJECTS
- количество симуляций, которые мы хотим провести
N_TRAJECTS = length(initial_conditions)
sim = solve(ensemble_prob, Tsit5(), EnsembleThreads(), trajectories=N_TRAJECTS, saveat=0.01);
Отрисуем фазовый портрет, пройдясь по каждому решению в наборе симуляции
plt = plot()
for sol in sim.u
plot!(plt, sol[1, :], sol[2, :], title="Фазовый портрет системы", xlabel="Жертва", ylabel="Хищник" ,legend=false)
end
display(plt)
Сравнение времени выполнения¶
Теперь сравним время выполнения кода в цикле и с использованием многопоточности.
Для этого подключим библиотеку BenchmarkTools
, зададим наши функции, которые мы хоти сравнить, и при помощи макроса @btime
наши функции прогонятся несколько раз и будет вычисленно их среднее время выполнения.
using BenchmarkTools
function for_loop_solution()
for y0 in range(10, stop=400, length=20)
u0 = [50; y0]
prob = ODEProblem(lotka, u0, tspan, p)
solve(prob, Tsit5(), saveat=0.01)
end
end
function ensemble_prob_solution()
prob = ODEProblem(lotka, u0, tspan, p)
initial_conditions = range(10, stop=400, length=20)
function prob_func(prob, i, repeat)
remake(prob, u0=[50, initial_conditions[i]])
end
ensemble_prob = EnsembleProblem(prob, prob_func=prob_func)
solve(ensemble_prob, Tsit5(), EnsembleThreads(), trajectories=20, saveat=0.01)
end
@btime for_loop_solution()
@btime ensemble_prob_solution()
[]; #заглушка чтобы не выводить EnsembleProblem в ячейку вывода
Заключение¶
Исходя из результатов оценки времени выполнения, ensemble_prob_solution
оказался быстрее в 4 раза, что соответствует количеству ядер, выделенных под задачу.
Из чего можно сделать вывод, что при решении системы дифференциальных уравнений с различными начальными условиями циклу for
стоит предпочитать EnsemleProblem
и распределённое вычисление, с использованием пакета Distributed
.