Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Combining ODE and PDE (part 2) #91

Closed
ctessum opened this issue Apr 30, 2022 · 3 comments · Fixed by #230
Closed

Combining ODE and PDE (part 2) #91

ctessum opened this issue Apr 30, 2022 · 3 comments · Fixed by #230

Comments

@ctessum
Copy link
Contributor

ctessum commented Apr 30, 2022

Hello,

Continuing from #83, and starting with this test:

@testset "Test 13: one linear diffusion with mixed BCs, one ODE" begin
# Method of Manufactured Solutions
u_exact = (x,t) -> exp.(-t) * sin.(x)
v_exact = (t) -> exp.(-t)
@parameters t x
@variables u(..) v(..)
Dt = Differential(t)
Dx = Differential(x)
Dxx = Dx^2
# 1D PDE and boundary conditions
eqs = [Dt(u(t,x)) ~ Dxx(u(t,x)),
Dt(v(t)) ~ -v(t)]
bcs = [u(0,x) ~ sin(x),
v(0) ~ 1,
u(t,0) ~ 0,
Dx(u(t,1)) ~ exp(-t) * cos(1)]
# Space and time domains
domains = [t Interval(0.0,1.0),
x Interval(0.0,1.0)]
# PDE system
@named pdesys = PDESystem(eqs,bcs,domains,[t,x],[u(t,x),v(t)])
# Method of lines discretization
l = 100
dx = range(0.0,1.0,length=l)
dx_ = dx[2]-dx[1]
order = 2
discretization = MOLFiniteDifference([x=>dx_],t)
# Convert the PDE problem into an ODE problem
prob = discretize(pdesys,discretization)
# Solve ODE problem
sol = solve(prob,Tsit5(),saveat=0.1)
x_sol = dx[2:end-1]
t_sol = sol.t
# Test against exact solution
for i in 1:length(sol)
@test all(isapprox.(u_exact(x_sol, t_sol[i]), sol.u[i][1:length(x_sol)], atol=0.01))
@test v_exact(t_sol[i]) sol.u[i][end] atol=0.01
end
end

I change the equation Dt(u(t,x)) ~ Dxx(u(t,x)) to be Dt(u(t,x)) ~ Dxx(u(t,x)) + v(t):

using ModelingToolkit, MethodOfLines, OrdinaryDiffEq, DomainSets

@parameters t x
@variables u(..) v(..)
Dt = Differential(t)
Dx = Differential(x)
Dxx = Dx^2

# 1D PDE and boundary conditions
eqs = [Dt(u(t,x)) ~ Dxx(u(t,x)) + v(t), # This is the only line that is significantly changed from the test.
        Dt(v(t)) ~ -v(t)]

bcs = [u(0,x) ~ sin(x),
        v(0) ~ 1,
        u(t,0) ~ 0,
        Dx(u(t,1)) ~ exp(-t) * cos(1)]
domains = [t  Interval(0.0,1.0),
            x  Interval(0.0,1.0)]
@named pdesys = PDESystem(eqs,bcs,domains,[t,x],[u(t,x),v(t)])
discretization = MOLFiniteDifference([x=>0.01],t)
prob = discretize(pdesys,discretization)

I was expecting that this would add a value to that equation that varied in t but was constant across all x. Instead, I get the following error: ERROR: ArgumentError: reducing over an empty collection is not allowed.

Is this a bug, or is it the expected behavior? Thanks for your help!

@xtalax
Copy link
Member

xtalax commented May 6, 2022

I know where this error comes from, I expected that this should already be supported. I will take a closer look at fixing this soon.

@ctessum
Copy link
Contributor Author

ctessum commented May 6, 2022

Thanks for checking!

@xtalax xtalax pinned this issue Jul 20, 2022
@xtalax xtalax mentioned this issue Jan 11, 2023
xtalax added a commit that referenced this issue Jan 11, 2023
@ctessum
Copy link
Contributor Author

ctessum commented Jan 11, 2023

Thanks!

@ChrisRackauckas ChrisRackauckas unpinned this issue May 3, 2023
Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment
Labels
None yet
Projects
None yet
Development

Successfully merging a pull request may close this issue.

2 participants