legacy_ar_models Module

Legacy Thermodynamic routines Module for a cubic eos system, made with the intention to keep compatiblity with legacy codes but with a better structure. this should be later adapted into a simple oop system where an eos object stores the relevant parameters (or some functional oriented approach)



Variables

Type Visibility Attributes Name Initial
real(kind=pr), public, allocatable :: ac(:)

Critical attractive parameter [bar (L/mol)^2]

real(kind=pr), public, allocatable :: b(:)

repulsive parameter [L]

real(kind=pr), public, allocatable :: bij(:,:)
real(kind=pr), public, allocatable :: dc(:)

Critical density [mol/L]

real(kind=pr), public, allocatable :: del1(:)

parameter

real(kind=pr), public, allocatable :: k(:)

Attractive parameter constant

real(kind=pr), public, allocatable :: kij(:,:)

Attractive BIP

real(kind=pr), public, allocatable :: kij0(:,:)
real(kind=pr), public, allocatable :: kinf(:,:)
real(kind=pr), public, allocatable :: lij(:,:)

Repulsive BIP

integer, public :: mixing_rule

What mixing rule to use

integer, public :: nc

Number of components

real(kind=pr), public, allocatable :: pc(:)

Critical pressure [bar]

real(kind=pr), public, allocatable :: tc(:)

Critical temperature [K]

integer, public :: tdep

Temperature dependance of kij

integer, public :: thermo_model

Which thermodynamic model to use

real(kind=pr), public, allocatable :: tstar(:,:)
real(kind=pr), public, allocatable :: w(:)

Acentric factor

real(kind=pr), public, allocatable :: z(:)

Mole fractions vector


Functions

public function cubic_v0(z, p, t)

Arguments

Type IntentOptional Attributes Name
real(kind=pr) :: z(:)
real(kind=pr) :: p
real(kind=pr) :: t

Return Value real(kind=pr)


Subroutines

public subroutine ArVnder(nc, nder, ntemp, z, V, T, ar, arv, artv, arv2, Arn, ArVn, ArTn, Arn2)

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: nc
integer, intent(in) :: nder
integer, intent(in) :: ntemp
real(kind=pr), intent(in) :: z(nc)
real(kind=pr), intent(in) :: V
real(kind=pr), intent(in) :: T
real(kind=pr), intent(out) :: ar
real(kind=pr), intent(out) :: arv
real(kind=pr), intent(out) :: artv
real(kind=pr), intent(out) :: arv2
real(kind=pr), intent(out), dimension(size(z)) :: Arn
real(kind=pr), intent(out), dimension(size(z)) :: ArVn
real(kind=pr), intent(out), dimension(size(z)) :: ArTn
real(kind=pr), intent(out) :: Arn2(size(z),size(z))

public subroutine Bnder(nc, rn, Bmix, dBi, dBij)

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: nc
real(kind=pr), intent(in) :: rn(nc)
real(kind=pr), intent(out) :: Bmix
real(kind=pr), intent(out) :: dBi(nc)
real(kind=pr), intent(out) :: dBij(nc,nc)

public subroutine DELTAnder(nc, rn, D1m, dD1i, dD1ij)

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: nc
real(kind=pr), intent(in) :: rn(nc)
real(kind=pr), intent(out) :: D1m
real(kind=pr), intent(out) :: dD1i(nc)
real(kind=pr), intent(out) :: dD1ij(nc,nc)

public subroutine DandTnder(ntd, nc, T, rn, D, dDi, dDiT, dDij, dDdT, dDdT2)

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: ntd
integer, intent(in) :: nc
real(kind=pr), intent(in) :: T
real(kind=pr), intent(in) :: rn(nc)
real(kind=pr), intent(out) :: D
real(kind=pr), intent(out) :: dDi(nc)
real(kind=pr), intent(out) :: dDiT(nc)
real(kind=pr), intent(out) :: dDij(nc,nc)
real(kind=pr), intent(out) :: dDdT
real(kind=pr), intent(out) :: dDdT2

public subroutine HelmRKPR(nco, NDE, NTD, rn, V, T, Ar, ArV, ArTV, ArV2, Arn, ArVn, ArTn, Arn2)

Calculate the reduced residual Helmholtz Energy and it’s derivatives with the RKPR EOS

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: nco
integer, intent(in) :: NDE
integer, intent(in) :: NTD
real(kind=pr), intent(in) :: rn(nco)
real(kind=pr), intent(in) :: V
real(kind=pr), intent(in) :: T
real(kind=pr), intent(out) :: Ar
real(kind=pr), intent(out) :: ArV
real(kind=pr), intent(out) :: ArTV
real(kind=pr), intent(out) :: ArV2
real(kind=pr), intent(out) :: Arn(nco)
real(kind=pr), intent(out) :: ArVn(nco)
real(kind=pr), intent(out) :: ArTn(nco)
real(kind=pr), intent(out) :: Arn2(nco,nco)

public subroutine HelmSRKPR(nc, nd, nt, rn, v, t, ar, arv, artv, arv2, Arn, ArVn, ArTn, Arn2)

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: nc

Number of components

integer, intent(in) :: nd

Compositional derivatives

integer, intent(in) :: nt

Temperature derivatives

real(kind=pr), intent(in) :: rn(nc)

Number of moles

real(kind=pr), intent(in) :: v

Volume [L]

real(kind=pr), intent(in) :: t

Temperature [K]

real(kind=pr), intent(out) :: ar

Residual Helmholtz

real(kind=pr), intent(out) :: arv

dAr/dV

real(kind=pr), intent(out) :: artv

dAr2/dTV

real(kind=pr), intent(out) :: arv2

dAr2/dV2

real(kind=pr), intent(out) :: Arn(nc)

dAr/dn

real(kind=pr), intent(out) :: ArVn(nc)

dAr2/dVn

real(kind=pr), intent(out) :: ArTn(nc)

dAr2/dTn

real(kind=pr), intent(out) :: Arn2(nc,nc)

dAr2/dn2

public subroutine PR76_factory(moles_in, ac_in, b_in, tc_in, pc_in, w_in, k_in)

PengRobinson 76 factory

Read more…

Arguments

Type IntentOptional Attributes Name
real(kind=pr), intent(in) :: moles_in(nc)
real(kind=pr), intent(in), optional :: ac_in(nc)
real(kind=pr), intent(in), optional :: b_in(nc)
real(kind=pr), intent(in), optional :: tc_in(nc)
real(kind=pr), intent(in), optional :: pc_in(nc)
real(kind=pr), intent(in), optional :: w_in(nc)
real(kind=pr), intent(in), optional :: k_in(nc)

public subroutine PR78_factory(moles_in, ac_in, b_in, tc_in, pc_in, w_in, k_in)

PengRobinson 78 factory

Read more…

Arguments

Type IntentOptional Attributes Name
real(kind=pr), intent(in) :: moles_in(nc)
real(kind=pr), intent(in), optional :: ac_in(nc)
real(kind=pr), intent(in), optional :: b_in(nc)
real(kind=pr), intent(in), optional :: tc_in(nc)
real(kind=pr), intent(in), optional :: pc_in(nc)
real(kind=pr), intent(in), optional :: w_in(nc)
real(kind=pr), intent(in), optional :: k_in(nc)

public subroutine SRK_factory(moles_in, ac_in, b_in, tc_in, pc_in, w_in, k_in)

SoaveRedlichKwong factory

Read more…

Arguments

Type IntentOptional Attributes Name
real(kind=pr), intent(in) :: moles_in(nc)
real(kind=pr), intent(in), optional :: ac_in(nc)
real(kind=pr), intent(in), optional :: b_in(nc)
real(kind=pr), intent(in), optional :: tc_in(nc)
real(kind=pr), intent(in), optional :: pc_in(nc)
real(kind=pr), intent(in), optional :: w_in(nc)
real(kind=pr), intent(in), optional :: k_in(nc)

public subroutine aTder(ac, Tc, k, T, a, dadT, dadT2)

Arguments

Type IntentOptional Attributes Name
real(kind=pr), intent(in) :: ac
real(kind=pr), intent(in) :: Tc
real(kind=pr), intent(in) :: k
real(kind=pr), intent(in) :: T
real(kind=pr), intent(out) :: a
real(kind=pr), intent(out) :: dadT
real(kind=pr), intent(out) :: dadT2

public subroutine aijTder(ntd, nc, T, aij, daijdT, daijdT2)

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: ntd
integer, intent(in) :: nc
real(kind=pr), intent(in) :: T
real(kind=pr), intent(out) :: aij(nc,nc)
real(kind=pr), intent(out) :: daijdT(nc,nc)
real(kind=pr), intent(out) :: daijdT2(nc,nc)

public subroutine ar_rkpr(z, v, t, ar, arv, artv, arv2, Arn, ArVn, ArTn, Arn2)

Arguments

Type IntentOptional Attributes Name
real(kind=pr), intent(in) :: z(:)

Number of moles

real(kind=pr), intent(in) :: v

Volume [L]

real(kind=pr), intent(in) :: t

Temperature [K]

real(kind=pr), intent(out) :: ar

Residual Helmholtz

real(kind=pr), intent(out) :: arv

dAr/dV

real(kind=pr), intent(out) :: artv

dAr2/dTV

real(kind=pr), intent(out) :: arv2

dAr2/dV2

real(kind=pr), intent(out) :: Arn(size(z))

dAr/dn

real(kind=pr), intent(out) :: ArVn(size(z))

dAr2/dVn

real(kind=pr), intent(out) :: ArTn(size(z))

dAr2/dTn

real(kind=pr), intent(out) :: Arn2(size(z),size(z))

dAr2/dn2

public subroutine ar_srkpr(z, v, t, ar, arv, artv, arv2, Arn, ArVn, ArTn, Arn2)

Wrapper subroutine to the SRK/PR Residula Helmholtz function to use the general interface

Arguments

Type IntentOptional Attributes Name
real(kind=pr), intent(in) :: z(:)

Number of moles

real(kind=pr), intent(in) :: v

Volume [L]

real(kind=pr), intent(in) :: t

Temperature [K]

real(kind=pr), intent(out) :: ar

Residual Helmholtz

real(kind=pr), intent(out) :: arv

dAr/dV

real(kind=pr), intent(out) :: artv

dAr2/dTV

real(kind=pr), intent(out) :: arv2

dAr2/dV2

real(kind=pr), intent(out) :: Arn(size(z))

dAr/dn

real(kind=pr), intent(out) :: ArVn(size(z))

dAr2/dVn

real(kind=pr), intent(out) :: ArTn(size(z))

dAr2/dTn

real(kind=pr), intent(out) :: Arn2(size(z),size(z))

dAr2/dn2

public subroutine get_Zc_OMa_OMb(del1, Zc, OMa, OMb)

Calculate Zc, OMa and OMb from the delta_1 parameter.

Arguments

Type IntentOptional Attributes Name
real(kind=pr), intent(in) :: del1(:)

delta_1 parameter

real(kind=pr), intent(out) :: Zc(:)

Critical compressibility factor

real(kind=pr), intent(out) :: OMa(:)

OMa

real(kind=pr), intent(out) :: OMb(:)

OMb

public subroutine setup(n, nmodel, ntdep, ncomb)

Setup the basics variables that describe the model.

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: n

Number of components

integer, intent(in) :: nmodel

Number of model

integer, intent(in) :: ntdep

Kij dependant of temperature

integer, intent(in) :: ncomb

Combining rule