co2_span_wagner_spline.F90       coverage:  0.00 %func     0.00 %block


     1)   module co2_span_wagner_spline_module
     2) 
     3)   use PFLOTRAN_Constants_module
     4) 
     5)   implicit none
     6) 
     7) #include "petsc/finclude/petscsys.h"
     8)   
     9)   save
    10) 
    11)   private
    12)   
    13)   PetscInt :: nptab,nttab0,ncrit_pts
    14)   PetscInt, allocatable :: nttab(:),ncrit(:)
    15)   PetscReal, allocatable :: p_tab(:),t_tab(:,:),r_tab(:,:), &
    16)                          h_tab(:,:),u_tab(:,:),s_tab(:,:),f_tab(:,:), &
    17)                          rr(:,:),hh(:,:),uu(:,:),ss(:,:),ff(:,:)
    18)   PetscReal, allocatable :: tcrit(:),pcrit(:),rhol(:),ul(:),hl(:),sl(:), &
    19)                          rhov(:),hv(:),uv(:),sv(:),fv(:)
    20)   
    21)   public sw_spline_read, sw_prop
    22)   
    23)   contains
    24) subroutine sw_spline_read
    25) 
    26)   use spline_module
    27) 
    28)   PetscInt :: i,ipx,j,n,iunit=9
    29)   
    30)   open (unit=iunit,file='co2_prop_TC.dat',status='old')
    31)   
    32)   read(iunit,*) nptab
    33)   
    34)   allocate(p_tab(nptab))
    35)   allocate(nttab(nptab))
    36)   allocate(ncrit(nptab))
    37) 
    38)   read(iunit,*) (p_tab(i),i=1,nptab)
    39)   
    40)   read(iunit,*) (nttab(i),i=1,nptab)
    41)   
    42) ! read: t, rho, h,u,f,s
    43)   nttab0 = 1
    44)   do i=1,nptab
    45)     nttab0 = max(nttab0,nttab(i))
    46)   enddo
    47)   print *,'nttab0=',nttab0
    48)   
    49)   allocate(t_tab(nptab,nttab0))
    50)   allocate(r_tab(nptab,nttab0))
    51)   allocate(h_tab(nptab,nttab0))
    52)   allocate(u_tab(nptab,nttab0))
    53)   allocate(s_tab(nptab,nttab0))
    54)   allocate(f_tab(nptab,nttab0))
    55)   allocate(rr(nptab,nttab0))
    56)   allocate(hh(nptab,nttab0))
    57)   allocate(uu(nptab,nttab0))
    58)   allocate(ss(nptab,nttab0))
    59)   allocate(ff(nptab,nttab0))
    60) 
    61)   do i=1,nptab
    62)     read(iunit,*) (t_tab(i,j),j=1,nttab(i))
    63)   enddo
    64) ! print *,'finished reading t'
    65)   
    66)   do i=1,nptab
    67)     read(iunit,*) (r_tab(i,j),j=1,nttab(i))
    68)   enddo
    69) ! print *,'finished reading r'
    70) 
    71)   do i=1,nptab
    72)     read(iunit,*) (h_tab(i,j),j=1,nttab(i))
    73)   enddo
    74) ! print *,'finished reading h'
    75) 
    76)   do i=1,nptab
    77)     read(iunit,*) (u_tab(i,j),j=1,nttab(i))
    78)   enddo
    79) ! print *,'finished reading u'
    80) 
    81)   do i=1,nptab
    82)     read(iunit,*) (f_tab(i,j),j=1,nttab(i))
    83)   enddo
    84) ! print *,'finished reading f'
    85)   
    86)   read(iunit,*) (ncrit(i),i=1,nptab)
    87)   read(iunit,*) ncrit_pts
    88) ! print *,'ncrit_pts=',ncrit_pts
    89)   
    90)   allocate(tcrit(ncrit_pts))
    91)   allocate(pcrit(ncrit_pts))
    92)   allocate(rhol(ncrit_pts))
    93)   allocate(ul(ncrit_pts))
    94)   allocate(hl(ncrit_pts))
    95)   allocate(sl(ncrit_pts))
    96)   allocate(rhov(ncrit_pts))
    97)   allocate(hv(ncrit_pts))
    98)   allocate(uv(ncrit_pts))
    99)   allocate(sv(ncrit_pts))
   100)   allocate(fv(ncrit_pts))
   101)   do i = 1, ncrit_pts
   102)     read(iunit,*) tcrit(i),pcrit(i),rhol(i),hl(i),ul(i),rhov(i),hv(i),uv(i),fv(i)
   103)   enddo
   104) 
   105) !---  Vapor temperature splines  ---
   106)   rr = 0.d0
   107)   do ipx = 1,nptab
   108)     n = nttab(ipx)-ncrit(ipx)+1
   109)     call spline(t_tab(ipx,ncrit(ipx):nttab(ipx)),r_tab(ipx,ncrit(ipx):nttab(ipx)),n,rr(ipx,ncrit(ipx):nttab(ipx)))
   110)     call spline(t_tab(ipx,ncrit(ipx):nttab(ipx)),h_tab(ipx,ncrit(ipx):nttab(ipx)),n,hh(ipx,ncrit(ipx):nttab(ipx)))
   111)     call spline(t_tab(ipx,ncrit(ipx):nttab(ipx)),u_tab(ipx,ncrit(ipx):nttab(ipx)),n,uu(ipx,ncrit(ipx):nttab(ipx)))
   112)     call spline(t_tab(ipx,ncrit(ipx):nttab(ipx)),f_tab(ipx,ncrit(ipx):nttab(ipx)),n,ff(ipx,ncrit(ipx):nttab(ipx)))
   113) !   call spline(t(ipx,ncrit(ipx)),s(ipx,ncrit(ipx)),n,ss(ipx,ncrit(ipx)))
   114) !   print *,'p= ',n,ipx,p_tab(ipx),rr(ipx,ncrit(ipx)),ncrit(ipx)
   115) !   if (ipx == 10 .or. ipx ==11) then
   116) !     print *,ncrit(ipx), nttab(ipx), rr(ipx,ncrit(ipx):nttab(ipx))
   117) !     print *,r_tab(ipx,ncrit(ipx):nttab(ipx))
   118) !   endif
   119)   enddo
   120) 
   121)   return
   122) end subroutine sw_spline_read
   123) 
   124) ! ************************************************************************** !
   125) 
   126) subroutine sw_prop(tx,px,rho,h,u,fg)
   127)       
   128)        use spline_module
   129) 
   130) !     density of liquid or vapor co2.
   131) 
   132) !     isrx liquid or vapor index: 1 - liquid 2 - vapor or supercritical
   133) 
   134) !     span, r., and w. wagner.  1996.  a new equation of state for
   135) !     carbon dioxide covering the fluid region from the triple-point
   136) !     to 1100 k at pressures up to 800 mpa.
   137) !     J. phys. chem. ref. data 25(6):1509-1588.
   138) 
   139)       implicit none
   140)       save
   141) 
   142)       PetscReal :: tx,px,rho, h,u,fg
   143)       PetscInt :: ipx,jpx,n
   144)       PetscReal :: tkx,pcx,ptx,tcx,ttx
   145)       PetscReal :: rtab(nptab+1) !,rtab2(nptab+1)
   146)       PetscReal :: htab(nptab+1) !,htab2(nptab+1)
   147)       PetscReal :: utab(nptab+1) !,utab2(nptab+1)
   148)       PetscReal :: fgtab(nptab+1) !,fgtab2(nptab+1)
   149) 
   150)       tkx = tx + 273.15d0
   151)       pcx = 7.3773d0
   152)       ptx = 0.51795d0
   153)       tcx = 304.1282d0
   154)       ttx = 216.592d0
   155) 
   156)           jpx = 0
   157)           do ipx = 1,nptab
   158)             n = nttab(ipx)-ncrit(ipx)+1
   159)             if (tkx.gt.t_tab(ipx,1)) then
   160)               jpx = jpx+1
   161)               call splint(t_tab(ipx,ncrit(ipx):nttab(ipx)),r_tab(ipx,ncrit(ipx):nttab(ipx)), &
   162)                 rr(ipx,ncrit(ipx):nttab(ipx)),n,tkx,rtab(jpx))
   163)  
   164)               call splint(t_tab(ipx,ncrit(ipx):nttab(ipx)),h_tab(ipx,ncrit(ipx):nttab(ipx)), &
   165)                 hh(ipx,ncrit(ipx):nttab(ipx)),n,tkx,htab(jpx))
   166) 
   167)               call splint(t_tab(ipx,ncrit(ipx):nttab(ipx)),u_tab(ipx,ncrit(ipx):nttab(ipx)), &
   168)                 uu(ipx,ncrit(ipx):nttab(ipx)),n,tkx,utab(jpx))
   169) 
   170)               call splint(t_tab(ipx,ncrit(ipx):nttab(ipx)),f_tab(ipx,ncrit(ipx):nttab(ipx)), &
   171)                 ff(ipx,ncrit(ipx):nttab(ipx)),n,tkx,fgtab(jpx))
   172) 
   173) !             print *,ipx,jpx,t_tab(ipx,ncrit(ipx)),r_tab(ipx,ncrit(ipx)),rtab(jpx)
   174)             endif
   175)           enddo
   176)           !call locate(p_tab,jpx,px,ipx)
   177)           !ipx = min(max(1,ipx),jpx-1)
   178) 
   179) #if 0
   180) ! Density
   181)           call spline(p_tab,rtab,nptab,rtab2)   
   182)           call splint(p_tab,rtab,rtab2,nptab,px,rho)
   183) ! H
   184)           call spline(p_tab,htab,nptab,htab2)   
   185)           call splint(p_tab,htab,htab2,nptab,px,h)
   186) ! U
   187)           call spline(p_tab,utab,nptab,utab2)   
   188)           call splint(p_tab,utab,utab2,nptab,px,u)
   189) ! fg
   190)           call spline(p_tab,fgtab,nptab,fgtab2)   
   191)           call splint(p_tab,fgtab,fgtab2,nptab,px,fg)
   192) #endif
   193) 
   194) ! ************** linear interpolation in pressure *******************
   195)           call locate(p_tab,jpx,px,ipx)
   196)           ipx = min(max(1,ipx),jpx-1)
   197)           rho = (rtab(ipx+1)-rtab(ipx))*(px-p_tab(ipx))/(p_tab(ipx+1)-p_tab(ipx)) + rtab(ipx)
   198)           h   = (htab(ipx+1)-htab(ipx))*(px-p_tab(ipx))/(p_tab(ipx+1)-p_tab(ipx)) + htab(ipx)
   199)           u   = (utab(ipx+1)-utab(ipx))*(px-p_tab(ipx))/(p_tab(ipx+1)-p_tab(ipx)) + utab(ipx)
   200)           fg  = (fgtab(ipx+1)-fgtab(ipx))*(px-p_tab(ipx))/(p_tab(ipx+1)-p_tab(ipx)) + fgtab(ipx)
   201) 
   202) !         print *,'density: ',tx,px,ipx,jpx,p_tab(ipx+1),p_tab(ipx),rtab(ipx+1),rtab(ipx),rho
   203) !************************************************************************************
   204)     return
   205) end subroutine sw_prop
   206)     
   207) end module co2_span_wagner_spline_module

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