UNIQUAC (universal quasichemical) Excess Gibbs free energy model.
With:
use yaeos, only: pr, setup_uniquac, UNIQUAC
integer, parameter :: nc = 3
real(pr) :: rs(nc), qs(nc)
real(pr) :: b(nc, nc)
real(pr) :: n(nc)
real(pr) :: ln_gammas(nc), T
type(UNIQUAC) :: model
rs = [0.92_pr, 2.1055_pr, 3.1878_pr]
qs = [1.4_pr, 1.972_pr, 2.4_pr]
T = 298.15_pr
! Calculate bij from DUij. We need -DU/R to get bij
b(1, :) = [0.0_pr, -526.02_pr, -309.64_pr]
b(2, :) = [318.06_pr, 0.0_pr, 91.532_pr]
b(3, :) = [-1325.1_pr, -302.57_pr, 0.0_pr]
model = setup_uniquac(qs, rs, bij=b)
n = [0.8_pr, 0.1_pr, 0.2_pr]
call model%ln_activity_coefficient(n, T, ln_gammas)
print *, exp(ln_gammas) ! [8.856, 0.860, 1.425]
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
real(kind=pr), | public, | allocatable | :: | aij(:,:) |
Interaction parameters matrix |
||
real(kind=pr), | public, | allocatable | :: | bij(:,:) |
Interaction parameters matrix |
||
real(kind=pr), | public, | allocatable | :: | cij(:,:) |
Interaction parameters matrix |
||
type(Substances), | public | :: | components |
Substances contained in the module |
|||
real(kind=pr), | public, | allocatable | :: | dij(:,:) |
Interaction parameters matrix |
||
real(kind=pr), | public, | allocatable | :: | eij(:,:) |
Interaction parameters matrix |
||
real(kind=pr), | public, | allocatable | :: | qs(:) |
Molecule’s relative areas |
||
real(kind=pr), | public, | allocatable | :: | rs(:) |
Molecule’s relative volumes |
||
real(kind=pr), | public | :: | z | = | 10.0_pr |
Model coordination number |
Calculate the excess Gibbs free energy and its derivatives of the UNIQUAC model.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(UNIQUAC), | intent(in) | :: | self |
UNIQUAC model |
||
real(kind=pr), | intent(in) | :: | n(:) |
Moles vector [mol] |
||
real(kind=pr), | intent(in) | :: | T |
Temperature [K] |
||
real(kind=pr), | intent(out), | optional | :: | Ge |
Excess Gibbs energy |
|
real(kind=pr), | intent(out), | optional | :: | GeT |
|
|
real(kind=pr), | intent(out), | optional | :: | GeT2 |
|
|
real(kind=pr), | intent(out), | optional | :: | Gen(size(n)) |
|
|
real(kind=pr), | intent(out), | optional | :: | GeTn(size(n)) |
|
|
real(kind=pr), | intent(out), | optional | :: | Gen2(size(n),size(n)) |
|
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(GeModel), | intent(in) | :: | self | |||
real(kind=pr), | intent(in) | :: | n(:) | |||
real(kind=pr), | intent(in) | :: | T | |||
real(kind=pr), | intent(out) | :: | lngamma(:) |
Calculate the temperature dependence term of the UNIQUAC model.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(UNIQUAC), | intent(in) | :: | self |
UNIQUAC model |
||
real(kind=pr), | intent(in) | :: | T |
Temperature [K] |
||
real(kind=pr), | intent(out), | optional | :: | tau(size(self%qs),size(self%qs)) |
UNIQUAC temperature dependence term |
|
real(kind=pr), | intent(out), | optional | :: | tauT(size(self%qs),size(self%qs)) |
|
|
real(kind=pr), | intent(out), | optional | :: | tauT2(size(self%qs),size(self%qs)) |
|