Skip to content

Commit

Permalink
michelsen tp-flash: if we don't know the default eq and vle fails, tr…
Browse files Browse the repository at this point in the history
…y lle
  • Loading branch information
longemen3000 committed Jun 10, 2024
1 parent f598922 commit 1027a63
Show file tree
Hide file tree
Showing 4 changed files with 38 additions and 74 deletions.
Original file line number Diff line number Diff line change
Expand Up @@ -45,7 +45,7 @@ end

numphases(::MichelsenTPFlash) = 2

function MichelsenTPFlash(;equilibrium = :vle,
function MichelsenTPFlash(;equilibrium = :unknown,
K0 = nothing,
x0 = nothing,
y0 = nothing,
Expand All @@ -56,7 +56,7 @@ function MichelsenTPFlash(;equilibrium = :vle,
second_order = false,
noncondensables = nothing,
nonvolatiles = nothing)
!(is_vle(equilibrium) | is_lle(equilibrium)) && throw(error("invalid equilibrium specification for MichelsenTPFlash"))
!(is_vle(equilibrium) | is_lle(equilibrium) | (equilibrium == :unknown)) && throw(error("invalid equilibrium specification for MichelsenTPFlash"))
if K0 == x0 == y0 === v0 == nothing #nothing specified
#is_lle(equilibrium)
T = Nothing
Expand Down Expand Up @@ -128,7 +128,7 @@ function tp_flash_michelsen(model::EoSModel, p, T, z; equilibrium=:vle, K0=nothi
z = z_full[z_nonzero]
end

if is_vle(equilibrium)
if is_vle(equilibrium) || (equilibrium == :unknown)
phasex = :liquid
phasey = :vapor
elseif is_lle(equilibrium)
Expand Down Expand Up @@ -180,6 +180,8 @@ function tp_flash_michelsen(model::EoSModel, p, T, z; equilibrium=:vle, K0=nothi
# Computing the initial guess for the K vector
x = similar(z)
y = similar(z)
x .= z
y .= z
if !isnothing(K0)
K = 1. * K0
lnK = log.(K)
Expand All @@ -189,9 +191,17 @@ function tp_flash_michelsen(model::EoSModel, p, T, z; equilibrium=:vle, K0=nothi
lnK = log.(x ./ y)
lnK,volx,voly,_ = update_K!(lnK,model,p,T,x,y,volx,voly,phasex,phasey,nothing,inx,iny)
K = exp.(lnK)
elseif is_vle(equilibrium)
elseif is_vle(equilibrium) || (equilibrium == :unknown)
# Wilson Correlation for K
K = tp_flash_K0(model,p,T)
#if we can't predict K, we use lle
if equilibrium == :unknown
Kmin,Kmax = extrema(K)

if Kmin >= 1 || Kmax <= 1
K = K0_lle_init(model,p,T,z)
end
end
lnK = log.(K)
# volx,voly = NaN*_1,NaN*_1
else
Expand Down Expand Up @@ -315,7 +325,6 @@ function tp_flash_michelsen(model::EoSModel, p, T, z; equilibrium=:vle, K0=nothi
x = index_expansion(x,z_nonzero)
y = index_expansion(y,z_nonzero)
end

return x, y, β, (vx,vy)
end

Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -44,7 +44,7 @@ function index_reduction(m::RRTPFlash,idx::AbstractVector)
return RRTPFlash{eltype(m)}(equilibrium,K0,x0,y0,v0,K_tol,max_iters,nacc,noncondensables,nonvolatiles)
end

function RRTPFlash(;equilibrium = :vle,
function RRTPFlash(;equilibrium = :unknown,
K0 = nothing,
x0 = nothing,
y0 = nothing,
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -13,7 +13,7 @@ function tp_flash_michelsen(model::ElectrolyteModel, p, T, z; equilibrium=:vle,
z = z_full[z_nonzero]
end

if is_vle(equilibrium)
if is_vle(equilibrium) || (equilibrium == :unknown)
phasex = :liquid
phasey = :vapor
append!(non_iny_list,ions)
Expand Down Expand Up @@ -75,9 +75,17 @@ function tp_flash_michelsen(model::ElectrolyteModel, p, T, z; equilibrium=:vle,
lnK = log.(x ./ y)
lnK,volx,voly,_ = update_K!(lnK,model,p,T,x,y,volx,voly,phasex,phasey,nothing,inx,iny)
K = exp.(lnK)
elseif is_vle(equilibrium)
elseif is_vle(equilibrium) || (equilibrium == :unknown)
# Wilson Correlation for K
K = tp_flash_K0(model,p,T)
#if we can't predict K, we use lle
if equilibrium == :unknown
Kmin,Kmax = extrema(K)

if Kmin > 1 || Kmax < 1
K = K0_lle_init(model,p,T,z)
end
end
lnK = log.(K)
# volx,voly = NaN*_1,NaN*_1
else
Expand Down
79 changes: 13 additions & 66 deletions src/methods/tpd.jl
Original file line number Diff line number Diff line change
Expand Up @@ -657,74 +657,21 @@ function K0_lle_act_coefficient(model::ActivityModel,p,T,z,μ_pure)
end

function K0_lle_init(model::EoSModel, p, T, z)
nc = length(model)
z_test = __z_test(z)
ntest = length(@view(z_test[1,:]))
μ_pure = K0_lle_init_cache(model::EoSModel,p,T)
γ1 = K0_lle_act_coefficient(model,p,T,@view(z_test[1,:]),μ_pure)
γ = fill(γ1,1)
for i in 2:ntest
push!(γ,K0_lle_act_coefficient(model,p,T,@view(z_test[i,:]),μ_pure))
end
γ2 = γ
_0 = zero(eltype(γ1))
err = Inf*one(_0)
i_min = 0
j_min = 0
ϕi = similar(γ1)
for i in 1:ntest
z_test_i = @view(z_test[i,:])
γi = γ[i]
for j in i+1:ntest
z_test_j = @view(z_test[j,:])
γj = γ[i]
ϕ = _0
ϕi_finite = 0
for k in 1:length(model)
ϕik = (z_test_i[k] - z[k])/(z_test_i[k] - z_test_j[k])
ϕi[k] = ϕik
if isfinite(ϕik)
ϕi_finite += 1
ϕ += ϕik
end
end
ϕ = ϕ/ϕi_finite

err_ij = _0
for k in 1:length(model)
zi,zj = z_test_i[k],z_test_j[k]
logγᵢzᵢ = log(γi[k] * zi)
logγⱼzⱼ = log(γj[k] * zj)
err_ij += logγᵢzᵢ
err_ij -= logγⱼzⱼ
err_ij += logγᵢzᵢ**zi + (1 - ϕ)*zj - z[k])
end
#err_ij = sum(log.(γ[i].*z_test[i,:]) .-log.(γ[j].*z_test[j,:])+log.(γ[i].*z_test[i,:]).*(ϕ*z_test[i,:]+(1-ϕ)*z_test[j,:].-z))
if err_ij < err
err = err_ij
i_min,j_min = i,j
end
comps,tpds,_,_ = tpd(model,p,T,z,lle = true, strategy = :pure, break_first = true)
if length(comps) == 1
w = comps[1]
β = one(eltype(w))
for i in 1:length(z)
β = min(β,z[i]/w[i])
end
β = 0.5*β
w2 = (z .- β .*w)/(1 .- β)
w2 ./= sum(w2)
K = w ./ w2
else
K = ones(eltype(eltype(comps)),length(z))
end

K0 = γ[i_min]./ γ[j_min]

#=
#extra step, reduce magnitudes if we aren't in a single phase
g0 = dot(z, K0) - 1.
g1 = 1. - sum(zi/Ki for (zi,Ki) in zip(z,K0))
if g0 < 0 && g1 <= 0 #we need to correct the value of g0
g0i = z .* K0
kz,idx = findmax(g0i)
#the maximum K value is too great
elseif g1 > 0 && g0 >= 0 #we need to correct the value of g1
g1i = z ./ K0
else #bail out here?
end =#
return K0
return K
end

export tpd

2 comments on commit 1027a63

@pw0908
Copy link
Member

@pw0908 pw0908 commented on 1027a63 Jun 10, 2024

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

@JuliaRegistrator
Copy link

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Registration pull request created: JuliaRegistries/General/108617

Tip: Release Notes

Did you know you can add release notes too? Just add markdown formatted text underneath the comment after the text
"Release notes:" and it will be added to the registry PR, and if TagBot is installed it will also be added to the
release that TagBot creates. i.e.

@JuliaRegistrator register

Release notes:

## Breaking changes

- blah

To add them here just re-invoke and the PR will be updated.

Tagging

After the above pull request is merged, it is recommended that a tag is created on this repository for the registered package version.

This will be done automatically if the Julia TagBot GitHub Action is installed, or can be done manually through the github interface, or via:

git tag -a v0.6.0 -m "<description of version>" 1027a63688edf36320c9c71b6909be5ba68d6199
git push origin v0.6.0

Please sign in to comment.