835 lines
24 KiB
Plaintext
835 lines
24 KiB
Plaintext
|
# Euler's method
|
|||
|
|
|||
|
This section uses these add-on packages:
|
|||
|
|
|||
|
```julia
|
|||
|
using CalculusWithJulia
|
|||
|
using Plots
|
|||
|
using SymPy
|
|||
|
using Roots
|
|||
|
```
|
|||
|
|
|||
|
```julia; echo=false; results="hidden"
|
|||
|
using CalculusWithJulia.WeaveSupport
|
|||
|
|
|||
|
|
|||
|
const frontmatter = (
|
|||
|
title = "Euler's method",
|
|||
|
description = "Calculus with Julia: Euler's method",
|
|||
|
tags = ["CalculusWithJulia", "odes", "euler's method"],
|
|||
|
);
|
|||
|
fig_size = (600, 400)
|
|||
|
nothing
|
|||
|
```
|
|||
|
|
|||
|
----
|
|||
|
|
|||
|
The following section takes up the task of numerically approximating solutions to differential equations. `Julia` has a huge set of state-of-the-art tools for this task starting with the [DifferentialEquations](https://github.com/SciML/DifferentialEquations.jl) package. We don't use that package in this section, focusing on simpler methods and implementations for pedagogical purposes, but any further exploration should utilize the tools provided therein. A brief introduction to the package follows in an upcoming [section](./differential_equations.html).
|
|||
|
|
|||
|
|
|||
|
----
|
|||
|
|
|||
|
|
|||
|
|
|||
|
Consider the differential equation:
|
|||
|
|
|||
|
```math
|
|||
|
y'(x) = y(x) \cdot x, \quad y(1)=1,
|
|||
|
```
|
|||
|
|
|||
|
which can be solved with `SymPy`:
|
|||
|
|
|||
|
```julia;
|
|||
|
@syms x, y, u()
|
|||
|
D = Differential(x)
|
|||
|
x0, y0 = 1, 1
|
|||
|
F(y,x) = y*x
|
|||
|
|
|||
|
dsolve(D(u)(x) - F(u(x), x))
|
|||
|
```
|
|||
|
|
|||
|
With the given initial condition, the solution becomes:
|
|||
|
|
|||
|
```julia;
|
|||
|
out = dsolve(D(u)(x) - F(u(x),x), u(x), ics=Dict(u(x0) => y0))
|
|||
|
```
|
|||
|
|
|||
|
|
|||
|
Plotting this solution over the slope field
|
|||
|
|
|||
|
```julia;
|
|||
|
p = plot(legend=false)
|
|||
|
vectorfieldplot!((x,y) -> [1, F(x,y)], xlims=(0, 2.5), ylims=(0, 10))
|
|||
|
plot!(rhs(out), linewidth=5)
|
|||
|
```
|
|||
|
|
|||
|
|
|||
|
we see that the vectors that are drawn seem to be tangent to the graph
|
|||
|
of the solution. This is no coincidence, the tangent lines to integral
|
|||
|
curves are in the direction of the slope field.
|
|||
|
|
|||
|
|
|||
|
What if the graph of the solution were not there, could we use this
|
|||
|
fact to *approximately* reconstruct the solution?
|
|||
|
|
|||
|
That is, if we stitched together pieces of the slope field, would we
|
|||
|
get a curve that was close to the actual answer?
|
|||
|
|
|||
|
```julia; hold=true; echo=false; cache=true
|
|||
|
## {{{euler_graph}}}
|
|||
|
function make_euler_graph(n)
|
|||
|
x, y = symbols("x, y")
|
|||
|
F(y,x) = y*x
|
|||
|
x0, y0 = 1, 1
|
|||
|
|
|||
|
h = (2-1)/5
|
|||
|
xs = zeros(n+1)
|
|||
|
ys = zeros(n+1)
|
|||
|
xs[1] = x0 # index is off by 1
|
|||
|
ys[1] = y0
|
|||
|
for i in 1:n
|
|||
|
xs[i + 1] = xs[i] + h
|
|||
|
ys[i + 1] = ys[i] + h * F(ys[i], xs[i])
|
|||
|
end
|
|||
|
|
|||
|
p = plot(legend=false)
|
|||
|
vectorfieldplot!((x,y) -> [1, F(y,x)], xlims=(1,2), ylims=(0,6))
|
|||
|
|
|||
|
## Add Euler soln
|
|||
|
plot!(p, xs, ys, linewidth=5)
|
|||
|
scatter!(p, xs, ys)
|
|||
|
|
|||
|
## add function
|
|||
|
out = dsolve(u'(x) - F(u(x), x), u(x), ics=(u, x0, y0))
|
|||
|
plot!(p, rhs(out), x0, xs[end], linewidth=5)
|
|||
|
|
|||
|
p
|
|||
|
end
|
|||
|
|
|||
|
|
|||
|
|
|||
|
|
|||
|
n = 5
|
|||
|
anim = @animate for i=1:n
|
|||
|
make_euler_graph(i)
|
|||
|
end
|
|||
|
|
|||
|
imgfile = tempname() * ".gif"
|
|||
|
gif(anim, imgfile, fps = 1)
|
|||
|
|
|||
|
|
|||
|
caption = """
|
|||
|
Illustration of a function stitching together slope field lines to
|
|||
|
approximate the answer to an initial-value problem. The other function drawn is the actual solution.
|
|||
|
"""
|
|||
|
|
|||
|
ImageFile(imgfile, caption)
|
|||
|
```
|
|||
|
|
|||
|
The illustration suggests the answer is yes, let's see. The solution
|
|||
|
is drawn over $x$ values $1$ to $2$. Let's try piecing together $5$
|
|||
|
pieces between $1$ and $2$ and see what we have.
|
|||
|
|
|||
|
The slope-field vectors are *scaled* versions of the vector `[1, F(y,x)]`. The `1`
|
|||
|
is the part in the direction of the $x$ axis, so here we would like
|
|||
|
that to be $0.2$ (which is $(2-1)/5$. So our vectors would be `0.2 *
|
|||
|
[1, F(y,x)]`. To allow for generality, we use `h` in place of the
|
|||
|
specific value $0.2$.
|
|||
|
|
|||
|
Then our first pieces would be the line connecting $(x_0,y_0)$ to
|
|||
|
|
|||
|
```math
|
|||
|
\langle x_0, y_0 \rangle + h \cdot \langle 1, F(y_0, x_0) \rangle.
|
|||
|
```
|
|||
|
|
|||
|
The above uses vector notation to add the piece scaled by $h$ to the
|
|||
|
starting point. Rather than continue with that notation, we will use
|
|||
|
subscripts. Let $x_1$, $y_1$ be the postion of the tip of the
|
|||
|
vector. Then we have:
|
|||
|
|
|||
|
```math
|
|||
|
x_1 = x_0 + h, \quad y_1 = y_0 + h F(y_0, x_0).
|
|||
|
```
|
|||
|
|
|||
|
With this notation, it is easy to see what comes next:
|
|||
|
|
|||
|
```math
|
|||
|
x_2 = x_1 + h, \quad y_2 = y_1 + h F(y_1, x_1).
|
|||
|
```
|
|||
|
|
|||
|
We just shifted the indices forward by $1$. But graphically what is
|
|||
|
this? It takes the tip of the first part of our "stitched" together
|
|||
|
solution, finds the slope filed there (`[1, F(y,x)]`) and then uses
|
|||
|
this direction to stitch together one more piece.
|
|||
|
|
|||
|
Clearly, we can repeat. The $n$th piece will end at:
|
|||
|
|
|||
|
```math
|
|||
|
x_{n+1} = x_n + h, \quad y_{n+1} = y_n + h F(y_n, x_n).
|
|||
|
```
|
|||
|
|
|||
|
For our example, we can do some numerics. We want $h=0.2$ and $5$
|
|||
|
pieces, so values of $y$ at $x_0=1, x_1=1.2, x_2=1.4, x_3=1.6,
|
|||
|
x_4=1.8,$ and $x_5=2$.
|
|||
|
|
|||
|
Below we do this in a loop. We have to be a bit careful, as in `Julia`
|
|||
|
the vector of zeros we create to store our answers begins indexing at
|
|||
|
$1$, and not $0$.
|
|||
|
|
|||
|
```julia;
|
|||
|
n=5
|
|||
|
h = (2-1)/n
|
|||
|
xs = zeros(n+1)
|
|||
|
ys = zeros(n+1)
|
|||
|
xs[1] = x0 # index is off by 1
|
|||
|
ys[1] = y0
|
|||
|
for i in 1:n
|
|||
|
xs[i + 1] = xs[i] + h
|
|||
|
ys[i + 1] = ys[i] + h * F(ys[i], xs[i])
|
|||
|
end
|
|||
|
```
|
|||
|
|
|||
|
So how did we do? Let's look graphically:
|
|||
|
|
|||
|
```julia;
|
|||
|
plot(exp(-1/2)*exp(x^2/2), x0, 2)
|
|||
|
plot!(xs, ys)
|
|||
|
```
|
|||
|
|
|||
|
Not bad. We wouldn't expect this to be exact - due to the concavity
|
|||
|
of the solution, each step is an underestimate. However, we see it is
|
|||
|
an okay approximation and would likely be better with a smaller $h$. A
|
|||
|
topic we pursue in just a bit.
|
|||
|
|
|||
|
Rather than type in the above command each time, we wrap it all up in
|
|||
|
a function. The inputs are $n$, $a=x_0$, $b=x_n$, $y_0$, and, most
|
|||
|
importantly, $F$. The output is massaged into a function through a
|
|||
|
call to `linterp`, rather than two vectors. The `linterp` function we define below just
|
|||
|
finds a function that linearly interpolates between the points and is
|
|||
|
`NaN` outside of the range of the $x$ values:
|
|||
|
|
|||
|
```julia;
|
|||
|
function linterp(xs, ys)
|
|||
|
function(x)
|
|||
|
((x < xs[1]) || (x > xs[end])) && return NaN
|
|||
|
for i in 1:(length(xs) - 1)
|
|||
|
if xs[i] <= x < xs[i+1]
|
|||
|
l = (x-xs[i]) / (xs[i+1] - xs[i])
|
|||
|
return (1-l) * ys[i] + l * ys[i+1]
|
|||
|
end
|
|||
|
end
|
|||
|
ys[end]
|
|||
|
end
|
|||
|
end
|
|||
|
```
|
|||
|
|
|||
|
With that, here is our function to find an approximate solution to $y'=F(y,x)$ with initial condition:
|
|||
|
|
|||
|
```julia;
|
|||
|
function euler(F, x0, xn, y0, n)
|
|||
|
h = (xn - x0)/n
|
|||
|
xs = zeros(n+1)
|
|||
|
ys = zeros(n+1)
|
|||
|
xs[1] = x0
|
|||
|
ys[1] = y0
|
|||
|
for i in 1:n
|
|||
|
xs[i + 1] = xs[i] + h
|
|||
|
ys[i + 1] = ys[i] + h * F(ys[i], xs[i])
|
|||
|
end
|
|||
|
linterp(xs, ys)
|
|||
|
end
|
|||
|
```
|
|||
|
|
|||
|
With `euler`, it becomes easy to explore different values.
|
|||
|
|
|||
|
For example, we thought the solution would look better with a smaller $h$ (or larger $n$). Instead of $n=5$, let's try $n=50$:
|
|||
|
|
|||
|
```julia;
|
|||
|
u₁₂ = euler(F, 1, 2, 1, 50)
|
|||
|
plot(exp(-1/2)*exp(x^2/2), x0, 2)
|
|||
|
plot!(u₁₂, x0, 2)
|
|||
|
```
|
|||
|
|
|||
|
It is more work for the computer, but not for us, and clearly a much better approximation to the actual answer is found.
|
|||
|
|
|||
|
|
|||
|
## The Euler method
|
|||
|
|
|||
|
|
|||
|
```julia; hold=true; echo=false
|
|||
|
imgfile ="figures/euler.png"
|
|||
|
caption = """
|
|||
|
Figure from first publication of Euler's method. From [Gander and Wanner](http://www.unige.ch/~gander/Preprints/Ritz.pdf).
|
|||
|
"""
|
|||
|
|
|||
|
ImageFile(:ODEs, imgfile, caption)
|
|||
|
```
|
|||
|
|
|||
|
|
|||
|
The name of our function reflects the [mathematician](https://en.wikipedia.org/wiki/Leonhard_Euler) associated with the iteration:
|
|||
|
|
|||
|
```math
|
|||
|
x_{n+1} = x_n + h, \quad y_{n+1} = y_n + h \cdot F(y_n, x_n),
|
|||
|
```
|
|||
|
|
|||
|
to approximate a solution to the first-order, ordinary differential
|
|||
|
equation with initial values: $y'(x) = F(y,x)$.
|
|||
|
|
|||
|
|
|||
|
[The Euler method](https://en.wikipedia.org/wiki/Euler_method) uses
|
|||
|
linearization. Each "step" is just an approximation of the function
|
|||
|
value $y(x_{n+1})$ with the value from the tangent line tangent to the
|
|||
|
point $(x_n, y_n)$.
|
|||
|
|
|||
|
|
|||
|
Each step introduces an error. The error in one step is known as the
|
|||
|
*local truncation error* and can be shown to be about equal to $1/2
|
|||
|
\cdot h^2 \cdot f''(x_{n})$ assuming $y$ has ``3`` or more derivatives.
|
|||
|
|
|||
|
The total error, or more commonly, *global truncation error*, is the
|
|||
|
error between the actual answer and the approximate answer at the end
|
|||
|
of the process. It reflects an accumulation of these local errors. This
|
|||
|
error is *bounded* by a constant times $h$. Since it gets smaller as
|
|||
|
$h$ gets smaller in direct proportion, the Euler method is called
|
|||
|
*first order*.
|
|||
|
|
|||
|
Other, somewhat more complicated, methods have global truncation errors that
|
|||
|
involve higher powers of $h$ - that is for the same size $h$, the
|
|||
|
error is smaller. In analogy is the fact that Riemann sums have
|
|||
|
error that depends on $h$, whereas other methods of approximating the
|
|||
|
integral have smaller errors. For example, Simpson's rule had error
|
|||
|
related to $h^4$. So, the Euler method may not be employed if there
|
|||
|
is concern about total resources (time, computer, ...), it is
|
|||
|
important for theoretical purposes in a manner similar to the role of the Riemann
|
|||
|
integral.
|
|||
|
|
|||
|
In the examples, we will see that for many problems the simple Euler
|
|||
|
method is satisfactory, but not always so. The task of numerically
|
|||
|
solving differential equations is not a one-size-fits-all one. In the
|
|||
|
following, a few different modifications are presented to the basic
|
|||
|
Euler method, but this just scratches the surface of the topic.
|
|||
|
|
|||
|
#### Examples
|
|||
|
|
|||
|
##### Example
|
|||
|
|
|||
|
|
|||
|
Consider the initial value problem $y'(x) = x + y(x)$ with initial
|
|||
|
condition $y(0)=1$. This problem can be solved exactly. Here we
|
|||
|
approximate over $[0,2]$ using Euler's method.
|
|||
|
|
|||
|
```julia;
|
|||
|
𝑭(y,x) = x + y
|
|||
|
𝒙0, 𝒙n, 𝒚0 = 0, 2, 1
|
|||
|
𝒇 = euler(𝑭, 𝒙0, 𝒙n, 𝒚0, 25)
|
|||
|
𝒇(𝒙n)
|
|||
|
```
|
|||
|
|
|||
|
We graphically compare our approximate answer with the exact one:
|
|||
|
|
|||
|
```julia;
|
|||
|
plot(𝒇, 𝒙0, 𝒙n)
|
|||
|
𝒐ut = dsolve(D(u)(x) - 𝑭(u(x),x), u(x), ics = Dict(u(𝒙0) => 𝒚0))
|
|||
|
plot(rhs(𝒐ut), 𝒙0, 𝒙n)
|
|||
|
plot!(𝒇, 𝒙0, 𝒙n)
|
|||
|
```
|
|||
|
|
|||
|
From the graph it appears our value for `f(xn)` will underestimate the
|
|||
|
actual value of the solution slightly.
|
|||
|
|
|||
|
##### Example
|
|||
|
|
|||
|
The equation $y'(x) = \sin(x \cdot y)$ is not separable, so need not have an
|
|||
|
easy solution. The default method will fail. Looking at the available methods with `sympy.classify_ode(𝐞qn, u(x))` shows a power series method which
|
|||
|
can return a power series *approximation* (a Taylor polynomial). Let's
|
|||
|
look at comparing an approximate answer given by the Euler method to
|
|||
|
that one returned by `SymPy`.
|
|||
|
|
|||
|
First, the `SymPy` solution:
|
|||
|
|
|||
|
```julia;
|
|||
|
𝐅(y,x) = sin(x*y)
|
|||
|
𝐞qn = D(u)(x) - 𝐅(u(x), x)
|
|||
|
𝐨ut = dsolve(𝐞qn, hint="1st_power_series")
|
|||
|
```
|
|||
|
|
|||
|
|
|||
|
|
|||
|
If we assume $y(0) = 1$, we can continue:
|
|||
|
|
|||
|
```julia;
|
|||
|
𝐨ut1 = dsolve(𝐞qn, u(x), ics=Dict(u(0) => 1), hint="1st_power_series")
|
|||
|
```
|
|||
|
|
|||
|
The approximate value given by the Euler method is
|
|||
|
|
|||
|
```julia;
|
|||
|
𝐱0, 𝐱n, 𝐲0 = 0, 2, 1
|
|||
|
|
|||
|
plot(legend=false)
|
|||
|
vectorfieldplot!((x,y) -> [1, 𝐅(y,x)], xlims=(𝐱0, 𝐱n), ylims=(0,5))
|
|||
|
plot!(rhs(𝐨ut1).removeO(), linewidth=5)
|
|||
|
|
|||
|
𝐮 = euler(𝐅, 𝐱0, 𝐱n, 𝐲0, 10)
|
|||
|
plot!(𝐮, linewidth=5)
|
|||
|
```
|
|||
|
|
|||
|
We see that the answer found from using a polynomial series matches that of Euler's method for a bit, but as time evolves, the approximate solution given by Euler's method more closely tracks the slope field.
|
|||
|
|
|||
|
##### Example
|
|||
|
|
|||
|
|
|||
|
The
|
|||
|
[Brachistochrone problem](http://www.unige.ch/~gander/Preprints/Ritz.pdf)
|
|||
|
was posed by Johann Bernoulli in 1696. It asked for the curve between
|
|||
|
two points for which an object will fall faster along that curve than
|
|||
|
any other. For an example, a bead sliding on a wire will take a certain amount of time to get from point $A$ to point $B$, the time depending on the shape of the wire. Which shape will take the least amount of time?
|
|||
|
|
|||
|
|
|||
|
```julia; hold=true; echo=false
|
|||
|
imgfile = "figures/bead-game.jpg"
|
|||
|
caption = """
|
|||
|
|
|||
|
A child's bead game. What shape wire will produce the shortest time for a bed to slide from a top to the bottom?
|
|||
|
|
|||
|
"""
|
|||
|
ImageFile(:ODEs, imgfile, caption)
|
|||
|
```
|
|||
|
|
|||
|
Restrict our attention to the $x$-$y$ plane, and consider a path,
|
|||
|
between the point $(0,A)$ and $(B,0)$. Let $y(x)$ be the distance from
|
|||
|
$A$, so $y(0)=0$ and at the end $y$ will be $A$.
|
|||
|
|
|||
|
|
|||
|
[Galileo](http://www-history.mcs.st-and.ac.uk/HistTopics/Brachistochrone.html)
|
|||
|
knew the straight line was not the curve, but incorrectly thought the
|
|||
|
answer was a part of a circle.
|
|||
|
|
|||
|
```julia; hold=true; echo=false
|
|||
|
imgfile = "figures/galileo.gif"
|
|||
|
caption = """
|
|||
|
As early as 1638, Galileo showed that an object falling along `AC` and then `CB` will fall faster than one traveling along `AB`, where `C` is on the arc of a circle.
|
|||
|
From the [History of Math Archive](http://www-history.mcs.st-and.ac.uk/HistTopics/Brachistochrone.html).
|
|||
|
"""
|
|||
|
ImageFile(:ODEs, imgfile, caption)
|
|||
|
```
|
|||
|
|
|||
|
|
|||
|
This simulation also suggests that a curved path is better than the shorter straight one:
|
|||
|
|
|||
|
```julia; hold=true; echo=false; cache=true
|
|||
|
##{{{brach_graph}}}
|
|||
|
|
|||
|
function brach(f, x0, vx0, y0, vy0, dt, n)
|
|||
|
m = 1
|
|||
|
g = 9.8
|
|||
|
|
|||
|
axs = Float64[0]
|
|||
|
ays = Float64[-g]
|
|||
|
vxs = Float64[vx0]
|
|||
|
vys = Float64[vy0]
|
|||
|
xs = Float64[x0]
|
|||
|
ys = Float64[y0]
|
|||
|
|
|||
|
for i in 1:n
|
|||
|
x = xs[end]
|
|||
|
vx = vxs[end]
|
|||
|
|
|||
|
ax = -f'(x) * (f''(x) * vx^2 + g) / (1 + f'(x)^2)
|
|||
|
ay = f''(x) * vx^2 + f'(x) * ax
|
|||
|
|
|||
|
push!(axs, ax)
|
|||
|
push!(ays, ay)
|
|||
|
|
|||
|
push!(vxs, vx + ax * dt)
|
|||
|
push!(vys, vys[end] + ay * dt)
|
|||
|
push!(xs, x + vxs[end] * dt)# + (1/2) * ax * dt^2)
|
|||
|
push!(ys, ys[end] + vys[end] * dt)# + (1/2) * ay * dt^2)
|
|||
|
end
|
|||
|
|
|||
|
[xs ys vxs vys axs ays]
|
|||
|
|
|||
|
end
|
|||
|
|
|||
|
|
|||
|
fs = [x -> 1 - x,
|
|||
|
x -> (x-1)^2,
|
|||
|
x -> 1 - sqrt(1 - (x-1)^2),
|
|||
|
x -> - (x-1)*(x+1),
|
|||
|
x -> 3*(x-1)*(x-1/3)
|
|||
|
]
|
|||
|
|
|||
|
|
|||
|
MS = [brach(f, 1/100, 0, 1, 0, 1/100, 100) for f in fs]
|
|||
|
|
|||
|
|
|||
|
function make_brach_graph(n)
|
|||
|
|
|||
|
p = plot(xlim=(0,1), ylim=(-1/3, 1), legend=false)
|
|||
|
for (i,f) in enumerate(fs)
|
|||
|
plot!(f, 0, 1)
|
|||
|
U = MS[i]
|
|||
|
x = min(1.0, U[n,1])
|
|||
|
scatter!(p, [x], [f(x)])
|
|||
|
end
|
|||
|
p
|
|||
|
|
|||
|
end
|
|||
|
|
|||
|
|
|||
|
|
|||
|
n = 4
|
|||
|
anim = @animate for i=[1,5,10,15,20,25,30,35,40,45,50,55,60]
|
|||
|
make_brach_graph(i)
|
|||
|
end
|
|||
|
|
|||
|
imgfile = tempname() * ".gif"
|
|||
|
gif(anim, imgfile, fps = 1)
|
|||
|
|
|||
|
|
|||
|
caption = """
|
|||
|
The race is on. An illustration of beads falling along a path, as can be seen, some paths are faster than others. The fastest path would follow a cycloid. See [Bensky and Moelter](https://pdfs.semanticscholar.org/66c1/4d8da6f2f5f2b93faf4deb77aafc7febb43a.pdf) for details on simulating a bead on a wire.
|
|||
|
"""
|
|||
|
|
|||
|
ImageFile(imgfile, caption)
|
|||
|
```
|
|||
|
|
|||
|
|
|||
|
|
|||
|
|
|||
|
Now, the natural question is which path is best? The solution can be
|
|||
|
[reduced](http://mathworld.wolfram.com/BrachistochroneProblem.html) to
|
|||
|
solving this equation for a positive $c$:
|
|||
|
|
|||
|
```math
|
|||
|
1 + (y'(x))^2 = \frac{c}{y}, \quad c > 0.
|
|||
|
```
|
|||
|
|
|||
|
Reexpressing, this becomes:
|
|||
|
|
|||
|
```math
|
|||
|
\frac{dy}{dx} = \sqrt{\frac{C-y}{y}}.
|
|||
|
```
|
|||
|
|
|||
|
This is a separable equation and can be solved, but even `SymPy` has
|
|||
|
trouble with this integral. However, the result has been known to be a piece of a cycloid since the insightful
|
|||
|
Jacob Bernoulli used an analogy from light bending to approach the problem. The answer is best described parametrically
|
|||
|
through:
|
|||
|
|
|||
|
```math
|
|||
|
x(u) = c\cdot u - \frac{c}{2}\sin(2u), \quad y(u) = \frac{c}{2}( 1- \cos(2u)), \quad 0 \leq u \leq U.
|
|||
|
```
|
|||
|
|
|||
|
The values of $U$ and $c$ must satisfy $(x(U), y(U)) = (B, A)$.
|
|||
|
|
|||
|
|
|||
|
Rather than pursue this, we will solve it numerically for a fixed
|
|||
|
value of $C$ over a fixed interval to see the shape.
|
|||
|
|
|||
|
|
|||
|
The equation can be written in terms of $y'=F(y,x)$, where
|
|||
|
|
|||
|
```math
|
|||
|
F(y,x) = \sqrt{\frac{c-y}{y}}.
|
|||
|
```
|
|||
|
|
|||
|
But as $y_0 = 0$, we immediately would have a problem with the first step, as there would be division by $0$.
|
|||
|
|
|||
|
This says that for the optimal solution, the bead picks up speed by first sliding straight down before heading off towards $B$. That's great for the physics, but runs roughshod over our Euler method, as the first step has an infinity.
|
|||
|
|
|||
|
For this, we can try the *backwards Euler* method which uses the slope at $(x_{n+1}, y_{n+1})$, rather than $(x_n, y_n)$. The update step becomes:
|
|||
|
|
|||
|
```math
|
|||
|
y_{n+1} = y_n + h \cdot F(y_{n+1}, x_{n+1}).
|
|||
|
```
|
|||
|
|
|||
|
Seems innocuous, but the value we are trying to find, $y_{n+1}$, is
|
|||
|
now on both sides of the equation, so is only *implicitly* defined. In
|
|||
|
this code, we use the `find_zero` function from the `Roots` package. The
|
|||
|
caveat is, this function needs a good initial guess, and the one we
|
|||
|
use below need not be widely applicable.
|
|||
|
|
|||
|
|
|||
|
```julia;
|
|||
|
function back_euler(F, x0, xn, y0, n)
|
|||
|
h = (xn - x0)/n
|
|||
|
xs = zeros(n+1)
|
|||
|
ys = zeros(n+1)
|
|||
|
xs[1] = x0
|
|||
|
ys[1] = y0
|
|||
|
for i in 1:n
|
|||
|
xs[i + 1] = xs[i] + h
|
|||
|
## solve y[i+1] = y[i] + h * F(y[i+1], x[i+1])
|
|||
|
ys[i + 1] = find_zero(y -> ys[i] + h * F(y, xs[i + 1]) - y, ys[i]+h)
|
|||
|
end
|
|||
|
linterp(xs, ys)
|
|||
|
end
|
|||
|
```
|
|||
|
|
|||
|
We then have with $C=1$ over the interval $[0,1.2]$ the following:
|
|||
|
|
|||
|
```julia;
|
|||
|
𝐹(y, x; C=1) = sqrt(C/y - 1)
|
|||
|
𝑥0, 𝑥n, 𝑦0 = 0, 1.2, 0
|
|||
|
cyc = back_euler(𝐹, 𝑥0, 𝑥n, 𝑦0, 50)
|
|||
|
plot(x -> 1 - cyc(x), 𝑥0, 𝑥n)
|
|||
|
```
|
|||
|
|
|||
|
Remember, $y$ is the displacement from the top, so it is
|
|||
|
non-negative. Above we flipped the graph to make it look more like
|
|||
|
expectation. In general, the trajectory may actually dip below the
|
|||
|
ending point and come back up. The above won't see this, for as
|
|||
|
written $dy/dx \geq 0$, which need not be the case, as the defining
|
|||
|
equation is in terms of $(dy/dx)^2$, so the derivative could have any
|
|||
|
sign.
|
|||
|
|
|||
|
|
|||
|
|
|||
|
##### Example: stiff equations
|
|||
|
|
|||
|
The Euler method is *convergent*, in that as $h$ goes to $0$, the
|
|||
|
approximate solution will converge to the actual answer. However, this
|
|||
|
does not say that for a fixed size $h$, the approximate value will be
|
|||
|
good. For example, consider the differential equation $y'(x) =
|
|||
|
-5y$. This has solution $y(x)=y_0 e^{-5x}$. However, if we try the
|
|||
|
Euler method to get an answer over $[0,2]$ with $h=0.5$ we don't see
|
|||
|
this:
|
|||
|
|
|||
|
```julia;
|
|||
|
ℱ(y,x) = -5y
|
|||
|
𝓍0, 𝓍n, 𝓎0 = 0, 2, 1
|
|||
|
𝓊 = euler(ℱ, 𝓍0, 𝓍n, 𝓎0, 4) # n =4 => h = 2/4
|
|||
|
vectorfieldplot((x,y) -> [1, ℱ(y,x)], xlims=(0, 2), ylims=(-5, 5))
|
|||
|
plot!(x -> y0 * exp(-5x), 0, 2, linewidth=5)
|
|||
|
plot!(𝓊, 0, 2, linewidth=5)
|
|||
|
```
|
|||
|
|
|||
|
What we see is that the value of $h$ is too big to capture the decay
|
|||
|
scale of the solution. A smaller $h$, can do much better:
|
|||
|
|
|||
|
```julia;
|
|||
|
𝓊₁ = euler(ℱ, 𝓍0, 𝓍n, 𝓎0, 50) # n=50 => h = 2/50
|
|||
|
plot(x -> y0 * exp(-5x), 0, 2)
|
|||
|
plot!(𝓊₁, 0, 2)
|
|||
|
```
|
|||
|
|
|||
|
This is an example of a
|
|||
|
[stiff equation](https://en.wikipedia.org/wiki/Stiff_equation). Such
|
|||
|
equations cause explicit methods like the Euler one problems, as small
|
|||
|
$h$s are needed to good results.
|
|||
|
|
|||
|
The implicit, backward Euler method does not have this issue, as we can see here:
|
|||
|
|
|||
|
```julia;
|
|||
|
𝓊₂ = back_euler(ℱ, 𝓍0, 𝓍n, 𝓎0, 4) # n =4 => h = 2/4
|
|||
|
vectorfieldplot((x,y) -> [1, ℱ(y,x)], xlims=(0, 2), ylims=(-1, 1))
|
|||
|
plot!(x -> y0 * exp(-5x), 0, 2, linewidth=5)
|
|||
|
plot!(𝓊₂, 0, 2, linewidth=5)
|
|||
|
```
|
|||
|
|
|||
|
|
|||
|
##### Example: The pendulum
|
|||
|
|
|||
|
|
|||
|
The differential equation describing the simple pendulum is
|
|||
|
|
|||
|
```math
|
|||
|
\theta''(t) = - \frac{g}{l}\sin(\theta(t)).
|
|||
|
```
|
|||
|
|
|||
|
The typical approach to solving for $\theta(t)$ is to use the small-angle approximation that $\sin(x) \approx x$, and then the differential equation simplifies to:
|
|||
|
$\theta''(t) = -g/l \cdot \theta(t)$, which is easily solved.
|
|||
|
|
|||
|
Here we try to get an answer numerically. However, the problem, as stated, is not a first order equation due to the $\theta''(t)$ term. If we let $u(t) = \theta(t)$ and $v(t) = \theta'(t)$, then we get *two* coupled first order equations:
|
|||
|
|
|||
|
```math
|
|||
|
v'(t) = -g/l \cdot \sin(u(t)), \quad u'(t) = v(t).
|
|||
|
```
|
|||
|
|
|||
|
We can try the Euler method here. A simple approach might be this iteration scheme:
|
|||
|
|
|||
|
```math
|
|||
|
\begin{align*}
|
|||
|
x_{n+1} &= x_n + h,\\
|
|||
|
u_{n+1} &= u_n + h v_n,\\
|
|||
|
v_{n+1} &= v_n - h \cdot g/l \cdot \sin(u_n).
|
|||
|
\end{align*}
|
|||
|
```
|
|||
|
|
|||
|
Here we need *two* initial conditions: one for the initial value
|
|||
|
$u(t_0)$ and the initial value of $u'(t_0)$. We have seen if we start at an angle $a$ and release the bob from rest, so $u'(0)=0$ we get a sinusoidal answer to the linearized model. What happens here? We let $a=1$, $L=5$ and $g=9.8$:
|
|||
|
|
|||
|
We write a function to solve this starting from $(x_0, y_0)$ and ending at $x_n$:
|
|||
|
|
|||
|
```julia;
|
|||
|
function euler2(x0, xn, y0, yp0, n; g=9.8, l = 5)
|
|||
|
xs, us, vs = zeros(n+1), zeros(n+1), zeros(n+1)
|
|||
|
xs[1], us[1], vs[1] = x0, y0, yp0
|
|||
|
h = (xn - x0)/n
|
|||
|
for i = 1:n
|
|||
|
xs[i+1] = xs[i] + h
|
|||
|
us[i+1] = us[i] + h * vs[i]
|
|||
|
vs[i+1] = vs[i] + h * (-g / l) * sin(us[i])
|
|||
|
end
|
|||
|
linterp(xs, us)
|
|||
|
end
|
|||
|
```
|
|||
|
|
|||
|
Let's take $a = \pi/4$ as the initial angle, then the approximate
|
|||
|
solution should be $\pi/4\cos(\sqrt{g/l}x)$ with period $T =
|
|||
|
2\pi\sqrt{l/g}$. We try first to plot then over 4 periods:
|
|||
|
|
|||
|
```julia;
|
|||
|
𝗅, 𝗀 = 5, 9.8
|
|||
|
𝖳 = 2pi * sqrt(𝗅/𝗀)
|
|||
|
𝗑0, 𝗑n, 𝗒0, 𝗒p0 = 0, 4𝖳, pi/4, 0
|
|||
|
plot(euler2(𝗑0, 𝗑n, 𝗒0, 𝗒p0, 20), 0, 4𝖳)
|
|||
|
```
|
|||
|
|
|||
|
Something looks terribly amiss. The issue is the step size, $h$, is
|
|||
|
too large to capture the oscillations. There are basically only $5$
|
|||
|
steps to capture a full up and down motion. Instead, we try to get $20$ steps per period
|
|||
|
so $n$ must be not $20$, but $4 \cdot 20 \cdot T \approx 360$. To this
|
|||
|
graph, we add the approximate one:
|
|||
|
|
|||
|
```julia;
|
|||
|
plot(euler2(𝗑0, 𝗑n, 𝗒0, 𝗒p0, 360), 0, 4𝖳)
|
|||
|
plot!(x -> pi/4*cos(sqrt(𝗀/𝗅)*x), 0, 4𝖳)
|
|||
|
```
|
|||
|
|
|||
|
Even now, we still see that something seems amiss, though the issue is
|
|||
|
not as dramatic as before. The oscillatory nature of the pendulum is
|
|||
|
seen, but in the Euler solution, the amplitude grows, which would
|
|||
|
necessarily mean energy is being put into the system. A familiar
|
|||
|
instance of a pendulum would be a child on a swing. Without pumping
|
|||
|
the legs - putting energy in the system - the height of the swing's
|
|||
|
arc will not grow. Though we now have oscillatory motion, this growth
|
|||
|
indicates the solution is still not quite right. The issue is likely
|
|||
|
due to each step mildly overcorrecting and resulting in an overall
|
|||
|
growth. One of the questions pursues this a bit further.
|
|||
|
|
|||
|
## Questions
|
|||
|
|
|||
|
##### Question
|
|||
|
|
|||
|
Use Euler's method with $n=5$ to approximate $u(1)$ where
|
|||
|
|
|||
|
```math
|
|||
|
u'(x) = x - u(x), \quad u(0) = 1
|
|||
|
```
|
|||
|
|
|||
|
```julia; hold=true; echo=false
|
|||
|
F(y,x) = x - y
|
|||
|
x0, xn, y0 = 0, 1, 1
|
|||
|
val = euler(F, x0, xn, y0, 5)(1)
|
|||
|
numericq(val)
|
|||
|
```
|
|||
|
|
|||
|
##### Question
|
|||
|
|
|||
|
Consider the equation
|
|||
|
|
|||
|
```math
|
|||
|
y' = x \cdot \sin(y), \quad y(0) = 1.
|
|||
|
```
|
|||
|
|
|||
|
Use Euler's method with $n=50$ to find the value of $y(5)$.
|
|||
|
|
|||
|
```julia; hold=true; echo=false
|
|||
|
F(y, x) = x * sin(y)
|
|||
|
x0, xn, y0 = 0, 5, 1
|
|||
|
n = 50
|
|||
|
u = euler(F, x0, xn, y0, n)
|
|||
|
numericq(u(xn))
|
|||
|
```
|
|||
|
|
|||
|
|
|||
|
##### Question
|
|||
|
|
|||
|
Consider the ordinary differential equation
|
|||
|
|
|||
|
```math
|
|||
|
\frac{dy}{dx} = 1 - 2\frac{y}{x}, \quad y(1) = 0.
|
|||
|
```
|
|||
|
|
|||
|
Use Euler's method to solve for $y(2)$ when $n=50$.
|
|||
|
|
|||
|
```julia; hold=true; echo=false
|
|||
|
F(y, x) = 1 - 2y/x
|
|||
|
x0, xn, y0 = 1, 2, 0
|
|||
|
n = 50
|
|||
|
u = euler(F, x0, xn, y0, n)
|
|||
|
numericq(u(xn))
|
|||
|
```
|
|||
|
|
|||
|
##### Question
|
|||
|
|
|||
|
|
|||
|
Consider the ordinary differential equation
|
|||
|
|
|||
|
```math
|
|||
|
\frac{dy}{dx} = \frac{y \cdot \log(y)}{x}, \quad y(2) = e.
|
|||
|
```
|
|||
|
|
|||
|
Use Euler's method to solve for $y(3)$ when $n=25$.
|
|||
|
|
|||
|
```julia; hold=true; echo=false
|
|||
|
F(y, x) = y*log(y)/x
|
|||
|
x0, xn, y0 = 2, 3, exp(1)
|
|||
|
n = 25
|
|||
|
u = euler(F, x0, xn, y0, n)
|
|||
|
numericq(u(xn))
|
|||
|
```
|
|||
|
|
|||
|
|
|||
|
##### Question
|
|||
|
|
|||
|
Consider the first-order non-linear ODE
|
|||
|
|
|||
|
```math
|
|||
|
y' = y \cdot (1-2x), \quad y(0) = 1.
|
|||
|
```
|
|||
|
|
|||
|
Use Euler's method with $n=50$ to approximate the solution $y$ over $[0,2]$.
|
|||
|
|
|||
|
What is the value at $x=1/2$?
|
|||
|
|
|||
|
```julia; hold=true; echo=false
|
|||
|
F(y, x) = y * (1-2x)
|
|||
|
x0, xn, y0 = 0, 2, 1
|
|||
|
n = 50
|
|||
|
u = euler(F, x0, xn, y0, n)
|
|||
|
numericq(u(1/2))
|
|||
|
```
|
|||
|
|
|||
|
What is the value at $x=3/2$?
|
|||
|
|
|||
|
```julia; hold=true; echo=false
|
|||
|
F(y, x) = y * (1-2x)
|
|||
|
x0, xn, y0 = 0, 2, 1
|
|||
|
n = 50
|
|||
|
u = euler(F, x0, xn, y0, n)
|
|||
|
numericq(u(3/2))
|
|||
|
```
|
|||
|
|
|||
|
##### Question: The pendulum revisited.
|
|||
|
|
|||
|
The issue with the pendulum's solution growing in amplitude can be
|
|||
|
addressed using a modification to the Euler method attributed to
|
|||
|
[Cromer](http://astro.physics.ncsu.edu/urca/course_files/Lesson14/index.html). The
|
|||
|
fix is to replace the term `sin(us[i])` in the line `vs[i+1] = vs[i] + h * (-g / l) *
|
|||
|
sin(us[i])` of the `euler2` function with `sin(us[i+1])`, which uses the updated angular
|
|||
|
velocity in the ``2``nd step in place of the value before the step.
|
|||
|
|
|||
|
Modify the `euler2` function to implement the Euler-Cromer method. What do you see?
|
|||
|
|
|||
|
```julia; hold=true; echo=false
|
|||
|
choices = [
|
|||
|
"The same as before - the amplitude grows",
|
|||
|
"The solution is identical to that of the approximation found by linearization of the sine term",
|
|||
|
"The solution has a constant amplitude, but its period is slightly *shorter* than that of the approximate solution found by linearization",
|
|||
|
"The solution has a constant amplitude, but its period is slightly *longer* than that of the approximate solution found by linearization"]
|
|||
|
ans = 4
|
|||
|
radioq(choices, ans, keep_order=true)
|
|||
|
```
|