Skip to content
This repository was archived by the owner on Jul 19, 2023. It is now read-only.

Commit 45f79f2

Browse files
Merge pull request #365 from tinosulzer/issue-352-order-choice
#352 allow higher-order central difference; add MOL_2D test with Dirichlet
2 parents f2b26cc + a1dcf61 commit 45f79f2

4 files changed

Lines changed: 97 additions & 23 deletions

File tree

src/MOLFiniteDifference/MOL_discretization.jl

Lines changed: 37 additions & 22 deletions
Original file line numberDiff line numberDiff line change
@@ -15,6 +15,7 @@ function SciMLBase.symbolic_discretize(pdesys::ModelingToolkit.PDESystem,discret
1515
pdeeqs = pdesys.eqs isa Vector ? pdesys.eqs : [pdesys.eqs]
1616
t = discretization.time
1717
nottime = filter(x->x.val != t.val,pdesys.indvars)
18+
nspace = length(nottime)
1819

1920
# Discretize space
2021
space = map(nottime) do x
@@ -33,13 +34,13 @@ function SciMLBase.symbolic_discretize(pdesys::ModelingToolkit.PDESystem,discret
3334
depvars = map(pdesys.depvars) do u
3435
[Num(Variable{Symbolics.FnType{Tuple{Any}, Real}}(Base.nameof(ModelingToolkit.operation(u.val)),II.I...))(t) for II in indices]
3536
end
36-
spacevals = map(y->[Pair(nottime[i],space[i][y.I[i]]) for i in 1:length(nottime)],indices)
37+
spacevals = map(y->[Pair(nottime[i],space[i][y.I[i]]) for i in 1:nspace],indices)
3738

3839

3940
### INITIAL AND BOUNDARY CONDITIONS ###
4041
# Build symbolic maps for boundaries
41-
edges = reduce(vcat,[[vcat([Colon() for j in 1:i-1],1,[Colon() for j in i+1:length(nottime)]),
42-
vcat([Colon() for j in 1:i-1],size(depvars[1],i),[Colon() for j in i+1:length(nottime)])] for i in 1:length(nottime)])
42+
edges = reduce(vcat,[[vcat([Colon() for j in 1:i-1],1,[Colon() for j in i+1:nspace]),
43+
vcat([Colon() for j in 1:i-1],size(depvars[1],i),[Colon() for j in i+1:nspace])] for i in 1:nspace])
4344

4445
#edgeindices = [indices[e...] for e in edges]
4546
edgevals = reduce(vcat,[[nottime[i]=>first(space[i]),nottime[i]=>last(space[i])] for i in 1:length(space)])
@@ -51,11 +52,11 @@ function SciMLBase.symbolic_discretize(pdesys::ModelingToolkit.PDESystem,discret
5152
subvar(depvar) = substitute.((depvar,),edgevals)
5253
depvarmaps = reduce(vcat,[subvar(depvar) .=> edgevars[i] for (i, depvar) in enumerate(pdesys.depvars)])
5354
# depvarderivmaps will dictate what to replace the Differential terms with
54-
if length(nottime) == 1
55+
if nspace == 1
5556
# 1D system
56-
left_weights(j) = DiffEqOperators.calculate_weights(discretization.upwind_order, 0.0, space[j][1:2])
57-
right_weights(j) = DiffEqOperators.calculate_weights(discretization.upwind_order, 0.0, space[j][end-1:end])
58-
central_neighbor_idxs(i,j) = [i+CartesianIndex([ifelse(l==j,-1,0) for l in 1:length(nottime)]...),i,i+CartesianIndex([ifelse(l==j,1,0) for l in 1:length(nottime)]...)]
57+
left_weights(j) = DiffEqOperators.calculate_weights(discretization.upwind_order, space[j][1], space[j][1:2])
58+
right_weights(j) = DiffEqOperators.calculate_weights(discretization.upwind_order, space[j][end], space[j][end-1:end])
59+
central_neighbor_idxs(i,j) = [i+CartesianIndex([ifelse(l==j,-1,0) for l in 1:nspace]...),i,i+CartesianIndex([ifelse(l==j,1,0) for l in 1:nspace]...)]
5960
left_idxs = central_neighbor_idxs(CartesianIndex(2),1)[1:2]
6061
right_idxs(j) = central_neighbor_idxs(CartesianIndex(length(space[j])-1),1)[end-1:end]
6162
# Constructs symbolic spatially discretized terms of the form e.g. au₂ - bu₁
@@ -86,7 +87,7 @@ function SciMLBase.symbolic_discretize(pdesys::ModelingToolkit.PDESystem,discret
8687

8788
# Replace Differential terms in the bc lhs with the symbolic spatially discretized terms
8889
# TODO: Fix Neumann and Robin on higher dimension
89-
lhs = length(nottime) == 1 ? substitute(bc.lhs,depvarderivmaps[i]) : bc.lhs
90+
lhs = nspace == 1 ? substitute(bc.lhs,depvarderivmaps[i]) : bc.lhs
9091

9192
# Replace symbol in the bc lhs with the spatial discretized term
9293
lhs = substitute(lhs,depvarmaps[i])
@@ -101,28 +102,42 @@ function SciMLBase.symbolic_discretize(pdesys::ModelingToolkit.PDESystem,discret
101102
### PDE EQUATIONS ###
102103
interior = indices[[2:length(s)-1 for s in space]...]
103104
eqs = vec(map(Base.product(interior,pdeeqs)) do p
104-
i,eq = p
105-
106-
# TODO: Number of points in the central_neighbor_idxs should be dependent
107-
# on discretization.centered_order
105+
II,eq = p
106+
107+
# Create a stencil in the required dimension centered around 0
108+
# e.g. (-1,0,1) for 2nd order, (-2,-1,0,1,2) for 4th order, etc
109+
order = discretization.centered_order
110+
stencil_indices(j) = [1:ifelse(l==j,order+1,1) for l in 1:nspace]
111+
stencil_shift(j) = [ifelse(l==j,order÷2+1,1) for l in 1:nspace]
112+
stencil(j) = CartesianIndices(Tuple(stencil_indices(j))) .- CartesianIndex(stencil_shift(j)...)
113+
108114
# TODO: Generalize central difference handling to allow higher even order derivatives
109-
central_neighbor_idxs(i,j) = [i+CartesianIndex([ifelse(l==j,-1,0) for l in 1:length(nottime)]...),i,i+CartesianIndex([ifelse(l==j,1,0) for l in 1:length(nottime)]...)]
110-
central_weights(i,j) = DiffEqOperators.calculate_weights(2, 0.0, space[j][i[j]-1:i[j]+1])
111-
central_deriv(i,j,k) = dot(central_weights(i,j),depvars[k][central_neighbor_idxs(i,j)])
112-
central_deriv_rules = [(Differential(s)^2)(pdesys.depvars[k]) => central_deriv(i,j,k) for (j,s) in enumerate(nottime), k in 1:length(pdesys.depvars)]
113-
valrules = vcat([pdesys.depvars[k] => depvars[k][i] for k in 1:length(pdesys.depvars)],
114-
[nottime[j] => space[j][i[j]] for j in 1:length(nottime)])
115-
115+
# The central neighbour indices should add the stencil to II, unless II is too close
116+
# to an edge in which case we need to shift away from the edge
117+
# Calculate buffers
118+
I1 = oneunit(first(indices))
119+
Imin = first(indices) + I1 * (order÷2)
120+
Imax = last(indices) - I1 * (order÷2)
121+
# Use max and min to apply buffers
122+
central_neighbor_idxs(II,j) = stencil(j) .+ max(Imin,min(II,Imax))
123+
central_neighbor_space(II,j) = vec(space[j][map(i->i[j],central_neighbor_idxs(II,j))])
124+
central_weights(II,j) = DiffEqOperators.calculate_weights(2, space[j][II[j]], central_neighbor_space(II,j))
125+
central_deriv(II,j,k) = dot(central_weights(II,j),depvars[k][central_neighbor_idxs(II,j)])
126+
127+
central_deriv_rules = [(Differential(s)^2)(u) => central_deriv(II,j,k) for (j,s) in enumerate(nottime), (k,u) in enumerate(pdesys.depvars)]
128+
valrules = vcat([pdesys.depvars[k] => depvars[k][II] for k in 1:length(pdesys.depvars)],
129+
[nottime[j] => space[j][II[j]] for j in 1:nspace])
130+
116131
# TODO: Use rule matching for nonlinear Laplacian
117-
132+
118133
# TODO: upwind rules needs interpolation into `@rule`
119134
#forward_weights(i,j) = DiffEqOperators.calculate_weights(discretization.upwind_order, 0.0, [space[j][i[j]],space[j][i[j]+1]])
120135
#reverse_weights(i,j) = DiffEqOperators.calculate_weights(discretization.upwind_order, 0.0, [space[j][i[j]-1],space[j][i[j]]])
121136
#upwinding_rules = [@rule(*(~~a,(Differential(nottime[j]))(u),~~b) => IfElse.ifelse(*(~~a..., ~~b...,)>0,
122137
# *(~~a..., ~~b..., dot(reverse_weights(i,j),depvars[k][central_neighbor_idxs(i,j)[1:2]])),
123138
# *(~~a..., ~~b..., dot(forward_weights(i,j),depvars[k][central_neighbor_idxs(i,j)[2:3]]))))
124-
# for j in 1:length(nottime), k in 1:length(pdesys.depvars)]
125-
139+
# for j in 1:nspace, k in 1:length(pdesys.depvars)]
140+
126141
rules = vcat(vec(central_deriv_rules),valrules)
127142
substitute(eq.lhs,rules) ~ substitute(eq.rhs,rules)
128143
end)
Lines changed: 7 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -36,14 +36,18 @@ using ModelingToolkit: Differential
3636
discretization = MOLFiniteDifference([x=>dx],t)
3737
# Explicitly specify order of centered difference
3838
discretization_centered = MOLFiniteDifference([x=>dx],t;centered_order=order)
39+
# Higher order centered difference
40+
discretization_centered_order4 = MOLFiniteDifference([x=>dx],t;centered_order=4)
3941

4042
# Convert the PDE problem into an ODE problem
4143
prob = discretize(pdesys,discretization)
4244
prob_centered = discretize(pdesys,discretization_centered)
45+
prob_centered_order4 = discretize(pdesys,discretization_centered_order4)
4346

4447
# Solve ODE problem
4548
sol = solve(prob,Tsit5(),saveat=0.1)
4649
sol_centered = solve(prob_centered,Tsit5(),saveat=0.1)
50+
sol_centered_order4 = solve(prob_centered_order4,Tsit5(),saveat=0.1)
4751

4852
x = dx[2:end-1]
4953
t = sol.t
@@ -53,11 +57,14 @@ using ModelingToolkit: Differential
5357
exact = u_exact(x, t[i])
5458
u_approx = sol.u[i]
5559
u_approx_centered = sol_centered.u[i]
60+
u_approx_centered_order4 =sol_centered_order4.u[i]
5661
@test u_approx exact atol=0.01
5762
@test u_approx_centered exact atol=0.01
63+
@test u_approx_centered_order4 exact atol=0.01
5864
end
5965
end
6066

67+
6168
@testset "Test 01: Dt(u(t,x)) ~ D*Dxx(u(t,x))" begin
6269
# Parameters, variables, and derivatives
6370
@parameters t x D

test/MOL/MOL_2D_Diffusion.jl

Lines changed: 51 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,51 @@
1+
# 2D diffusion problem
2+
3+
# Packages and inclusions
4+
using ModelingToolkit,DiffEqOperators,LinearAlgebra,Test,OrdinaryDiffEq
5+
using ModelingToolkit: Differential
6+
7+
# Tests
8+
@testset "Test 00: Dt(u(t,x)) ~ Dxx(u(t,x,y)) + Dyy(u(t,x,y))" begin
9+
@parameters t x y
10+
@variables u(..)
11+
Dxx = Differential(x)^2
12+
Dyy = Differential(y)^2
13+
Dt = Differential(t)
14+
t_min= 0.
15+
t_max = 2.0
16+
x_min = 0.
17+
x_max = 2.
18+
y_min = 0.
19+
y_max = 2.
20+
21+
# 3D PDE
22+
eq = Dt(u(t,x,y)) ~ Dxx(u(t,x,y)) + Dyy(u(t,x,y))
23+
24+
analytic_sol_func(t,x,y) = exp(x+y)*cos(x+y+4t)
25+
# Initial and boundary conditions
26+
bcs = [u(t_min,x,y) ~ analytic_sol_func(t_min,x,y),
27+
u(t,x_min,y) ~ analytic_sol_func(t,x_min,y),
28+
u(t,x_max,y) ~ analytic_sol_func(t,x_max,y),
29+
u(t,x,y_min) ~ analytic_sol_func(t,x,y_min),
30+
u(t,x,y_max) ~ analytic_sol_func(t,x,y_max)]
31+
32+
# Space and time domains
33+
domains = [t IntervalDomain(t_min,t_max),
34+
x IntervalDomain(x_min,x_max),
35+
y IntervalDomain(y_min,y_max)]
36+
pdesys = PDESystem([eq],bcs,domains,[t,x,y],[u(t,x,y)])
37+
38+
# Method of lines discretization
39+
dx = 0.1; dy = 0.2
40+
for order in [2,4,6]
41+
discretization = MOLFiniteDifference([x=>dx,y=>dy],t;centered_order=order)
42+
prob = ModelingToolkit.discretize(pdesys,discretization)
43+
sol = solve(prob,Tsit5())
44+
45+
# Test against exact solution
46+
# TODO: do this properly when sol[u] with reshape etc works
47+
@test sol.u[1][1] analytic_sol_func(sol.t[1],0.1,0.2)
48+
@test sol.u[1][2] analytic_sol_func(sol.t[1],0.2,0.2)
49+
50+
end
51+
end

test/runtests.jl

Lines changed: 2 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -38,7 +38,8 @@ end
3838
if GROUP == "All" || GROUP == "MOLFiniteDifference"
3939
@time @safetestset "MOLFiniteDifference Interface" begin include("MOL/MOLtest.jl") end
4040
#@time @safetestset "MOLFiniteDifference Interface: Linear Convection" begin include("MOL/MOL_1D_Linear_Convection.jl") end
41-
@time @safetestset "MOLFiniteDifference Interface: Linear Diffusion" begin include("MOL/MOL_1D_Linear_Diffusion.jl") end
41+
@time @safetestset "MOLFiniteDifference Interface: 1D Diffusion" begin include("MOL/MOL_1D_Diffusion.jl") end
42+
@time @safetestset "MOLFiniteDifference Interface: 2D Diffusion" begin include("MOL/MOL_2D_Diffusion.jl") end
4243
end
4344

4445
if GROUP == "All" || GROUP == "Misc"

0 commit comments

Comments
 (0)