eos_water.F90       coverage:  75.00 %func     49.50 %block


     1) module EOS_Water_module
     2)  
     3)   use PFLOTRAN_Constants_module
     4) 
     5)   implicit none
     6) 
     7)   private
     8)   
     9) #include "petsc/finclude/petscsys.h"
    10) 
    11)   ! module variables
    12)   PetscReal :: constant_density
    13)   PetscReal :: constant_enthalpy
    14)   PetscReal :: constant_viscosity
    15)   PetscReal :: constant_steam_density
    16)   PetscReal :: constant_steam_enthalpy
    17) 
    18)   ! exponential
    19)   PetscReal :: exponent_reference_density
    20)   PetscReal :: exponent_reference_pressure
    21)   PetscReal :: exponent_water_compressibility
    22) 
    23)   ! In order to support generic EOS subroutines, we need the following:
    24)   ! 1. An interface declaration that defines the argument list (best to have 
    25)   !    "Dummy" appended.
    26)   ! 2. A procedure pointer that is initially set to null.  This pointer is
    27)   !    pointed to the appropriate subroutine later on (e.g. EOSWaterInit())
    28)   ! 3. An interface for derivative/non-derivative versions
    29) 
    30)   ! procedure pointer declarations
    31)   ! standard versions
    32)   procedure(EOSWaterViscosityDummy), pointer :: EOSWaterViscosityPtr => null()
    33)   procedure(EOSWaterSatPressDummy), pointer :: &
    34)     EOSWaterSaturationPressurePtr => null()
    35)   procedure(EOSWaterDensityDummy), pointer :: EOSWaterDensityPtr => null()
    36)   procedure(EOSWaterEnthalpyDummy), pointer :: EOSWaterEnthalpyPtr => null()
    37)   procedure(EOSWaterSteamDenEnthDummy), pointer :: &
    38)     EOSWaterSteamDensityEnthalpyPtr => null()
    39)   procedure(EOSWaterDensityIceDummy), pointer :: &
    40)     EOSWaterDensityIcePtr => null()
    41)   ! extended versions
    42)   procedure(EOSWaterViscosityExtDummy), pointer :: &
    43)     EOSWaterViscosityExtPtr => null()
    44)   procedure(EOSWaterDensityExtDummy), pointer :: &
    45)     EOSWaterDensityExtPtr => null()
    46)   procedure(EOSWaterEnthalpyExtDummy), pointer :: &
    47)     EOSWaterEnthalpyExtPtr => null()
    48)     
    49)   ! interface blocks
    50)   interface
    51)   ! standard versions
    52)     subroutine EOSWaterViscosityDummy(T, P, PS, dPS_dT, &
    53)                                       calculate_derivatives, VW, &
    54)                                       dVW_dT, dVW_dP, ierr)
    55)       implicit none
    56)       PetscReal, intent(in) :: T, P, PS, dPS_dT
    57)       PetscBool, intent(in) :: calculate_derivatives
    58)       PetscReal, intent(out) :: VW
    59)       PetscReal, intent(out) :: dVW_dT, dVW_dP
    60)       PetscErrorCode, intent(out) :: ierr
    61)     end subroutine EOSWaterViscosityDummy
    62)     subroutine EOSWaterSatPressDummy(T, calculate_derivatives, &
    63)                                      PS, dPS_dT, ierr)
    64)       implicit none
    65)       PetscReal, intent(in) :: T
    66)       PetscBool, intent(in) :: calculate_derivatives
    67)       PetscReal, intent(out) :: PS, dPS_dT
    68)       PetscErrorCode, intent(out) :: ierr
    69)     end subroutine EOSWaterSatPressDummy
    70)     subroutine EOSWaterDensityDummy(t,p,calculate_derivatives, &
    71)                                     dw,dwmol,dwp,dwt,ierr)
    72)       implicit none
    73)       PetscReal, intent(in) :: t
    74)       PetscReal, intent(in) :: p
    75)       PetscBool, intent(in) :: calculate_derivatives
    76)       PetscReal, intent(out) :: dw,dwmol,dwp,dwt
    77)       PetscErrorCode, intent(out) :: ierr
    78)     end subroutine EOSWaterDensityDummy
    79)     subroutine EOSWaterEnthalpyDummy(t,p,calculate_derivatives, &
    80)                                      hw,hwp,hwt,ierr)
    81)       implicit none
    82)       PetscReal, intent(in) :: t
    83)       PetscReal, intent(in) :: p
    84)       PetscBool, intent(in) :: calculate_derivatives
    85)       PetscReal, intent(out) :: hw,hwp,hwt
    86)       PetscErrorCode, intent(out) :: ierr
    87)     end subroutine EOSWaterEnthalpyDummy
    88)     subroutine EOSWaterSteamDenEnthDummy(t,p,calculate_derivatives, &
    89)                                          dg,dgmol,hg,dgp,dgt,hgp,hgt,ierr)
    90)       implicit none
    91)       PetscReal, intent(in) :: t
    92)       PetscReal, intent(in) :: p
    93)       PetscBool, intent(in) :: calculate_derivatives
    94)       PetscReal, intent(out) :: dg,dgmol,dgp,dgt
    95)       PetscReal, intent(out) :: hg,hgp,hgt
    96)       PetscErrorCode, intent(out) :: ierr 
    97)     end subroutine EOSWaterSteamDenEnthDummy
    98)     subroutine EOSWaterDensityIceDummy(t,p,calculate_derivatives, &
    99)                                        dw,dwp,dwt,ierr)
   100)       implicit none
   101)       PetscReal, intent(in) :: t
   102)       PetscReal, intent(in) :: p
   103)       PetscBool, intent(in) :: calculate_derivatives
   104)       PetscReal, intent(out) :: dw,dwp,dwt
   105)       PetscErrorCode, intent(out) :: ierr
   106)     end subroutine EOSWaterDensityIceDummy
   107)     ! Extended versions
   108)     subroutine EOSWaterViscosityExtDummy(T, P, PS, dPS_dT, aux, &
   109)                                          calculate_derivatives, VW, &
   110)                                          dVW_dT, dVW_dP, ierr)
   111)       implicit none
   112)       PetscReal, intent(in) :: T, P, PS, dPS_dT, aux(*)
   113)       PetscBool, intent(in) :: calculate_derivatives
   114)       PetscReal, intent(out) :: VW
   115)       PetscReal, intent(out) :: dVW_dT, dVW_dP
   116)       PetscErrorCode, intent(out) :: ierr
   117)     end subroutine EOSWaterViscosityExtDummy
   118)     subroutine EOSWaterDensityExtDummy(t,p,aux,calculate_derivatives, &
   119)                                        dw,dwmol,dwp,dwt,ierr)
   120)       implicit none
   121)       PetscReal, intent(in) :: t
   122)       PetscReal, intent(in) :: p
   123)       PetscReal, intent(in) :: aux(*)
   124)       PetscBool, intent(in) :: calculate_derivatives
   125)       PetscReal, intent(out) :: dw,dwmol,dwp,dwt
   126)       PetscErrorCode, intent(out) :: ierr
   127)     end subroutine EOSWaterDensityExtDummy
   128)     subroutine EOSWaterEnthalpyExtDummy(t,p,aux,calculate_derivatives, &
   129)                                         hw,hwp,hwt,ierr)
   130)       implicit none
   131)       PetscReal, intent(in) :: t
   132)       PetscReal, intent(in) :: p
   133)       PetscReal, intent(in) :: aux(*)
   134)       PetscBool, intent(in) :: calculate_derivatives
   135)       PetscReal, intent(out) :: hw,hwp,hwt
   136)       PetscErrorCode, intent(out) :: ierr
   137)     end subroutine EOSWaterEnthalpyExtDummy
   138)   end interface
   139)   
   140)   ! interfaces for derivative/non-derivative versions that are visible outside
   141)   ! the module.
   142)   ! standard versions
   143)   interface EOSWaterViscosity
   144)     procedure EOSWaterViscosityNoDerive
   145)     procedure EOSWaterViscosityDerive
   146)   end interface
   147)   interface EOSWaterSaturationPressure
   148)     procedure EOSWaterSatPresNoDerive
   149)     procedure EOSWaterSatPresDerive
   150)   end interface
   151)   interface EOSWaterDensity
   152)     procedure EOSWaterDensityNoDerive
   153)     procedure EOSWaterDensityDerive
   154)   end interface
   155)   interface EOSWaterEnthalpy
   156)     procedure EOSWaterEnthalpyNoDerive
   157)     procedure EOSWaterEnthalpyDerive
   158)   end interface
   159)   interface EOSWaterSteamDensityEnthalpy
   160)     procedure EOSWaterSteamDenEnthNoDerive
   161)     procedure EOSWaterSteamDenEnthDerive
   162)   end interface
   163)   interface EOSWaterDensityIce
   164)     procedure EOSWaterDensityIceNoDerive
   165)     procedure EOSWaterDensityIceDerive
   166)   end interface
   167)   ! Extended versions
   168)   interface EOSWaterViscosityExt
   169)     procedure EOSWaterViscosityExtNoDerive
   170)     procedure EOSWaterViscosityExtDerive
   171)   end interface
   172)   interface EOSWaterDensityExt
   173)     procedure EOSWaterDensityExtNoDerive
   174)     procedure EOSWaterDensityExtDerive
   175) !geh: very useful for debuggin
   176) !    procedure EOSWaterDensityExtNumericalDerive
   177)   end interface
   178)   interface EOSWaterEnthalpyExt
   179)     procedure EOSWaterEnthalpyExtNoDerive
   180)     procedure EOSWaterEnthalpyExtDerive
   181)   end interface
   182) 
   183)   ! the "public" definition that makes subroutines visible outside.
   184)   public :: EOSWaterInit, &
   185)             EOSWaterVerify, &
   186)             EOSWaterViscosity, &
   187)             EOSWaterSaturationPressure, &
   188)             EOSWaterDensity, &
   189)             EOSWaterEnthalpy, &
   190)             EOSWaterSteamDensityEnthalpy, &
   191)             EOSWaterDuanMixture, &
   192)             EOSWaterViscosityNaCl, &
   193)             EOSWaterInternalEnergyIce, &
   194)             EOSWaterDensityIcePainter, &
   195)             EOSWaterSaturationTemperature, &
   196)             EOSWaterDensityIce, &
   197)             EOSWaterDensityTGDPB01, &
   198)             EOSWaterViscosityExt, &
   199)             EOSWaterDensityExt, &
   200)             EOSWaterEnthalpyExt, &
   201)             EOSWaterInputRecord
   202) 
   203)   public :: EOSWaterSetDensity, &
   204)             EOSWaterSetEnthalpy, &
   205)             EOSWaterSetViscosity, &
   206)             EOSWaterSetSteamDensity, &
   207)             EOSWaterSetSteamEnthalpy
   208)             
   209)             
   210)   public :: TestEOSWaterBatzleAndWang, &
   211)             EOSWaterTest
   212)  
   213)   contains
   214) 
   215) ! ************************************************************************** !
   216) 
   217) subroutine EOSWaterInit()
   218) 
   219)   implicit none
   220)   
   221)   constant_density = UNINITIALIZED_DOUBLE
   222)   constant_viscosity = UNINITIALIZED_DOUBLE
   223)   constant_enthalpy = UNINITIALIZED_DOUBLE
   224)   constant_steam_density = UNINITIALIZED_DOUBLE
   225)   constant_steam_enthalpy = UNINITIALIZED_DOUBLE
   226)   exponent_reference_density = UNINITIALIZED_DOUBLE
   227)   exponent_reference_pressure = UNINITIALIZED_DOUBLE
   228)   exponent_water_compressibility = UNINITIALIZED_DOUBLE
   229)   
   230)   ! standard versions  
   231)   EOSWaterDensityPtr => EOSWaterDensityIFC67
   232)   EOSWaterEnthalpyPtr => EOSWaterEnthalpyIFC67
   233)   EOSWaterViscosityPtr => EOSWaterViscosity1
   234)   EOSWaterSaturationPressurePtr => EOSWaterSaturationPressureIFC67
   235)   EOSWaterSteamDensityEnthalpyPtr => EOSWaterSteamDensityEnthalpyIFC67
   236)   EOSWaterDensityIcePtr => EOSWaterDensityIcePainter
   237)   
   238)   ! extended versions
   239)   EOSWaterViscosityExtPtr => EOSWaterViscosityKestinExt
   240)   EOSWaterDensityExtPtr => EOSWaterDensityBatzleAndWangExt
   241)   
   242) end subroutine EOSWaterInit
   243) 
   244) ! ************************************************************************** !
   245) 
   246) subroutine EOSWaterVerify(ierr,error_string)
   247) 
   248)   implicit none
   249)   
   250)   PetscErrorCode, intent(out) :: ierr
   251)   character(len=MAXSTRINGLENGTH), intent(out) :: error_string
   252)   
   253)   ierr = 0
   254)   
   255)   error_string = ''
   256)   if ((associated(EOSWaterDensityPtr,EOSWaterDensityIFC67) .and. &
   257)         Initialized(constant_density)) .or. &
   258)       (associated(EOSWaterEnthalpyPtr,EOSWaterEnthalpyIFC67) .and. &
   259)         Initialized(constant_enthalpy)) &
   260)      ) then
   261)     ierr = 1
   262)   endif
   263) 
   264)   if (associated(EOSWaterDensityPtr,EOSWaterDensityConstant) .and. &
   265)       Uninitialized(constant_density)) then
   266)     error_string = trim(error_string) // &
   267)       ' CONSTANT density not set.'
   268)     ierr = 1
   269)   endif
   270)   
   271)   if (associated(EOSWaterEnthalpyPtr,EOSWaterEnthalpyConstant) .and. &
   272)       Uninitialized(constant_enthalpy)) then
   273)     error_string = trim(error_string) // &
   274)       ' CONSTANT enthalpy not set.'
   275)     ierr = 1
   276)   endif
   277)   
   278)   if (associated(EOSWaterDensityPtr,EOSWaterDensityExponential) .and. &
   279)       (Uninitialized(exponent_reference_density) .or. & 
   280)        Uninitialized(exponent_reference_pressure) .or. &
   281)        Uninitialized(exponent_water_compressibility))) then
   282)     error_string = trim(error_string) // &
   283)       ' Exponential parameters incorrect.'
   284)     ierr = 1
   285)   endif
   286) 
   287)   if ((associated(EOSWaterViscosityPtr, &
   288)                   EOSWaterViscosityConstant) .and. &
   289)        Uninitialized(constant_viscosity)) .or. &
   290)       (associated(EOSWaterViscosityPtr, &
   291)                   EOSWaterViscosity1) .and. &
   292)        Initialized(constant_viscosity))) then
   293)     ierr = 1
   294)   endif
   295)   
   296)   if (associated(EOSWaterSteamDensityEnthalpyPtr, &
   297)                  EOSWaterSteamDenEnthConstant) .and. &
   298)       (Uninitialized(constant_steam_density) .or. &
   299)        Uninitialized(constant_steam_enthalpy))) then
   300)     if (Uninitialized(constant_steam_density)) then
   301)       error_string = trim(error_string) // &
   302)         ' CONSTANT steam density not set.'
   303)     endif
   304)     if (Uninitialized(constant_steam_enthalpy)) then
   305)       error_string = trim(error_string) // &
   306)         ' CONSTANT steam enthalpy not set.'
   307)     endif
   308)     ierr = 1
   309)   endif
   310)   
   311) end subroutine EOSWaterVerify
   312) 
   313) ! ************************************************************************** !
   314) 
   315) subroutine EOSWaterSetDensity(keyword,aux)
   316) 
   317)   implicit none
   318)   
   319)   character(len=*) :: keyword
   320)   PetscReal, optional :: aux(*)
   321)   
   322)   select case(keyword)
   323)     case('CONSTANT')
   324)       constant_density = aux(1)  
   325)       EOSWaterDensityPtr => EOSWaterDensityConstant
   326)     case('DEFAULT','IFC67')
   327)       EOSWaterDensityPtr => EOSWaterDensityIFC67
   328)       EOSWaterDensityExtPtr => EOSWaterDensityBatzleAndWangExt
   329)     case('EXPONENTIAL')
   330)       exponent_reference_density = aux(1)
   331)       exponent_reference_pressure = aux(2)
   332)       exponent_water_compressibility = aux(3)  
   333)       EOSWaterDensityPtr => EOSWaterDensityExponential    
   334)     case('TGDPB01')
   335)       EOSWaterDensityPtr => EOSWaterDensityTGDPB01
   336)     case('PAINTER')
   337)       EOSWaterDensityPtr => EOSWaterDensityPainter
   338)     case('BATZLE_AND_WANG')
   339)       EOSWaterDensityPtr => EOSWaterDensityBatzleAndWang
   340)       EOSWaterDensityExtPtr => EOSWaterDensityBatzleAndWangExt
   341)     case default
   342)       print *, 'Unknown pointer type "' // trim(keyword) // &
   343)         '" in EOSWaterSetDensity().'      
   344)       stop
   345)   end select
   346)   
   347) end subroutine EOSWaterSetDensity
   348) 
   349) ! ************************************************************************** !
   350) 
   351) subroutine EOSWaterSetEnthalpy(keyword,aux)
   352) 
   353)   implicit none
   354)   
   355)   character(len=*) :: keyword
   356)   PetscReal, optional :: aux(*)
   357)   
   358)   select case(keyword)
   359)     case('CONSTANT')
   360)       constant_enthalpy = aux(1)  
   361)       EOSWaterEnthalpyPtr => EOSWaterEnthalpyConstant
   362)     case('DEFAULT','IFC67')
   363)       EOSWaterEnthalpyPtr => EOSWaterEnthalpyIFC67
   364)     case('PAINTER')
   365)       EOSWaterEnthalpyPtr => EOSWaterEnthalpyPainter      
   366)     case default
   367)       print *, 'Unknown pointer type "' // trim(keyword) // &
   368)         '" in EOSWaterSetEnthalpy().'
   369)       stop
   370)   end select
   371)   
   372)   
   373) end subroutine EOSWaterSetEnthalpy
   374) 
   375) ! ************************************************************************** !
   376) 
   377) subroutine EOSWaterSetViscosity(keyword,aux)
   378) 
   379)   implicit none
   380)   
   381)   character(len=*) :: keyword
   382)   PetscReal, optional :: aux(*)
   383)   
   384)   select case(keyword)
   385)     case('CONSTANT')
   386)       constant_viscosity = aux(1)  
   387)       EOSWaterViscosityPtr => EOSWaterViscosityConstant
   388)     case('DEFAULT')
   389)       EOSWaterViscosityPtr => EOSWaterViscosity1
   390)       EOSWaterViscosityExtPtr => EOSWaterViscosityKestinExt
   391)     case('BATZLE_AND_WANG')
   392)       EOSWaterViscosityPtr => EOSWaterViscosityBatzleAndWang
   393)       EOSWaterViscosityExtPtr => EOSWaterViscosityBatzleAndWangExt
   394)     case default
   395)       print *, 'Unknown pointer type "' // trim(keyword) // &
   396)         '" in EOSWaterSetViscosity().'  
   397)       stop 
   398)   end select
   399)   
   400) end subroutine EOSWaterSetViscosity
   401) 
   402) 
   403) ! ************************************************************************** !
   404) 
   405) subroutine EOSWaterSetSteamDensity(keyword,aux)
   406) 
   407)   implicit none
   408)   
   409)   character(len=*) :: keyword
   410)   PetscReal, optional :: aux(*)
   411)   
   412)   select case(keyword)
   413)     case('CONSTANT')
   414)       constant_steam_density = aux(1)  
   415)       EOSWaterSteamDensityEnthalpyPtr => EOSWaterSteamDenEnthConstant
   416)     case('IFC67')
   417)       EOSWaterSteamDensityEnthalpyPtr => EOSWaterSteamDensityEnthalpyIFC67
   418)     case default
   419)       print *, 'Unknown pointer type "' // trim(keyword) // &
   420)         '" in EOSWaterSetSteamDensity().'
   421)       stop
   422)   end select
   423)   
   424) end subroutine EOSWaterSetSteamDensity
   425) 
   426) ! ************************************************************************** !
   427) 
   428) subroutine EOSWaterSetSteamEnthalpy(keyword,aux)
   429) 
   430)   implicit none
   431)   
   432)   character(len=*) :: keyword
   433)   PetscReal, optional :: aux(*)
   434)   
   435)   select case(keyword)
   436)     case('CONSTANT')
   437)       constant_steam_enthalpy = aux(1)  
   438)       EOSWaterSteamDensityEnthalpyPtr => EOSWaterSteamDenEnthConstant
   439)     case('DEFAULT','IFC67')
   440)       EOSWaterSteamDensityEnthalpyPtr => EOSWaterSteamDensityEnthalpyIFC67
   441)     case default
   442)       print *, 'Unknown pointer type "' // trim(keyword) // &
   443)         '" in EOSWaterSetSteamEnthalpy().'
   444)       stop
   445)   end select
   446)   
   447) end subroutine EOSWaterSetSteamEnthalpy
   448) 
   449) ! ************************************************************************** !
   450) 
   451) subroutine EOSWaterViscosityNoDerive(T, P, PS, VW, ierr)
   452) 
   453)   implicit none
   454) 
   455)   PetscReal, intent(in) :: T, P, PS ! temperature, pressure, saturation_press
   456)   PetscReal, intent(out) :: VW ! water viscosity
   457)   PetscErrorCode, intent(out) :: ierr
   458)   
   459)   PetscReal :: dPS_dT ! derivative of PS with respect to temp
   460)   PetscReal :: dum1, dum2
   461)   
   462)   dPS_dT = 0.d0
   463)   call EOSWaterViscosityPtr(T, P, PS, dPS_dT, PETSC_FALSE, VW, &
   464)                             dum1, dum2, ierr)
   465)   
   466) end subroutine EOSWaterViscosityNoDerive
   467) 
   468) ! ************************************************************************** !
   469) 
   470) subroutine EOSWaterViscosityDerive(T, P, PS, dPS_dT, VW, dVW_dT, &
   471)                                    dVW_dP, ierr)
   472) 
   473)   implicit none
   474) 
   475)   PetscReal, intent(in) :: T, P, PS ! temperature, pressure, saturation_press
   476)   PetscReal, intent(in) :: dPS_dT ! derivative of PS with respect to temp
   477)   PetscReal, intent(out) :: VW ! water viscosity
   478)   PetscReal, intent(out) :: dVW_dT, dVW_dP ! derivatives
   479)   PetscErrorCode, intent(out) :: ierr
   480)   
   481)   call EOSWaterViscosityPtr(T, P, PS, dPS_dT, PETSC_TRUE, VW, &
   482)                             dVW_dT, dVW_dP, ierr)
   483)   
   484) end subroutine EOSWaterViscosityDerive
   485) 
   486) ! ************************************************************************** !
   487) 
   488) subroutine EOSWaterSatPresNoDerive(T, PS, ierr)
   489) 
   490)   implicit none
   491) 
   492)   PetscReal, intent(in) :: T ! temperature
   493)   PetscReal, intent(out) :: PS ! Saturation pres. and derivative
   494)   PetscErrorCode, intent(out) :: ierr
   495)   
   496)   PetscReal :: dummy
   497)   
   498)   call EOSWaterSaturationPressurePtr(T, PETSC_FALSE, PS, dummy, ierr)
   499)   
   500) end subroutine EOSWaterSatPresNoDerive
   501) 
   502) ! ************************************************************************** !
   503) 
   504) subroutine EOSWaterSatPresDerive(T, PS, dPS_dT, ierr)
   505) 
   506)   implicit none
   507) 
   508)   PetscReal, intent(in) :: T ! temperature
   509)   PetscReal, intent(out) :: PS, dPS_dT ! Saturation pres. and derivative
   510)   PetscErrorCode, intent(out) :: ierr
   511)   
   512)   call EOSWaterSaturationPressurePtr(T, PETSC_TRUE, PS, dPS_dT, ierr)
   513)   
   514) end subroutine EOSWaterSatPresDerive
   515) 
   516) ! ************************************************************************** !
   517) 
   518) subroutine EOSWaterDensityNoDerive(t,p,dw,dwmol,ierr)
   519) 
   520)   implicit none
   521) 
   522)   PetscReal, intent(in) :: t
   523)   PetscReal, intent(in) :: p
   524)   PetscReal, intent(out) :: dw,dwmol
   525)   PetscErrorCode, intent(out) :: ierr
   526)   
   527)   PetscReal :: dum1, dum2
   528)   
   529)   call EOSWaterDensityPtr(t,p,PETSC_FALSE,dw,dwmol,dum1,dum2,ierr)
   530)   
   531) end subroutine EOSWaterDensityNoDerive
   532) 
   533) ! ************************************************************************** !
   534) 
   535) subroutine EOSWaterDensityDerive(t,p,dw,dwmol,dwp,dwt,ierr)
   536) 
   537)   implicit none
   538) 
   539)   PetscReal, intent(in) :: t
   540)   PetscReal, intent(in) :: p
   541)   PetscReal, intent(out) :: dw,dwmol,dwp,dwt
   542)   PetscErrorCode, intent(out) :: ierr
   543)   
   544)   call EOSWaterDensityPtr(t,p,PETSC_TRUE,dw,dwmol,dwp,dwt,ierr)
   545)   
   546) end subroutine EOSWaterDensityDerive
   547) 
   548) ! ************************************************************************** !
   549) 
   550) subroutine EOSWaterEnthalpyNoDerive(t,p,hw,ierr)
   551) 
   552)   implicit none
   553) 
   554)   PetscReal, intent(in) :: t
   555)   PetscReal, intent(in) :: p
   556)   PetscReal, intent(out) :: hw
   557)   PetscErrorCode, intent(out) :: ierr
   558)   
   559)   PetscReal :: dum1, dum2
   560)   
   561)   call EOSWaterEnthalpyPtr(t,p,PETSC_FALSE,hw,dum1,dum2,ierr)
   562)   
   563) end subroutine EOSWaterEnthalpyNoDerive
   564) 
   565) ! ************************************************************************** !
   566) 
   567) subroutine EOSWaterEnthalpyDerive(t,p,hw,hwp,hwt,ierr)
   568)   implicit none
   569) 
   570)   PetscReal, intent(in) :: t
   571)   PetscReal, intent(in) :: p
   572)   PetscReal, intent(out) :: hw,hwp,hwt
   573)   PetscErrorCode, intent(out) :: ierr
   574)   
   575)   call EOSWaterEnthalpyPtr(t,p,PETSC_TRUE,hw,hwp,hwt,ierr)
   576)   
   577) end subroutine EOSWaterEnthalpyDerive
   578) 
   579) ! ************************************************************************** !
   580) 
   581) subroutine EOSWaterViscosityExtNoDerive(T, P, PS, aux, VW, ierr)
   582) 
   583)   implicit none
   584) 
   585)   PetscReal, intent(in) :: T, P, PS ! temperature, pressure, saturation_press
   586)   PetscReal, intent(in) :: aux(*)
   587)   PetscReal, intent(out) :: VW ! water viscosity
   588)   PetscErrorCode, intent(out) :: ierr
   589)   
   590)   PetscReal :: dPS_dT ! derivative of PS with respect to temp
   591)   PetscReal :: dum1, dum2
   592)   
   593)   dPS_dT = 0.d0
   594)   call EOSWaterViscosityExtPtr(T, P, PS, dPS_dT, aux, PETSC_FALSE, VW, &
   595)                                dum1, dum2, ierr)
   596)   
   597) end subroutine EOSWaterViscosityExtNoDerive
   598) 
   599) ! ************************************************************************** !
   600) 
   601) subroutine EOSWaterViscosityExtDerive(T, P, PS, dPS_dT, aux, VW, dVW_dT, &
   602)                                       dVW_dP, ierr)
   603) 
   604)   implicit none
   605) 
   606)   PetscReal, intent(in) :: T, P, PS ! temperature, pressure, saturation_press
   607)   PetscReal, intent(in) :: dPS_dT ! derivative of PS with respect to temp
   608)   PetscReal, intent(in) :: aux(*)
   609)   PetscReal, intent(out) :: VW ! water viscosity
   610)   PetscReal, intent(out) :: dVW_dT, dVW_dP ! derivatives
   611)   PetscErrorCode, intent(out) :: ierr
   612)   
   613)   call EOSWaterViscosityExtPtr(T, P, PS, dPS_dT, aux, PETSC_TRUE, VW, &
   614)                                dVW_dT, dVW_dP, ierr)
   615)   
   616) end subroutine EOSWaterViscosityExtDerive
   617) 
   618) ! ************************************************************************** !
   619) 
   620) subroutine EOSWaterDensityExtNoDerive(t,p,aux,dw,dwmol,ierr)
   621) 
   622)   implicit none
   623) 
   624)   PetscReal, intent(in) :: t
   625)   PetscReal, intent(in) :: p
   626)   PetscReal, intent(in) :: aux(*)
   627)   PetscReal, intent(out) :: dw,dwmol
   628)   PetscErrorCode, intent(out) :: ierr
   629)   
   630)   PetscReal :: dum1, dum2
   631)   
   632)   call EOSWaterDensityExtPtr(t,p,aux,PETSC_FALSE,dw,dwmol,dum1,dum2,ierr)
   633)   
   634) end subroutine EOSWaterDensityExtNoDerive
   635) 
   636) ! ************************************************************************** !
   637) 
   638) subroutine EOSWaterDensityExtDerive(t,p,aux,dw,dwmol,dwp,dwt,ierr)
   639)   implicit none
   640) 
   641)   PetscReal, intent(in) :: t
   642)   PetscReal, intent(in) :: p
   643)   PetscReal, intent(in) :: aux(*)
   644)   PetscReal, intent(out) :: dw,dwmol,dwp,dwt
   645)   PetscErrorCode, intent(out) :: ierr
   646)   
   647)   call EOSWaterDensityExtPtr(t,p,aux,PETSC_TRUE,dw,dwmol,dwp,dwt,ierr)
   648)   
   649) end subroutine EOSWaterDensityExtDerive
   650) 
   651) ! ************************************************************************** !
   652) 
   653) subroutine EOSWaterEnthalpyExtNoDerive(t,p,aux,hw,ierr)
   654) 
   655)   implicit none
   656) 
   657)   PetscReal, intent(in) :: t
   658)   PetscReal, intent(in) :: p
   659)   PetscReal, intent(in) :: aux(*)  
   660)   PetscReal, intent(out) :: hw
   661)   PetscErrorCode, intent(out) :: ierr
   662)   
   663)   PetscReal :: dum1, dum2
   664)   
   665)   call EOSWaterEnthalpyExtPtr(t,p,aux,PETSC_FALSE,hw,dum1,dum2,ierr)
   666)   
   667) end subroutine EOSWaterEnthalpyExtNoDerive
   668) 
   669) ! ************************************************************************** !
   670) 
   671) subroutine EOSWaterEnthalpyExtDerive(t,p,aux,hw,hwp,hwt,ierr)
   672)   implicit none
   673) 
   674)   PetscReal, intent(in) :: t
   675)   PetscReal, intent(in) :: p
   676)   PetscReal, intent(in) :: aux(*)  
   677)   PetscReal, intent(out) :: hw,hwp,hwt
   678)   PetscErrorCode, intent(out) :: ierr
   679)   
   680)   call EOSWaterEnthalpyExtPtr(t,p,aux,PETSC_TRUE,hw,hwp,hwt,ierr)
   681)   
   682) end subroutine EOSWaterEnthalpyExtDerive
   683) 
   684) ! ************************************************************************** !
   685) 
   686) subroutine EOSWaterViscosity1(T, P, PS, dPS_dT, calculate_derivatives, &
   687)                               VW, dVW_dT, dVW_dP, ierr)
   688) 
   689) ! Calculates the viscosity of water and derivatives as a function of 
   690) ! temperature, pressure, and saturation pressure.
   691) 
   692)   implicit none
   693)     
   694)   PetscReal, intent(in) :: T, P, PS ! temperature, pressure, saturation_press
   695)   PetscReal, intent(in) :: dPS_dT ! derivative of PS with respect to temp
   696)   PetscBool, intent(in) :: calculate_derivatives
   697)   PetscReal, intent(out) :: VW ! water viscosity
   698)   PetscReal, intent(out) :: dVW_dT, dVW_dP ! derivatives
   699)   PetscErrorCode, intent(out) :: ierr
   700)   
   701)   PetscReal :: EX, PHI, AM, pwr, aln10
   702)   
   703)   EX  = 247.8d0/(T+133.15d0)
   704)   PHI = 1.0467d0*(T-31.85d0)
   705)   !geh: added max(P,PS)-PS in place of P-PS
   706)   !geh: here P should be the maximum of Pl and Pg
   707)   AM  = 1.d0+PHI*(max(P,PS)-PS)*1.d-11
   708)   pwr = 10.d0**EX
   709)   VW = 1.d-7*AM*241.4d0*pwr
   710)     
   711)   if (calculate_derivatives) then
   712)     aln10 = log(10.d0)
   713)     dVW_dT = VW/AM*1.d-11* &
   714)             ! dAM_PHI_dT       dAM_PS_dT
   715)             (1.0467d0*(P-PS) - PHI*dPS_dT) - &
   716)             ! dpwr_EX_dT
   717)             VW*aln10*247.8d0/(T+133.15d0)**2
   718)     dVW_dP = VW/AM*PHI*1.d-11
   719)   else
   720)     dVW_dT = UNINITIALIZED_DOUBLE
   721)     dVW_dP = UNINITIALIZED_DOUBLE
   722)   endif
   723)   ierr = 0
   724)  
   725) end subroutine EOSWaterViscosity1
   726) 
   727) ! ************************************************************************** !
   728) 
   729) subroutine EOSWaterViscosityConstant(T, P, PS, dPS_dT, &
   730)                                      calculate_derivatives, &
   731)                                       VW, dVW_dT, dVW_dP, ierr)
   732) 
   733) ! Calculates the viscosity of water and derivatives as a function of 
   734) ! temperature, pressure, and saturation pressure.
   735) 
   736)   implicit none
   737)     
   738)   PetscReal, intent(in) :: T, P, PS ! temperature, pressure, saturation_press
   739)   PetscReal, intent(in) :: dPS_dT ! derivative of PS with respect to temp
   740)   PetscBool, intent(in) :: calculate_derivatives
   741)   PetscReal, intent(out) :: VW ! water viscosity
   742)   PetscReal, intent(out) :: dVW_dT, dVW_dP ! derivatives
   743)   PetscErrorCode, intent(out) :: ierr
   744)   
   745)   VW = constant_viscosity
   746)   
   747)   dVW_dT = 0.d0
   748)   dVW_dP = 0.d0
   749)   
   750) end subroutine EOSWaterViscosityConstant
   751) 
   752) ! ************************************************************************** !
   753) 
   754) subroutine EOSWaterSaturationPressureIFC67(T, calculate_derivatives, &
   755)                                            PS, dPS_dT, ierr)
   756) 
   757)   implicit none
   758) 
   759)   PetscReal, intent(in) :: T ! temperature
   760)   PetscBool, intent(in) :: calculate_derivatives
   761)   PetscReal, intent(out) :: PS, dPS_dT ! Saturation pres. and derivative
   762)   PetscErrorCode, intent(out) :: ierr
   763)   
   764)   PetscReal, save, dimension(9) :: A(9)
   765)   PetscReal :: TC, SC, PCAP, E1, E2
   766)   PetscReal :: one_m_tc, one_m_tc_sq, E2_bottom
   767)   PetscReal :: dTC_dT, dSC_dTC, dE1_dTC, dE2_dTC, dPC_dSC, dPC_dTC
   768)     
   769)   DATA A/ &
   770)     -7.691234564d0,-2.608023696d1,-1.681706546d2,6.423285504d1, &
   771)     -1.189646225d2,4.167117320d0,2.097506760E1,1.d9,6.d0/
   772) 
   773)   if (T .GT. 500.d0) then
   774)     ierr = 1
   775)     return
   776)   end if
   777)   TC = (T+273.15d0)/H2O_CRITICAL_TEMPERATURE
   778)   one_m_tc = 1.d0-TC
   779)   one_m_tc_sq = one_m_tc*one_m_tc
   780)   SC = A(1)*one_m_tc+A(2)*one_m_tc_sq+A(3)*one_m_tc**3.d0+ &
   781)         A(4)*one_m_tc**4.d0+A(5)*one_m_tc**5.d0
   782)   E1 = TC*(1.d0+A(6)*one_m_tc+A(7)*one_m_tc_sq)
   783)   E2_bottom = A(8)*one_m_tc_sq+A(9)
   784)   E2 = one_m_tc/E2_bottom
   785)   PCAP = EXP(SC/E1-E2)
   786)    
   787)   PS = PCAP*H2O_CRITICAL_PRESSURE
   788)     
   789)   if (calculate_derivatives) then
   790)     dTC_dT = 1.d0/H2O_CRITICAL_TEMPERATURE
   791)     dSC_dTC = -A(1)-2.d0*A(2)*one_m_tc-3.d0*A(3)*one_m_tc_sq- &
   792)               4.d0*A(4)*one_m_tc**3.-5.d0*A(5)*one_m_tc**4.
   793)     dE1_dTC = (1.d0+A(6)*one_m_tc+A(7)*one_m_tc_sq)+ &
   794)               TC*(-A(6)-2.d0*A(7)*one_m_tc)
   795)     dE2_dTC = -1.d0/E2_bottom+one_m_tc/(E2_bottom*E2_bottom)*2.d0*one_m_tc
   796)     dPC_dTC = (-SC/(E1*E1)*dE1_dTC-dE2_dTC)*PCAP
   797)     dPC_dSC = 1.d0/E1*PCAP
   798)     dPS_dT = (dPC_dSC*dSC_dTC+dPC_dTC)*dTC_dT*H2O_CRITICAL_PRESSURE
   799)   else
   800)     dPS_dT = UNINITIALIZED_DOUBLE
   801)   endif
   802)   ierr = 0
   803) 
   804) end subroutine EOSWaterSaturationPressureIFC67
   805) 
   806) ! ************************************************************************** !
   807) 
   808) subroutine EOSWaterDensityIFC67(t,p,calculate_derivatives,dw,dwmol, &
   809)                                 dwp,dwt,ierr)
   810) 
   811) !  This subroutine calculates water and steam-gas mixture properties.
   812) !  The water and steam properties are valid in the range of:
   813) !
   814) !            0 < p < 165.4 * 10^5 pascals (165.4 bars)
   815) !            0 < t < 350 centigrade (623.15 Kelvin)
   816) !
   817) !  The properties cover densities, enthalpies, internal energies,
   818) !  and partial derivatives of these quanties with respect to
   819) !  pressure and temperature.
   820) !
   821) !  For saturated fluid, it will also calculate water saturation
   822) !  temperature for the specified pressure using Newton-Raphson and
   823) !  the derivative dts/dp (=tsp) or Ps for a given temperature.
   824) !
   825) !  Ref.: International Formulation Committee of the Sixth International
   826) !       Conference on Properties of Steam (1967).
   827) 
   828)   implicit none
   829)   
   830)   PetscReal, intent(in) :: t   ! Temperature in centigrade
   831)   PetscReal, intent(in) :: p   ! Pressure in Pascals
   832)   PetscBool, intent(in) :: calculate_derivatives
   833)   PetscReal, intent(out) :: dw,dwmol,dwp,dwt
   834)   PetscErrorCode, intent(out) :: ierr
   835)   
   836)   PetscInt :: i
   837)     
   838)   PetscReal, save :: aa(0:22)
   839)   PetscReal, save :: a1,a2,a3,a4,a5,a6,a7,a8,a9,a10,a11,a12
   840)   
   841)   PetscReal :: beta,beta2x,beta4,theta,utheta,theta2x,theta18,theta20
   842)   PetscReal :: xx,yy,zz
   843)   PetscReal :: u0,u1,u2,u3,u4,u5,u6,u7,u8,u9
   844)   PetscReal :: tempreal
   845) !  PetscReal :: v0_1, v1_1, v2_1, v3_1, v4_1
   846) !  PetscReal :: v1_2, v2_2, v3_2, v4_2, v20_2, v40_2
   847) !  PetscReal :: v1_3, v2_3, v3_3, v4_3
   848) !  PetscReal :: v1_4, v2_4, v3_4
   849) !  PetscReal :: v1_5, v2_5
   850) !  PetscReal :: v1_6
   851) !  PetscReal :: term1,term2,term2t,term3,term3t,term3p,term4,term4t,term4p, &
   852) !               term5,term5t,term5p,term6,term6t,term6p,term7,term7t,term7p
   853) !  PetscReal :: dv2t,dv2p,dv3t
   854)   PetscReal :: vr,ypt,yptt,zpt,zpp,vrpt,vrpp,cnv
   855)   PetscReal :: tc1,pc1,vc1,utc1,upc1,vc1mol
   856)   PetscReal :: d2z_dp2    ! 2nd derivative of z w.r.t. pressure
   857)   PetscReal :: d2vr_dp2   ! 2nd derivative of vr w.r.t. pressure
   858)   PetscReal :: d2wmol_dp2 ! 2nd derivative of density w.r.t. pressure
   859)   PetscReal, parameter :: zero = 0.d0
   860)   PetscReal, parameter :: one = 1.d0
   861)   PetscReal, parameter :: two = 2.d0
   862)   PetscReal, parameter :: three = 3.d0
   863)   PetscReal, parameter :: four = 4.d0
   864)   PetscReal, parameter :: five = 5.d0
   865)   PetscReal, parameter :: six = 6.d0
   866)   PetscReal, parameter :: seven = 7.d0
   867)   PetscReal, parameter :: eight = 8.d0
   868)   PetscReal, parameter :: nine = 9.d0
   869)   PetscReal, parameter :: ten = 10.d0
   870)   
   871)   data aa/ &
   872) !-----data aa0,aa1,aa2,aa3/
   873)         6.824687741d03,-5.422063673d02,-2.096666205d04, 3.941286787d04, &
   874) !-----data aa4,aa5,aa6,aa7/
   875)         -6.733277739d04, 9.902381028d04,-1.093911774d05, 8.590841667d04, &
   876) !-----data aa8,aa9,aa10,aa11/
   877)         -4.511168742d04, 1.418138926d04,-2.017271113d03, 7.982692717d00, &
   878) !-----data aa12,aa13,aa14,aa15/
   879)         -2.616571843d-2, 1.522411790d-3, 2.284279054d-2, 2.421647003d02, &
   880) !-----data aa16,aa17,aa18,aa19/
   881)         1.269716088d-10,2.074838328d-7, 2.174020350d-8, 1.105710498d-9, &
   882) !-----data aa20,aa21,aa22/    
   883)         1.293441934d01, 1.308119072d-5, 6.047626338d-14/
   884) 
   885)   data a1,a2,a3,a4/ &
   886)   8.438375405d-1, 5.362162162d-4, 1.720000000d00, 7.342278489d-2/
   887)   data a5,a6,a7,a8/ &
   888)   4.975858870d-2, 6.537154300d-1, 1.150000000d-6, 1.510800000d-5/
   889)   data a9,a10,a11,a12/ &
   890)   1.418800000d-1, 7.002753165d00, 2.995284926d-4, 2.040000000d-1/
   891)     
   892)   ierr = 0
   893) 
   894)   tc1 = H2O_CRITICAL_TEMPERATURE    ! K
   895)   pc1 = H2O_CRITICAL_PRESSURE     ! Pa 
   896)   vc1 = 0.00317d0  ! m^3/kg
   897)   utc1 = one/tc1   ! 1/C
   898)   upc1 = one/pc1   ! 1/Pa
   899)   vc1mol = vc1*FMWH2O ! m^3/kmol
   900)     
   901)   theta = (t+273.15d0)*utc1
   902)   theta2x = theta*theta
   903)   theta18 = theta**18.d0
   904)   theta20 = theta18*theta2x
   905)     
   906)   beta = p*upc1
   907)   beta2x = beta*beta
   908)   beta4  = beta2x*beta2x
   909)     
   910)   yy = one-a1*theta2x-a2*theta**(-6.d0)
   911)   xx = a3*yy*yy-two*(a4*theta-a5*beta)
   912)     
   913) !   Note: xx may become negative near the critical point-pcl.
   914)   if (xx.gt.zero) then
   915)     xx = sqrt(xx)
   916)   else
   917)     write(*,*) 'Warning: negative term in density (eos_water.F90:&
   918)       &EOSWaterDensityIFC67):'
   919)     write(*,*) 't= ',t,' p= ',p,' xx= ',xx
   920)     ierr = 1
   921)     xx = 1.d-6               !set arbitrarily
   922)   end if
   923)   zz = yy + xx                                     
   924)   u0 = -five/17.d0
   925)   u1 = aa(11)*a5*zz**u0
   926)   u2 = one/(a8+theta**11.d0)
   927)   u3 = aa(17)+(two*aa(18)+three*aa(19)*beta)*beta
   928)   u4 = one/(a7+theta18*theta)
   929)   u5 = (a10+beta)**(-4.d0)
   930)   u6 = a11-three*u5
   931)   u7 = aa(20)*theta18*(a9+theta2x)
   932)   u8 = aa(15)*(a6-theta)**9.d0
   933)     
   934)   vr = u1+aa(12)+theta*(aa(13)+aa(14)*theta)+u8*(a6-theta) &
   935)         +aa(16)*u4-u2*u3-u6*u7+(three*aa(21)*(a12-theta) &
   936)         +four*aa(22)*beta/theta20)*beta2x
   937)     
   938)   dwmol = one/(vr*vc1mol) ! kmol/m^3
   939)   dw = one/(vr*vc1) ! kg/m^3
   940) 
   941)   ! ypt used for enthalpy even if derivative not calculated
   942)   ypt = six*a2*theta**(-7.d0)-two*a1*theta
   943)   
   944)   !---calculate derivatives for water density
   945)   if (calculate_derivatives) then
   946)     zpt = ypt+(a3*yy*ypt-a4)/xx
   947)     zpp = a5/xx
   948)     u9 = u0*u1/zz
   949)     vrpt = u9*zpt+aa(13)+two*aa(14)*theta-ten*u8 &
   950)         -19.d0*aa(16)*u4*u4*theta18+11.d0*u2*u2*u3*theta**10.d0 &
   951)         -aa(20)*u6*(18.d0*a9*theta18+20.d0*theta20)/theta &
   952)         -(three*aa(21)+80.d0*aa(22)*beta/(theta20*theta))*beta2x
   953)     
   954)     vrpp = u9*zpp-u2*(two*aa(18)+six*aa(19)*beta)-12.d0*u7*u5/ &
   955)         (a10+beta)+(six*aa(21)*(a12-theta)+12.d0*aa(22)*beta/ &
   956)         theta20)*beta
   957)     
   958)     cnv = -one/(vc1mol*vr*vr)
   959)     dwt = cnv*vrpt*utc1 ! kmol/m^3/C
   960)     dwp = cnv*vrpp*upc1 ! kmol/m^3/Pa
   961) 
   962)     ! 2nd derivative w.r.t pressure
   963)     d2z_dp2 = -a5*a5/xx/xx/xx
   964) 
   965)     d2vr_dp2 = u9*(u0-one)/zz*zpp*zpp + &
   966)                u0*u1/zz*d2z_dp2 + &
   967)                six*u2*aa(19) + &
   968)                60.d0*u7*u5/(a10+beta)/(a10+beta) + &
   969)                six*(a12 - theta) + &
   970)                24.d0*aa(22)*beta/theta20
   971) 
   972)     d2wmol_dp2 = -cnv*upc1*upc1*(2.d0/vr*vrpp*vrpp -  d2vr_dp2) ! kmol/m^3/Pa^2
   973) 
   974)   else
   975)     dwt = UNINITIALIZED_DOUBLE
   976)     dwp = UNINITIALIZED_DOUBLE
   977)     d2wmol_dp2 = UNINITIALIZED_DOUBLE
   978)   endif
   979) 
   980) end subroutine EOSWaterDensityIFC67
   981) 
   982) ! ************************************************************************** !
   983) 
   984) subroutine EOSWaterEnthalpyIFC67(t,p,calculate_derivatives,hw, &
   985)                                  hwp,hwt,ierr)
   986) 
   987) !  This subroutine calculates water and steam-gas mixture properties.
   988) !  The water and steam properties are valid in the range of:
   989) !
   990) !            0 < p < 165.4 * 10^5 pascals (165.4 bars)
   991) !            0 < t < 350 centigrade (623.15 Kelvin)
   992) !
   993) !  The properties cover densities, enthalpies, internal energies,
   994) !  and partial derivatives of these quanties with respect to
   995) !  pressure and temperature.
   996) !
   997) !  For saturated fluid, it will also calculate water saturation
   998) !  temperature for the specified pressure using Newton-Raphson and
   999) !  the derivative dts/dp (=tsp) or Ps for a given temperature.
  1000) !
  1001) !  Ref.: International Formulation Committee of the Sixth International
  1002) !       Conference on Properties of Steam (1967).
  1003) 
  1004)   implicit none
  1005)   
  1006)   PetscReal, intent(in) :: t   ! Temperature in centigrade
  1007)   PetscReal, intent(in) :: p   ! Pressure in Pascals
  1008)   PetscBool, intent(in) :: calculate_derivatives
  1009)   PetscReal, intent(out) :: hw,hwp,hwt
  1010)   PetscErrorCode, intent(out) :: ierr
  1011)   
  1012)   PetscInt :: i
  1013)     
  1014)   PetscReal, save :: aa(0:22)
  1015)   PetscReal, save :: a1,a2,a3,a4,a5,a6,a7,a8,a9,a10,a11,a12
  1016)   
  1017)   PetscReal :: beta,beta2x,beta4,theta,utheta,theta2x,theta18,theta20
  1018)   PetscReal :: xx,yy,zz
  1019)   PetscReal :: u0,u1!,u2,u3,u4,u5,u6,u7,u8,u9
  1020)   PetscReal :: tempreal
  1021)   PetscReal :: v0_1, v1_1, v2_1, v3_1, v4_1
  1022)   PetscReal :: v1_2, v2_2, v3_2, v4_2, v20_2, v40_2
  1023)   PetscReal :: v1_3, v2_3, v3_3, v4_3
  1024)   PetscReal :: v1_4, v2_4, v3_4
  1025)   PetscReal :: v1_5, v2_5
  1026)   PetscReal :: v1_6
  1027)   PetscReal :: term1,term2,term2t,term3,term3t,term3p,term4,term4t,term4p, &
  1028)                term5,term5t,term5p,term6,term6t,term6p,term7,term7t,term7p
  1029)   PetscReal :: dv2t,dv2p,dv3t
  1030) !  PetscReal :: vr,ypt,yptt,zpt,zpp,vrpt,vrpp,cnv
  1031)   PetscReal :: ypt,yptt,zpt,zpp
  1032)   PetscReal :: tc1,pc1,vc1,utc1,upc1,vc1mol
  1033)   PetscReal, parameter :: zero = 0.d0
  1034)   PetscReal, parameter :: one = 1.d0
  1035)   PetscReal, parameter :: two = 2.d0
  1036)   PetscReal, parameter :: three = 3.d0
  1037)   PetscReal, parameter :: four = 4.d0
  1038)   PetscReal, parameter :: five = 5.d0
  1039)   PetscReal, parameter :: six = 6.d0
  1040)   PetscReal, parameter :: seven = 7.d0
  1041)   PetscReal, parameter :: eight = 8.d0
  1042)   PetscReal, parameter :: nine = 9.d0
  1043)   PetscReal, parameter :: ten = 10.d0
  1044)   
  1045)   data aa/ &
  1046) !-----data aa0,aa1,aa2,aa3/
  1047)         6.824687741d03,-5.422063673d02,-2.096666205d04, 3.941286787d04, &
  1048) !-----data aa4,aa5,aa6,aa7/
  1049)         -6.733277739d04, 9.902381028d04,-1.093911774d05, 8.590841667d04, &
  1050) !-----data aa8,aa9,aa10,aa11/
  1051)         -4.511168742d04, 1.418138926d04,-2.017271113d03, 7.982692717d00, &
  1052) !-----data aa12,aa13,aa14,aa15/
  1053)         -2.616571843d-2, 1.522411790d-3, 2.284279054d-2, 2.421647003d02, &
  1054) !-----data aa16,aa17,aa18,aa19/
  1055)         1.269716088d-10,2.074838328d-7, 2.174020350d-8, 1.105710498d-9, &
  1056) !-----data aa20,aa21,aa22/    
  1057)         1.293441934d01, 1.308119072d-5, 6.047626338d-14/
  1058) 
  1059)   data a1,a2,a3,a4/ &
  1060)   8.438375405d-1, 5.362162162d-4, 1.720000000d00, 7.342278489d-2/
  1061)   data a5,a6,a7,a8/ &
  1062)   4.975858870d-2, 6.537154300d-1, 1.150000000d-6, 1.510800000d-5/
  1063)   data a9,a10,a11,a12/ &
  1064)   1.418800000d-1, 7.002753165d00, 2.995284926d-4, 2.040000000d-1/
  1065)     
  1066)   ierr = 0
  1067) 
  1068)   tc1 = H2O_CRITICAL_TEMPERATURE    ! K
  1069)   pc1 = H2O_CRITICAL_PRESSURE     ! Pa 
  1070)   vc1 = 0.00317d0  ! m^3/kg
  1071)   utc1 = one/tc1   ! 1/C
  1072)   upc1 = one/pc1   ! 1/Pa
  1073)   vc1mol = vc1*FMWH2O ! m^3/kmol
  1074)     
  1075)   theta = (t+273.15d0)*utc1
  1076)   theta2x = theta*theta
  1077)   theta18 = theta**18.d0
  1078)   theta20 = theta18*theta2x
  1079)     
  1080)   beta = p*upc1
  1081)   beta2x = beta*beta
  1082)   beta4  = beta2x*beta2x
  1083)     
  1084)   yy = one-a1*theta2x-a2*theta**(-6.d0)
  1085)   xx = a3*yy*yy-two*(a4*theta-a5*beta)
  1086) !   Note: xx may become negative near the critical point-pcl.
  1087)   if (xx.gt.zero) then
  1088)     xx = sqrt(xx)
  1089)   else
  1090)     write(*,*) 'Warning: negative term in density (eos_water.F90:&
  1091)       &EOSWaterEnthalpyIFC67):'
  1092)     write(*,*) 't= ',t,' p= ',p,' xx= ',xx
  1093)     ierr = 1
  1094)     xx = 1.d-6               !set arbitrarily
  1095)   end if
  1096)   zz = yy + xx                                     
  1097)   u0 = -five/17.d0
  1098)   u1 = aa(11)*a5*zz**u0
  1099) #if 0  
  1100)   u2 = one/(a8+theta**11)
  1101)   u3 = aa(17)+(two*aa(18)+three*aa(19)*beta)*beta
  1102)   u4 = one/(a7+theta18*theta)
  1103)   u5 = (a10+beta)**(-4)
  1104)   u6 = a11-three*u5
  1105)   u7 = aa(20)*theta18*(a9+theta2x)
  1106)   u8 = aa(15)*(a6-theta)**9
  1107)     
  1108)   vr = u1+aa(12)+theta*(aa(13)+aa(14)*theta)+u8*(a6-theta) &
  1109)         +aa(16)*u4-u2*u3-u6*u7+(three*aa(21)*(a12-theta) &
  1110)         +four*aa(22)*beta/theta20)*beta2x
  1111)     
  1112)   dwmol = one/(vr*vc1mol) ! kmol/m^3
  1113)   dw = one/(vr*vc1) ! kg/m^3
  1114) 
  1115)   !---calculate derivatives for water density
  1116)   if (calculate_derivatives) then
  1117)     u9 = u0*u1/zz
  1118)     vrpt = u9*zpt+aa(13)+two*aa(14)*theta-ten*u8 &
  1119)         -19.d0*aa(16)*u4*u4*theta18+11.d0*u2*u2*u3*theta**10.d0 &
  1120)         -aa(20)*u6*(18.d0*a9*theta18+20.d0*theta20)/theta &
  1121)         -(three*aa(21)+80.d0*aa(22)*beta/(theta20*theta))*beta2x
  1122)     
  1123)     vrpp = u9*zpp-u2*(two*aa(18)+six*aa(19)*beta)-12.d0*u7*u5/ &
  1124)         (a10+beta)+(six*aa(21)*(a12-theta)+12.d0*aa(22)*beta/ &
  1125)         theta20)*beta
  1126)     
  1127)     cnv = -one/(vc1mol*vr*vr)
  1128)     dwt = cnv*vrpt*utc1 ! kmol/m^3/C
  1129)     dwp = cnv*vrpp*upc1 ! kmol/m^3/Pa
  1130)   else
  1131)     dwt = UNINITIALIZED_DOUBLE
  1132)     dwp = UNINITIALIZED_DOUBLE
  1133)   endif
  1134) #endif  
  1135) 
  1136)   ! ypt used for enthalpy even if derivative not calculated
  1137)   ypt = six*a2*theta**(-7.d0)-two*a1*theta
  1138)   
  1139) 
  1140) !---compute enthalpy internal energy and derivatives for water
  1141)   utheta = one/theta
  1142)   term1 = aa(0)*theta
  1143)   term2 = -aa(1)
  1144)   ! term2t is part of the derivative calc., but left here to avoid
  1145)   ! recomputing the expensive do loop
  1146)   term2t = zero
  1147)   do i = 3,10
  1148)     tempreal = dfloat(i-2)*aa(i)*theta**(i-1)
  1149)     term2t = term2t+tempreal*utheta*dfloat(i-1)
  1150)     term2 = term2+tempreal                            
  1151)   end do
  1152)     
  1153)   ! "v" section 1
  1154)   v0_1 = u1/a5
  1155)   v2_1 = 17.d0*(zz/29.d0-yy/12.d0)+five*theta*ypt/12.d0
  1156)   v3_1 = a4*theta-(a3-one)*theta*yy*ypt
  1157)   v1_1 = zz*v2_1+v3_1
  1158)   term3 = v0_1*v1_1
  1159)   
  1160)   ! block 1 removed from here
  1161) 
  1162)   ! "v" section 2
  1163)   v1_2 = nine*theta+a6
  1164)   v20_2 = (a6-theta)
  1165)   v2_2 = v20_2**9.d0
  1166)   v3_2 = a7+20.d0*theta**19.d0
  1167)   v40_2 = a7+theta**19.d0
  1168)   v4_2 = one/(v40_2*v40_2)
  1169)   ! term4p is a derivative, but left due to dependency in term4
  1170)   term4p = aa(12)-aa(14)*theta2x+aa(15)*v1_2*v2_2+aa(16)*v3_2*v4_2
  1171)   term4 = term4p*beta
  1172)   
  1173)   ! block 2 removed from here
  1174)     
  1175)   ! "v" section 3
  1176)   v1_3 = beta*(aa(17)+aa(18)*beta+aa(19)*beta2x)
  1177)   v2_3 = 12.d0*theta**11.d0+a8
  1178)   v4_3 = one/(a8+theta**11.d0)
  1179)   v3_3 = v4_3*v4_3
  1180)   term5 = v1_3*v2_3*v3_3
  1181)   
  1182)   ! block 3 removed from here
  1183)     
  1184)   ! "v" section 4
  1185)   v1_4 = (a10+beta)**(-3.d0)+a11*beta
  1186)   v3_4 = (17.d0*a9+19.d0*theta2x)
  1187)   v2_4 = aa(20)*theta18*v3_4       
  1188)   term6 = v1_4*v2_4
  1189) 
  1190)   ! block 4 removed from here
  1191)     
  1192)   ! "v" section 5
  1193)   v1_5 = 21.d0*aa(22)/theta20*beta4
  1194)   v2_5 = aa(21)*a12*beta2x*beta
  1195)   term7 = v1_5+v2_5  
  1196)     
  1197)   ! "v" section 6
  1198)   v1_6 = pc1*vc1mol
  1199)   hw = (term1-term2+term3+term4-term5+term6+term7)*v1_6
  1200)     
  1201)   if (calculate_derivatives) then
  1202) 
  1203)     zpt = ypt+(a3*yy*ypt-a4)/xx
  1204)     zpp = a5/xx
  1205)   
  1206)     ! block 1
  1207)     yptt = -two*a1-42.d0*a2/theta**8.d0
  1208)     dv2t = 17.d0*(zpt/29.d0-ypt/12.d0)+five/12.d0*(ypt+theta*yptt) 
  1209)     dv3t = a4-(a3-one)*(theta*yy*yptt+yy*ypt+theta*ypt*ypt)
  1210)     dv2p = 17.d0*zpp/29.d0
  1211)     v4_1 = five*v1_1/(17.d0*zz)       
  1212)     term3t = v0_1*(zz*dv2t+(v2_1-v4_1)*zpt+dv3t)
  1213)     term3p = v0_1*(zz*dv2p+(v2_1-v4_1)*zpp)
  1214)   
  1215)     ! block 2
  1216)     term4t = (-two*aa(14)*theta+nine*aa(15)*(v2_2-v1_2*v2_2/v20_2) &
  1217)              +38.d0*theta18*aa(16)*(ten*v4_2-v3_2*v4_2/v40_2))*beta
  1218) 
  1219)     ! block 3
  1220)     term5p = v3_3*v2_3*(aa(17)+two*aa(18)*beta+three*aa(19)*beta2x)
  1221)     term5t = v1_3*(132.d0*v3_3*theta**10.d0-22.d0*v2_3*v3_3*v4_3*theta**10.d0)
  1222)   
  1223)     ! block 4
  1224)     term6p = v2_4*(a11-three*(a10+beta)**(-4.d0))
  1225)     term6t = v1_4*aa(20)*theta18*(18.d0*v3_4*utheta+38.d0*theta)
  1226)     
  1227)     ! block 5
  1228)     term7p = beta2x*(three*aa(21)*a12+84.d0*aa(22)*beta/theta20)
  1229)     term7t = -420.d0*aa(22)*beta4/(theta20*theta)
  1230) 
  1231)     hwp = (term3p+term4p-term5p+term6p+term7p)*vc1mol
  1232)     hwt = (aa(0)-term2t+term3t+term4t-term5t+term6t+term7t)*v1_6*utc1
  1233)   else
  1234)     hwp = UNINITIALIZED_DOUBLE
  1235)     hwt = UNINITIALIZED_DOUBLE
  1236)   endif
  1237)     
  1238) end subroutine EOSWaterEnthalpyIFC67
  1239) 
  1240) ! ************************************************************************** !
  1241) 
  1242) subroutine EOSWaterDensityConstant(t,p,calculate_derivatives,dw,dwmol, &
  1243)                                    dwp,dwt,ierr)
  1244)   implicit none
  1245)   
  1246)   PetscReal, intent(in) :: t
  1247)   PetscReal, intent(in) :: p
  1248)   PetscBool, intent(in) :: calculate_derivatives
  1249)   PetscReal, intent(out) :: dw,dwmol,dwp,dwt
  1250)   PetscErrorCode, intent(out) :: ierr
  1251)   
  1252)   dw = constant_density ! kg/m^3
  1253)   dwmol = dw/FMWH2O ! kmol/m^3
  1254)   
  1255)   dwp = 0.d0
  1256)   dwt = 0.d0
  1257) 
  1258) end subroutine EOSWaterDensityConstant
  1259) 
  1260) ! ************************************************************************** !
  1261) 
  1262) subroutine EOSWaterEnthalpyConstant(t,p,calculate_derivatives, &
  1263)                                     hw,hwp,hwt,ierr)
  1264)   implicit none
  1265)   
  1266)   PetscReal, intent(in) :: t   ! Temperature in centigrade
  1267)   PetscReal, intent(in) :: p   ! Pressure in Pascals
  1268)   PetscBool, intent(in) :: calculate_derivatives
  1269)   PetscReal, intent(out) :: hw,hwp,hwt
  1270)   PetscErrorCode, intent(out) :: ierr
  1271)   
  1272)   hw = constant_enthalpy ! J/kmol
  1273)   
  1274)   hwp = 0.d0
  1275)   hwt = 0.d0
  1276)   
  1277) end subroutine EOSWaterEnthalpyConstant
  1278) 
  1279) ! ************************************************************************** !
  1280) 
  1281) subroutine EOSWaterDensityExponential(t,p,calculate_derivatives, &
  1282)                                       dw,dwmol,dwp,dwt,ierr)
  1283)   implicit none
  1284)   
  1285)   PetscReal, intent(in) :: t   ! Temperature in centigrade
  1286)   PetscReal, intent(in) :: p   ! Pressure in Pascals
  1287)   PetscBool, intent(in) :: calculate_derivatives
  1288)   PetscReal, intent(out) :: dw,dwmol,dwp,dwt
  1289)   PetscErrorCode, intent(out) :: ierr
  1290)   
  1291)   ! kg/m^3
  1292)   dw = exponent_reference_density*exp(exponent_water_compressibility* &
  1293)                                      (p-exponent_reference_pressure))
  1294)   dwmol = dw/FMWH2O ! kmol/m^3
  1295)   
  1296)   if (calculate_derivatives) then
  1297)     dwp = dwmol*exponent_water_compressibility !kmol/m^3/Pa
  1298)   else
  1299)     dwp = 0.d0
  1300)   endif
  1301)   dwt = 0.d0
  1302) 
  1303) end subroutine EOSWaterDensityExponential
  1304) 
  1305) ! ************************************************************************** !
  1306) 
  1307) subroutine EOSWaterSteamDenEnthNoDerive(t,p,dg,dgmol,hg,ierr)
  1308) 
  1309)   implicit none
  1310) 
  1311)   PetscReal, intent(in) :: t   ! Temperature in centigrade.
  1312)   PetscReal, intent(in) :: p   ! Vapor Pressure in Pascals.
  1313)   PetscReal, intent(out) :: dg,dgmol
  1314)   PetscReal, intent(out) :: hg
  1315)   PetscErrorCode, intent(out) :: ierr
  1316)   
  1317)   PetscReal :: dum1, dum2, dum3, dum4
  1318)   
  1319)   call EOSWaterSteamDensityEnthalpyPtr(t,p,PETSC_FALSE,dg,dgmol,hg, &
  1320)                                        dum1,dum2,dum3,dum4,ierr)
  1321)   
  1322) end subroutine EOSWaterSteamDenEnthNoDerive
  1323) 
  1324) ! ************************************************************************** !
  1325) 
  1326) subroutine EOSWaterSteamDenEnthDerive(t,pv,dg,dgmol,hg, &
  1327)                                       dgp,dgt,hgp,hgt,ierr)
  1328)   implicit none
  1329) 
  1330)   PetscReal, intent(in) :: t   ! Temperature in centigrade.
  1331)   PetscReal, intent(in) :: pv  ! Vapor Pressure in Pascals.
  1332)   PetscReal, intent(out) :: dg,dgmol,dgp,dgt
  1333)   PetscReal, intent(out) :: hg,hgp,hgt
  1334)   PetscErrorCode, intent(out) :: ierr
  1335)   
  1336)   call EOSWaterSteamDensityEnthalpyPtr(t,pv,PETSC_TRUE,dg,dgmol,hg,dgp, &
  1337)                                        dgt,hgp,hgt,ierr)
  1338)   
  1339) end subroutine EOSWaterSteamDenEnthDerive
  1340) 
  1341) ! ************************************************************************** !
  1342) 
  1343) subroutine EOSWaterSteamDensityEnthalpyIFC67(t,pv,calculate_derivatives, &
  1344)                                              dg,dgmol,hg, &
  1345)                                              dgp,dgt,hgp,hgt,ierr)
  1346) ! t/C  p/Pa dgmol/(mol/m^3)  h/J/kmol
  1347)   implicit none
  1348)   
  1349)   PetscReal, intent(in) :: t   ! Temperature in centigrade.
  1350)   PetscReal, intent(in) :: pv  ! Vapor Pressure in Pascals.
  1351)   PetscBool, intent(in) :: calculate_derivatives
  1352)   PetscReal, intent(out) :: dg,dgmol,dgp,dgt
  1353)   PetscReal, intent(out) :: hg,hgp,hgt
  1354)   PetscErrorCode, intent(out) :: ierr
  1355)   
  1356)   PetscInt, save :: n(8),ll(8),x(8,2),z(8,3)
  1357)   PetscReal, save :: xi1,xl0,xl1,xl2
  1358)   PetscReal, save :: b(8,2),bb(0:9,0:6)
  1359)   PetscReal :: sumbx(8),sumbxt(8)
  1360)   
  1361)   PetscInt :: i,j
  1362)   PetscReal, save :: delp,delt
  1363)   PetscReal :: beta,betap,betal,betalp,betalp1,betal1,ubeta,theta, &
  1364)             thetap,utheta
  1365)   PetscReal :: xx,xxp
  1366)   PetscReal :: f1,f2,fim1,fim2,sum,sumt,sum1,sum1t,sum1p, &
  1367)             sum2,sum2t
  1368)   PetscReal :: u1,u1t,u1p,u2,u2p,u2t,u3,u3p,u3t,v1,v1t
  1369)   PetscReal :: term,term1,term1t,term1p,term2,term2t,term2p, &
  1370)             term3,term3t,term3p,term4,term4t,term4p,term5,term5t,term5p
  1371)   PetscReal :: hr,hrp,hrt,hrpt,hrpp
  1372)   PetscReal :: vr,vrpt,vrpp
  1373)   PetscReal :: tc1,pc1,vc1,utc1,upc1,vc1mol
  1374)   PetscReal, parameter :: zero = 0.d0
  1375)   PetscReal, parameter :: one = 1.d0
  1376)   PetscReal, parameter :: two = 2.d0
  1377)   PetscReal, parameter :: three = 3.d0
  1378)   PetscReal, parameter :: four = 4.d0
  1379)   PetscReal, parameter :: five = 5.d0
  1380)   PetscReal, parameter :: six = 6.d0
  1381)   PetscReal, parameter :: seven = 7.d0
  1382)   PetscReal, parameter :: eight = 8.d0
  1383)   PetscReal, parameter :: nine = 9.d0
  1384)   PetscReal, parameter :: ten = 10.d0
  1385)  
  1386)   data delt,delp/1.d-6,1.d-6/
  1387)     
  1388)   data n/2,3,2,2,3,2,2,2/
  1389)   data ll/0,0,0,0,0,1,1,2/
  1390)   data x/0,0,0,0,0,14,19,54, &
  1391)           0,0,0,0,0, 0, 0,27/
  1392)   data z/13,18,18,25,32,12,24,24, &
  1393)           3, 2,10,14,28,11,18,14, &
  1394)           0, 1, 0, 0,24, 0, 0, 0/
  1395)     
  1396)   data b/7.6333333333d-1,0.d0,0.d0,0.d0,0.d0,4.006073948d-1, &
  1397)           8.636081627d-2,-8.532322921d-1, &
  1398)           0.d0,0.d0,0.d0,0.d0,0.d0,0.d0,0.d0,3.460208861d-1/
  1399) 
  1400) !---sub-region 2
  1401)   data bb/ &
  1402) !---bb(0-9,0)
  1403)           1.683599274d1,  8*0.d0        , 1.936587558d2, &
  1404) !---bb(0-9,1)
  1405)           2.856067796d01, 6.670375918d-2, 8.390104328d-2, &
  1406)           4.520918904d-1,-5.975336707d-1, 5.958051609d-1, &
  1407)           1.190610271d-1, 1.683998803d-1, 6.552390126d-3, &
  1408)           -1.388522425d+3, &
  1409) !---bb(0-9,2)
  1410)           -5.438923329d01, 1.388983801d0 , 2.614670893d-2, &
  1411)           1.069036614d-1,-8.847535804d-2,-5.159303373d-1, &
  1412)           -9.867174132d-2,-5.809438001d-2, 5.710218649d-4, &
  1413)           4.126607219d03, &
  1414) !---bb(0-9,3)
  1415)           4.330662834d-1, 0.d0          ,-3.373439453d-2, &
  1416)           0.d0          , 0.d0          , 2.075021122d-1, &
  1417)           0.d0          , 0.d0          , 0.d0          , &
  1418)           -6.508211677d03, &
  1419) !---bb(0-9,4)
  1420)           -6.547711697d-1, 8*0.d0        , 5.745984054d03, &
  1421) !---bb(0-9,5)
  1422)           8.565182058d-2, 8*0.d0        ,-2.693088365d03, &
  1423) !---bb(0-9,6)
  1424)                           9*0.d0        , 5.235718623d02/
  1425)     
  1426)   data xi1/4.260321148d0/
  1427)   data xl0,xl1,xl2/15.74373327d0,-34.17061978d0,19.31380707d0/
  1428)     
  1429)   ierr = 0
  1430)   
  1431)   tc1 = H2O_CRITICAL_TEMPERATURE       ! K
  1432)   pc1 = H2O_CRITICAL_PRESSURE        ! Pa
  1433)   vc1 = 0.00317d0     ! m^3/kg
  1434)   utc1 = one/tc1
  1435)   upc1 = one/pc1
  1436)   vc1mol = vc1*FMWH2O ! m^3/kmol
  1437) 
  1438)   theta  = (t+273.15d0)*utc1
  1439)   beta   = pv*upc1
  1440)   ubeta  = one/beta
  1441)   utheta = one/theta
  1442)   xx = exp(b(1,1)*(one-theta))
  1443)     
  1444) !---compute steam density and derivatives
  1445)     
  1446)   term1 = xi1*theta*ubeta
  1447)   term1t = xi1*ubeta
  1448)   term1p = -term1*ubeta        
  1449)     
  1450)   do i = 1,8
  1451)     sum = zero
  1452)     sumt = zero
  1453)     do j = 1,n(i)
  1454)       u1 = bb(i,j)*xx**z(i,j)     
  1455)       sumt = sumt-b(1,1)*z(i,j)*u1                     
  1456)       sum = sum + u1                     
  1457)     end do
  1458)     sumbx(i) = sum
  1459)     sumbxt(i) = sumt
  1460)   end do
  1461)     
  1462)   term2  = sumbx(1)+beta*(two*sumbx(2)+beta*(three*sumbx(3) &
  1463)             +beta*(four*sumbx(4)+five*beta*sumbx(5))))
  1464)   term2t = sumbxt(1)+beta*(two*sumbxt(2)+beta*(three*sumbxt(3) &
  1465)             +beta*(four*sumbxt(4)+five*beta*sumbxt(5))))
  1466)   term2p = two*sumbx(2)+beta*(six*sumbx(3)+beta*(12.d0*sumbx(4) &
  1467)             +20.d0*beta*sumbx(5)))
  1468)     
  1469)   term3  = zero
  1470)   term3p = zero
  1471)   term3t = zero
  1472)     
  1473)   do i = 6,8
  1474)     fim1 = dfloat(i-1)
  1475)     fim2 = fim1-one    
  1476)     
  1477)     sum2 = zero
  1478)     sum2t = zero
  1479)     do j = 1,ll(i)
  1480)       u1 = b(i,j)*xx**x(i,j)
  1481)       sum2t = sum2t-b(1,1)*x(i,j)*u1
  1482)       sum2 = sum2 + u1
  1483)     end do
  1484)     
  1485)     f1 = fim2*beta**(1-i)
  1486)     f2 = beta**(2-i)
  1487)     u1 = one/(f2+sum2)**2
  1488)     term = u1*f1*sumbx(i)
  1489)     term3 = term3+term      
  1490)     u2 = two*term/(f2+sum2)
  1491)     term3t = term3t+u1*(f1*sumbxt(i)+u2*sum2t)    
  1492)     term3p = term3p-u1*(f1*fim1*sumbx(i)+u2*fim2*f2)*ubeta
  1493)   end do
  1494)     
  1495)   term4 = bb(9,0)
  1496)   term4t = zero  
  1497)   do i = 1,6
  1498)     u1 = bb(9,i)*xx**i 
  1499)     term4t = term4t-b(1,1)*dfloat(i)*u1
  1500)     term4 = term4 + u1             
  1501)   end do
  1502)     
  1503)   betal = xl0+theta*(xl1+xl2*theta)
  1504)   betalp = xl1+ two*xl2*theta
  1505)   u1 = 11.d0*(beta/betal)**10
  1506)   term4t = u1*(term4t-10.d0*betalp*term4/betal)
  1507)   term4 = term4*u1                  
  1508)   term4p = 10.d0*term4*ubeta             
  1509)     
  1510)   vr = term1-term2-term3+term4
  1511)   vrpt = term1t-term2t-term3t+term4t
  1512)   vrpp = term1p-term2p-term3p+term4p
  1513)     
  1514)   u1 = -one/(vc1mol*vr)
  1515)   dgmol = -u1
  1516)   dg = one/(vc1*vr)
  1517)   
  1518)   if (calculate_derivatives) then
  1519)     u1 = u1/vr
  1520)     dgt = u1*vrpt*utc1
  1521)     dgp = u1*vrpp*upc1
  1522)   else
  1523)     dgt = UNINITIALIZED_DOUBLE
  1524)     dgp = UNINITIALIZED_DOUBLE
  1525)   endif
  1526)     
  1527) !---compute steam enthalpy, internal energy, and derivatives
  1528)   thetap = theta+delt 
  1529)   betap  = beta +delp 
  1530)   xxp    = exp(b(1,1)*(one-thetap))
  1531)     
  1532)   term1  = bb(0,0)*theta
  1533)   term1t = bb(0,0)*thetap
  1534)     
  1535)   term2  = zero
  1536)   term2t = zero
  1537)   do i = 1,5
  1538)     u1 = bb(0,i)*(dfloat(i)-two)
  1539)     term2t = term2t+u1*thetap**(i-1)
  1540)     term2  = term2+u1*theta**(i-1)
  1541)   end do
  1542)     
  1543)   term3  = zero
  1544)   term3t = zero
  1545)   term3p = zero
  1546)   u1 = one
  1547)   u1p = one
  1548)   do i = 1,5
  1549)     u1 = u1*beta
  1550)     u1p = u1p*betap
  1551)     
  1552)     sum = zero
  1553)     sumt = zero
  1554)     do j = 1,n(i)
  1555)       sumt = sumt+bb(i,j)*(one+z(i,j)*b(1,1)*thetap)*xxp**z(i,j)
  1556)       sum  = sum+bb(i,j)*(one+z(i,j)*b(1,1)*theta)*xx**z(i,j)
  1557)     end do
  1558)     
  1559)     term3t = term3t+u1*sumt
  1560)     term3p = term3p+u1p*sum
  1561)     term3  = term3+u1*sum
  1562)   end do
  1563)     
  1564)   term4  = zero
  1565)   term4t = zero
  1566)   term4p = zero
  1567)   do i = 6,8
  1568)     
  1569)     sum1  = zero
  1570)     sum2  = zero
  1571)     sum1t = zero
  1572)     sum2t = zero
  1573)     
  1574)     do j = 1,ll(i)
  1575)       u1 = b(i,j)*xxp**x(i,j)
  1576)       sum1t = sum1t+x(i,j)*u1               
  1577)       sum2t = sum2t+u1                 
  1578)       u1 = b(i,j)*xx**x(i,j)
  1579)       sum1 = sum1+x(i,j)*u1               
  1580)       sum2 = sum2+u1                 
  1581)     end do
  1582)     
  1583)     u1 = one/(beta**(2-i)+sum2)
  1584)     u2 = one-b(1,1)*theta*sum1*u1                  
  1585)     u3 = b(1,1)*theta
  1586)     
  1587)     u1t = one/(beta**(2-i)+sum2t)
  1588)     u2t = one-b(1,1)*thetap*sum1t*u1t                 
  1589)     u3t = b(1,1)*thetap
  1590)     
  1591)     u1p = one/(betap**(2-i)+sum2)
  1592)     u2p = one-b(1,1)*theta*sum1*u1p                   
  1593)     u3p = u3          
  1594)     
  1595)     sum1 = zero
  1596)     sum1t = zero
  1597)     sum1p = zero
  1598)     do j = 1,n(i)
  1599)       sum1t = sum1t + bb(i,j)*xxp**z(i,j)*(u2t+z(i,j)*u3t)
  1600)       sum1p = sum1p + bb(i,j)*xx **z(i,j)*(u2p+z(i,j)*u3p)
  1601)       sum1  = sum1  + bb(i,j)*xx **z(i,j)*(u2 +z(i,j)*u3)
  1602)     end do
  1603)     
  1604)     term4t = term4t+sum1t*u1t
  1605)     term4p = term4p+sum1p*u1p
  1606)     term4  = term4+sum1*u1
  1607)   end do
  1608)     
  1609)   u1 = ten*betalp/betal
  1610)   term5 = (one+theta*u1)*bb(9,0)
  1611)     
  1612)   betal1  = xl0+thetap*(xl1+xl2*thetap)
  1613)   betalp1 = xl1+two*xl2*thetap
  1614)   u1t     = ten*betalp1/betal1
  1615)   term5t  = (one+thetap*u1t)*bb(9,0)
  1616)     
  1617)   do i = 1,6
  1618)     v1     = one+theta*(u1+dfloat(i)*b(1,1))
  1619)     v1t    = one+thetap*(u1t+dfloat(i)*b(1,1))
  1620)     term5t = v1t*bb(9,i)*xxp**i+term5t                        
  1621)     term5  =  v1*bb(9,i)*xx **i+term5                         
  1622)   end do
  1623)     
  1624)   term5  = term5*beta*(beta/betal)**10
  1625)   term5t = term5t*beta*(beta/betal1)**10
  1626)   term5p = term5*(betap/beta)**11        
  1627)     
  1628)   hr   = term1 -term2 -term3 -term4 +term5
  1629)   hrt  = term1t-term2t-term3t-term4t+term5t
  1630)   hrp  = term1 -term2 -term3p-term4p+term5p
  1631)   hrpt = (hrt-hr)/delt
  1632)   hrpp = (hrp-hr)/delp
  1633)     
  1634)   v1 = pc1*vc1mol  ! Pa = (nRT/V) = J/m^3 : J/m^3 * m^3/kmol = J/kmol
  1635)   hg = hr*v1       ! J/kmol
  1636)   
  1637)   if (calculate_derivatives) then
  1638)     hgt = hrpt*v1*utc1
  1639)     hgp = hrpp*vc1mol
  1640)   else
  1641)     hgt = UNINITIALIZED_DOUBLE
  1642)     hgp = UNINITIALIZED_DOUBLE
  1643)   endif
  1644) 
  1645) end subroutine EOSWaterSteamDensityEnthalpyIFC67
  1646) 
  1647) ! ************************************************************************** !
  1648) 
  1649) subroutine EOSWaterSteamDenEnthConstant(t,pv,calculate_derivatives, &
  1650)                                         dg,dgmol,hg,dgp,dgt,hgp,hgt,ierr)
  1651) ! t/C  p/Pa dgmol/(mol/m^3)  h/J/kmol
  1652)   implicit none
  1653)   
  1654)   PetscReal, intent(in) :: t   ! Temperature in centigrade.
  1655)   PetscReal, intent(in) :: pv  ! Vapor Pressure in Pascals.
  1656)   PetscBool, intent(in) :: calculate_derivatives
  1657)   PetscReal, intent(out) :: dg,dgmol,dgp,dgt
  1658)   PetscReal, intent(out) :: hg,hgp,hgt
  1659)   PetscErrorCode, intent(out) :: ierr
  1660)   
  1661)   ierr = 0
  1662)   
  1663)   dg = constant_steam_density
  1664)   dgmol = dg/FMWH2O
  1665)   hg = constant_steam_enthalpy
  1666)   
  1667)   if (calculate_derivatives) then
  1668)     dgt = 0.d0
  1669)     dgp = 0.d0
  1670)     hgt = 0.d0
  1671)     hgp = 0.d0
  1672)   else
  1673)     dgt = UNINITIALIZED_DOUBLE
  1674)     dgp = UNINITIALIZED_DOUBLE
  1675)     hgt = UNINITIALIZED_DOUBLE
  1676)     hgp = UNINITIALIZED_DOUBLE
  1677)   endif
  1678)   
  1679) end subroutine EOSWaterSteamDenEnthConstant
  1680) 
  1681) ! ************************************************************************** !
  1682) 
  1683) subroutine EOSWaterDuanMixture(t,p,xmol,y_nacl,avgmw,dw_kg,denmix)
  1684) 
  1685) ! Duan et al. (2008) Energy and Fuels, v 22, 1666-1674.
  1686) 
  1687)   implicit none
  1688) 
  1689)   PetscReal :: t,tk,p,xco2,xmol,x1,y_nacl,vphi_a1,vphi_a2,vphi,denmix,pw_kg,dw_kg,avgmw
  1690) 
  1691)   PetscReal :: fmwh2o = 18.01534d0
  1692)   PetscReal :: fmwco2 = 44.0098d0
  1693)   PetscReal :: fmwnacl = 58.44277d0
  1694)   PetscReal :: dummy
  1695)   PetscErrorCode :: ierr
  1696) 
  1697)   !duan mixing **************************
  1698)   tk = t + 273.15D0; xco2 = xmol;
  1699)   call EOSWaterDensity(t,p,pw_kg,dummy,ierr)
  1700)   x1 = 1.D0-xco2;
  1701)   vphi_a1 = (0.3838402D-3*tk - 0.5595385D0)*tk + 0.30429268D3 + &
  1702)             (-0.72044305D5 + 0.63003388D7/tk)/tk;  
  1703)   vphi_a2 = (-0.57709332D-5*tk + 0.82764653D-2)*tk - 0.43813556D1 + &
  1704)             (0.10144907D4 - 0.86777045D5/tk)/tk;  
  1705)   vphi = (1.D0 + vphi_a1 + vphi_a2*p*1.D-6)*(fmwh2o*1.D-3/pw_kg); 
  1706)   vphi = x1*((1.D0 - y_nacl)*fmwh2o + y_nacl*fmwnacl)*1.D-3/dw_kg + xco2*vphi;
  1707)   denmix = (x1*((1.D0 - y_nacl)*fmwh2o + y_nacl*fmwnacl) + xco2*fmwco2)*1.D-3/vphi;
  1708)   denmix = denmix/avgmw
  1709)   
  1710) end subroutine EOSWaterDuanMixture
  1711) 
  1712) ! ************************************************************************** !
  1713) 
  1714) subroutine EOSWaterViscosityNaCl (t,p_Pa,xnacl,visnacl)
  1715) 
  1716)   !viscosity: Kestin et al. (1981)
  1717) 
  1718)   implicit none
  1719) 
  1720)   PetscReal, intent(in) :: t        ! [C]
  1721)   PetscReal, intent(in) :: p_Pa     ! [Pa]
  1722)   PetscReal, intent(in) :: xnacl    ! [-]
  1723)   PetscReal, intent(out) :: visnacl ! [Pa-s]
  1724) 
  1725) 
  1726)   PetscReal, save :: a1,a2,a3,b1,b2,b3,c1,c2,c3,c4,wnacl
  1727)   PetscReal :: ak,bk,ck
  1728)   PetscReal :: beta,betap,betas,betaw
  1729)   PetscReal :: tt,mnacl,fac,mu0,ms
  1730)   PetscReal :: p_GPa
  1731) 
  1732)   data a1,a2,a3 / 3.324d-2, 3.624d-3, -1.879d-4 /
  1733)   data b1,b2,b3 / -3.96d-2, 1.02d-2, -7.02d-4 /
  1734)   data c1,c2,c3,c4 / 1.2378d0, -1.303d-3, 3.06d-6, 2.55d-8 /
  1735) 
  1736)   data wnacl / 58.44277d-3 / ! (kg/mol NaCl)
  1737) 
  1738)   !convert pressure to GPa
  1739)   p_GPa = p_Pa*1.d-9
  1740) 
  1741)   mnacl = xnacl/(1.d0-xnacl)/wnacl
  1742) 
  1743)   tt = 20.d0-t
  1744)   ck = (c1 + (c2 + (c3+c4*tt)*tt)*tt)*tt/(96.d0+t)
  1745)   ak = (a1 + (a2 + a3*mnacl)*mnacl)*mnacl
  1746)   bk = (b1 + (b2 + b3*mnacl)*mnacl)*mnacl
  1747) 
  1748)   ms = 6.044d0 + (2.8d-3 + 3.6d-5*t)*t
  1749)   fac = mnacl/ms
  1750)   betaw = -1.297d0 + (5.74d-2 + (-6.97d-4 + (4.47d-6 - 1.05d-8*t)*t)*t)*t
  1751)   betas = 0.545d0 + 2.8d-3 * t - betaw
  1752)   betap = (2.5d0 + (-2.d0 + 0.5d0*fac)*fac)*fac
  1753)   beta = betas*betap + betaw
  1754) 
  1755)   mu0 = 1001.74d-6 * 10.d0**(ak + ck*(bk + 1.d0))
  1756) 
  1757)   visnacl = mu0*(1.d0 + beta*p_GPa)
  1758) 
  1759) end subroutine EOSWaterViscosityNaCl
  1760) 
  1761) ! ************************************************************************** !
  1762) 
  1763) subroutine EOSWaterViscosityKestinExt(T, P, PS, dPS_dT, aux, &
  1764)                                       calculate_derivatives, VW, &
  1765)                                       dVW_dT, dVW_dP, ierr)
  1766) 
  1767)   !viscosity: Kestin et al. (1981)
  1768) 
  1769)   implicit none
  1770) 
  1771)   PetscReal, intent(in) :: T   ! C
  1772)   PetscReal, intent(in) :: P   ! Pa
  1773)   PetscReal, intent(in) :: PS  ! Pa
  1774)   PetscReal, intent(in) :: dPS_dT
  1775)   PetscReal, intent(in) :: aux(*)
  1776)   PetscBool, intent(in) :: calculate_derivatives
  1777)   PetscReal, intent(out) :: VW ! Pa-s
  1778)   PetscReal, intent(out) :: dVW_dT, dVW_dP
  1779)   PetscErrorCode, intent(out) :: ierr
  1780) 
  1781)   PetscReal, save :: a1,a2,a3,b1,b2,b3,c1,c2,c3,c4,wnacl
  1782)   PetscReal :: xnacl    ! [-]
  1783)   PetscReal :: ak,bk,ck
  1784)   PetscReal :: beta,betap,betas,betaw
  1785)   PetscReal :: tt,mnacl,fac,mu0,ms
  1786)   PetscReal :: p_GPa
  1787) 
  1788)   data a1,a2,a3 / 3.324d-2, 3.624d-3, -1.879d-4 /
  1789)   data b1,b2,b3 / -3.96d-2, 1.02d-2, -7.02d-4 /
  1790)   data c1,c2,c3,c4 / 1.2378d0, -1.303d-3, 3.06d-6, 2.55d-8 /
  1791) 
  1792)   data wnacl / 58.44277d-3 / ! (kg/mol NaCl)
  1793) 
  1794)   if (calculate_derivatives) then
  1795)     print *, 'Derivatives not set up in EOSWaterViscosityKestinExt().'
  1796)     stop
  1797)   endif
  1798) 
  1799)   !convert pressure to GPa
  1800)   p_GPa = P*1.d-9
  1801)   xnacl = aux(1)
  1802) 
  1803)   mnacl = xnacl/(1.d0-xnacl)/wnacl
  1804) 
  1805)   tt = 20.d0-t
  1806)   ck = (c1 + (c2 + (c3+c4*tt)*tt)*tt)*tt/(96.d0+t)
  1807)   ak = (a1 + (a2 + a3*mnacl)*mnacl)*mnacl
  1808)   bk = (b1 + (b2 + b3*mnacl)*mnacl)*mnacl
  1809) 
  1810)   ms = 6.044d0 + (2.8d-3 + 3.6d-5*t)*t
  1811)   fac = mnacl/ms
  1812)   betaw = -1.297d0 + (5.74d-2 + (-6.97d-4 + (4.47d-6 - 1.05d-8*t)*t)*t)*t
  1813)   betas = 0.545d0 + 2.8d-3 * t - betaw
  1814)   betap = (2.5d0 + (-2.d0 + 0.5d0*fac)*fac)*fac
  1815)   beta = betas*betap + betaw
  1816) 
  1817)   mu0 = 1001.74d-6 * 10.d0**(ak + ck*(bk + 1.d0))
  1818) 
  1819)   VW = mu0*(1.d0 + beta*p_GPa)
  1820) 
  1821) end subroutine EOSWaterViscosityKestinExt
  1822) 
  1823) ! ************************************************************************** !
  1824) 
  1825) subroutine EOSWaterSaturationTemperature(ps,ts_guess,ts,t_ps,ierr)
  1826) 
  1827)   !  This function calculates saturation temperature for a given Ps c  
  1828)   !  Ref.: International Formulation Committee of the Sixth International
  1829)   !       Conference on Properties of Steam (1967).
  1830) 
  1831)   !    ps  = saturation pressure (pascals)
  1832)   !    ts  = saturation temperature (deg. C)
  1833)   !    tsp = estimated ts on entry and dT/dps on return
  1834) 
  1835)   implicit none
  1836) 
  1837)   PetscReal, intent(in) :: ps
  1838)   PetscReal, intent(in) :: ts_guess
  1839)   PetscReal, intent(out) :: ts
  1840)   PetscReal, intent(out) :: t_ps
  1841)   PetscErrorCode :: ierr
  1842)   
  1843) 
  1844)   PetscReal, parameter :: epsilon = 1.d-10
  1845)   PetscReal, parameter :: tc1 = H2O_CRITICAL_TEMPERATURE
  1846)   PetscReal, parameter :: pc1 = H2O_CRITICAL_PRESSURE
  1847)       
  1848)   PetscReal :: theta, beta, u1, err
  1849)   PetscReal :: t1num, t1nump
  1850)   PetscReal :: t1, t1den, t1denp, term1, term1p, t2, term2, term2p
  1851)   PetscReal :: f, fp
  1852)   PetscReal :: kn(9)
  1853)   PetscInt :: itr
  1854) 
  1855)   data kn / -7.691234564d0,-2.608023696d1,-1.681706546d2, &
  1856)             6.423285504d1,-1.189646225d2, 4.167117320d0, &
  1857)             2.097506760d1, 1.d9         , 6.d0/
  1858) 
  1859)   ierr = 0
  1860) 
  1861) !geh  if (ipvtab.eq.0 .or. tsp.gt.369.d0) then
  1862) 
  1863) !-------newton-raphson iteration for calculating ts by analytical funcs
  1864) 
  1865) !-------on entry, ts_guess = estimated ts
  1866)   theta = (ts_guess+273.15d0)/tc1
  1867)   beta  = ps/pc1
  1868) 
  1869)   u1  = 1.d0-theta
  1870)   itr = 0
  1871) 
  1872)   do         
  1873)     itr = itr + 1
  1874) 
  1875)     t1num  = u1*(kn(1)+u1*(kn(2)+u1*(kn(3)+u1*(kn(4)+u1*kn(5)))))
  1876)     t1nump = u1*(2.d0*kn(2)+u1*(3.d0*kn(3)+u1*(4.d0*kn(4)+5.d0* &
  1877)              u1*kn(5))))+kn(1)
  1878)     t1     = 1.d0+u1*(kn(6)+kn(7)*u1)
  1879)     t1den  = 1.d0/(theta*t1)
  1880)     t1denp = theta*(kn(6)+2.d0*kn(7)*u1)-t1
  1881)     term1  = t1num*t1den
  1882)     term1p = (t1nump-term1*t1denp)*t1den
  1883) 
  1884)     t2     = 1.d0/(kn(8)*u1*u1+kn(9))
  1885)     term2  = u1*t2
  1886)     term2p = t2*(1.d0-2.d0*kn(8)*u1*u1*t2)
  1887)     f      = exp(term1-term2)-beta
  1888)     fp     = (f+beta)*(term1p-term2p)
  1889)     err    = f/fp
  1890)     u1     = u1-err
  1891)     theta  = 1.d0-u1
  1892)     
  1893)     if (dabs(err) <= epsilon .or. itr >= 20) exit
  1894) 
  1895)   enddo
  1896) 
  1897)   ts = theta*tc1-273.15d0
  1898) 
  1899) !-------Note-(dbeta/dtheta) = -fp ; tsp = dT/dps
  1900)   t_ps = -tc1/(pc1*fp)
  1901) 
  1902) end subroutine EOSWaterSaturationTemperature
  1903) 
  1904) ! ************************************************************************** !
  1905) 
  1906) subroutine EOSWaterDensityIcePainter(T, P, calculate_derivatives, &
  1907)                                      den_ice, dden_ice_dT, dden_ice_dP, ierr)
  1908)   ! Subroutine to calculate the density of ice at given temperature
  1909)   ! and pressure
  1910)   ! T is in deg C, P is in Pa, density is in kmol/m3
  1911)   implicit none
  1912)   
  1913)   PetscReal, intent(in) :: T
  1914)   PetscReal, intent(in) :: P
  1915)   PetscBool, intent(in) :: calculate_derivatives
  1916) 
  1917)   PetscReal, intent(out) :: den_ice
  1918)   PetscReal, intent(out) :: dden_ice_dT
  1919)   PetscReal, intent(out) :: dden_ice_dP
  1920)   PetscErrorCode, intent(out) :: ierr
  1921) 
  1922)   PetscReal, parameter :: P_ref = 1.d5
  1923)   PetscReal, parameter :: alpha = 3.3d-10
  1924)   PetscReal, parameter :: beta = 1.53d-4
  1925) 
  1926)   den_ice = 5.09424d1*(1.d0 + alpha*(P - P_ref) - beta*(T)) !in Kmol/m3
  1927)   if (calculate_derivatives) then
  1928)     dden_ice_dT = 5.09424d1*(-beta)
  1929)     dden_ice_dP = 5.09424d1*alpha
  1930)   else
  1931)     dden_ice_dT = UNINITIALIZED_DOUBLE
  1932)     dden_ice_dP = UNINITIALIZED_DOUBLE
  1933)   endif
  1934)   
  1935) end subroutine EOSWaterDensityIcePainter
  1936) 
  1937) ! ************************************************************************** !
  1938) 
  1939) subroutine EOSWaterInternalEnergyIce(T, u_ice, du_ice_dT)
  1940)   ! Subroutine to calculate the internal energy of ice at given temperature and 
  1941)   ! pressure
  1942)   ! T is in deg C, internal energy is in J/mol
  1943)   implicit none
  1944) 
  1945)   PetscReal, intent(in) :: T
  1946)   PetscReal, intent(out) :: u_ice
  1947)   PetscReal, intent(out) :: du_ice_dT
  1948)   PetscErrorCode :: ierr
  1949)   
  1950)   PetscReal, parameter :: a = -10.6644d0
  1951)   PetscReal, parameter :: b = 0.1698d0
  1952)   PetscReal, parameter :: c = 198148.d0
  1953)   PetscReal, parameter :: T_ref = 273.15d0
  1954) 
  1955)   ! from Maier-Kelly type fit (integrated tref to t)
  1956)   ! in J/mol
  1957) 
  1958)   u_ice = a*(T) + b/2.d0*((T + T_ref)**(2.d0) - T_ref**(2.d0)) + &
  1959)           c*(1.d0/T_ref - 1.d0/(T + T_ref))
  1960)   u_ice = u_ice - HEAT_OF_FUSION*FMWH2O*1.d-3   ! kJ/kmol
  1961)   du_ice_dT = a + b*(T + T_ref) + c/((T + T_ref)**(2.d0)) !kJ/kmol/K
  1962)   
  1963) end subroutine EOSWaterInternalEnergyIce
  1964) 
  1965) ! ************************************************************************** !
  1966) 
  1967) subroutine EOSWaterDensityPainter(t,p,calculate_derivatives,dw,dwmol, &
  1968)                                   dwp,dwt,ierr)
  1969) 
  1970) ! wateos_simple: Simple water equation of state from Scott Painter
  1971) ! Author: Satish Karra, LANL
  1972) ! Date: 02/1/12
  1973) ! T in C, P in Pa
  1974)   implicit none
  1975) 
  1976)   PetscReal, intent(in) :: t   ! Temperature in centigrade
  1977)   PetscReal, intent(in) :: p   ! Pressure in Pascals
  1978)   PetscBool, intent(in) :: calculate_derivatives
  1979)   PetscReal, intent(out) :: dw,dwmol,dwp,dwt
  1980)   PetscErrorCode, intent(out) :: ierr
  1981) 
  1982)   PetscReal, parameter :: a = 999.915d0
  1983)   PetscReal, parameter :: b = 0.0416516d0
  1984)   PetscReal, parameter :: c = -0.0100836d0
  1985)   PetscReal, parameter :: d = 0.000206355
  1986)   PetscReal, parameter :: alpha = 5.0d-10     ! in Pa^(-1)
  1987)   PetscReal, parameter :: T_ref = 273.15d0    ! in K
  1988)   PetscReal, parameter :: P_ref = 1.0d5       ! in Pa
  1989) 
  1990)   PetscReal :: den_w_one_bar, T_K
  1991)   PetscReal :: u_J_kg, h_J_kg
  1992)   PetscReal :: du_dt
  1993) 
  1994)   ierr = 0
  1995) 
  1996)   ! Density of water
  1997)   T_K = T + T_ref    ! convert to Kelvin
  1998)   den_w_one_bar = a + b*(T_K - T_ref) + c*(T_K - T_ref)**(2.d0) + &
  1999)                   d*(T_K - T_ref)**(3.d0)
  2000)   dw = den_w_one_bar*(1 + alpha*(P - P_ref))
  2001)   dwmol = dw/FMWH2O     ! in mol
  2002) 
  2003)   ! Internal energy
  2004)   u_J_kg = 4.217*1.0d3*(T_K - T_ref)    ! in J/kg
  2005)   h_J_kg = u_J_kg + P/dw    ! in J/kg
  2006) 
  2007)   if (calculate_derivatives) then
  2008)     ! Derivatives of density
  2009)     dwp = 1/FMWH2O*den_w_one_bar*alpha    ! in Kmol/Pa
  2010)     dwt = 1/FMWH2O*(1 + alpha*(P - P_ref))*(b + 2.d0*c*(T_K - T_ref) + &
  2011)                               3.d0*d*(T_K - T_ref)**(2.d0))      ! in Kmol/K
  2012)   else
  2013)     dwp = UNINITIALIZED_DOUBLE
  2014)     dwt = UNINITIALIZED_DOUBLE
  2015)   endif
  2016) 
  2017) end subroutine EOSWaterDensityPainter
  2018) 
  2019) ! ************************************************************************** !
  2020) 
  2021) subroutine EOSWaterEnthalpyPainter(T, P, calculate_derivatives, &
  2022)                                    h_J_kmol, dh_dp, dh_dt, ierr)
  2023) 
  2024) ! wateos_simple: Simple water equation of state from Scott Painter
  2025) ! Author: Satish Karra, LANL
  2026) ! Date: 02/1/12
  2027) ! T in C, P in Pa
  2028)   implicit none
  2029) 
  2030)   PetscReal, intent(in) :: T
  2031)   PetscReal, intent(in) :: P
  2032)   PetscBool, intent(in) :: calculate_derivatives
  2033)   PetscReal, intent(out) :: h_J_kmol
  2034)   PetscReal :: den_water_kg, den_water_kmol
  2035)   PetscReal :: dden_water_dp, dden_water_dt
  2036)   PetscReal, intent(out) :: dh_dp, dh_dt
  2037) 
  2038)   PetscErrorCode, intent(out) :: ierr
  2039) 
  2040)   PetscReal, parameter :: a = 999.915d0
  2041)   PetscReal, parameter :: b = 0.0416516d0
  2042)   PetscReal, parameter :: c = -0.0100836d0
  2043)   PetscReal, parameter :: d = 0.000206355
  2044)   PetscReal, parameter :: alpha = 5.0d-10     ! in Pa^(-1)
  2045)   PetscReal, parameter :: T_ref = 273.15d0    ! in K
  2046)   PetscReal, parameter :: P_ref = 1.0d5       ! in Pa
  2047) 
  2048)   PetscReal :: den_w_one_bar, T_K
  2049)   PetscReal :: u_J_kg, h_J_kg
  2050)   PetscReal :: du_dt
  2051) 
  2052)   ierr = 0
  2053) 
  2054)   ! Density of water
  2055)   T_K = T + T_ref    ! convert to Kelvin
  2056)   den_w_one_bar = a + b*(T_K - T_ref) + c*(T_K - T_ref)**(2.d0) + &
  2057)                   d*(T_K - T_ref)**(3.d0)
  2058)   den_water_kg = den_w_one_bar*(1 + alpha*(P - P_ref))
  2059)   den_water_kmol = den_water_kg/FMWH2O     ! in mol
  2060) 
  2061)   ! Internal energy
  2062)   u_J_kg = 4.217*1.0d3*(T_K - T_ref)    ! in J/kg
  2063)   h_J_kg = u_J_kg + P/den_water_kg    ! in J/kg
  2064)   h_J_kmol = h_J_kg*FMWH2O     ! in J/kmol
  2065) 
  2066)   if (calculate_derivatives) then
  2067)     ! Derivatives of density
  2068)     dden_water_dp = 1/FMWH2O*den_w_one_bar*alpha    ! in Kmol/Pa
  2069)     dden_water_dt = 1/FMWH2O*(1 + alpha*(P - P_ref))*(b + 2.d0*c*(T_K - T_ref) + &
  2070)                               3.d0*d*(T_K - T_ref)**(2.d0))      ! in Kmol/K
  2071) 
  2072)     ! Derivatives of enthalpy
  2073)     dh_dp = FMWH2O/den_water_kg   ! in J/kmol/Pa
  2074)     du_dt = 4.217*1.d3                  ! in J/kg/K
  2075)     dh_dt = FMWH2O*(du_dt + P*(-1.d0/den_water_kg**(2.d0))* &
  2076)                     dden_water_dt*FMWH2O)    ! in MJ/kmol/K
  2077)   else
  2078)     dden_water_dp = UNINITIALIZED_DOUBLE
  2079)     dden_water_dp = UNINITIALIZED_DOUBLE
  2080)     dh_dp = UNINITIALIZED_DOUBLE
  2081)     du_dt = UNINITIALIZED_DOUBLE
  2082)     dh_dt = UNINITIALIZED_DOUBLE
  2083)   endif
  2084) 
  2085) end subroutine EOSWaterEnthalpyPainter
  2086) 
  2087) ! ************************************************************************** !
  2088) 
  2089) subroutine EOSWaterDensityIceNoDerive(t,p,dw,ierr)
  2090) 
  2091)   implicit none
  2092) 
  2093)   PetscReal, intent(in) :: t
  2094)   PetscReal, intent(in) :: p
  2095)   PetscReal, intent(out) :: dw
  2096)   PetscErrorCode, intent(out) :: ierr
  2097) 
  2098)   PetscReal :: dum1, dum2
  2099) 
  2100)   call EOSWaterDensityIcePtr(t,p,PETSC_FALSE,dw,dum1,dum2,ierr)
  2101) 
  2102) end subroutine EOSWaterDensityIceNoDerive
  2103) 
  2104) ! ************************************************************************** !
  2105) 
  2106) subroutine EOSWaterDensityIceDerive(t,p,dw,dwp,dwt,ierr)
  2107)   implicit none
  2108) 
  2109)   PetscReal, intent(in) :: t
  2110)   PetscReal, intent(in) :: p
  2111)   PetscReal, intent(out) :: dw,dwp,dwt
  2112)   PetscErrorCode, intent(out) :: ierr
  2113) 
  2114)   call EOSWaterDensityIcePtr(t,p,PETSC_TRUE,dw,dwp,dwt,ierr)
  2115) 
  2116) end subroutine EOSWaterDensityIceDerive
  2117) 
  2118) ! ************************************************************************** !
  2119) 
  2120) subroutine EOSWaterDensityTGDPB01(t, p, calculate_derivatives, &
  2121)                                   dw, dwmol, dwp, dwt, ierr)
  2122) 
  2123)   ! 
  2124)   ! Tanaka M. , G. Girard, R. Davis, A. Peuto, and N. Bignell. 2001.
  2125)   ! Recommended table for the density of water between 0 °C
  2126)   ! and 40 °C based on recent experimental reports. Metrologia,
  2127)   ! 38:301-309 [doi:10.1088/0026-1394/38/4/3].
  2128)   ! 
  2129)   ! Author: Gautam Bisht, LBNL
  2130)   ! Date: 24/06/15
  2131)   ! 
  2132)   implicit none
  2133) 
  2134)   PetscReal, intent(in) :: t   ! Temperature in centigrade
  2135)   PetscReal, intent(in) :: p   ! Pressure in Pascals
  2136)   PetscBool, intent(in) :: calculate_derivatives
  2137)   PetscReal, intent(out) :: dw,dwmol,dwp,dwt
  2138)   PetscErrorCode, intent(out) :: ierr
  2139) 
  2140)   PetscReal,parameter :: a1 = -3.983035d0     ! [degC]
  2141)   PetscReal,parameter :: a2 = 301.797d0       ! [degC]
  2142)   PetscReal,parameter :: a3 = 522528.9d0      ! [degC^{2}]
  2143)   PetscReal,parameter :: a4 = 69.34881d0      ! [degC]
  2144)   PetscReal,parameter :: a5 = 999.974950d0    ! [kg m^{-3}]
  2145)   PetscReal,parameter :: k0 = 50.74d-11       ! [Pa^{-1}]
  2146)   PetscReal,parameter :: k1 = -0.326d-11      ! [Pa^{-1} degC^{-1}]
  2147)   PetscReal,parameter :: k2 = 0.00416d-11     ! [Pa^{-1} degC^{-2}]
  2148)   PetscReal,parameter :: p0 = 101325.d0       ! [Pa]
  2149)   PetscReal :: t_c
  2150)   PetscReal :: dent
  2151)   PetscReal :: kappa
  2152)   PetscReal :: ddent_dt
  2153)   PetscReal :: ddent_dt_1
  2154)   PetscReal :: ddent_dt_2
  2155)   PetscReal :: ddent_dt_3
  2156)   PetscReal :: ddent_dp
  2157)   PetscReal :: dkappa_dp
  2158)   PetscReal :: dkappa_dt
  2159)   PetscReal :: dden_dt
  2160) 
  2161)   ! Density of water as function of temperature
  2162)   dent = a5*(1.d0 - ((t + a1)**2.d0)*(t + a2)/a3/(t + a4))
  2163) 
  2164)   ! Compressibility of water
  2165)   kappa = (1.d0 + (k0 + k1*t + k2*t**2.d0)*(p - p0))
  2166) 
  2167)   ! Density of water
  2168)   dw    = dent*kappa ! [kg m^{-3}]
  2169)   dwmol = dw/FMWH2O  ! [kmol m^{-3}]
  2170) 
  2171)   if (calculate_derivatives) then
  2172)     ! Derivative
  2173)     ddent_dp = 0.d0
  2174)     ddent_dt_1 = -((t + a1)**2.d0)/a3/(t + a4)
  2175)     ddent_dt_2 = -2.d0*(t + a1)*(t + a2)/a3/(t + a4)
  2176)     ddent_dt_3 =  ((t + a1)**2.d0)*(t + a2)/a3/((t + a4)**2.d0)
  2177)     ddent_dt   = a5*(ddent_dt_1 + ddent_dt_2 + ddent_dt_3)
  2178) 
  2179)     dkappa_dp = (k0 + k1*t + k2*t**2.d0)
  2180)     dkappa_dt = (k1 + 2.d0*k2*t)*(p - p0)
  2181) 
  2182)     dwt = (ddent_dt*kappa + dent*dkappa_dt)/FMWH2O ! [kmol m^{-3} degC^{-1}]
  2183)     dwp = (ddent_dp*kappa + dent*dkappa_dp)/FMWH2O ! [kmol m^{-3} Pa^{-1}]
  2184)   else
  2185)     dwt = UNINITIALIZED_DOUBLE
  2186)     dwp = UNINITIALIZED_DOUBLE
  2187)   endif
  2188) 
  2189) end subroutine EOSWaterDensityTGDPB01
  2190) 
  2191) ! ************************************************************************** !
  2192) 
  2193) subroutine EOSWaterDensityBatzleAndWang(tin, pin, calculate_derivatives, &
  2194)                                         dw, dwmol, dwp, dwt, ierr)
  2195) 
  2196)   ! 
  2197)   ! From Batlze M. and Z. Wang (1992) Seismic properties of fluids, Geophysics,
  2198)   ! Vol. 57, No. 11, Pg. 1396-1408.
  2199)   !
  2200)   ! Equation 27a
  2201)   !
  2202)   ! Author: Glenn Hammond
  2203)   ! Date: 02/08/16
  2204)   ! 
  2205)   implicit none
  2206) 
  2207)   PetscReal, intent(in) :: tin   ! Temperature in centigrade
  2208)   PetscReal, intent(in) :: pin   ! Pressure in Pascal
  2209)   PetscBool, intent(in) :: calculate_derivatives
  2210)   PetscReal, intent(out) :: dw ! kg/m^3
  2211)   PetscReal, intent(out) :: dwmol ! kmol/m^3
  2212)   PetscReal, intent(out) :: dwp ! kmol/m^3-Pa
  2213)   PetscReal, intent(out) :: dwt ! kmol/m^3-C
  2214)   PetscErrorCode, intent(out) :: ierr
  2215)   
  2216)   PetscReal, parameter :: g_cm3_to_kg_m3 = 1.d3
  2217)   PetscReal, parameter :: Pa_to_MPa = 1.d-6
  2218)   PetscReal :: t ! temperature in Celcius
  2219)   PetscReal :: t_sq, t_cub
  2220)   PetscReal :: p_MPa, p_MPa_sq
  2221)   
  2222)   t = tin
  2223)   t_sq = t*t
  2224)   t_cub = t*t_sq
  2225)   p_MPa = pin*Pa_to_MPa
  2226)   p_MPa_sq = p_MPa*p_MPa
  2227)   
  2228)   ! temperature is in C and pressure in MPa
  2229)   ! g/cm^3
  2230)   ! Eq. 27a
  2231)   dw = 1.d0 + 1.d-6*(-80.d0*t - 3.3d0*t_sq + 1.75d-3*t_cub + &
  2232)                      489.d0*p_MPa - 2.d0*t*p_MPa + 1.6d-2*t_sq*p_MPa - &
  2233)                      1.3d-5*t_cub*p_MPa - 3.33d-1*p_MPa_sq - 2.d-3*t*p_MPa_sq)
  2234)   ! convert from g/cm^3 to kg/m^3
  2235)   dw = dw * g_cm3_to_kg_m3
  2236)   dwmol = dw/FMWH2O ! kmol/m^3 
  2237)   
  2238)   if (calculate_derivatives) then
  2239)     dwp = 1.d-6*(489.d0 - 2.d0*t + 1.6d-2*t_sq - 1.3d-5*t_cub - &
  2240)                  6.66d-1*p_MPa - 4.d-3*t*p_MPa) *  &
  2241)           g_cm3_to_kg_m3/FMWH2O
  2242)     ! convert from kmol/m^3-MPa to kmol/m^3-Pa
  2243)     dwp = dwp*Pa_to_MPa
  2244)     dwt = 1.d-6*(-80.d0 - 6.6d0*t + 5.25d-3*t_sq - 2.d0*p_MPa + &
  2245)                  3.2d-2*t*p_MPa - 3.9d-5*t_sq*p_MPa - 2.d-3*p_MPa_sq) * &
  2246)           g_cm3_to_kg_m3/FMWH2O
  2247)   else
  2248)     dwp = 0.d0
  2249)     dwt = 0.d0
  2250)   endif
  2251) 
  2252) end subroutine EOSWaterDensityBatzleAndWang
  2253) 
  2254) ! ************************************************************************** !
  2255) 
  2256) subroutine EOSWaterDensityBatzleAndWangExt(tin, pin, aux, &
  2257)                                            calculate_derivatives, &
  2258)                                            dw, dwmol, dwp, dwt, ierr)
  2259) 
  2260)   ! 
  2261)   ! From Batlze M. and Z. Wang (1992) Seismic properties of fluids, Geophysics,
  2262)   ! Vol. 57, No. 11, Pg. 1396-1408.
  2263)   !
  2264)   ! Equation 27b
  2265)   !
  2266)   ! Author: Glenn Hammond
  2267)   ! Date: 02/08/16
  2268)   ! 
  2269)   implicit none
  2270) 
  2271)   PetscReal, intent(in) :: tin   ! Temperature in centigrade
  2272)   PetscReal, intent(in) :: pin   ! Pressure in Pascal
  2273)   PetscReal, intent(in) :: aux(*)
  2274)   PetscBool, intent(in) :: calculate_derivatives
  2275)   PetscReal, intent(out) :: dw ! kg/m^3
  2276)   PetscReal, intent(out) :: dwmol ! kmol/m^3
  2277)   PetscReal, intent(out) :: dwp ! kmol/m^3-Pa
  2278)   PetscReal, intent(out) :: dwt ! kmol/m^3-C
  2279)   PetscErrorCode, intent(out) :: ierr
  2280)   
  2281)   PetscReal, parameter :: g_cm3_to_kg_m3 = 1.d3
  2282)   PetscReal, parameter :: Pa_to_MPa = 1.d-6
  2283)   PetscReal :: t_C ! temperature in Celcius
  2284)   PetscReal :: p_MPa
  2285)   PetscReal :: s
  2286)   
  2287)   t_C = tin
  2288)   p_MPa = pin*Pa_to_MPa
  2289)   s = aux(1)
  2290)   
  2291)   call EOSWaterDensityPtr(tin, pin, calculate_derivatives, &
  2292)                           dw, dwmol, dwp, dwt, ierr)
  2293)   
  2294)   ! temperature is in C and pressure in MPa
  2295)   ! kg/m^3
  2296)   ! Eq. 27b
  2297)   dw = dw + &
  2298)        s*(0.668d0 + 0.44d0*s + &
  2299)           1.d-6*(300.d0*p_MPa - 2400.d0*p_MPa*s + &
  2300)                  t_C*(80.d0 + 3.d0*t_C - 3300.d0*s - 13.d0*p_MPa + &
  2301)                       47.d0*p_Mpa*s))) * &
  2302)        g_cm3_to_kg_m3
  2303)   dwmol = dw/FMWH2O ! kmol/m^3 
  2304)   
  2305)   if (calculate_derivatives) then
  2306)         ! v - this dwp is in the correct units of kmol/m^3-Pa
  2307)     dwp = dwp + &
  2308)           s*(1.d-6*(300.d0 - 2400.d0*s + t_C*(-13.d0 + 47.d0*s))) * &
  2309)                                  ! v - convert from kmol/m^3-MPa to kmol/m^3-Pa
  2310)           g_cm3_to_kg_m3/FMWH2O*Pa_to_MPa
  2311)     dwt = dwt + &
  2312)           s*(1.d-6*(80.d0 + 6.d0*t_C - 3300.d0*s - 13.d0*p_MPa + &
  2313)                     47.d0*p_Mpa*s)) * &
  2314)           g_cm3_to_kg_m3/FMWH2O
  2315)   else
  2316)     dwp = 0.d0
  2317)     dwt = 0.d0
  2318)   endif
  2319) 
  2320) end subroutine EOSWaterDensityBatzleAndWangExt
  2321) 
  2322) ! ************************************************************************** !
  2323) 
  2324) subroutine EOSWaterViscosityBatzleAndWang(T, P, PS, dPS_dT, &
  2325)                                           calculate_derivatives, VW, &
  2326)                                           dVW_dT, dVW_dP, ierr)
  2327)   ! 
  2328)   ! From Batlze M. and Z. Wang (1992) Seismic properties of fluids, Geophysics,
  2329)   ! Vol. 57, No. 11, Pg. 1396-1408.
  2330)   !
  2331)   ! Equation 32
  2332)   !
  2333)   ! Author: Glenn Hammond
  2334)   ! Date: 02/08/16
  2335)   ! 
  2336)   implicit none
  2337)   PetscReal, intent(in) :: T       ! C
  2338)   PetscReal, intent(in) :: P       ! Pa
  2339)   PetscReal, intent(in) :: PS      ! Pa
  2340)   PetscReal, intent(in) :: dPS_dT  ! Pa/C
  2341)   PetscBool, intent(in) :: calculate_derivatives
  2342)   PetscReal, intent(out) :: VW     ! Pa-s
  2343)   PetscReal, intent(out) :: dVW_dT, dVW_dP
  2344)   PetscErrorCode, intent(out) :: ierr
  2345)   
  2346)   ! convert from centipoise to Pa-s (1 cP = 1.d-3 Pa-s)
  2347)   PetscReal, parameter :: centipoise_to_Pa_s = 1.d-3
  2348)   PetscReal :: t_C
  2349)   PetscReal :: exponential_term
  2350)   PetscReal :: temperature_term
  2351)   
  2352)   t_C = T
  2353)   
  2354)   ! this is Eq. 32 without all the salt terms.
  2355)   ! -0.057138d0 = -1.d0*0.42d0*(-0.17d0)**2.d0+0.045d0  
  2356)   exponential_term = -0.057138d0*t_C**0.8d0
  2357)   temperature_term = 1.65d0*exp(exponential_term)
  2358)   VW = 0.1d0 + temperature_term
  2359)   VW = VW * centipoise_to_Pa_s
  2360)        
  2361)   if (calculate_derivatives) then
  2362)     dVW_dP = 0.d0
  2363)     dVW_dT = 0.8d0*temperature_term*exponential_term/t_C*centipoise_to_Pa_s
  2364)   else
  2365)     dVW_dP = 0.d0
  2366)     dVW_dT = 0.d0
  2367)   endif
  2368) 
  2369) end subroutine EOSWaterViscosityBatzleAndWang
  2370) 
  2371) ! ************************************************************************** !
  2372) 
  2373) subroutine EOSWaterViscosityBatzleAndWangExt(T, P, PS, dPS_dT, aux, &
  2374)                                              calculate_derivatives, VW, &
  2375)                                              dVW_dT, dVW_dP, ierr)
  2376)   ! 
  2377)   ! From Batlze M. and Z. Wang (1992) Seismic properties of fluids, Geophysics,
  2378)   ! Vol. 57, No. 11, Pg. 1396-1408.
  2379)   !
  2380)   ! Equation 32
  2381)   !
  2382)   ! Author: Glenn Hammond
  2383)   ! Date: 02/08/16
  2384)   ! 
  2385)   implicit none
  2386) 
  2387)   PetscReal, intent(in) :: T, P, PS, dPS_dT
  2388)   PetscReal, intent(in) :: aux(*)
  2389)   PetscBool, intent(in) :: calculate_derivatives
  2390)   PetscReal, intent(out) :: VW
  2391)   PetscReal, intent(out) :: dVW_dT, dVW_dP
  2392)   PetscErrorCode, intent(out) :: ierr
  2393)   
  2394)   ! convert from centipoise to Pa-s (1 cP = 1.d-3 Pa-s)
  2395)   PetscReal, parameter :: centipoise_to_Pa_s = 1.d-3
  2396)   PetscReal :: t_C
  2397)   PetscReal :: s
  2398)   PetscReal :: exponential_term
  2399)   PetscReal :: temperature_term
  2400)   
  2401)   s = aux(1)
  2402)   t_C = T
  2403)   
  2404)   exponential_term = -1.d0*(0.42d0*(s**0.8d0-0.17d0)**2.d0 + 0.045d0)* &
  2405)                      t_C**0.8d0
  2406)   temperature_term = (1.65d0 + 91.9d0*s**3.d0)*exp(exponential_term)
  2407)   VW = 0.1d0 + 0.333d0*s + temperature_term
  2408)   VW = VW * centipoise_to_Pa_s
  2409)        
  2410)   if (calculate_derivatives) then
  2411)     dVW_dP = 0.d0
  2412)     dVW_dT = 0.8d0*temperature_term*exponential_term/t_C*centipoise_to_Pa_s
  2413)   else
  2414)     dVW_dP = 0.d0
  2415)     dVW_dT = 0.d0
  2416)   endif
  2417)   
  2418) end subroutine EOSWaterViscosityBatzleAndWangExt
  2419) 
  2420) ! ************************************************************************** !
  2421) 
  2422) subroutine TestEOSWaterBatzleAndWang()
  2423)   ! 
  2424)   ! From Batlze M. and Z. Wang (1992) Seismic properties of fluids, Geophysics,
  2425)   ! Vol. 57, No. 11, Pg. 1396-1408.
  2426)   !
  2427)   ! Code for testing
  2428)   !
  2429)   ! Author: Glenn Hammond
  2430)   ! Date: 02/15/16
  2431)   ! 
  2432)   PetscReal :: p, t, dw, dwmol, dwp, dwt, ps, dps_dt, vw, dvw_dt, dvw_dp
  2433)   PetscReal :: t2, dw2, p2, vw2
  2434)   PetscReal :: aux(1)
  2435)   PetscErrorCode :: ierr
  2436)   
  2437)   p = 0.101325d6
  2438)   t = 20.d0
  2439)   ps = -999.d0
  2440)   dps_dt = -999.d0
  2441)   aux(1) = 0.14d0
  2442)   
  2443)   call EOSWaterSetDensity('BATZLE_AND_WANG')
  2444)   call EOSWaterSetViscosity('BATZLE_AND_WANG')
  2445)   
  2446)   call EOSWaterDensityBatzleAndWang(t,p, PETSC_TRUE, &
  2447)                                     dw, dwmol, dwp, dwt, ierr)
  2448)   print *, 'dw:    ', dw
  2449)   print *, 'dwmol: ', dwmol
  2450)   print *, 'dwp:   ', dwp
  2451)   print *, 'dwt:   ', dwt
  2452) 
  2453)   t2 = t + 1.d-6*t
  2454)   call EOSWaterDensityBatzleAndWang(t2,p, PETSC_TRUE, &
  2455)                                     dw2, dwmol, dwp, dwt, ierr)
  2456)   print *
  2457)   print *, 'dw2:   ', dw2
  2458)   print *, 'dw(t): ', (dw2-dw)/(t2-t)
  2459)   
  2460)   p2 = p + 1.d-6*p
  2461)   call EOSWaterDensityBatzleAndWang(t,p2, PETSC_TRUE, &
  2462)                                     dw2, dwmol, dwp, dwt, ierr)
  2463)   print *
  2464)   print *, 'dw2:   ', dw2
  2465)   print *, 'dw(p): ', (dw2-dw)/(p2-p)
  2466)   
  2467)   
  2468)   call EOSWaterDensityBatzleAndWangExt(t,p,aux, PETSC_TRUE, &
  2469)                                     dw, dwmol, dwp, dwt, ierr)
  2470)   print *, 'Density-Ext'
  2471)   print *, 'dw:    ', dw
  2472)   print *, 'dwmol: ', dwmol
  2473)   print *, 'dwp:   ', dwp
  2474)   print *, 'dwt:   ', dwt
  2475) 
  2476)   call EOSWaterViscosityBatzleAndWang(t, p, PS, dPS_dT, &
  2477)                                       PETSC_TRUE, vw, &
  2478)                                       dvw_dt, dvw_dp, ierr)  
  2479)   print *
  2480)   print *, 'vw:      ', vw
  2481)   print *, 'dvw_dp:  ', dvw_dp
  2482)   print *, 'dvw_dt:  ', dvw_dt
  2483)   
  2484)   
  2485)   call EOSWaterViscosityBatzleAndWangExt(t, p, PS, dPS_dT, aux, &
  2486)                                          PETSC_TRUE, vw, &
  2487)                                          dvw_dt, dvw_dp, ierr) 
  2488)   print *, 'Ext-'
  2489)   print *, 'vw:      ', vw
  2490)   print *, 'dvw_dp:  ', dvw_dp
  2491)   print *, 'dvw(t)t:  ', dvw_dt
  2492)   
  2493)   call EOSWaterViscosityBatzleAndWangExt(t2, p, PS, dPS_dT, aux, &
  2494)                                          PETSC_TRUE, vw2, &
  2495)                                          dvw_dt, dvw_dp, ierr) 
  2496)   
  2497)   print *, 'Ext-numerical'
  2498)   print *, 'vw:      ', vw2
  2499)   print *, 'dvw(t)t:  ', (vw2-vw)/(t2-t)
  2500)   
  2501)   call EOSWaterViscosityBatzleAndWangExt(t, p, PS, dPS_dT, aux, &
  2502)                                          PETSC_TRUE, vw, &
  2503)                                          dvw_dt, dvw_dp, ierr) 
  2504)   print *, 'Ext-S'
  2505)   print *, 'vw:      ', vw
  2506)   print *, 'dvw_dp:  ', dvw_dp
  2507)   print *, 'dvw(t)t:  ', dvw_dt
  2508)   
  2509) end subroutine TestEOSWaterBatzleAndWang
  2510) 
  2511) ! ************************************************************************** !
  2512) 
  2513) subroutine EOSWaterDensityExtNumericalDerive(t,p,aux,dw,dwmol,dwp,dwt,ierr)
  2514) 
  2515)   implicit none
  2516) 
  2517)   PetscReal, intent(in) :: t     ! Temperature in centigrade
  2518)   PetscReal, intent(in) :: p     ! Pressure in Pascal
  2519)   PetscReal, intent(in) :: aux(*)
  2520)   PetscReal, intent(out) :: dw ! kg/m^3
  2521)   PetscReal, intent(out) :: dwmol ! kmol/m^3
  2522)   PetscReal, intent(out) :: dwp ! kmol/m^3-Pa
  2523)   PetscReal, intent(out) :: dwt ! kmol/m^3-C
  2524)   PetscErrorCode, intent(out) :: ierr
  2525) 
  2526)   PetscReal :: dwp_analytical, dwt_analytical
  2527)   PetscReal :: dwp_numerical, dwt_numerical
  2528)   PetscReal :: dw_no_salt, dwp_no_salt, dwt_no_salt
  2529)   PetscReal :: dum1, dum2, dum3
  2530)   PetscReal :: dwmol_tpert, dwmol_ppert
  2531)   PetscReal :: tpert, ppert, t_plus_tpert, p_plus_ppert
  2532)   PetscReal :: salinity(1)
  2533)   PetscReal :: pert_tol = 1.d-6
  2534)   
  2535)   tpert = t*pert_tol
  2536)   ppert = p*pert_tol
  2537)   t_plus_tpert = t + tpert
  2538)   p_plus_ppert = p + ppert
  2539)   salinity(1) = aux(1)
  2540) 
  2541) #if 0 
  2542)   ! test against non-extended version
  2543)   call EOSWaterDensityPtr(t,p,PETSC_TRUE,dw,dwmol,dwp_analytical, &
  2544)                           dwt_analytical,ierr)
  2545)   dwp = dwp_analytical
  2546)   dwt = dwt_analytical
  2547) #else
  2548)   call EOSWaterDensityExtPtr(t,p,salinity,PETSC_TRUE, &
  2549)                              dw,dwmol,dwp_analytical,dwt_analytical,ierr)
  2550)   dwp = dwp_analytical
  2551)   dwt = dwt_analytical
  2552)   call EOSWaterDensityExtPtr(t_plus_tpert,p,salinity,PETSC_FALSE, &
  2553)                              dum1,dwmol_tpert,dum2,dum3,ierr)
  2554)   call EOSWaterDensityExtPtr(t,p_plus_ppert,salinity,PETSC_FALSE, &
  2555)                              dum1,dwmol_ppert,dum2,dum3,ierr)
  2556) 
  2557)   dwp_numerical = (dwmol_ppert-dwmol)/ppert
  2558)   dwt_numerical = (dwmol_tpert-dwmol)/tpert
  2559) 
  2560)   if (.not.PETSC_TRUE) then
  2561)     dwp = dwp_numerical
  2562)     dwt = dwt_numerical
  2563)   else
  2564)     dwp = dwp_analytical
  2565)     dwt = dwt_analytical
  2566)   endif
  2567) 
  2568)   if (dabs((dwp_numerical-dwp_analytical)/dwp_numerical) > 1.d-4) then
  2569)     print *, p, t, salinity(1), dw, dwmol, dwp, dwp_analytical, dwp_numerical
  2570)   endif
  2571) #endif
  2572)   
  2573) end subroutine EOSWaterDensityExtNumericalDerive
  2574) 
  2575) ! **************************************************************************** !
  2576) 
  2577) subroutine EOSWaterInputRecord()
  2578)   ! 
  2579)   ! Prints ingested equation of state information to the input record file.
  2580)   ! 
  2581)   ! Author: Jenn Frederick
  2582)   ! Date: 05/04/2016
  2583)   ! 
  2584)   
  2585)   implicit none
  2586)   
  2587)   character(len=MAXWORDLENGTH) :: word1, word2
  2588)   character(len=MAXSTRINGLENGTH) :: string
  2589)   PetscInt :: id = INPUT_RECORD_UNIT
  2590) 
  2591)   write(id,'(a)') '---------------------------------------------------------&
  2592)                   &-----------------------'
  2593)   write(id,'(a29)',advance='no') '---------------------------: '
  2594)   write(id,'(a)') 'WATER'
  2595)   
  2596)   ! water density [kg/m^3]
  2597)   if (associated(EOSWaterDensityPtr,EOSWaterDensityConstant)) then
  2598)     write(id,'(a29)',advance='no') 'water density: '
  2599)     write(word1,*) constant_density
  2600)     write(id,'(a)') 'constant, ' // adjustl(trim(word1)) // ' kg/m^3'
  2601)   endif
  2602)   if (associated(EOSWaterDensityPtr,EOSWaterDensityExponential)) then
  2603)     write(id,'(a29)',advance='no') 'water density: '
  2604)     write(id,'(a)') 'exponential'
  2605)     write(id,'(a29)',advance='no') 'exp. ref. density: '
  2606)     write(word1,*) exponent_reference_density
  2607)     write(id,'(a)') adjustl(trim(word1)) // ' kg/m^3'
  2608)     write(id,'(a29)',advance='no') 'exp. ref. pressure: '
  2609)     write(word1,*) exponent_reference_pressure
  2610)     write(id,'(a)') adjustl(trim(word1)) // ' Pa'
  2611)     write(id,'(a29)',advance='no') 'exp. water compressibility: '
  2612)     write(word1,*) exponent_water_compressibility
  2613)     write(id,'(a)') adjustl(trim(word1)) // ' 1/Pa'
  2614)   endif
  2615)   if (associated(EOSWaterDensityPtr,EOSWaterDensityIFC67)) then
  2616)     write(id,'(a29)',advance='no') 'water density: '
  2617)     write(id,'(a)') 'default, IFC67'
  2618)   endif
  2619)   if (associated(EOSWaterDensityPtr,EOSWaterDensityTGDPB01)) then
  2620)     write(id,'(a29)',advance='no') 'water density: '
  2621)     write(id,'(a)') 'TGDPB01'
  2622)   endif
  2623)   if (associated(EOSWaterDensityPtr,EOSWaterDensityPainter)) then
  2624)     write(id,'(a29)',advance='no') 'water density: '
  2625)     write(id,'(a)') 'PAINTER'
  2626)   endif
  2627)   if (associated(EOSWaterDensityPtr,EOSWaterDensityBatzleAndWang) .and. &
  2628)       associated(EOSWaterDensityExtPtr,EOSWaterDensityBatzleAndWangExt)) then
  2629)     write(id,'(a29)',advance='no') 'water density: '
  2630)     write(id,'(a)') 'Batzle and Wang'
  2631)   endif
  2632)   
  2633)   ! water viscosity [Pa-s]
  2634)   if (associated(EOSWaterViscosityPtr,EOSWaterViscosityConstant)) then
  2635)     write(id,'(a29)',advance='no') 'water viscosity: '
  2636)     write(word1,*) constant_viscosity
  2637)     write(id,'(a)') 'constant, ' // trim(word1) // ' Pa-sec'
  2638)   endif
  2639)   if (associated(EOSWaterViscosityPtr,EOSWaterViscosity1)) then
  2640)     write(id,'(a29)',advance='no') 'water viscosity: '
  2641)     write(id,'(a)') 'default'
  2642)   endif
  2643)   if (associated(EOSWaterViscosityPtr,EOSWaterViscosityBatzleAndWang)) then
  2644)     write(id,'(a29)',advance='no') 'water viscosity: '
  2645)     write(id,'(a)') 'Batzle and Wang'
  2646)   endif
  2647)   
  2648)   ! water enthalpy [J/kmol]
  2649)   if (associated(EOSWaterViscosityPtr,EOSWaterEnthalpyConstant)) then
  2650)     write(id,'(a29)',advance='no') 'water enthalpy: '
  2651)     write(word1,*) constant_enthalpy
  2652)     write(id,'(a)') 'constant, ' // trim(word1) // ' J/kmol'
  2653)   endif
  2654)   if (associated(EOSWaterViscosityPtr,EOSWaterEnthalpyIFC67)) then
  2655)     write(id,'(a29)',advance='no') 'water enthalpy: '
  2656)     write(id,'(a)') 'default, IFC67'
  2657)   endif
  2658)   if (associated(EOSWaterViscosityPtr,EOSWaterEnthalpyPainter)) then
  2659)     write(id,'(a29)',advance='no') 'water enthalpy: '
  2660)     write(id,'(a)') 'PAINTER'
  2661)   endif
  2662)   
  2663)   ! steam density
  2664)   if (associated(EOSWaterSteamDensityEnthalpyPtr, &
  2665)                  EOSWaterSteamDenEnthConstant)) then
  2666)     write(id,'(a29)',advance='no') 'steam density: '
  2667)     write(word1,*) constant_steam_density
  2668)     write(id,'(a)') 'constant, ' // trim(word1) // ' kg/m^3'
  2669)   endif
  2670)   if (associated(EOSWaterSteamDensityEnthalpyPtr, &
  2671)                  EOSWaterSteamDensityEnthalpyIFC67)) then
  2672)     write(id,'(a29)',advance='no') 'steam density: '
  2673)     write(id,'(a)') 'default, IFC67'
  2674)   endif
  2675)   
  2676)   ! steam enthalpy [J/kmol]
  2677)   if (associated(EOSWaterSteamDensityEnthalpyPtr, &
  2678)                  EOSWaterSteamDenEnthConstant)) then
  2679)     write(id,'(a29)',advance='no') 'steam enthalpy: '
  2680)     write(word1,*) constant_steam_enthalpy
  2681)     write(id,'(a)') 'constant, ' // trim(word1) // ' J/kmol'
  2682)   endif
  2683)   if (associated(EOSWaterSteamDensityEnthalpyPtr, &
  2684)                  EOSWaterSteamDensityEnthalpyIFC67)) then
  2685)     write(id,'(a29)',advance='no') 'steam enthalpy: '
  2686)     write(id,'(a)') 'default, IFC67'
  2687)   endif
  2688)   
  2689)   write(id,'(a)') '---------------------------------------------------------&
  2690)                   &-----------------------'
  2691)   
  2692) end subroutine EOSWaterInputRecord
  2693) 
  2694) ! ************************************************************************** !
  2695) 
  2696) subroutine EOSWaterTest(temp_low,temp_high,pres_low,pres_high, &
  2697)                         ntemp,npres,uniform_temp,uniform_pres,filename)
  2698) 
  2699)   implicit none
  2700) 
  2701)   PetscReal :: temp_low
  2702)   PetscReal :: temp_high
  2703)   PetscReal :: pres_low
  2704)   PetscReal :: pres_high
  2705)   PetscInt :: npres
  2706)   PetscInt :: ntemp
  2707)   PetscBool :: uniform_temp
  2708)   PetscBool :: uniform_pres
  2709)   character(len=MAXWORDLENGTH) :: filename
  2710) 
  2711)   PetscReal, allocatable :: temp(:)
  2712)   PetscReal, allocatable :: pres(:)
  2713)   PetscReal, allocatable :: density_kg(:,:)
  2714)   PetscReal, allocatable :: enthalpy(:,:)
  2715)   PetscReal, allocatable :: viscosity(:,:)
  2716)   PetscReal, allocatable :: saturation_pressure_array(:)
  2717)   PetscReal :: dum1, dum2, dum3, dum4
  2718)   PetscInt :: itemp, ipres
  2719)   PetscReal :: ln_low, ln_high
  2720)   PetscReal :: saturation_pressure
  2721)   PetscReal :: NaN
  2722)   character(len=MAXWORDLENGTH) :: eos_density_name
  2723)   character(len=MAXWORDLENGTH) :: eos_enthalpy_name
  2724)   character(len=MAXWORDLENGTH) :: eos_viscosity_name
  2725)   character(len=MAXWORDLENGTH) :: eos_saturation_pressure_name
  2726)   character(len=MAXSTRINGLENGTH) :: header
  2727) 
  2728)   PetscErrorCode :: ierr
  2729) 
  2730)   NaN = 0.d0
  2731)   NaN = 1.d0/NaN
  2732)   NaN = 0.d0*NaN
  2733) 
  2734)   allocate(temp(ntemp))
  2735)   temp = UNINITIALIZED_DOUBLE
  2736)   allocate(pres(ntemp))
  2737)   pres = UNINITIALIZED_DOUBLE
  2738)   allocate(density_kg(npres,ntemp))
  2739)   density_kg = UNINITIALIZED_DOUBLE
  2740)   allocate(viscosity(npres,ntemp))
  2741)   viscosity = UNINITIALIZED_DOUBLE
  2742)   allocate(enthalpy(npres,ntemp))
  2743)   enthalpy = UNINITIALIZED_DOUBLE
  2744)   allocate(saturation_pressure_array(ntemp))
  2745)   saturation_pressure_array = UNINITIALIZED_DOUBLE
  2746) 
  2747)   if (uniform_pres) then
  2748)     do ipres = 1, npres
  2749)       pres(ipres) = (pres_high-pres_low)/dble(npres-1) * (ipres-1) + pres_low
  2750)     enddo
  2751)   else
  2752)     ln_high = log(pres_high)
  2753)     ln_low = log(pres_low)
  2754)     do ipres = 1, npres
  2755)       pres(ipres) = exp((ln_high-ln_low)/dble(npres-1) * (ipres-1) + ln_low)
  2756)     enddo
  2757)   endif
  2758) 
  2759)   if (uniform_temp) then
  2760)     do itemp = 1, ntemp
  2761)       temp(itemp) = (temp_high-temp_low)/dble(ntemp-1) * (itemp-1) + temp_low
  2762)     enddo
  2763)   else
  2764)     ln_high = log(temp_high)
  2765)     ln_low = log(temp_low)
  2766)     do itemp = 1, ntemp
  2767)       temp(itemp) = exp((ln_high-ln_low)/dble(ntemp-1) * (itemp-1) + ln_low)
  2768)     enddo
  2769)   endif
  2770) 
  2771)   ! density
  2772)   if (associated(EOSWaterDensityPtr,EOSWaterDensityConstant)) then
  2773)     eos_density_name = 'Constant'
  2774)   else if (associated(EOSWaterDensityPtr,EOSWaterDensityExponential)) then
  2775)     eos_density_name = 'Exponential'
  2776)   else if (associated(EOSWaterDensityPtr,EOSWaterDensityIFC67)) then
  2777)     eos_density_name = 'IFC67'
  2778)   else if (associated(EOSWaterDensityPtr,EOSWaterDensityTGDPB01)) then
  2779)     eos_density_name = 'TGDPB01'
  2780)   else if (associated(EOSWaterDensityPtr,EOSWaterDensityPainter)) then
  2781)     eos_density_name = 'Painter'
  2782)   else if (associated(EOSWaterDensityPtr,EOSWaterDensityBatzleAndWang)) then
  2783)     eos_density_name = 'Batzle and Wang'
  2784)   else 
  2785)     eos_density_name = 'Unknown'
  2786)   endif
  2787) 
  2788)   ! enthalpy
  2789)   if (associated(EOSWaterEnthalpyPtr,EOSWaterEnthalpyConstant)) then
  2790)     eos_enthalpy_name = 'Constant'
  2791)   else if (associated(EOSWaterEnthalpyPtr,EOSWaterEnthalpyIFC67)) then
  2792)     eos_enthalpy_name = 'IFC67'
  2793)   else if (associated(EOSWaterEnthalpyPtr,EOSWaterEnthalpyPainter)) then
  2794)     eos_enthalpy_name = 'Painter'
  2795)   else
  2796)     eos_enthalpy_name = 'Unknown'
  2797)   endif
  2798) 
  2799)   ! viscosity
  2800)   if (associated(EOSWaterViscosityPtr,EOSWaterViscosityConstant)) then
  2801)     eos_viscosity_name = 'Constant'
  2802)   else if (associated(EOSWaterViscosityPtr,EOSWaterViscosity1)) then
  2803)     eos_viscosity_name = 'Default'
  2804)   else if (associated(EOSWaterViscosityPtr,EOSWaterViscosityBatzleAndWang)) then
  2805)     eos_viscosity_name = 'Batzle and Wang'
  2806)   else
  2807)     eos_viscosity_name = 'Unknown'
  2808)   endif
  2809) 
  2810)   ! saturation pressure
  2811)   if (associated(EOSWaterSaturationPressurePtr, &
  2812)                  EOSWaterSaturationPressureIFC67)) then
  2813)     eos_saturation_pressure_name = 'IFC67'
  2814)   else
  2815)     eos_saturation_pressure_name = 'Unknown'
  2816)   endif
  2817) 
  2818)   do itemp = 1, ntemp
  2819)     do ipres = 1, npres
  2820)       ! EOSWaterSaturationPressurePtr() must come before call to 
  2821)       ! EOSWaterViscosityPtr()
  2822)       call EOSWaterSaturationPressurePtr(temp(itemp),PETSC_FALSE, &
  2823)                                          saturation_pressure,dum1,ierr)
  2824)       if (ipres == 1) &
  2825)         saturation_pressure_array(itemp) = saturation_pressure
  2826)       call EOSWaterDensityPtr(temp(itemp),pres(ipres),PETSC_FALSE, &
  2827)                               density_kg(ipres,itemp), &
  2828)                               dum1,dum2,dum3,ierr)
  2829)       call EOSWaterEnthalpyPtr(temp(itemp),pres(ipres),PETSC_FALSE, &
  2830)                                enthalpy(ipres,itemp),dum1,dum2,ierr)
  2831)       call EOSWaterViscosityPtr(temp(itemp),pres(ipres),saturation_pressure, &
  2832)                                 dum1,PETSC_FALSE,viscosity(ipres,itemp), &
  2833)                                 dum2,dum3,ierr)
  2834)     enddo
  2835)   enddo
  2836) 
  2837) 100 format(100es16.8)
  2838)   if (len_trim(filename) == 0) filename = 'eos_water_test.txt'
  2839)   open(unit=IUNIT_TEMP,file=filename)
  2840)   header = 'T[C], P[Pa], &
  2841)     &Density (' // trim(eos_density_name) // ') [kg/m^3], &
  2842)     &Enthalpy (' // trim(eos_enthalpy_name) // ') [J/kmol], &
  2843)     &Viscosity (' // trim(eos_viscosity_name) // ') [Pa-s], &
  2844)     &Saturation Pressure (' // trim(eos_saturation_pressure_name) // ') [Pa]'
  2845)   write(IUNIT_TEMP,'(a)') trim(header)
  2846)   write(IUNIT_TEMP,'(100i9)') ntemp, npres
  2847)   do itemp = 1, ntemp
  2848)     do ipres = 1, npres
  2849)       write(IUNIT_TEMP,100) temp(itemp), pres(ipres), &
  2850)             density_kg(ipres,itemp), enthalpy(ipres,itemp), &
  2851)             viscosity(ipres,itemp), saturation_pressure_array(itemp)
  2852)     enddo
  2853)   enddo
  2854)   close(IUNIT_TEMP)
  2855) 
  2856)   deallocate(temp)
  2857)   deallocate(pres)
  2858)   deallocate(density_kg)
  2859)   deallocate(enthalpy)
  2860)   deallocate(viscosity)
  2861)   deallocate(saturation_pressure_array)
  2862) 
  2863) end subroutine EOSWaterTest
  2864) 
  2865) ! ************************************************************************** !
  2866) 
  2867) end module EOS_Water_module

generated by
Intel(R) C++/Fortran Compiler code-coverage tool
Web-Page Owner: Nobody