熱方程式の数値解:有限差分法
📂数値解析熱方程式の数値解:有限差分法
数値的な熱方程式の解
以下のような一次元の熱方程式が与えられたとしよう。
∂t∂u=∂x2∂2u,0≤x≤1,t≥0(1)
目標は、有限の点でこのソリューションを近似することだ。 次のように時空間ドメインを分割しよう。
{(ℓΔx,nΔt):ℓ=0,1,…,d+1, n≥0} where Δx=d+11
u(ℓΔx,nΔt)の近似を uℓnと表記する。 上付き文字の nが冪数ではなく、時間のインデックスであることに注意しよう。 (1)の右辺と左辺は、それぞれ次のように近似される。
∂x2∂2u(x,t)=(Δx)21[u(x−Δx,t)−2u(x,t)+u(x+Δx,t)]+O((Δx)2)
∂t∂u(x,t)=Δt1[u(x,t+Δt)−u(x,t)]+O(Δt)
これで、μ=(Δx)2Δtとしよう。 この二つの式を (1)に代入して、Δt=μ(Δx)2を掛けると次を得る。
Δtu(x,t+Δt)−u(x,t)=(Δx)21[u(x−Δx,t)−2u(x,t)+u(x+Δx,t)]+O((Δx)2)
⟹u(x,t+Δt)=u(x,t)+(Δx)2μ(Δx)2[u(x−Δx,t)−2u(x,t)+u(x+Δx,t)]+O((Δx)4)
⟹u(x,t+Δt)=u(x,t)+μ[u(x−Δx,t)−2u(x,t)+u(x+Δx,t)]+O((Δx)4)
⟹uℓn+1=uℓn+μ(uℓ−1n−2uℓn+uℓ+1n),ℓ=1,2,…,d, n=0,1,…(2)
ここで、空間間隔と時間間隔の比 μ=(Δx)2Δtは、クーラン数Courant numberと呼ばれる重要な定数であり、(2)の収束を決定する。
説明
以下の定理から、(2)は μ≤21の時に収束する。これは、空間間隔に応じて、一度に進められるタイムステップの上限が存在することを意味する。つまり、空間ドメインを細かく分割するほど、一度に進められるタイムステップもそれだけ小さくなる。単純な例では問題ないかもしれないが、時間ドメインの長さが長い場合、計算時間が長くなる原因になり得る。例えば、空間ドメインを大きさ100のベクトルで分割した場合、時間ドメインはなんと1733個に分割される。
収束性
μ≤21であれば、数値解(2)は収束する。
証明
t∗>0を任意の定数としよう。各点での近似値と実際のソリューションのエラーを次のように定義しよう。
eℓn:=uℓn−u(ℓΔx,nΔt),ℓ=0,1,…,d+1,n=0,1,…,nΔt
ここで、nΔt=⌊t∗/Δt⌋=⌊t∗/(μ(Δx)2)⌋はnの終点である。 方法が収束するとは、エラーeℓnに対して次の式が成立することとして表現できる。
Δx→0lim[n=0,1,…,nΔtmax(ℓ=0,1,…,d+1max∣eℓn∣)]=0
括弧内の値を次のようにしよう。
ηn:=ℓ=0,1,…,d+1max∣eℓn∣,n=0,1,…,nΔt(3)
その場合、収束に対する表現は、
Δx→0lim(n=0,1,…,nΔtmaxηn)=0(4)
今、u~ℓn=u(ℓΔx,nΔt)をuℓnの真値としよう。
uℓn+1u~ℓn+1=uℓn+μ(uℓ−1n−2uℓn+uℓ+1n)=u~ℓn+μ(u~ℓ−1n−2u~ℓn+u~ℓ+1n)+O((Δx)4)ℓ=0,1,…,d+1n=0,1,…,nΔt
これを引くと、次を得る。
eℓn+1=eℓn+μ(eℓ−1n−2eℓn+eℓ+1n)+O((Δx)4)≤eℓn+μ(eℓ−1n−2eℓn+eℓ+1n)+c(Δx)4,c>0
不等式は漸近記法の定義によって成立する。三角不等式により、次を得る。
eℓn+1≤eℓn+μ(eℓ−1n−2eℓn+eℓ+1n)+c(Δx)4≤eℓn+μ(eℓ−1n−2eℓn+eℓ+1n)+c(Δx)4≤μeℓ−1n+∣1−2μ∣∣eℓn∣+μeℓ+1n+c(Δx)4≤(2μ+∣1−2μ∣)ηn+c(Δx)4,n=0,1,…,nΔt−1
最後の不等式は(3)によって成立する。μ≤1/2であるので、(2μ+∣1−2μ∣)=1であり、再び(3)から次を得る。
ηn+1=ℓ=0,1,…,d+1maxeℓn+1≤ηn+c(Δx)4,n=0,1,…,nΔt−1
これを帰納的に繰り返すと、
ηn+1≤ηn+c(Δx)4≤ηn−1+2c(Δx)4≤ηn−2+3c(Δx)4≤⋯
⟹ηn≤η0+nc(Δx)4,n=0,1,…,nΔt
初期値がη0=0であり、n(Δx)2=nΔt/μ=t∗/μであるので、
ηn≤μct∗(Δx)2,n=0,1,…,nΔt
したがって、すべてのnに対して、
ηn→0 as Δx→0
これと(4)を組み合わせると証明終わり。
■
例
これをJuliaで次のように実装できる。
1次元
今、次のような1次元の熱方程式を考えよう。
Domain:Heat equation:Exact solution:Ω=[−1,1]×[0,0.35]ut−uxx=0u(x,0)=sin(πx)u(−1,t)=0,u(1,t)=0u(x,t)=sin(πx)e−π2t
ここで、Δx=991−(−1)そしてΔt=μ(Δ)2=0.5(Δx)2として設定し、(2)で解くと、次のような結果が得られる。

一方、Δt=μ(Δ)2=0.51(Δx)2として解くと、解が正しく得られないことが分かる。

using Plots
using LinearAlgebra
function Solver_Heat_1d(Δt, T, x, f, g)
t = Vector(0:Δt:T)
n = length(t)
ℓ = length(x)
Δx = x[2]-x[1]
μ = Δt/((Δx)^2)
μ = Δt/((Δx)^2)
U = zeros(ℓ,n)
U[1,:] .= g[1]
U[end,:] .= g[end]
U[:,1] = f
for i ∈ 2:n
U[2:end-1, i] = (1-2μ)*U[2:end-1, i-1] + μ*(U[1:end-2, i-1] + U[3:end, i-1])
end
return U
end
function ref_sol(x, t, exact_sol)
xt =[kron(ones(length(t)),x);; kron(t,ones(length(x)))]
ref_U = reshape(exact_sol(xt),length(x),length(t))
return ref_U
end
function results(t, x, U, ref_U, μ)
h1 = heatmap(t, x, U, colormap=:viridis, clim=(-1,1))
title!("FDM solution")
h2 = heatmap(t, x, ref_U, colormap=:viridis, clim=(-1,1))
title!("exact solution")
h3 = heatmap(t, x, U - ref_U, colormap=:viridis)
title!("error")
h4 = heatmap(framestyle=:none)
annotate!([0.5],[0.5], text("μ = $(round(μ, digits=3))", 30))
plot(h1, h2, h3, h4, size=(1500,750))
end
for μ ∈ [0.5, 0.51]
x = Vector(LinRange(-1., 1, 100))
Δx = x[2]-x[1]
initial_f = sin.(π*x)
Δt = μ*((Δx)^2)
T = 0.35
bdry_g = [0.0, 0.0]
t = Vector(0:Δt:T)
exact_sol(xt) = sin.(π*xt[:,1]) .* exp.(- (π)^2 * xt[:,2])
ref_u = ref_sol(x, t, exact_sol)
U = Solver_Heat_1d(Δt, T, x, initial_f, bdry_g)
results(t, x, U, ref_u, μ)
savefig("mu=" * string(μ) * ".png")
end
2次元
Domain:Heat equation:Ω=[−2,2]×[−2,2]×[0,1]ut−uxx=0u(x,y,0)=41e−x2−y2u(x,y,t)=0 on ∂([−2,2]×[−2,2])
同じ方法で上の微分方程式を解くと、

using Plots
cd = @__DIR__
function Solver_Heat_2d(Δt, T, x, f)
t = Vector(0:Δt:T)
n = length(t)
ℓ = length(x)
Δx = x[2]-x[1]
μ = Δt/((Δx)^2)
U = zeros(ℓ,ℓ,n)
U[:,:,1] = f
for i ∈ 2:n
U[2:end-1, 2:end-1, i] = (1-4μ)*U[2:end-1, 2:end-1, i-1] + μ*(U[1:end-2, 2:end-1, i-1] + U[3:end, 2:end-1, i-1] + U[2:end-1, 1:end-2, i-1] + U[2:end-1, 3:end, i-1])
end
return U
end
x = range(-2., stop=2, length=100)'
y = range(-2., stop=2, length=100)
initial_f = @. (1/4)exp(-x^2) * exp(-y^2)
Δx = x[2]-x[1]
Δt = 0.49*((Δx)^2)/2
T = 1.0
U = Solver_Heat_2d(Δt, T, x, initial_f)
anim = @animate for i ∈ range(1, stop=size(U)[3], step=100)
surface(U[:,:,i], zlims=(0,0.5), camera=(30,30), colorbar=:none, showaxis = false, grid=false, size=(500,500))
end
gif(anim, cd*"/2d.gif", fps=10)
環境
- OS: Windows11
- Version: Julia v1.8.3, Plots v1.38.6