nrtl_hv.f90 Source File


This file depends on

sourcefile~~nrtl_hv.f90~~EfferentGraph sourcefile~nrtl_hv.f90 nrtl_hv.f90 sourcefile~autodiff.f90 autodiff.f90 sourcefile~nrtl_hv.f90->sourcefile~autodiff.f90 sourcefile~base.f90 base.f90 sourcefile~nrtl_hv.f90->sourcefile~base.f90 sourcefile~constants.f90 constants.f90 sourcefile~nrtl_hv.f90->sourcefile~constants.f90 sourcefile~ge_models.f90 ge_models.f90 sourcefile~nrtl_hv.f90->sourcefile~ge_models.f90 sourcefile~armodel_adiff_api.f90 armodel_adiff_api.f90 sourcefile~autodiff.f90->sourcefile~armodel_adiff_api.f90 sourcefile~hyperdual.f90 hyperdual.f90 sourcefile~autodiff.f90->sourcefile~hyperdual.f90 sourcefile~base.f90->sourcefile~constants.f90 sourcefile~ge_models.f90->sourcefile~constants.f90 sourcefile~base.f90~2 base.f90 sourcefile~ge_models.f90->sourcefile~base.f90~2 sourcefile~armodel_adiff_api.f90->sourcefile~constants.f90 sourcefile~armodel_adiff_api.f90->sourcefile~hyperdual.f90 sourcefile~ar_models.f90 ar_models.f90 sourcefile~armodel_adiff_api.f90->sourcefile~ar_models.f90 sourcefile~substance.f90 substance.f90 sourcefile~base.f90~2->sourcefile~substance.f90 sourcefile~hyperdual.f90->sourcefile~constants.f90 sourcefile~ar_models.f90->sourcefile~constants.f90 sourcefile~ar_models.f90->sourcefile~base.f90~2 sourcefile~math.f90 math.f90 sourcefile~ar_models.f90->sourcefile~math.f90 sourcefile~substance.f90->sourcefile~constants.f90 sourcefile~math.f90->sourcefile~constants.f90 sourcefile~auxiliar.f90 auxiliar.f90 sourcefile~math.f90->sourcefile~auxiliar.f90 sourcefile~continuation.f90 continuation.f90 sourcefile~math.f90->sourcefile~continuation.f90 sourcefile~linalg.f90 linalg.f90 sourcefile~math.f90->sourcefile~linalg.f90 sourcefile~auxiliar.f90->sourcefile~constants.f90 sourcefile~continuation.f90->sourcefile~constants.f90 sourcefile~continuation.f90->sourcefile~auxiliar.f90 sourcefile~continuation.f90->sourcefile~linalg.f90 sourcefile~linalg.f90->sourcefile~constants.f90 sourcefile~linalg.f90->sourcefile~auxiliar.f90

Files dependent on this one

sourcefile~~nrtl_hv.f90~~AfferentGraph sourcefile~nrtl_hv.f90 nrtl_hv.f90 sourcefile~huron_vidal.f90 huron_vidal.f90 sourcefile~huron_vidal.f90->sourcefile~nrtl_hv.f90 sourcefile~implementations.f90 implementations.f90 sourcefile~implementations.f90->sourcefile~nrtl_hv.f90 sourcefile~implementations.f90~2 implementations.f90 sourcefile~implementations.f90~2->sourcefile~huron_vidal.f90 sourcefile~implementations.f90~2->sourcefile~implementations.f90 sourcefile~models.f90 models.f90 sourcefile~models.f90->sourcefile~huron_vidal.f90 sourcefile~models.f90->sourcefile~implementations.f90 sourcefile~models.f90->sourcefile~implementations.f90~2 sourcefile~gerg2008.f90 gerg2008.f90 sourcefile~models.f90->sourcefile~gerg2008.f90 sourcefile~gerg2008.f90->sourcefile~implementations.f90~2 sourcefile~yaeos.f90 yaeos.f90 sourcefile~yaeos.f90->sourcefile~models.f90

Source Code

module yaeos__models_ge_nrtlhv
   !! NRTL-HV model for excess Gibbs energy
   use yaeos__constants, only: pr, R
   use yaeos__autodiff
   use yaeos__models_ge, only: GeModel

   implicit none

   type, extends(GeModel) :: NRTLHV
      !! # NRTLHV
      !! Huron-Vidal modification of the NRTL model for excess Gibbs energy.
      !!
      !! # Description
      !! Huron-Vidal modified the NRTL model, including the covolume parameter
      !! into the model. This helps when running the model coupled with a cubic
      !! equation of state. Because it can be reduced to the classical
      !! quadratic mixing rules when some parameters are defined accordingly
      !!
      !! It uses the temperature function:
      !!
      !! \[ \tau_{ij}(T) =  \frac{g_{ij}^0}{RT} + \frac{g_{ij}^T T}{RT} \]
      !! # Examples
      !!
      !! # References
      real(pr), allocatable :: b(:) 
         !! Covolume parameter
      real(pr), allocatable :: alpha(:, :) 
         !! \( \alpha \) matrix
      real(pr), allocatable :: gji0(:, :) 
         !! \( g_{ij}^0 \) matrix
      real(pr), allocatable :: gjiT(:, :) 
         !! \( g_{ij}^T \) matrix
   contains
      procedure :: excess_gibbs => excess_gibbs
   end type NRTLHV

contains

   subroutine excess_gibbs(self, n, T, Ge, GeT, GeT2, Gen, GeTn, Gen2)
      !! Calculate Excess Gibbs and its derivatives.
      use yaeos__models_ge_base, only: nrtl_hv_ge, nrtl_hv_tdep, nrtl_hv_tdep_linear
      class(NRTLHV), intent(in) :: self !! Model
      real(pr), intent(in) ::n(:) !! Moles vector
      real(pr), intent(in) :: T !! Temperature [K]
      real(pr), optional, intent(out) :: Ge !! Excess Gibbs free energy
      real(pr), optional, intent(out) :: GeT !! \(\frac{dG^E}{dT}\)
      real(pr), optional, intent(out) :: GeT2 !! \(\frac{d^2G^E}{dT^2}\)
      real(pr), optional, intent(out) :: Gen(size(n)) !! \(\frac{dG^E}{dn}\)
      real(pr), optional, intent(out) :: GeTn(size(n))
      real(pr), optional, intent(out) :: Gen2(size(n), size(n))
      real(pr) :: tau(size(n), size(n)), dtaudt(size(n), size(n)), dtaudt2(size(n), size(n))

      call nrtl_hv_tdep_linear(T, self%gjiT/R, self%gji0/R, tau, dtaudt, dtaudt2)
      call nrtl_hv_ge(n=n, T=T,&
         b=self%b, alpha=self%alpha, &
         tau=tau, dtaudt=dtaudt, dtaudt2=dtaudt2, &
         Ge=Ge, Gen=Gen, GeT=GeT, GeT2=GeT2, GeTn=GeTn, Gen2=Gen2)
   end subroutine excess_gibbs

end module yaeos__models_ge_nrtlhv