diff --git a/docs/src/finite.md b/docs/src/finite.md index c9f3c49d..88bb37e6 100644 --- a/docs/src/finite.md +++ b/docs/src/finite.md @@ -8,6 +8,7 @@ CurrentModule = Fronts FiniteProblem FiniteDirichletProblem FiniteFluxProblem +FiniteReservoirProblem solve(prob::FiniteProblem{<:DiffusionEquation{1}}, tstop) FiniteSolution ``` diff --git a/src/Fronts.jl b/src/Fronts.jl index 89396358..4c8f7de9 100644 --- a/src/Fronts.jl +++ b/src/Fronts.jl @@ -78,7 +78,7 @@ export solve export Solution, rb, flux, sorptivity export SolvingError export inverse -export FiniteProblem, FiniteDirichletProblem, FiniteFluxProblem, FiniteSolution +export FiniteProblem, FiniteDirichletProblem, FiniteFluxProblem, FiniteReservoirProblem, FiniteSolution include("ParamEstim.jl") diff --git a/src/finite.jl b/src/finite.jl index fdb59f73..934636e9 100644 --- a/src/finite.jl +++ b/src/finite.jl @@ -65,6 +65,32 @@ function Base.show(io::IO, prob::FiniteFluxProblem) print(io, "⎩ q(0,t) = ", prob.qb, ", t>0") end +""" + FiniteReservoirProblem(eq, rstop; i, b, capacity) <: FiniteProblem + +Models `eq` in the domain `0 ≤ r ≤ rstop` with initial condition `i` and a reservoir boundary condition with capacity `capacity`. +""" +struct FiniteReservoirProblem{Teq, _Trstop, _Ti, _Tb, _Tcapacity} <: FiniteProblem{Teq} + eq::Teq + rstop::_Trstop + i::_Ti + b::_Tb + capacity::_Tcapacity + + function FiniteReservoirProblem(eq::DiffusionEquation{1}, rstop; i, b, capacity) + @argcheck rstop > 0 + new{typeof(eq),typeof(rstop),typeof(i),typeof(b),typeof(capacity)}(eq, rstop, i, b, capacity) + end +end + +""" + FiniteReservoirProblem(D, rstop; i, b, capacity) <: FiniteProblem + +Shortcut for `FiniteReservoirProblem(DiffusionEquation(D), rstop, i=i, b=b, capacity=capacity)`. +""" + +FiniteReservoirProblem(D, rstop; i, b, capacity) = FiniteReservoirProblem(DiffusionEquation(D), rstop, i=i, b=b, capacity=capacity) + """ solve(prob::FiniteProblem{<:DiffusionEquation{1}}, tstop[; N, tol, Δt]) -> FiniteSolution @@ -83,13 +109,21 @@ function solve(prob::FiniteProblem{<:DiffusionEquation{1}}, tstop; N=500, tol=1e Δr = step(r) Δr² = Δr^2 + if prob isa FiniteReservoirProblem + used = zero(prob.capacity) + end + θ = similar(r) t = 0 isol = nothing # Solve with Fronts as much as possible - if prob.i isa Number && prob isa FiniteDirichletProblem + if prob.i isa Number && (prob isa FiniteDirichletProblem || prob isa FiniteReservoirProblem) isol = solve(DirichletProblem(prob.eq, i=prob.i, b=prob.b), itol=tol) t = min(tstop, (prob.rstop/isol.ϕi)^2) + if prob isa FiniteReservoirProblem + t = min(t, (prob.capacity/sorptivity(isol))^2) + used = sorptivity(isol)*√t + end θ .= isol.(r, t) else θ .= prob.i @@ -135,14 +169,18 @@ function solve(prob::FiniteProblem{<:DiffusionEquation{1}}, tstop; N=500, tol=1e B .= θ # Apply boundary conditions - if prob isa FiniteDirichletProblem + if prob isa FiniteReservoirProblem + influx = min(-Df[begin]*(θ[begin+1] - prob.b)/Δr*Δt, prob.capacity - used) + end + + if prob isa FiniteDirichletProblem || (prob isa FiniteReservoirProblem && influx < prob.capacity - used) A[begin,begin] = 1 A[begin,begin+1] = 0 B[begin] = prob.b - elseif (prob isa FiniteFluxProblem && prob.qb != zero(prob.qb)) + elseif (prob isa FiniteFluxProblem && prob.qb != zero(prob.qb)) || (prob isa FiniteReservoirProblem && influx > zero(influx)) A[begin,begin] = Df[begin]/Δr A[begin,begin+1] = -Df[begin]/Δr - B[begin] = prob isa FiniteFluxProblem ? prob.qb : influx + B[begin] = prob isa FiniteReservoirProblem ? influx/Δt : prob.qb end θ_prev_sweep .= θ @@ -156,6 +194,11 @@ function solve(prob::FiniteProblem{<:DiffusionEquation{1}}, tstop; N=500, tol=1e push!(ts, t) push!(θs, copy(θ)) + if prob isa FiniteReservoirProblem + influx = -Df[begin]*(θ[begin+1] - θ[begin])/Δr*Δt + used += influx + end + if sweeps < 3 Δt *= 1.3 end diff --git a/test/test_finite.jl b/test/test_finite.jl index 8bd15378..309d3e43 100644 --- a/test/test_finite.jl +++ b/test/test_finite.jl @@ -32,4 +32,25 @@ @test NumericalIntegration.integrate(r, θ.(r, t) - θ.(r, 0)) ≈ prob.qb*t atol=1e-3 end end + + @testset "FiniteReservoirProblem" begin + # Wetting of a Whatman No. 1 paper strip, LETd model + # Reference: Gerlero et al. (2022) + # https://doi.org/10.1007/s11242-021-01724-w + θs = 0.7 + θi = 0.025 + ϵ = 1e-7 + + pm = LETd(L=0.004569, E=12930, T=1.505, Dwt=4.660e-4, θr=0.019852, θs=θs) + + r = range(0, 0.05, length=500) + + prob = FiniteReservoirProblem(pm, r[end], i=θi, b=θs-ϵ, capacity=1e-2) + + θ = solve(prob, 200, N=length(r)) + + for t in [100, 150, 200] + @test NumericalIntegration.integrate(r, θ.(r, t) .- θi) ≈ prob.capacity atol=1e-4 + end + end end