Module that contains the automatic differentiation logic for an Ar model.
All that is needed to define an Ar model that uses automatic differentiation with hyperdual numbers is to define a new derived type that overloads the method to the Ar function that you want to use. A minimal example follows:
module newmodel
use yaeos__adiff_hyperdual_ar_api, only: ArModelAdiff
type, extends(ArModelAdiff) :: YourNewModel
type(Substances) :: composition
real(8) :: parameters(:)
contains
procedure :: Ar => arfun
procedure :: get_v0 => v0
end type
contains
subroutine arfun(self, n, v, t, Ar)
class(YourNewModel), intent(in) :: self
type(hyperdual), intent(in) :: n(:) ! Number of moles
type(hyperdual), intent(in) :: v ! Volume [L]
type(hyperdual), intent(in) :: t ! Temperature [K]
type(hyperdual), intent(out) :: ar_value ! Residual Helmholtz Energy
! A very complicated residual helmholtz function of a mixture
Ar = sum(n) * v * t
end subroutine
function v0(self, n, p, t)
class(YourNewModel), intent(in) :: self
real(pr), intent(in) :: n(:) ! Number of moles
real(pr), intent(in) :: p ! Pressure [bar]
real(pr), intent(in) :: t ! Temperature [K]
real(pr) :: v0
v0 = self%parameters(3)
end function
A complete implementation of the PR76 Equation of State can me found in
example/adiff/adiff_pr76.f90
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ArModelAdiff) | :: | self | ||||
type(hyperdual), | intent(in) | :: | n(:) | |||
type(hyperdual), | intent(in) | :: | v | |||
type(hyperdual), | intent(in) | :: | t |
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
type(Substances), | public | :: | components |
Substances contained in the module |
|||
character(len=:), | public, | allocatable | :: | name |
Name of the model |
procedure(hyperdual_Ar), public, deferred :: Ar | |
procedure, public :: Cp_residual_vt | |
procedure, public :: Cv_residual_vt | |
procedure, public :: enthalpy_residual_vt | |
procedure, public :: entropy_residual_vt | |
procedure(abs_volume_initializer), public, deferred :: get_v0 | |
procedure, public :: gibbs_residual_vt => gibbs_residual_VT | |
procedure, public :: lnphi_pt => fugacity_pt | |
procedure, public :: lnphi_vt => fugacity_vt | |
procedure, public :: pressure | |
procedure, public :: residual_helmholtz | |
procedure, public :: volume |
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ArModelAdiff), | intent(in) | :: | self | |||
real(kind=pr), | intent(in) | :: | n(:) | |||
real(kind=pr), | intent(in) | :: | v | |||
real(kind=pr), | intent(in) | :: | t | |||
real(kind=pr), | intent(out), | optional | :: | Ar | ||
real(kind=pr), | intent(out), | optional | :: | ArV | ||
real(kind=pr), | intent(out), | optional | :: | ArT | ||
real(kind=pr), | intent(out), | optional | :: | ArTV | ||
real(kind=pr), | intent(out), | optional | :: | ArV2 | ||
real(kind=pr), | intent(out), | optional | :: | ArT2 | ||
real(kind=pr), | intent(out), | optional, | dimension(size(n)) | :: | Arn | |
real(kind=pr), | intent(out), | optional, | dimension(size(n)) | :: | ArVn | |
real(kind=pr), | intent(out), | optional, | dimension(size(n)) | :: | ArTn | |
real(kind=pr), | intent(out), | optional | :: | Arn2(size(n),size(n)) |