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

Commit 5cd91bf

Browse files
#300 add more error checking
1 parent 964c593 commit 5cd91bf

2 files changed

Lines changed: 52 additions & 17 deletions

File tree

src/MOL_discretization.jl

Lines changed: 29 additions & 13 deletions
Original file line numberDiff line numberDiff line change
@@ -8,7 +8,7 @@ end
88

99
function throw_bc_err(bc)
1010
throw(BoundaryConditionError(
11-
"Could not recognize boundary condition '$(bc.lhs) ~ $(bc.rhs)'"
11+
"Could not read boundary condition '$(bc.lhs) ~ $(bc.rhs)'"
1212
))
1313
end
1414
# Get boundary conditions from an array
@@ -28,13 +28,16 @@ function get_bcs(bcs,tdomain,domain)
2828
# Extract the variable from the lhs
2929
if operation(lhs) isa Sym
3030
# Dirichlet boundary condition
31-
var = operation(lhs).name
31+
var = nameof(operation(lhs))
3232
α = 1.0
3333
β = 0.0
3434
bc_args = lhs.args
3535
elseif operation(lhs) isa Differential
3636
# Neumann boundary condition
37-
var = operation(lhs.args[1]).name
37+
# Check that we don't have a second-order derivative in the
38+
# boundary condition, by checking that the argument is a Sym
39+
@assert operation(lhs.args[1]) isa Sym throw_bc_err(bcs[i])
40+
var = nameof(operation(lhs.args[1]))
3841
α = 0.0
3942
β = 1.0
4043
bc_args = lhs.args[1].args
@@ -44,38 +47,50 @@ function get_bcs(bcs,tdomain,domain)
4447
# Left side of the expression should be Sym or α * Sym
4548
if operation(lhs_l) isa Sym
4649
α = 1.0
47-
var_l = operation(lhs_l).name
50+
var_l = nameof(operation(lhs_l))
4851
bc_args_l = lhs_l.args
4952
elseif operation(lhs_l) isa typeof(*)
5053
α = lhs_l.args[1]
5154
# Convert α to a Float64 if it is an Int, leave unchanged otherwise
5255
α = α isa Int ? Float64(α) : α
53-
@assert operation(lhs_l.args[2]) isa Sym
54-
var_l = operation(lhs_l.args[2]).name
56+
@assert operation(lhs_l.args[2]) isa Sym throw_bc_err(bcs[i])
57+
var_l = nameof(operation(lhs_l.args[2]))
5558
bc_args_l = lhs_l.args[2].args
5659
else
5760
throw_bc_err(bcs[i])
5861
end
5962
# Right side of the expression should be Differential or β * Differential
6063
if operation(lhs_r) isa Differential
64+
# Check that we don't have a second-order derivative in the
65+
# boundary condition
66+
@assert operation(lhs_r.args[1]) isa Sym throw_bc_err(bcs[i])
6167
β = 1.0
62-
var_r = operation(lhs_r.args[1]).name
68+
var_r = nameof(operation(lhs_r.args[1]))
6369
bc_args_r = lhs_r.args[1].args
6470
elseif operation(lhs_r) isa typeof(*)
6571
β = lhs_r.args[1]
6672
# Convert β to a Float64 if it is an Int, leave unchanged otherwise
6773
β = β isa Int ? Float64(β) : β
68-
@assert operation(lhs_r.args[2]) isa Differential
69-
var_r = operation(lhs_r.args[2].args[1]).name
74+
# Check that the bc is a derivative
75+
@assert operation(lhs_r.args[2]) isa Differential throw_bc_err(bcs[i])
76+
# But not second order (argument should be a Sym)
77+
@assert operation(lhs_r.args[2].args[1]) isa Sym throw_bc_err(bcs[i])
78+
var_r = nameof(operation(lhs_r.args[2].args[1]))
7079
bc_args_r = lhs_r.args[2].args[1].args
7180
else
7281
throw_bc_err(bcs[i])
7382
end
7483
# Check var and bc_args are the same in lhs and rhs, and if so assign
7584
# the unique value
76-
@assert var_l == var_r
85+
@assert var_l == var_r throw(BoundaryConditionError(
86+
"mismatched variables '$var_l' and '$var_r' "
87+
* "in Robin BC '$(bcs[i].lhs) ~ $(bcs[i].rhs)'"
88+
))
7789
var = var_l
78-
@assert bc_args_l == bc_args_r
90+
@assert bc_args_l == bc_args_r throw(BoundaryConditionError(
91+
"mismatched args $bc_args_l and $bc_args_r "
92+
* "in Robin BC '$(bcs[i].lhs) ~ $(bcs[i].rhs)'"
93+
))
7994
bc_args = bc_args_l
8095
else
8196
throw_bc_err(bcs[i])
@@ -93,8 +108,9 @@ function get_bcs(bcs,tdomain,domain)
93108
key = "right bc"
94109
else
95110
throw(BoundaryConditionError(
96-
"Boundary condition not recognized. "
97-
* "Should be applied at t=$(tdomain.lower), x=$(domain.lower), or x=$(domain.upper)"
111+
"Boundary condition '$(bcs[i].lhs) ~ $(bcs[i].rhs)' could not be read. "
112+
* "BCs should be applied at t=$(tdomain.lower), "
113+
* "x=$(domain.lower), or x=$(domain.upper)"
98114
))
99115
end
100116
# Create value

test/MOL_1D_Linear_Diffusion.jl

Lines changed: 23 additions & 4 deletions
Original file line numberDiff line numberDiff line change
@@ -359,6 +359,13 @@ end
359359
Dx(u(t,0.5)) ~ 0.0]
360360
pdesys = PDESystem(eq,bcs,domains,[t,x],[u])
361361
@test_throws BoundaryConditionError discretize(pdesys,discretization)
362+
363+
# Second-order derivative in BC
364+
bcs = [u(0,x) ~ 0.5 + sin(2pi*x),
365+
Dxx(u(t,0)) ~ 0.0,
366+
Dx(u(t,0.5)) ~ 0.0]
367+
pdesys = PDESystem(eq,bcs,domains,[t,x],[u])
368+
@test_throws BoundaryConditionError discretize(pdesys,discretization)
362369

363370
# Wrong format for Robin BCs
364371
bcs = [u(0,x) ~ 0.5 + sin(2pi*x),
@@ -367,6 +374,18 @@ end
367374
pdesys = PDESystem(eq,bcs,domains,[t,x],[u])
368375
@test_throws BoundaryConditionError discretize(pdesys,discretization)
369376

377+
bcs = [u(0,x) ~ 0.5 + sin(2pi*x),
378+
Dx(u(t,0)) ~ 0.0,
379+
u(t,1) + Dxx(u(t,1)) ~ 0.0]
380+
pdesys = PDESystem(eq,bcs,domains,[t,x],[u])
381+
@test_throws BoundaryConditionError discretize(pdesys,discretization)
382+
383+
bcs = [u(0,x) ~ 0.5 + sin(2pi*x),
384+
Dx(u(t,0)) ~ 0.0,
385+
u(t,1) + 2Dxx(u(t,1)) ~ 0.0]
386+
pdesys = PDESystem(eq,bcs,domains,[t,x],[u])
387+
@test_throws BoundaryConditionError discretize(pdesys,discretization)
388+
370389
bcs = [u(0,x) ~ 0.5 + sin(2pi*x),
371390
Dx(u(t,0)) ~ 0.0,
372391
Dx(u(t,1)) + u(t,1) ~ 0.0]
@@ -390,13 +409,13 @@ end
390409
Dx(u(t,0)) ~ 0.0,
391410
u(t,0) + Dx(u(t,1)) ~ 0.0]
392411
pdesys = PDESystem(eq,bcs,domains,[t,x],[u])
393-
@test_throws AssertionError discretize(pdesys,discretization)
412+
@test_throws BoundaryConditionError discretize(pdesys,discretization)
394413

395414
# Mismatching variables
396415
bcs = [u(0,x) ~ 0.5 + sin(2pi*x),
397-
Dx(u(t,0)) ~ 0.0,
398-
u(t,1) + Dx(v(t,1)) ~ 0.0]
416+
Dx(u(t,0)) ~ 0.0,
417+
u(t,1) + Dx(v(t,1)) ~ 0.0]
399418
pdesys = PDESystem(eq,bcs,domains,[t,x],[u])
400-
@test_throws AssertionError discretize(pdesys,discretization)
419+
@test_throws BoundaryConditionError discretize(pdesys,discretization)
401420

402421
end

0 commit comments

Comments
 (0)