co2_span_wagner.F90       coverage:  6.67 %func     13.38 %block


     1)   module co2_span_wagner_module
     2) 
     3) 
     4)   ! module contains only 1 interface with other part, read as:
     5) 
     6)   ! co2_span_wagner(p,t,rho,dddt,dddp,fg,dfgdp,dfgdt,
     7)   !                 eng,ent,dhdt,dhdp,visc,dvdt,dvdp)
     8) 
     9)   !  where units are assigned as:   
    10)   !     P [MPa]       T [K]         
    11)   !     rho [kg/m^3] Energy [MJ/kg] Enthalpy [MJ/kg] Vis [Pa s]
    12) 
    13)   use PFLOTRAN_Constants_module
    14) 
    15)       implicit none
    16) 
    17) #include "petsc/finclude/petscsys.h"
    18) 
    19)       save
    20)       
    21) !     table lookup parameters t[k], p[MPa]
    22) 
    23) !     t = 35 - 410 C, p = 0.01 - 250 bars
    24) !     PetscInt,   public :: ntab_t = 150, ntab_p = 500
    25) !     PetscReal,  public :: t0_tab = 35.d0+273.15D0, p0_tab = 0.01d0
    26) !     PetscReal,  public :: dt_tab = 2.5d0, dp_tab = 0.5d0
    27) 
    28) !     t = 0 - 375 C, p = 0.01 - 1250.01 bars
    29)       PetscInt,   public :: ntab_t = 150, ntab_p = 500
    30)       PetscReal,  public :: t0_tab = 0.d0+273.15D0, p0_tab = 0.01d0
    31)       PetscReal,  public :: dt_tab = 2.5d0, dp_tab = 2.5d0
    32) 
    33)       PetscReal, private :: n(42),ti(40),gamma(5),phic(8),c(40),d(40),a(8)
    34)       PetscReal, private :: alpha(5),beta(8),delta(4),epsilon(5)
    35)       PetscReal, private :: aco2(4),bco2(4),capa(5),capb(5),capc(5),capd(5)
    36)       PetscReal, private :: av(0:4)
    37) 
    38)       PetscReal, private :: denc,tc,rg,pc
    39)       PetscReal ,private, allocatable :: co2_prop_spwag(:,:,:)
    40)       PetscReal, private :: p,t,rhosav
    41) 
    42)       public initialize_span_wagner, co2_span_wagner, vappr
    43) 
    44)       private
    45) 
    46)       contains
    47)  
    48) 
    49) ! ************************************************************************** !
    50) 
    51) subroutine initialize_span_wagner(itable,myrank,option)
    52) 
    53)       use Option_module
    54) 
    55)       implicit none
    56)       PetscInt, optional :: itable
    57)       
    58)       PetscReal :: pl,tl,tmp,tmp2,dtmp,dtemp,dpres,dddt,dddp
    59)       
    60)       PetscReal :: rhodp,rhodt,fgdp,fgdt,engdp,engdt,entdp,entdt,vdp,vdt
    61)        
    62)       PetscInt :: iitable,i,j
    63)       PetscMPIInt :: myrank
    64)       PetscInt :: iflag = 1
    65)       
    66)       character*3 :: q
    67)       character*1 :: tab
    68)       
    69)       PetscReal :: temparray(15)
    70)       PetscInt :: status
    71) 
    72)       character(len=MAXSTRINGLENGTH) :: co2_database_filename
    73) 
    74)       type(option_type) :: option
    75)       
    76)       tab = char(9)
    77)       q = '","'
    78) 
    79)       if (.not. allocated(co2_prop_spwag)) &
    80)         allocate(co2_prop_spwag(0:ntab_p,0:ntab_t,1:15))
    81)       
    82)       iitable = 0      
    83)       if (present(itable)) iitable=itable
    84) 
    85)       denc = 467.6d0
    86)       tc = 304.1282d0
    87)       rg = 0.1889241d0
    88)       pc = 7.3773d0
    89) 
    90)       av(0) = 0.235156d0
    91)       av(1) = -0.491266d0
    92)       av(2) = 5.211155d-2
    93)       av(3) = 5.347906d-2
    94)       av(4) = -1.537102d-2
    95) 
    96)       a(1) = 8.37304456d0
    97)       a(2) = -3.70454304d0
    98)       a(3) = 2.5d0
    99)       a(4) = 1.99427042d0
   100)       a(5) = 0.62105248d0
   101)       a(6) = 0.41195293d0
   102)       a(7) = 1.04028922d0
   103)       a(8) = 0.08327678d0
   104) 
   105)       phic(1) = 0.0d0
   106)       phic(2) = 0.0d0
   107)       phic(3) = 0.0d0
   108)       phic(4) = 3.15163d0
   109)       phic(5) = 6.1119d0
   110)       phic(6) = 6.77708d0
   111)       phic(7) = 11.32384d0
   112)       phic(8) = 27.08792d0
   113)       n(1) = 0.38856823203161d0
   114)       n(2) = 0.2938547594274d1
   115)       n(3) = -0.55867188534934d1
   116)       n(4) = -0.76753199592477d0
   117)       n(5) = 0.31729005580416d0
   118)       n(6) = 0.54803315897767d0
   119)       n(7) = 0.12279411220335d0
   120)       n(8) = 0.2165896154322d1
   121)       n(9) = 0.15841735109724d1
   122)       n(10) = -0.23132705405503d0
   123)       n(11) = 0.58116916431436d-1
   124)       n(12) = -0.55369137205382d0
   125)       n(13) = 0.48946615909422d0
   126)       n(14) = -0.24275739843501d-1
   127)       n(15) = 0.62494790501678d-1
   128)       n(16) = -0.12175860225246d0
   129)       n(17) = -0.37055685270086d0
   130)       n(18) = -0.16775879700426d-1
   131)       n(19) = -0.11960736637987d0
   132)       n(20) = -0.45619362508778d-1
   133)       n(21) = 0.35612789270346d-1
   134)       n(22) = -0.74427727132052d-2
   135)       n(23) = -0.17395704902432d-2
   136)       n(24) = -0.21810121289527d-1
   137)       n(25) = 0.24332166559236d-1
   138)       n(26) = -0.37440133423463d-1
   139)       n(27) = 0.14338715756878d0
   140)       n(28) = -0.13491969083286d0
   141)       n(29) = -0.2315122505348d-1
   142)       n(30) = 0.12363125492901d-1
   143)       n(31) = 0.2105832197294d-2
   144)       n(32) = -0.33958519026368d-3
   145)       n(33) = 0.55993651771592d-2
   146)       n(34) = -0.30335118055646d-3
   147)       n(35) = -0.2136548868832d3
   148)       n(36) = 0.26641569149272d5
   149)       n(37) = -0.24027212204557d5
   150)       n(38) = -0.28341603423999d3
   151)       n(39) = 0.21247284400179d3
   152)       n(40) = -0.66642276540751d0
   153)       n(41) = 0.72608632349897d0
   154)       n(42) = 0.55068668612842d-1
   155)       c(1) = 0.d0
   156)       c(2) = 0.d0
   157)       c(3) = 0.d0
   158)       c(4) = 0.d0
   159)       c(5) = 0.d0
   160)       c(6) = 0.d0
   161)       c(7) = 0.d0
   162)       c(8) = 1.d0
   163)       c(9) = 1.d0
   164)       c(10) = 1.d0
   165)       c(11) = 1.d0
   166)       c(12) = 1.d0
   167)       c(13) = 1.d0
   168)       c(14) = 1.d0
   169)       c(15) = 1.d0
   170)       c(16) = 1.d0
   171)       c(17) = 2.d0
   172)       c(18) = 2.d0
   173)       c(19) = 2.d0
   174)       c(20) = 2.d0
   175)       c(21) = 2.d0
   176)       c(22) = 2.d0
   177)       c(23) = 2.d0
   178)       c(24) = 3.d0
   179)       c(25) = 3.d0
   180)       c(26) = 3.d0
   181)       c(27) = 4.d0
   182)       c(28) = 4.d0
   183)       c(29) = 4.d0
   184)       c(30) = 4.d0
   185)       c(31) = 4.d0
   186)       c(32) = 4.d0
   187)       c(33) = 5.d0
   188)       c(34) = 6.d0
   189)       d(1) = 1.d0
   190)       d(2) = 1.d0
   191)       d(3) = 1.d0
   192)       d(4) = 1.d0
   193)       d(5) = 2.d0
   194)       d(6) = 2.d0
   195)       d(7) = 3.d0
   196)       d(8) = 1.d0
   197)       d(9) = 2.d0
   198)       d(10) = 4.d0
   199)       d(11) = 5.d0
   200)       d(12) = 5.d0
   201)       d(13) = 5.d0
   202)       d(14) = 6.d0
   203)       d(15) = 6.d0
   204)       d(16) = 6.d0
   205)       d(17) = 1.d0
   206)       d(18) = 1.d0
   207)       d(19) = 4.d0
   208)       d(20) = 4.d0
   209)       d(21) = 4.d0
   210)       d(22) = 7.d0
   211)       d(23) = 8.d0
   212)       d(24) = 2.d0
   213)       d(25) = 3.d0
   214)       d(26) = 3.d0
   215)       d(27) = 5.d0
   216)       d(28) = 5.d0
   217)       d(29) = 6.d0
   218)       d(30) = 7.d0
   219)       d(31) = 8.d0
   220)       d(32) = 10.d0
   221)       d(33) = 4.d0
   222)       d(34) = 8.d0
   223)       d(35) = 2.d0
   224)       d(36) = 2.d0
   225)       d(37) = 2.d0
   226)       d(38) = 3.d0
   227)       d(39) = 3.d0
   228)       ti(1) = 0.d0
   229)       ti(2) = 0.75d0
   230)       ti(3) = 1.0d0
   231)       ti(4) = 2.0d0
   232)       ti(5) = 0.75d0
   233)       ti(6) = 2.0d0
   234)       ti(7) = 0.75d0
   235)       ti(8) = 1.5d0
   236)       ti(9) = 1.5d0
   237)       ti(10) = 2.5d0
   238)       ti(11) = 0.0d0
   239)       ti(12) = 1.5d0
   240)       ti(13) = 2.0d0
   241)       ti(14) = 0.0d0
   242)       ti(15) = 1.0d0
   243)       ti(16) = 2.0d0
   244)       ti(17) = 3.0d0
   245)       ti(18) = 6.0d0
   246)       ti(19) = 3.0d0
   247)       ti(20) = 6.0d0
   248)       ti(21) = 8.0d0
   249)       ti(22) = 6.0d0
   250)       ti(23) = 0.0d0
   251)       ti(24) = 7.0d0
   252)       ti(25) = 12.0d0
   253)       ti(26) = 16.0d0
   254)       ti(27) = 22.0d0
   255)       ti(28) = 24.0d0
   256)       ti(29) = 16.0d0
   257)       ti(30) = 24.0d0
   258)       ti(31) = 8.0d0
   259)       ti(32) = 2.0d0
   260)       ti(33) = 28.0d0
   261)       ti(34) = 14.0d0
   262)       ti(35) = 1.0d0
   263)       ti(36) = 0.0d0
   264)       ti(37) = 1.0d0
   265)       ti(38) = 3.0d0
   266)       ti(39) = 3.0d0
   267)       alpha(1) = 25.d0
   268)       alpha(2) = 25.d0
   269)       alpha(3) = 25.d0
   270)       alpha(4) = 15.d0
   271)       alpha(5) = 20.d0
   272)       beta(1) = 325.d0
   273)       beta(2) = 300.d0
   274)       beta(3) = 300.d0
   275)       beta(4) = 275.d0
   276)       beta(5) = 275.d0
   277)       beta(6) = 0.3d0
   278)       beta(7) = 0.3d0
   279)       beta(8) = 0.3d0
   280)       gamma(1) = 1.16d0
   281)       gamma(2) = 1.19d0
   282)       gamma(3) = 1.19d0
   283)       gamma(4) = 1.25d0
   284)       gamma(5) = 1.22d0
   285)       epsilon(1) = 1.d0
   286)       epsilon(2) = 1.d0
   287)       epsilon(3) = 1.d0
   288)       epsilon(4) = 1.d0
   289)       epsilon(5) = 1.d0
   290)       aco2(1)=3.5d0
   291)       aco2(2)=3.5d0
   292)       aco2(3)=3.0d0
   293)       bco2(1)=0.875d0
   294)       bco2(2)=0.925d0
   295)       bco2(3)=0.875d0
   296)       capa(1) = 0.7d0
   297)       capa(2) = 0.7d0
   298)       capa(3) = 0.7d0
   299)       capb(1) = 0.3d0
   300)       capb(2) = 0.3d0
   301)       capb(3) = 1.0d0
   302)       capc(1) = 10.d0
   303)       capc(2) = 10.d0
   304)       capc(3) = 12.5d0
   305)       capd(1) = 275.d0
   306)       capd(2) = 275.d0
   307)       capd(3) = 275.d0
   308)   
   309)   if (myrank==0) print *,'Preparing Table...',iitable
   310)   if (iitable == 1) then
   311)       ! Generate the table: 3 dimensional array
   312)          !  1 p,     2  T
   313)          !  3 rho    4 dddt,  5 dddp,
   314)          !  6 fg,    7 dfgdp, 8 dfgdt
   315)          !  9 eng,  
   316)          ! 10 ent,  11 dhdt, 12 dhdp,
   317)          ! 13 visc, 14 dvdt, 15 dvdp
   318)           
   319)     tmp2=0.D0    
   320)     do i = 0, ntab_p
   321)       tmp=tmp2
   322)       pl = p0_tab + dp_tab * real(i)
   323)       
   324)       do j = 0, ntab_t
   325)         tl = t0_tab + dt_tab * real(j)
   326)         
   327)         co2_prop_spwag(i,j,1) = pl
   328)         co2_prop_spwag(i,j,2) = tl
   329)         co2_prop_spwag(i,j,3) = tmp 
   330)         call co2_span_wagner(pl,tl,co2_prop_spwag(i,j,3),co2_prop_spwag(i,j,4),&
   331)         co2_prop_spwag(i,j,5),co2_prop_spwag(i,j,6),co2_prop_spwag(i,j,7),&
   332)         co2_prop_spwag(i,j,8),co2_prop_spwag(i,j,9),co2_prop_spwag(i,j,10),&
   333)         co2_prop_spwag(i,j,11),co2_prop_spwag(i,j,12),co2_prop_spwag(i,j,13),&
   334)         co2_prop_spwag(i,j,14),co2_prop_spwag(i,j,15),iflag)
   335) 
   336) !       p = p + dp
   337)         dpres = 1.e-3
   338)         call co2_span_wagner(pl+dpres,tl,rhodp,co2_prop_spwag(i,j,4),&
   339)         co2_prop_spwag(i,j,5),fgdp,co2_prop_spwag(i,j,7),&
   340)         co2_prop_spwag(i,j,8),engdp,entdp,&
   341)         co2_prop_spwag(i,j,11),co2_prop_spwag(i,j,12),vdp,&
   342)         co2_prop_spwag(i,j,14),co2_prop_spwag(i,j,15),iflag)
   343) 
   344) !       t = t + dt
   345)         dtemp = 1.e-4
   346)         call co2_span_wagner(pl,tl+dtemp,rhodt,co2_prop_spwag(i,j,4),&
   347)         co2_prop_spwag(i,j,5),fgdt,co2_prop_spwag(i,j,7),&
   348)         co2_prop_spwag(i,j,8),engdt,entdt,&
   349)         co2_prop_spwag(i,j,11),co2_prop_spwag(i,j,12),vdt,&
   350)         co2_prop_spwag(i,j,14),co2_prop_spwag(i,j,15),iflag)
   351) 
   352) ! compute derivatives numerically: 1-p,2-T,3-d,4-dddt,5-dddp,6-fg,7-dfgdp,8-dfgdt,
   353) !                 9-energy,10-enthalpy,11-dhdt,12-dhdp,13-visc,14-dvdt,15-dvdp
   354) !#if 0
   355)           
   356)           !dddt
   357)           co2_prop_spwag(i,j,4) = (rhodt-co2_prop_spwag(i,j,3))/dtemp
   358)           
   359)           !dfgdt
   360)           co2_prop_spwag(i,j,8) = (fgdt-co2_prop_spwag(i,j,6))/dtemp
   361)           
   362)           !dhdt
   363)           co2_prop_spwag(i,j,11) = (entdt-co2_prop_spwag(i,j,10))/dtemp
   364)           
   365)           !dvdt
   366)           co2_prop_spwag(i,j,14) = (vdt-co2_prop_spwag(i,j,13))/dtemp
   367) 
   368)           !dddp
   369)           co2_prop_spwag(i,j,5) = (rhodp-co2_prop_spwag(i,j,3))/dpres
   370)           
   371)           !dfgdp
   372)           co2_prop_spwag(i,j,7) = (fgdp-co2_prop_spwag(i,j,6))/dpres
   373)           
   374)           !dhdp
   375)           co2_prop_spwag(i,j,12) = (entdp-co2_prop_spwag(i,j,10))/dpres
   376)           
   377)           !dvdp
   378)           co2_prop_spwag(i,j,15) = (vdp-co2_prop_spwag(i,j,13))/dpres
   379) !#endif
   380) #if 0
   381)         if (j>0) then
   382)           dtemp = tl-co2_prop_spwag(i,j-1,2)
   383)           
   384)           !dddt
   385)           co2_prop_spwag(i,j,4) = (co2_prop_spwag(i,j,3)-co2_prop_spwag(i,j-1,3))/dtemp
   386)           
   387)           !dfgdt
   388)           co2_prop_spwag(i,j,8) = (co2_prop_spwag(i,j,6)-co2_prop_spwag(i,j-1,6))/dtemp
   389)           
   390)           !dhdt
   391)           co2_prop_spwag(i,j,11) = (co2_prop_spwag(i,j,10)-co2_prop_spwag(i,j-1,10))/dtemp
   392)           
   393)           !dvdt
   394)           co2_prop_spwag(i,j,14) = (co2_prop_spwag(i,j,13)-co2_prop_spwag(i,j-1,13))/dtemp
   395)         endif
   396)         if (i>0) then
   397)           dpres = pl-co2_prop_spwag(i-1,j,1)
   398)           !dddp
   399)           co2_prop_spwag(i,j,5) = (co2_prop_spwag(i,j,3)-co2_prop_spwag(i-1,j,3))/dpres
   400)           
   401)           !dfgdp
   402)           co2_prop_spwag(i,j,7) = (co2_prop_spwag(i,j,6)-co2_prop_spwag(i-1,j,6))/dpres
   403)           
   404)           !dhdp
   405)           co2_prop_spwag(i,j,12) = (co2_prop_spwag(i,j,10)-co2_prop_spwag(i-1,j,10))/dpres
   406)           
   407)           !dvdp
   408)           co2_prop_spwag(i,j,15) = (co2_prop_spwag(i,j,13)-co2_prop_spwag(i-1,j,13))/dpres
   409)         endif
   410) #endif
   411)         tmp = co2_prop_spwag(i,j,3)
   412)         if (j==0) tmp2 = co2_prop_spwag(i,j,3)
   413)       ! print *, co2_prop_spwag(i,j,:)
   414)       enddo
   415)     enddo
   416)   endif
   417)   
   418)   if (len_trim(option%co2_database_filename) < 2) then
   419)     option%io_buffer = 'CO2 database filename not included in input deck.'
   420)     call printErrMsg(option)
   421)   endif
   422)   
   423)   if (myrank == 0) then
   424)     if (iitable == 1) then
   425)       co2_database_filename = 'co2data.dat'
   426)       print *,'Writing Table lookup file in working directory...'
   427)       if (myrank==0) print *,'--> open CO2 database file: ', &
   428)                              trim(co2_database_filename)
   429)       open(unit=IUNIT_TEMP,file=co2_database_filename,status='unknown',iostat=status)
   430)       if (status /= 0) then
   431)         print *, 'file: ', trim(co2_database_filename), ' not found.'
   432)         stop
   433)       endif
   434)       write(IUNIT_TEMP,'(''TITLE= "'',''co2data.dat'',''"'')')
   435)       write(IUNIT_TEMP,'(''VARIABLES= "'',a6,100(a3,a6))') &
   436)           'p',q,'T',q,'d',q,'dddT',q,'dddp',q,'fg',q,'dfgdp',q,'dfgdT',q, &
   437)           'u',q,'h',q,'dhdT',q,'dhdp',q,'vis',q,'dvdT',q,'dvdp','"'
   438)       write(IUNIT_TEMP,'(''ZONE T= "'',''",'','' I='',i4,'' , J='',i4)') ntab_t+1,ntab_p+1
   439)       do i = 0, ntab_p
   440)         tmp=tmp2
   441)         pl = p0_tab + dp_tab * real(i)
   442)         do j = 0, ntab_t
   443)           tl = t0_tab + dt_tab * real(j)
   444)           write(IUNIT_TEMP,'(1p15e14.6)') co2_prop_spwag(i,j,1:15)
   445)         enddo
   446)       enddo
   447)       close (IUNIT_TEMP)
   448)     endif
   449)   endif
   450)   
   451)   if (iitable == 2) then
   452)     if (myrank == 0) print *,'Reading Table ...'
   453)     if (myrank == 0) print *,'--> CO2 database file: ', &
   454)                              trim(option%co2_database_filename)
   455)     open(unit = IUNIT_TEMP,file=option%co2_database_filename,status='old',iostat=status)
   456)     if (status /= 0) then
   457)       print *, 'file: ', trim(option%co2_database_filename), ' not found.'
   458)       stop
   459)     endif
   460) !   open(unit=IUNIT_TEMP,file='co2data0.dat',status='old')
   461)     read(IUNIT_TEMP,*)
   462)     read(IUNIT_TEMP,*)
   463)     read(IUNIT_TEMP,*)
   464)     do i = 0, ntab_p
   465)       do j = 0, ntab_t
   466) #ifdef PC_BUG
   467)         read(IUNIT_TEMP,'(1p15e14.6)') temparray
   468)         co2_prop_spwag(i,j,1:15) = temparray(:)
   469) #else
   470)         read(IUNIT_TEMP,'(1p15e14.6)') co2_prop_spwag(i,j,1:15)
   471) #endif
   472)       enddo
   473)     enddo
   474)     close (IUNIT_TEMP)
   475)   endif
   476) 
   477) end subroutine initialize_span_wagner
   478) 
   479) ! ************************************************************************** !
   480) 
   481) subroutine co2_span_wagner(pl,tl,rho,dddt,dddp,fg,dfgdp,dfgdt, &
   482)       eng,ent,dhdt,dhdp,visc,dvdt,dvdp,iflag,itable)
   483)       
   484)   use co2_sw_rtsafe_module
   485) 
   486)   implicit none
   487)       
   488)       PetscReal :: pl,tl,rho,eng,ent,dhdt,dhdp,dddt,dddp,visc,dvdt,dvdp
   489)       PetscReal :: fiot,fiott,fird,firt,firdt,firdd,firtt,ftau,fdel,del,tau,fir
   490) !     PetscReal :: vpartial, rho_h2o, Cs, xco2, rho_wco2
   491)       PetscReal :: fg,dfgdp,fg1,dfgdt
   492)       PetscReal :: rho1,rho2
   493)       PetscInt :: iflag
   494)       
   495) !     PetscInt :: it
   496)       PetscInt, optional :: itable
   497)        
   498)       PetscInt :: iitable,i1,i2,j1,j2,i,isucc
   499)       PetscReal :: iindex,jindex,tmp,factor(1:4)
   500) 
   501) 
   502)       p=pl; t=tl; iitable = 0
   503)       if (present(itable)) iitable = itable
   504) 
   505) !     co2data0.dat - units: 
   506) !      1-P      : MPa
   507) !      2-T      : K
   508) !      3-d (den): kg/m^3
   509) !      4-dddT   : kg/m3/C
   510) !      5-dddp   : kg/m3/MPa
   511) !      6-fg     : -
   512) !      7-dfgdp  : 1/MPa
   513) !      8-dfgdT  : 1/C
   514) !      9-u(eng) : MJ/Kg
   515) !     10-h(ent) : MJ/Kg
   516) !     11-dhdT   : MJ/kg/C
   517) !     12-dhdp   : MJ/Kg/MPa
   518) !     13-vis    : Pa.s
   519) !     14-dvdT   : Pa.s/C
   520) !     15-dvdp   : Pa.s/Mpa
   521) 
   522) !     print *,'span_wag: ',p,t,tc,pc,rg,iitable
   523) 
   524) !     call co2_density (rho,it)
   525) !     call co2_density0 (p,t,rho)
   526) 
   527) !*************************** Table Lookup *********************
   528) 
   529)   tablelookup: if (iitable >= 1) then
   530) 
   531)   isucc=1
   532)   tmp = (p - p0_tab) / dp_tab; i1 = floor(tmp); i2 = ceiling(tmp); iindex=tmp 
   533)   tmp = (t - t0_tab) / dt_tab; j1 = floor(tmp); j2 = ceiling(tmp); jindex=tmp 
   534) 
   535)   if (iindex > ntab_p .or. iindex < 0.d0 .or. jindex < 0.d0 .or. jindex > ntab_t) then
   536)     print  *,' Out of Table Bounds (Span-Wagner): ', 'p [MPa] =',p, &
   537)     ' t [C] =',t-273.15,' i=',iindex,' j=',jindex
   538) !geh    isucc=0
   539)     iflag = -1
   540)     return
   541)   endif
   542) 
   543)   if (isucc > 0) then
   544) 
   545)     if (i1==i2 .and. j1==j2) then
   546)       factor(1) = 1.D0
   547)       factor(2:4) = 0.D0
   548)     elseif (i2==i1) then
   549)       factor(1) = -(jindex-j2)
   550)       factor(2) = 0.D0
   551)       factor(3) = (jindex-j1)
   552)       factor(4) = 0.D0
   553)     elseif (j2==j1) then
   554)       factor(1) = -(iindex-i2)
   555)       factor(2) = (iindex-i1)
   556)       factor(3) = 0.0
   557)       factor(4) = 0.D0
   558)     else
   559)       factor(1) =  (iindex-i2) * (jindex-j2)
   560)       factor(2) = -(iindex-i1) * (jindex-j2)
   561)       factor(3) = -(iindex-i2) * (jindex-j1)
   562)       factor(4) =  (iindex-i1) * (jindex-j1)
   563)     endif
   564) 
   565)     i=1
   566)     tmp = factor(1)*co2_prop_spwag(i1,j1,i) + factor(2)*co2_prop_spwag(i2,j1,i) &
   567)          + factor(3)*co2_prop_spwag(i1,j2,i) + factor(4)*co2_prop_spwag(i2,j2,i)
   568)     if (dabs(tmp-p) > 1.D-10 ) then
   569)       print *,' Error in interpolation: P ',tmp,p,iindex,factor;isucc=0
   570)     endif
   571)    !print *, 'Table: P ',iindex,jindex, factor,i  
   572)     i=i+1
   573)     tmp = factor(1)*co2_prop_spwag(i1,j1,i) + factor(2)*co2_prop_spwag(i2,j1,i) &
   574)          + factor(3)*co2_prop_spwag(i1,j2,i) + factor(4)*co2_prop_spwag(i2,j2,i)
   575)     if (dabs(tmp-t) > 1.D-10 ) then
   576)       print *,' Error in interpolation: T', tmp,t,jindex,factor; isucc=0
   577)     endif
   578)   endif
   579) 
   580) ! if (i <= 3) &
   581) ! print *, 'Table: T',i1,i2,j1,j2,iindex,jindex,factor,i,isucc,itable
   582) 
   583)   if (isucc == 1) then
   584)     i=i+1
   585)     rho = factor(1)*co2_prop_spwag(i1,j1,i) + factor(2)*co2_prop_spwag(i2,j1,i) &
   586)          + factor(3)*co2_prop_spwag(i1,j2,i) + factor(4)*co2_prop_spwag(i2,j2,i)
   587)     i=i+1 
   588)     dddt = factor(1)*co2_prop_spwag(i1,j1,i) + factor(2)*co2_prop_spwag(i2,j1,i) &
   589)          + factor(3)*co2_prop_spwag(i1,j2,i) + factor(4)*co2_prop_spwag(i2,j2,i)
   590) 
   591)     i=i+1
   592)     dddp= factor(1)*co2_prop_spwag(i1,j1,i) + factor(2)*co2_prop_spwag(i2,j1,i) &
   593)          + factor(3)*co2_prop_spwag(i1,j2,i) + factor(4)*co2_prop_spwag(i2,j2,i)
   594) 
   595)     i=i+1
   596)     fg= factor(1)* co2_prop_spwag(i1,j1,i) + factor(2)*co2_prop_spwag(i2,j1,i) &
   597)          + factor(3)* co2_prop_spwag(i1,j2,i) + factor(4)*co2_prop_spwag(i2,j2,i)
   598) 
   599)     i=i+1
   600)     dfgdp= factor(1)* co2_prop_spwag(i1,j1,i) + factor(2)* co2_prop_spwag(i2,j1,i) &
   601)          + factor(3)* co2_prop_spwag(i1,j2,i) + factor(4)* co2_prop_spwag(i2,j2,i)
   602) 
   603)     i=i+1
   604)     dfgdt= factor(1)* co2_prop_spwag(i1,j1,i) + factor(2)* co2_prop_spwag(i2,j1,i) &
   605)          + factor(3)* co2_prop_spwag(i1,j2,i) + factor(4)* co2_prop_spwag(i2,j2,i)
   606) 
   607)     i=i+1
   608)     eng= factor(1)* co2_prop_spwag(i1,j1,i) + factor(2)* co2_prop_spwag(i2,j1,i) &
   609)          + factor(3)* co2_prop_spwag(i1,j2,i) + factor(4)* co2_prop_spwag(i2,j2,i)
   610) 
   611)     i=i+1
   612)     ent= factor(1)* co2_prop_spwag(i1,j1,i) + factor(2)* co2_prop_spwag(i2,j1,i) &
   613)          + factor(3)* co2_prop_spwag(i1,j2,i) + factor(4)* co2_prop_spwag(i2,j2,i)
   614) 
   615)     i=i+1
   616)     dhdt= factor(1)* co2_prop_spwag(i1,j1,i) + factor(2)* co2_prop_spwag(i2,j1,i) &
   617)          + factor(3)* co2_prop_spwag(i1,j2,i) + factor(4)* co2_prop_spwag(i2,j2,i)
   618) 
   619)     i=i+1
   620)     dhdp= factor(1)* co2_prop_spwag(i1,j1,i) + factor(2)* co2_prop_spwag(i2,j1,i) &
   621)          + factor(3)* co2_prop_spwag(i1,j2,i) + factor(4)* co2_prop_spwag(i2,j2,i)
   622) 
   623)     i=i+1
   624)     visc= factor(1)* co2_prop_spwag(i1,j1,i) + factor(2)* co2_prop_spwag(i2,j1,i) &
   625)          + factor(3)* co2_prop_spwag(i1,j2,i) + factor(4)* co2_prop_spwag(i2,j2,i)
   626) 
   627)     i=i+1
   628)     dvdt= factor(1)* co2_prop_spwag(i1,j1,i) + factor(2)* co2_prop_spwag(i2,j1,i) &
   629)          + factor(3)* co2_prop_spwag(i1,j2,i) + factor(4)* co2_prop_spwag(i2,j2,i)
   630) 
   631)     i=i+1
   632)     dvdp= factor(1)* co2_prop_spwag(i1,j1,i) + factor(2)* co2_prop_spwag(i2,j1,i) &
   633)          + factor(3)* co2_prop_spwag(i1,j2,i) + factor(4)* co2_prop_spwag(i2,j2,i)
   634) ! print *,'sp:tab: ', p,t,rho,ent, factor
   635)     return 
   636)     endif     
   637)     
   638)     
   639)   endif tablelookup
   640) 
   641) 
   642) !*************************** Solving EOS *********************
   643)       call guess(rho1,rho2)
   644) !     print *,'spanwag-guess: ',p,t,rho1,rho2,iitable
   645) 
   646)       call bracket(co2den,rho1,rho2)
   647) 
   648)       rho = rtsafe(co2den,rho1,rho2,1.d-10)
   649) 
   650)       del = rho/denc
   651)       tau = tc/t
   652)       call phir(fir,del,tau)
   653)       call dphiodtau(fiot,del,tau)
   654)       call dphirddel(fird,del,tau)
   655)       call dphirdtau(firt,del,tau)
   656)       call dphirdddel(firdd,del,tau)
   657)       call dphirdtautau(firtt,del,tau)
   658)       call dphirddeldtau(firdt,del,tau)
   659)       call dphiodtautau(fiott,del,tau)
   660)       
   661)       ent = rg*t*(1.0d0+tau*(fiot+firt)+del*fird) + 506.78d0
   662)       eng = ent - p/rho*1.d3
   663) 
   664) !     dhdt = -rg*tau*tau*((-1.d0/(tau*tau))+fiott+firtt+ &
   665) !         ((del*firdt)/tau)-((del*fird)/(tau*tau)))
   666) 
   667)       dhdt = rg*((1.d0+del*(fird-tau*firdt))**2/(1.d0+del*(2.d0*fird+del*firdd))- &
   668)              tau*tau*(fiott+firtt))
   669) 
   670)       dhdp = (1000.d0/denc)*((tau*firdt)+(del*firdd)+fird)/ &
   671)           (1.d0+(del*del*firdd)+(2.d0*del*fird))
   672) 
   673) !     change units to MJ/kg from KJ/Kg
   674)       ent = ent*1.d-3
   675)       eng = eng*1.d-3
   676)       dhdt = dhdt*1.d-3
   677)       dhdp = dhdp*1.d-3
   678) 
   679)       ftau = (del*del*firdt)-(1000.d0*p/(denc*rg*tc))
   680)       fdel = (2.d0*del*fird)+(del*del*firdd)+1.d0
   681)       dddt = (denc*tc/(t*t))*(ftau/fdel)
   682)       dddp = 1000.d0/(rg*t*(1.d0+(del*del*firdd)+(2.d0*del*fird)))
   683)       
   684) !     print *,'co2_span: ','p=',p,'t=',t,'rho=',rho,'dddt=',dddt,'dddp=',dddp,'del=',del, &
   685) !     'ftau=',ftau,'fdel=',fdel,'firdt=',firdt,'firdd=',firdd,'fird=',fird
   686) 
   687) !     fugacity
   688)       fg1 = exp(fir+(del*fird)-log(1.d0+(del*fird)))
   689)       fg = fg1*p
   690)       dfgdp = (fird*(1.d0+(2.d0*del*fird)+(del*del*firdd))/ &
   691)           (1.d0+(del*fird)))*dddp/denc
   692)       dfgdp = fg1+(fg*dfgdp)
   693)       
   694)       dfgdt = (firt+(del*firt*fird)+(del*del*firdt*fird))/ &
   695)           (1.d0+(del*fird))
   696)       dfgdt = -fg*dfgdt*tc/(t*t)
   697)       
   698) !      write(*,*) "Molality?"
   699) !      read(*,*) mol
   700) !     call dissco2(p,t,mco2,fg1,mol)
   701) 
   702) !     density of aqueous CO2 solution
   703) !     vpartial = 37.51d0-(9.585d-2*(t-273.15d0))+ &
   704) !     (8.74d-4*((t-273.15d0)**2.d0))-(5.044d-7*((t-273.15d0)**3))
   705) 
   706) !      write(*,*) "Density of water?"
   707) !      read(*,*) rho_h2o
   708) 
   709) !     Cs = mco2*rho_h2o
   710) !     Xco2 = Cs/(Cs+(rho_h2o/18.d-3))
   711) !     rho_wco2 = rho_h2o + (44.d-3*Cs) - (vpartial*1d-6*rho_h2o*Cs)
   712) 
   713) !      write(*,*) mco2, rho_wco2, xco2
   714) !     Xco2 = .0292d0
   715) !     rho_wco2 = rho_h2o + (1.96d2*Xco2) + (1.54d4*Xco2*Xco2)
   716) !      write(*,*) rho_wco2
   717) 
   718)       call viscosity(pl,tl,rho,dddp,dddt,visc,dvdt,dvdp)
   719) 
   720) !     change units to Pa.s from microPa.s
   721)       visc = visc*1d-6
   722)       dvdt = dvdt*1d-6
   723)       dvdp = dvdp*1d-6
   724)       
   725)       if (allocated(co2_prop_spwag)) deallocate(co2_prop_spwag)
   726) 
   727)       
   728) end subroutine co2_span_wagner
   729) 
   730) ! ************************************************************************** !
   731) 
   732) subroutine guess(lguess,uguess)
   733) 
   734)       implicit none
   735) 
   736)       PetscReal :: ts,dum1,dum2,uguess,lguess,rhov,rhol
   737)       PetscReal :: a1,a2,a3,a4,a5,t1,t2,t3,t4,t5,tr
   738)       PetscInt, parameter :: TWELVE_INTEGER = 12
   739) 
   740)       if (p.le.pc .and. t.le.tc) then
   741)          call vappr (ts,p,dum1,dum2,TWELVE_INTEGER)
   742) 
   743) !vapor density
   744)             a1 = -1.7074879d0
   745)             a2 = -0.8227467d0
   746)             a3 = -4.6008549d0
   747)             a4 = -10.111178d0
   748)             a5 = -29.742252d0
   749)             t1 = 0.34d0
   750)             t2 = 0.5d0
   751)             t3 = 1.d0
   752)             t4 = 7.d0/3.d0
   753)             t5 = 14.d0/3.d0
   754) !           tr = 1.d0-(ts+273.15d0)/tc
   755)             tr = 1.d0-ts/tc
   756)             rhov =a1*(tr**t1)+a2*(tr**t2)+a3*(tr**t3)+ &
   757)                     a4*(tr**t4)+a5*(tr**t5)
   758)             rhov =denc*exp(rhov)
   759) 
   760) !liquid density
   761)             a1 = 1.9245108d0
   762)             a2 = -0.62385555d0
   763)             a3 = -0.32731127d0
   764)             a4 = -0.39245142d0
   765)             t1 = 0.34d0
   766)             t2 = 0.5d0
   767)             t3 = 5.d0/3.d0
   768)             t4 = 11.d0/6.d0
   769)             rhol =a1*(tr**t1)+a2*(tr**t2)+a3*(tr**t3)+ &
   770)                     a4*(tr**t4)
   771)             rhol =denc*exp(rhol)
   772) 
   773) !           print *,'guess: ',p,t,ts,tc,denc,tr,rhov,rhol
   774)  
   775)         if (t.le.ts) then ! sub-critical liquid region
   776)           if (p.lt.6.d0) then
   777)             lguess = rhol+2.d0*(ts-t)
   778)           else if (p.lt.7.d0) then
   779)             lguess = rhol+5.d0*(ts-t)
   780)           else if (p.lt.8.d0) then
   781)             lguess = rhol+9.d0*(ts-t)
   782)           else
   783)             lguess = rhol+10.d0*(ts-t)
   784)           endif
   785) !         iflag = 1
   786)           uguess = 2000.d0
   787)         else             ! vapor region
   788)           uguess = 1.d-5
   789)           lguess = 1.2*rhov
   790) !         iflag = 2
   791)         endif
   792)       else if (t.le.tc .and. p.gt.pc) then
   793)         if (p.le.8.d0) then
   794)           uguess = 2000.d0
   795)           if (t.lt.275.d0) then
   796)             lguess = 950.d0
   797)           else if (t.le.290.d0) then
   798)             lguess = 800.d0
   799)           else
   800)             lguess = 700.d0
   801)           endif
   802) !         iflag = 3
   803)         else if (p.le.9.5d0) then
   804)           uguess = 2000.d0
   805)           if (t.lt.275.d0) then
   806)             lguess = 950.d0
   807)           else if (t.le.290.d0) then
   808)             lguess = 800.d0
   809)           else
   810)             lguess = 750.d0
   811)           endif
   812) !         iflag = 4
   813)         else if (p.le.100.d0) then
   814)           uguess = 2000.d0
   815)           if (t.lt.275.d0) then
   816)             lguess = 950.d0
   817)           else if (t.le.290.d0) then
   818)             lguess = 800.d0
   819)           else
   820)             lguess = 700.d0
   821)           endif
   822) !         iflag = 5
   823)         else
   824)           uguess = 2000.d0
   825)           lguess = 900.d0
   826) !         iflag = 6
   827)         endif
   828)       else if (t.gt.tc .and. p.le.pc) then
   829)         uguess = 2000.d0
   830)         lguess = 1.d-5
   831) !       iflag = 7
   832)       else if (p.gt.pc) then ! supercritical
   833)         uguess = 2000.d0
   834)         lguess = 25.d0
   835) !       iflag = 8
   836)       endif
   837) 
   838) end subroutine guess
   839) 
   840) ! ************************************************************************** !
   841) 
   842) subroutine co2den(den,f,df)
   843)      
   844)     IMPLICIT NONE
   845)       PetscReal :: den,tau1,del1
   846)       PetscReal :: f1,df1,f,df
   847) 
   848)       tau1 = tc/t
   849)       del1 = den/denc
   850) 
   851)       call dphirddel(f1,del1,tau1)
   852)       call dphirdddel(df1,del1,tau1)
   853)       f = del1*(1.d0+del1*f1)-(p/(denc*0.001d0*rg*t))
   854)       df = 1.d0+(2.d0*del1*f1)+(del1*del1*df1)
   855) 
   856)       return
   857) end subroutine co2den
   858) 
   859) ! ************************************************************************** !
   860) 
   861)       PetscReal function psi(i,del2,tau2)
   862)       implicit none
   863) !     PetscReal :: psi
   864) !     PetscReal aco2(4),bco2(4),capa(5),capb(5),capc(5),capd(5)
   865)       PetscReal :: del2,tau2
   866)       PetscInt :: i
   867) !     common/params3/aco2,bco2,capa,capb,capc,capd
   868) 
   869)       psi=-(capc(i)*((del2-1.d0)**2.d0))-(capd(i)*((tau2-1.d0)**2.d0))
   870)       psi=exp(psi)
   871) 
   872)       end function psi
   873) 
   874) ! ************************************************************************** !
   875) 
   876) subroutine dphiodtau(dr,del2,tau2)
   877)       implicit none
   878)       PetscInt :: i
   879)       PetscReal :: del2,tau2,dr
   880)       PetscReal :: ideal_helm, derti_helm
   881) 
   882)       ideal_helm = log(del2)+a(1)+a(2)*tau2+a(3)*log(tau2)
   883)       do i = 4, 8
   884)         ideal_helm = ideal_helm+(a(i)*log(1.d0-exp(-phic(i)*tau2)))
   885)       enddo
   886) 
   887) !     Table 34 derivative of ideal helmholtz wrt tau
   888) 
   889)       derti_helm = a(2) + (a(3)/tau2)
   890)       do i = 4, 8
   891)        derti_helm = derti_helm+a(i)*phic(i)*(1.d0/(1.d0-exp(-phic(i)*tau2))-1.d0)
   892)       enddo
   893)       dr = derti_helm
   894) 
   895) end subroutine dphiodtau
   896) 
   897) ! ************************************************************************** !
   898) 
   899) subroutine dphiodtautau(dr,del2,tau2)
   900) 
   901)       implicit none
   902)       
   903) !     Span & Wagner (1996) Table 32
   904) 
   905)       PetscInt :: i
   906)       PetscReal :: del2,tau2,dr,dihelm_dtautau
   907) 
   908) !     Table 34 double derivative of ideal helmholtz wrt tau
   909) 
   910)       dihelm_dtautau = -a(3)/(tau2*tau2)
   911)       do i = 4, 8
   912)         dihelm_dtautau = dihelm_dtautau-a(i)*(phic(i)**2)* &
   913)         (1.d0-exp(-phic(i)*tau2))**(-2)*exp(-phic(i)*tau2)
   914)       enddo
   915)       dr = dihelm_dtautau
   916) 
   917) end subroutine dphiodtautau
   918) 
   919) ! ************************************************************************** !
   920) 
   921) subroutine phir(r_helm,del2,tau2)
   922) 
   923) !     residual helmholtz energy: Span & Wagner (1996), p. 1544, eq. (6.5)
   924) !     Table 32
   925) 
   926)       implicit none
   927)       
   928)       PetscReal :: del2,tau2,r_helm,psi1
   929) !     PetscReal :: x1,x2,x3,x4,x5,x6,x7,x8,x10
   930) !     PetscReal :: e1,e2,e3,e4,e5,e6
   931) !     PetscReal :: tsqr,t1,t2,t3,t6,t7,t8,t12,t14,t16,t22,t24,t28
   932)       PetscReal :: capdel1
   933)       PetscInt :: i
   934) 
   935) !     equation 5.3 for residual part of Helmholtz function
   936) 
   937)       r_helm = 0.0d0
   938) 
   939)       do i = 1, 7
   940)         r_helm = r_helm+(n(i)*(del2**d(i))*(tau2**ti(i)))
   941)       enddo
   942) 
   943)       do i = 8, 34
   944)         r_helm = r_helm+(n(i)*(del2**d(i))*(tau2**ti(i))*exp(-del2**c(i))) 
   945)       enddo
   946) 
   947)       do i = 35, 39
   948)         r_helm = r_helm+(n(i)*(del2**d(i))*(tau2**ti(i))* &
   949)         exp((-alpha(i-34)*((del2-epsilon(i-34))**2))-(beta(i-34)* &
   950)         ((tau2-gamma(i-34))**2))))
   951)       enddo
   952) 
   953)       do i = 40, 42
   954)         psi1=psi(i-39,del2,tau2)
   955)         capdel1=capdel(i-39,del2,tau2)
   956)         r_helm=r_helm+(n(i)*(capdel1**bco2(i-39))*del2*psi1)
   957)       enddo
   958) 
   959)       return
   960) end subroutine phir
   961) 
   962) ! ************************************************************************** !
   963) 
   964) subroutine dphirddel(dr,del2,tau2)
   965) 
   966) !     Span & Wagner (1996) Table 32
   967) 
   968)       implicit none
   969)       PetscInt :: i
   970)       PetscReal :: del2,tau2,dr,derdr_helm
   971)       PetscReal :: psi1,capdel1,dsidd,ddelbdd
   972) !     PetscReal :: d1,d2,d3,d4,d5,d6,d7,d9
   973) !     PetscReal :: e1,e2,e3,e4,e5,e6
   974) !     PetscReal :: tsqr,t1,t2,t3,t6,t7,t8,t12,t14,t16,t22,t24,t28
   975) 
   976) !     table 32, derivative of residual helmholtz wrt delta
   977) 
   978)       derdr_helm = 0.0d0
   979)       do i = 1, 7
   980)         derdr_helm=derdr_helm+(n(i)*d(i)*(del2**(d(i)-1.d0))*(tau2**ti(i)))
   981)       enddo
   982)       do i = 8, 34
   983)         derdr_helm=derdr_helm+(n(i)*(del2**(d(i)-1.d0))*(tau2**ti(i)) &
   984)         *exp(-del2**c(i))*(d(i)-(c(i)*(del2**c(i))))) 
   985)       enddo
   986)       do i = 35, 39
   987)         derdr_helm = derdr_helm+(n(i)*(del2**d(i))*(tau2**ti(i))* &
   988)         exp((-alpha(i-34)*((del2-epsilon(i-34))**2.0d0))-(beta(i-34)* &
   989)         ((tau2-gamma(i-34))**2.0d0)))* &
   990)         ((d(i)/del2)-(2.0d0*alpha(i-34)*(del2-epsilon(i-34)))))
   991)       enddo
   992)       do i = 40, 42
   993)         psi1    = psi(i-39,del2,tau2)
   994)         capdel1 = capdel(i-39,del2,tau2)
   995)         dsidd   = dpsiddel(i-39,del2,tau2)
   996)         ddelbdd = ddelbiddel(i-39,del2,tau2)
   997)         derdr_helm = derdr_helm+(n(i)*(((capdel1**bco2(i-39))*(psi1+ &
   998)         (del2*dsidd)))+(del2*psi1*ddelbdd)))
   999)       enddo
  1000)       dr = derdr_helm
  1001)       
  1002) end subroutine dphirddel
  1003) 
  1004) ! ************************************************************************** !
  1005) 
  1006) subroutine dphirdddel(dpdd,del2,tau2)
  1007)       implicit none
  1008)       PetscInt :: i
  1009)       PetscReal :: del2,tau2,derdr_helm
  1010)       PetscReal :: dpdd,psi1,capdel1
  1011)       PetscReal :: dsidd,d2sidd,ddelbdd
  1012)       PetscReal :: d2delbdd
  1013) !     PetscReal :: d1,d2,d3,d4,d5,d6,d7,d8,d9
  1014) !     PetscReal :: e1,e2,e3,e4,e5,e6
  1015) !     PetscReal :: tsqr,t1,t2,t3,t6,t7,t8,t12,t14,t16,t22,t24,t28
  1016) 
  1017) !     table 32, derivative of residual helmholtz wrt delta
  1018) 
  1019)       derdr_helm = 0.0d0
  1020)       do i = 1, 7
  1021)         derdr_helm=derdr_helm+(n(i)*d(i)*(d(i)-1.d0)* &
  1022)         (del2**(d(i)-2.d0))*(tau2**ti(i)))
  1023)       enddo
  1024)       
  1025)       do i = 8, 34
  1026)         derdr_helm=derdr_helm+(n(i)*(del2**(d(i)-2.d0))*(tau2**ti(i))* &
  1027)         exp(-del2**c(i))*(((d(i)-(c(i)*(del2**c(i))))* &
  1028)         (d(i)-1.d0-(c(i)*(del2**c(i)))))-((c(i)**2.d0)* &
  1029)         (del2**c(i)))))
  1030)       enddo
  1031) 
  1032)       do i = 35, 39
  1033)         derdr_helm = derdr_helm+(n(i)*(tau2**ti(i))* &
  1034)         exp((-alpha(i-34)*((del2-epsilon(i-34))**2.0d0))-(beta(i-34)* &
  1035)         ((tau2-gamma(i-34))**2.0d0)))* &
  1036)         ((-2.d0*alpha(i-34)*(del2**d(i)))+(4.d0*(alpha(i-34)**2.d0)* &
  1037)         (del2**d(i))*((del2-epsilon(i-34))**2.d0))-(4.d0*d(i)* &
  1038)         alpha(i-34)*(del2**(d(i)-1.d0))*(del2-epsilon(i-34)))+ &
  1039)         (d(i)*(d(i)-1.d0)*(del2**(d(i)-2.d0)))))
  1040)       enddo
  1041) 
  1042)       do i = 40, 42
  1043)         psi1     = psi(i-39,del2,tau2)
  1044)         capdel1  = capdel(i-39,del2,tau2)
  1045)         dsidd    = dpsiddel(i-39,del2,tau2)
  1046)         d2sidd   = d2psiddel2(i-39,del2,tau2)
  1047)         ddelbdd  = ddelbiddel(i-39,del2,tau2)
  1048)         d2delbdd = d2delbiddel2(i-39,del2,tau2)
  1049)         derdr_helm = derdr_helm+(n(i)*(((capdel1**bco2(i-39))*((2.d0* &
  1050)         dsidd)+(del2*d2sidd)))+(2.d0*ddelbdd*(psi1+(del2*dsidd))) &
  1051)         +(d2delbdd*del2*psi1)))
  1052)         
  1053) !       print *,'dphirdddel: ',i,psi1,capdel1,dsidd,d2sidd,ddelbdd,derdr_helm
  1054)       enddo
  1055)       
  1056)       dpdd = derdr_helm
  1057)       
  1058) end subroutine dphirdddel
  1059) 
  1060) ! ************************************************************************** !
  1061) 
  1062) subroutine dphirdtau(dpdtau,del2,tau2)
  1063) 
  1064)       implicit none
  1065)       PetscInt :: i
  1066)       PetscReal :: del2,tau2,derdr_helm
  1067)       PetscReal :: dpdtau
  1068)       PetscReal :: psi1,dbidt,capdel1,dsidt
  1069) 
  1070) !     table 32, derivative of residual helmholtz wrt delta
  1071) 
  1072)       derdr_helm = 0.0d0
  1073)       do i = 1, 7
  1074)         derdr_helm=derdr_helm+(n(i)*ti(i)* &
  1075)         (tau2**(ti(i)-1))*(del2**d(i)))
  1076)       enddo
  1077)       
  1078)       do i = 8, 34
  1079)         derdr_helm=derdr_helm+(n(i)*ti(i)*(del2**d(i))* &
  1080)         (tau2**(ti(i)-1))*exp(-del2**c(i)))
  1081)       enddo
  1082)       
  1083)       do i = 35, 39
  1084)         derdr_helm = derdr_helm+(n(i)*(del2**d(i))*(tau2**ti(i))* &
  1085)         exp((-alpha(i-34)*((del2-epsilon(i-34))**2.0d0))-(beta(i-34)* &
  1086)         ((tau2-gamma(i-34))**2.0d0)))* &
  1087)         ((ti(i)/tau2)-(2.0d0*beta(i-34)*(tau2-gamma(i-34)))))
  1088)       enddo
  1089)       do i = 40, 42
  1090)         psi1=psi(i-39,del2,tau2)
  1091)         dbidt=ddelbidtau(i-39,del2,tau2)
  1092)         capdel1=capdel(i-39,del2,tau2)
  1093)         dsidt=dpsidtau(i-39,del2,tau2)
  1094)         derdr_helm=derdr_helm+(n(i)*del2*((dbidt*psi1)+((capdel1** &
  1095)         bco2(i-39))*dsidt)))
  1096)       enddo
  1097)       dpdtau = derdr_helm
  1098)       
  1099) end subroutine dphirdtau
  1100) 
  1101) ! ************************************************************************** !
  1102) 
  1103) subroutine dphirdtautau(dpdtt,del2,tau2)
  1104) 
  1105)       implicit none
  1106)       PetscInt :: i
  1107)       PetscReal :: del2,tau2,derdr_helm
  1108)       PetscReal :: dpdtt
  1109)       PetscReal :: psi1
  1110)       PetscReal :: d2bidtt,capdel1,dbidt,dsidt,d2sidtt
  1111) 
  1112) !     table 36, derivative of residual helmholtz wrt delta
  1113) 
  1114)       derdr_helm = 0.0d0
  1115)       do i = 1, 7
  1116)         derdr_helm=derdr_helm+(n(i)*ti(i)*(ti(i)-1.d0)* &
  1117)         (tau2**(ti(i)-2.d0))*(del2**d(i)))
  1118)       enddo
  1119)       
  1120)       do i = 8, 34
  1121)         derdr_helm=derdr_helm+(n(i)*ti(i)*(ti(i)-1.d0)*(del2**d(i))* &
  1122)         (tau2**(ti(i)-2))*exp(-del2**c(i)))
  1123)       enddo
  1124)       
  1125)       do i = 35, 39
  1126)         derdr_helm = derdr_helm+(n(i)*(del2**d(i))*(tau2**ti(i))* &
  1127)         exp((-alpha(i-34)*((del2-epsilon(i-34))**2.0d0))-(beta(i-34)* &
  1128)         ((tau2-gamma(i-34))**2.0d0)))*( &
  1129)         (((ti(i)/tau2)-(2.0d0*beta(i-34)*(tau2-gamma(i-34))))**2.d0) &
  1130)         -(ti(i)/(tau2*tau2))-(2.d0*beta(i-34))))
  1131)       enddo
  1132)       
  1133)       do i = 40, 42
  1134)         psi1    = psi(i-39,del2,tau2)
  1135)         capdel1 = capdel(i-39,del2,tau2)
  1136)         d2bidtt = d2delbidtau2(i-39,del2,tau2)
  1137)         dbidt   = ddelbidtau(i-39,del2,tau2)
  1138)         dsidt   = dpsidtau(i-39,del2,tau2)
  1139)         d2sidtt = d2psidtau2(i-39,del2,tau2)
  1140)         derdr_helm = derdr_helm+(n(i)*del2*((d2bidtt*psi1)+ &
  1141)         (2.d0*dbidt*dsidt)+((capdel1**bco2(i-39))*d2sidtt)))
  1142)       enddo
  1143)       
  1144)       dpdtt = derdr_helm
  1145)       
  1146) end subroutine dphirdtautau
  1147) 
  1148) ! ************************************************************************** !
  1149) 
  1150) subroutine dphirddeldtau(dpddt,del2,tau2)
  1151) 
  1152)       implicit none
  1153)       PetscInt :: i
  1154)       PetscReal :: del2,tau2,derdr_helm
  1155)       PetscReal :: dpddt
  1156)       PetscReal :: psi1,capdel1,dsidt,dsiddt,dbidd,dbidt,dsidd,d2biddt
  1157) 
  1158) !     table 36, derivative of residual helmholtz wrt delta
  1159) 
  1160)       derdr_helm = 0.0d0
  1161)       do i = 1, 7
  1162)         derdr_helm=derdr_helm+(n(i)*ti(i)*d(i)* &
  1163)         (tau2**(ti(i)-1))*(del2**(d(i)-1)))
  1164)       enddo
  1165)       
  1166)       do i = 8, 34
  1167)         derdr_helm=derdr_helm+(n(i)*ti(i)*(del2**(d(i)-1))* &
  1168)         (tau2**(ti(i)-1))*exp(-del2**c(i))* &
  1169)         (d(i)-(c(i)*(del2**c(i)))))
  1170)       enddo
  1171)       
  1172)       do i = 35, 39
  1173)         derdr_helm = derdr_helm+(n(i)*(del2**d(i))*(tau2**ti(i))* &
  1174)         exp((-alpha(i-34)*((del2-epsilon(i-34))**2.0d0))-(beta(i-34)* &
  1175)         ((tau2-gamma(i-34))**2.0d0)))* &
  1176)         ((ti(i)/tau2)-(2.0d0*beta(i-34)*(tau2-gamma(i-34))))* &
  1177)         ((d(i)/del2)-(2.d0*alpha(i-34)*(del2-epsilon(i-34)))))
  1178)       enddo
  1179)       
  1180)       do i = 40, 42
  1181)         psi1    = psi(i-39,del2,tau2)
  1182)         capdel1 = capdel(i-39,del2,tau2)
  1183)         dsidt   = dpsidtau(i-39,del2,tau2)
  1184)         dsiddt  = d2psiddeltau(i-39,del2,tau2)
  1185)         dbidd   = ddelbiddel(i-39,del2,tau2)
  1186)         dbidt   = ddelbidtau(i-39,del2,tau2)
  1187)         dsidd   = dpsiddel(i-39,del2,tau2)
  1188)         d2biddt = d2delbiddeltau(i-39,del2,tau2)
  1189)         derdr_helm = derdr_helm+(n(i)*(((capdel1**bco2(i-39))*(dsidt+ &
  1190)         (del2*dsiddt)))+(del2*dbidd*dsidt)+(dbidt* &
  1191)         (psi1+(del2*dsidd)))+(d2biddt*del2*psi1)))
  1192)       enddo
  1193)       
  1194)       dpddt = derdr_helm
  1195)       
  1196) end subroutine dphirddeldtau
  1197) 
  1198) ! ************************************************************************** !
  1199) 
  1200) function dpsiddel(i,del2,tau2)
  1201)       implicit none
  1202)       PetscReal :: dpsiddel
  1203)       PetscReal :: del2,tau2,psi1
  1204)       PetscInt :: i
  1205) 
  1206)       psi1=psi(i,del2,tau2)
  1207)       dpsiddel = -2.d0*capc(i)*(del2-1.d0)*psi1
  1208)       
  1209) !     print *,'dpsiddel: ',i,dpsiddel,psi1,del2,tau2,capc(i)
  1210) 
  1211)       end function dpsiddel
  1212) 
  1213) ! ************************************************************************** !
  1214) 
  1215)       function d2psiddel2(i,del2,tau2)
  1216)       implicit none
  1217)       PetscReal :: d2psiddel2
  1218)       PetscReal :: del2,tau2,psi1
  1219)       PetscInt :: i
  1220) 
  1221)       psi1=psi(i,del2,tau2)
  1222)       d2psiddel2 = (2.d0*capc(i)*(del2-1.d0)**2-1.d0)*2.d0*capc(i)*psi1
  1223) 
  1224)       end function d2psiddel2
  1225) 
  1226) ! ************************************************************************** !
  1227) 
  1228)       function dpsidtau(i,del2,tau2)
  1229)       implicit none
  1230)       PetscReal :: dpsidtau
  1231)       PetscReal :: del2,tau2,psi1
  1232)       PetscInt :: i
  1233) 
  1234)       psi1=psi(i,del2,tau2)
  1235)       dpsidtau = -2.d0*capd(i)*(tau2-1.d0)*psi1
  1236) 
  1237)       end function dpsidtau
  1238) 
  1239) ! ************************************************************************** !
  1240) 
  1241)       function d2psidtau2(i,del2,tau2)
  1242)       implicit none
  1243)       PetscReal :: d2psidtau2
  1244)       PetscReal :: del2,tau2,psi1
  1245)       PetscInt :: i
  1246) 
  1247)       psi1=psi(i,del2,tau2)
  1248)       d2psidtau2 = ((2.d0*capd(i)*(tau2-1.d0)**2)-1.d0)*2.d0* &
  1249)       capd(i)*psi1
  1250) 
  1251)       end function d2psidtau2
  1252) 
  1253) ! ************************************************************************** !
  1254) 
  1255)       function d2psiddeltau(i,del2,tau2)
  1256)       implicit none
  1257)       PetscReal :: d2psiddeltau
  1258)       PetscReal :: del2,tau2,psi1
  1259)       PetscInt :: i
  1260) 
  1261)       psi1=psi(i,del2,tau2)
  1262)       d2psiddeltau = 4.d0*capc(i)*capd(i)*(del2-1.d0)*(tau2-1.d0)*psi1
  1263) 
  1264)       end function d2psiddeltau
  1265) 
  1266) ! ************************************************************************** !
  1267) 
  1268)       function theta(i,del2,tau2)
  1269)       implicit none
  1270)       PetscReal :: theta
  1271)       PetscReal :: del2,tau2
  1272)       PetscInt :: i
  1273) 
  1274)       theta=1.d0-tau2+capa(i)*((del2-1.d0)**2)**(1.d0/(2.d0*beta(i+5)))
  1275) 
  1276)       end function theta
  1277) 
  1278) ! ************************************************************************** !
  1279) 
  1280)       function capdel(i,del2,tau2)
  1281)       implicit none
  1282)       PetscReal :: capdel
  1283)       PetscReal :: theta1,del2,tau2
  1284)       PetscInt :: i
  1285) 
  1286)       theta1=theta(i,del2,tau2)
  1287)       capdel=theta1*theta1+capb(i)*(((del2-1.d0)**2)**aco2(i))
  1288)       end function capdel
  1289) 
  1290) ! ************************************************************************** !
  1291) 
  1292)       function dcapdelddel(i,del2,tau2)
  1293)       implicit none
  1294)       PetscReal :: dcapdelddel
  1295)       PetscReal :: del2,tau2,theta1
  1296)       PetscInt :: i
  1297) 
  1298)       theta1=theta(i,del2,tau2)
  1299) !     dcapdelddel=(del2-1.d0)*((capa(i)*theta1*(2.d0/beta(i+5)) &
  1300)       dcapdelddel=((capa(i)*theta1*(2.d0/beta(i+5)) & ! remove factor del2-1
  1301)       *((del2-1.d0)**2 &
  1302)       )**((1.d0/(2.d0*beta(i+5)))-1.d0))+(2.d0*capb(i)*aco2(i)* &
  1303)       (((del2-1.d0)**2)**(aco2(i)-1.d0))))
  1304) 
  1305)       end function dcapdelddel
  1306) 
  1307) ! ************************************************************************** !
  1308) 
  1309)       function d2capdelddel2(i,del2,tau2) ! d^2 Delta/d delta^2
  1310)       implicit none
  1311)       PetscReal :: d2capdelddel2
  1312)       PetscReal :: tmp1,del2,tau2,theta1
  1313)       PetscReal :: ddd
  1314)       PetscInt :: i
  1315) 
  1316)       theta1=theta(i,del2,tau2)
  1317)       tmp1=4.d0*capb(i)*aco2(i)*(aco2(i)-1.d0)*(((del2-1.d0) &
  1318)       **2)**(aco2(i)-2.d0))
  1319)  
  1320)       tmp1=tmp1+(2.d0*capa(i)*capa(i)*((1.d0/beta(i))**2)* &
  1321)       ((((del2-1.d0)**2)**((1.d0/(2.d0*beta(i)))-1.d0))**2))
  1322) 
  1323)       tmp1=tmp1+capa(i)*theta1*(4.d0/beta(i))*((1.d0/(2.d0*beta(i))) &
  1324)       -1.d0)*(((del2-1.d0)**2)**((1.d0/(2.d0*beta(i)))-2.d0))
  1325) 
  1326)       tmp1=tmp1*((del2-1.d0)**2.d0)
  1327) 
  1328)       ddd=dcapdelddel(i,del2,tau2)
  1329) !     d2capdelddel2=tmp1+((1.d0/(del2-1.d0))*ddd)
  1330)       d2capdelddel2=tmp1+ddd
  1331) 
  1332)       end function d2capdelddel2
  1333) 
  1334) ! ************************************************************************** !
  1335) 
  1336)       function ddelbiddel(i,del2,tau2)
  1337)       implicit none
  1338)       PetscReal :: ddelbiddel
  1339)       PetscInt :: i
  1340)       PetscReal :: del2,tau2,capdel1,ddd
  1341)       
  1342)       capdel1=capdel(i,del2,tau2)
  1343)       ddd=dcapdelddel(i,del2,tau2)*(del2-1.d0)
  1344)       ddelbiddel=bco2(i)*(capdel1**(bco2(i)-1.d0))*ddd
  1345) 
  1346)       end function ddelbiddel
  1347) 
  1348) ! ************************************************************************** !
  1349) 
  1350)       function d2delbiddel2(i,del2,tau2)
  1351)       implicit none
  1352)       PetscReal :: d2delbiddel2
  1353)       PetscInt :: i
  1354)       PetscReal :: del2,tau2,ddd1,ddd2,capdel1
  1355) 
  1356)       ddd1=dcapdelddel(i,del2,tau2)*(del2-1.d0)
  1357)       ddd2=d2capdelddel2(i,del2,tau2)
  1358)       capdel1=capdel(i,del2,tau2)
  1359)       d2delbiddel2=bco2(i)*ddd2*capdel1**(bco2(i)-1.d0) &
  1360)       +(bco2(i)-1.d0)*(capdel1**(bco2(i)-2.d0))*ddd1**2
  1361) 
  1362)       end function d2delbiddel2
  1363) 
  1364) ! ************************************************************************** !
  1365) 
  1366)       function ddelbidtau(i,del2,tau2)
  1367)       implicit none
  1368)       PetscReal :: ddelbidtau
  1369)       PetscInt :: i
  1370)       PetscReal :: del2,tau2,theta1,capdel1
  1371) 
  1372)       theta1=theta(i,del2,tau2)
  1373)       capdel1=capdel(i,del2,tau2)
  1374) 
  1375)       ddelbidtau=-2.d0*theta1*bco2(i)*capdel1**(bco2(i)-1.d0)
  1376) 
  1377)       end function ddelbidtau
  1378) 
  1379) ! ************************************************************************** !
  1380) 
  1381)       function d2delbidtau2(i,del2,tau2)
  1382)       implicit none
  1383)       PetscReal :: d2delbidtau2
  1384)       PetscReal :: del2,tau2
  1385)       PetscInt :: i
  1386)       PetscReal :: capdel1,theta1
  1387) 
  1388)       capdel1=capdel(i,del2,tau2)
  1389)       theta1=theta(i,del2,tau2)
  1390)       d2delbidtau2=(2.d0*bco2(i)*(capdel1**(bco2(i)-1.d0)))+(4.d0* &
  1391)       (theta1**2)*bco2(i)*(bco2(i)-1.d0)*(capdel1**(bco2(i)-2.d0)))
  1392) 
  1393)       end function d2delbidtau2
  1394) 
  1395) ! ************************************************************************** !
  1396) 
  1397)       function d2delbiddeltau(i,del2,tau2)
  1398)       implicit none
  1399)       PetscReal :: d2delbiddeltau
  1400)       PetscReal :: tmp3,del2,tau2
  1401)       PetscInt :: i
  1402)       PetscReal :: capdel1,ddd,theta1
  1403) 
  1404)       capdel1=capdel(i,del2,tau2)
  1405)       theta1=theta(i,del2,tau2)
  1406)       ddd=dcapdelddel(i,del2,tau2)*(del2-1.d0)
  1407)       tmp3=-capa(i)*bco2(i)*(2.d0/beta(i))*(capdel1**(bco2(i)-1.d0))* &
  1408)       (del2-1.d0)*(((del2-1.d0)**2)**((1.d0/(2.d0*beta(i)))-1.d0))
  1409) 
  1410)       tmp3=tmp3-(2.d0*theta1*bco2(i)*(bco2(i)-1.d0)*(capdel1** &
  1411)       (bco2(i)-2.d0))*ddd)
  1412)     
  1413)       d2delbiddeltau=tmp3
  1414) 
  1415) end function d2delbiddeltau
  1416) 
  1417) ! ************************************************************************** !
  1418) 
  1419) subroutine vappr(tm,ps,dertp,derpt,ifl1)
  1420) 
  1421)       implicit none
  1422)       
  1423) !     co2 vapor pressure curve (Span & Wagner, 1996, p. 1524, eq. (3.13)
  1424) 
  1425)       PetscInt :: j,maxit
  1426)       PetscReal :: tm,ps,nu,dertp,derpt,ps1,uguess,lguess,xacc
  1427)       PetscReal :: a1,a2,a3,a4
  1428)       PetscInt :: ifl1
  1429)       PetscReal :: dnu,ps2,f,df
  1430) 
  1431)       parameter(maxit=1000)
  1432)       parameter(xacc=1.d-7)
  1433)       a1 = -7.0602087d0
  1434)       a2 =  1.9391218d0
  1435)       a3 = -1.6463597d0
  1436)       a4 = -3.2995634
  1437)       if (ifl1 .eq. 11) then
  1438)         nu = 1.d0-(tm/tc)
  1439) 
  1440)         ps1 = a1*nu+a2*nu**1.5d0+a3*nu**2+a4*nu**4.d0
  1441)         ps2 = a1+1.5d0*a2*nu**0.5d0+2.0d0*a3*nu+4.d0*a4*nu**3
  1442) 
  1443)         ps = ps1*tc/tm
  1444) 
  1445)         ps = exp(ps)*pc
  1446) 
  1447)         dertp = (ps/tm)*(-ps2-((ps1*tc)/tm))
  1448)       else
  1449)         lguess = 50.
  1450)         uguess = 400.
  1451)         tm = (uguess-lguess)/2.d0
  1452)         nu = 1.d0-(tm/tc)            
  1453)         do j = 1, maxit
  1454)           ps1 = a1*nu+(a2*(nu**1.5d0))+(a3*(nu**2))+ &
  1455)                 (a4*(nu**4))
  1456)           ps2 = a1+(1.5d0*a2*(nu**0.5d0))+(2.0d0*a3*nu)+ &
  1457)                 (4.d0*a4*(nu**3))
  1458)           f = log(ps/pc)-(ps1/(1.d0-nu))
  1459)           df=-(ps1+(ps2*(1.d0-nu)))/((1.d0-nu)**2)
  1460)           dnu=f/df
  1461)           nu=nu-dnu
  1462)           if (dabs(dnu).lt.xacc) goto 10
  1463)         enddo
  1464)  10     continue
  1465)         tm = tc*(1.d0-nu)
  1466)         ps1 = a1*nu+(a2*(nu**1.5d0))+(a3*(nu**2))+ &
  1467)              (a4*(nu**4.d0))
  1468)         ps2 = a1+(1.5d0*a2*(nu**0.5d0))+(2.0d0*a3*nu)+ &
  1469)              (4.d0*a4*(nu**3))
  1470)         derpt = (ps/tm)*(-ps2-((ps1*tc)/tm))
  1471)         derpt=1/derpt
  1472)       endif
  1473)       
  1474) end subroutine vappr
  1475) 
  1476) ! ************************************************************************** !
  1477) 
  1478) subroutine viscosity(p,t,rho,drhodp,drhodt,mu,dmudt,dmudp)
  1479) 
  1480) ! Fenghour, A., W. A. Wakeham, and V. Vesovic,
  1481) ! The viscosity of carbon dioxide,
  1482) ! J. Phys. Chem. Ref. Data, 27(1), 31−44, 1998.
  1483) 
  1484)       implicit none
  1485)       PetscReal :: p, t, rho
  1486)       PetscReal :: xsection1, xsection2, drhodt, drhodp,dmudp,dmudt
  1487)       PetscReal :: dtxsection,dt_zerodenmu,dp_zerodenmu,drho_excessmu
  1488)       PetscReal :: dp_excessmu, dt_excessmu
  1489)       PetscReal :: t_star, xsection, zeroden_mu, excess_mu, mu
  1490)       PetscReal :: lnstr,r2,r4,r5,r6,r7,r8
  1491) 
  1492) !     zero density viscosity
  1493) !     av(0) = 0.235156d0
  1494) !     av(1) = -.491266d0
  1495) !     av(2) = 5.211155d-2
  1496) !     av(3) = 5.347906d-2
  1497) !     av(4) = -1.537102d-2
  1498) 
  1499)       t_star = t/251.196
  1500)       
  1501)       lnstr = log(t_star)
  1502) 
  1503)       xsection = av(0) + (av(1) + (av(2) + (av(3) + &
  1504)           av(4)*lnstr)*lnstr)*lnstr)*lnstr
  1505)       
  1506)       xsection1 = xsection
  1507)       xsection = exp(xsection)
  1508) 
  1509)       zeroden_mu = 1.00697d0*(t**0.5d0)/xsection
  1510)       
  1511)       r2 = rho*rho
  1512)       r4 = r2*r2
  1513)       r5 = rho*r4
  1514)       r6 = rho*r5
  1515)       r7 = rho*r6
  1516)       r8 = r2*r6
  1517) 
  1518)       excess_mu = 0.4071119d-2*rho+0.7198037d-4*r2+ &
  1519)           0.2411697d-16*r6/t_star**3+ &
  1520)           0.2971072d-22*r8+ &
  1521)           (-0.1627888d-22*r8/t_star)
  1522) 
  1523)       mu = zeroden_mu+excess_mu
  1524) 
  1525)       xsection2 = av(1) + (2.d0*av(2) + &
  1526)           (3.d0*av(3) + 4.d0*av(4)*lnstr)*lnstr)*lnstr
  1527) 
  1528)       dtxsection = xsection2/(251.196d0*t_star)
  1529)       dt_zerodenmu =  exp(xsection1)*dtxsection
  1530) 
  1531)       dp_zerodenmu = 0.d0
  1532) 
  1533)       drho_excessmu = 0.4071119d-2+2.d0*0.7198037d-4*rho+ &
  1534)           6.d0*0.2411697d-16*r5/(t_star**3)+ &
  1535)           8.d0*0.2971072d-22*r7+ &
  1536)           8.d0*(-0.1627888d-22)*r7/t_star
  1537) 
  1538)       dp_excessmu = drho_excessmu
  1539) 
  1540)       dt_excessmu = 0.4071119d-2*drhodt+ &
  1541)           2.d0*0.7198037d-4*rho*drhodt+ &
  1542)           (0.2411697d-16*r5*((6.d0*t_star*drhodt) &
  1543)           -(3.d0*rho/251.196d0))/ &
  1544)           (t_star**4))+ &
  1545)           (8.d0*0.2971072d-22*r7*drhodt)+ &
  1546)           (-0.1627888d-22*r7*((8.d0*t_star*drhodt)- &
  1547)           (rho/251.196d0))/(t_star**2))
  1548) 
  1549)       dmudt = dt_excessmu + dt_zerodenmu
  1550) 
  1551)       dmudp = dp_zerodenmu + (dp_excessmu*drhodp)
  1552) 
  1553) end subroutine viscosity
  1554) 
  1555) ! ************************************************************************** !
  1556) 
  1557) subroutine dissco2(p,t,mco2,fg,mol)
  1558) 
  1559)       implicit none
  1560)       PetscReal :: p, t, mco2, fg, mol
  1561)       PetscReal :: pc_h2o, tc_h2o, t1,ph2o,liq_cp, lambdaco2_na
  1562)       PetscReal :: tauco2_na_cl, rhs
  1563) 
  1564) !     calculate solubility of CO2 based on Duan model
  1565) 
  1566) !     partial pressure of water based on empirical model provided
  1567) !     in Duan paper, Eqn. B1
  1568) 
  1569) !     convert pressure in bar and temperature in K
  1570) 
  1571)       p = p*10.d0
  1572)       pc_h2o = 220.85d0
  1573)       tc_h2o = 647.29d0
  1574) 
  1575)       t1 = (t/tc_h2o)-1.d0
  1576) 
  1577)       ph2o = (pc_h2o*t/tc_h2o)*(1.d0+(-38.640844d0*((-t1)**1.9d0)) &
  1578)           +(5.894842d0*t1)+(59.876516*t1*t1)+(26.654627*(t1**3.d0)) &
  1579)           +(10.637097*(t1**4.d0)))
  1580) 
  1581)       liq_cp = 28.9447706d0 + (-0.0354581768d0*t) + (-4770.67077d0/t) &
  1582)       +(1.02782768e-5*t*t) + (33.8126098/(630-t)) + (9.0403714d-3*p) &
  1583)       +(-1.14934031d-3*p*log(t)) + (-0.307405726*p/t) &
  1584)       + (-0.0907301486*p/(630.d0-t)) &
  1585)       + (9.32713393d-4*p*p/((630.d0-t)**2.d0))
  1586) 
  1587)       lambdaco2_na = -0.411370585d0 + (6.07632013d-4*t) &
  1588)       + (97.5347708d0/t) + (-0.0237622469d0*p/t) &
  1589)       + (0.0170656236*p/(630.d0-t)) + (1.41335834d-5*t*(log(p)))
  1590) 
  1591)       tauco2_na_cl = 3.36389723d-4 + (-1.9829898d-5*t) &
  1592)       + (2.1222083d-3*p/t) + (-5.24873303d-3*p/(630.d0-t))
  1593) 
  1594)       rhs = liq_cp - log(fg) + (2.d0*lambdaco2_na*mol) &
  1595)       + (tauco2_na_cl*mol*mol)
  1596) 
  1597)       mco2 = (p-ph2o)/exp(rhs)
  1598) 
  1599)       p = p*0.1d0
  1600) 
  1601) end subroutine dissco2
  1602) 
  1603) end module co2_span_wagner_module

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