module modgps05refstruct 6,4
#if defined (DOC)
  !
  ! Structure containing the refractivity and derivatives
  ! of a model profile.
  !
  ! Functions:
  !      Construct the cached profiles.
  !
  ! Josep M. Aparicio
  ! Meteorological Service of Canada, 2006.
  !
#endif
  use modgps00base   , only : i4, dp, ngpssize
  use modgps03diff
  use modgps04profile, only : gpsprofile
  use modgps06gravity, only : gpsgravitysurf, gpsgravityabove
  
  implicit none
  
contains


  subroutine gpspre(prf) 8
    type(gpsprofile) :: prf

    integer(i4)      :: i, ngpslev
    real(dp)         :: rpt, rpr, rcf, pr1, termA, termB

    if (prf%bpst .EQV. .false.) then
       ngpslev = prf%ngpslev
       rpt = prf%rPT
       rpr = prf%rPR
       rcf = prf%rCF
       if ( abs(rcf-1._dp) .lt. 0.01_dp ) then
          do i = 1, ngpslev
             prf%pst(i)%Var               = 0.01_dp * (rpt+prf%rHYB(i)*(prf%rP0-rpt))
             prf%pst(i)%DVar              = 0._dp
             prf%pst(i)%DVar(2*ngpslev+1) = 0.01_dp * prf%rHYB(i)
          enddo
       else
          prf%pst(1)%Var               = 0.01_dp * rpt
          prf%pst(1)%DVar              = 0._dp
          prf%pst(1)%DVar(2*ngpslev+1) = 0._dp
          pr1 = 1._dp/(1._dp - rpt/rpr)
          do i = 2, ngpslev
             termB = ( (prf%rHYB(i) - rpt/rpr)*pr1 ) ** rcf
             termA = rpr * ( prf%rHYB(i) - termB )
             prf%pst(i)%Var               = 0.01_dp * (termA + termB * prf%rP0)
             prf%pst(i)%DVar              = 0._dp
             prf%pst(i)%DVar(2*ngpslev+1) = 0.01_dp * termB
          enddo
       endif
       prf%bpst=.true.
    endif
  end subroutine gpspre


  subroutine gpstem(prf) 6
    type(gpsprofile) :: prf

    integer(i4)      :: i, ngpslev

    if (prf%btst .EQV. .false.) then
       ngpslev = prf%ngpslev
       do i = 1, ngpslev
          prf%tst(i)%Var               = prf%rTT(i)+273.15_dp
          prf%tst(i)%DVar              = 0._dp
          prf%tst(i)%DVar(i)           = 1._dp
       enddo
       prf%btst=.true.
    endif
  end subroutine gpstem


  subroutine gpsQ(prf) 7
    type(gpsprofile) :: prf

    integer(i4)      :: i, ngpslev

    if (prf%bqst .EQV. .false.) then
       ngpslev = prf%ngpslev
       do i = 1, ngpslev
          prf%qst(i)%Var               = exp(prf%rLQ(i))
          prf%qst(i)%DVar              = 0._dp
          prf%qst(i)%DVar(ngpslev+i)   = exp(prf%rLQ(i))
       enddo
       prf%bqst=.true.
    endif
  end subroutine gpsQ


  subroutine gpsref(prf) 2,3
    type(gpsprofile) :: prf

    integer(i4)      :: i, ngpslev

    type(gpsdiff)               :: p, t, q
    type(gpsdiff)               :: pd
    type(gpsdiff)               :: pw
    type(gpsdiff)               :: tc

    type(gpsdiff)               :: nd
    type(gpsdiff)               :: nw
    type(gpsdiff)               :: n

    call gpspre(prf)
    call gpstem(prf)
    call gpsQ  (prf)
    if (prf%brst .EQV. .false.) then
       ngpslev = prf%ngpslev
       do i = 1, ngpslev
          p  = prf%pst(i)
          t  = prf%tst(i)
          q  = prf%qst(i)
          pw = 1.6069_dp*q*p/(1._dp+0.6069_dp*q)
          pd = p - pw
          nd = (77.6890_dp*pd/t                   )
          nw = (71.2952_dp*pw/t+375463._dp*pw/t**2)
          prf%rst(i)=nd+nw
       enddo
       prf%brst=.true.
    endif
  end subroutine gpsref


  subroutine gpsprew(prf, prew),3
    type(gpsprofile) :: prf
    type(gpsdiff)   , intent(out):: prew(:)

    integer(i4)      :: i, ngpslev

    type(gpsdiff)               :: p, q
    type(gpsdiff)               :: pw

    call gpspre(prf)
    call gpstem(prf)
    call gpsQ  (prf)
    ngpslev = prf%ngpslev
    do i = 1, ngpslev
       p  = prf%pst(i)
       q  = prf%qst(i)
       pw = 1.6069_dp*q*p/(1._dp+0.6069_dp*q)
       prew(i) = pw
    enddo
  end subroutine gpsprew


  subroutine gpsiprew(prf, iprew),3
    type(gpsprofile) :: prf
    type(gpsdiff)   , intent(out):: iprew

    integer(i4)      :: i, ngpslev

    type(gpsdiff)    :: pm, qm, pp, qp
    real(dp)         :: g

    call gpspre(prf)
    call gpsQ  (prf)
    ngpslev = prf%ngpslev
    iprew=0._dp
    do i = 1, ngpslev-1
       pm = prf%pst(i)
       qm = prf%qst(i)
       pp = prf%pst(i+1)
       qp = prf%qst(i+1)
       iprew=iprew+(qm+qp)*(pp-pm)
    enddo
    g  = gpsgravitysurf(prf%rLat)
    iprew = (100._dp * 0.5_dp/g) * iprew
  end subroutine gpsiprew


  subroutine gpsrefw(prf, refw) 1,3
    type(gpsprofile) :: prf
    type(gpsdiff)   , intent(out):: refw(:)

    integer(i4)      :: i, ngpslev

    type(gpsdiff)               :: p, t, q
    type(gpsdiff)               :: pd
    type(gpsdiff)               :: pw
    type(gpsdiff)               :: tc

    type(gpsdiff)               :: nd
    type(gpsdiff)               :: nw
    type(gpsdiff)               :: n

    call gpspre(prf)
    call gpstem(prf)
    call gpsQ  (prf)
!    if (prf%brst .EQV. .false.) then
    ngpslev = prf%ngpslev
    do i = 1, ngpslev
       p  = prf%pst(i)
       t  = prf%tst(i)
       q  = prf%qst(i)
       pw = 1.6069_dp*q*p/(1._dp+0.6069_dp*q)
       pd = p - pw
       nd = (77.6890_dp*pd/t                   )
       nw = (71.2952_dp*pw/t+375463._dp*pw/t**2)
       refw(i)=nw
!          prf%rst(i)=nd+nw
    enddo
!   prf%brst=.true.
!   endif
  end subroutine gpsrefw


  subroutine gpsZTD(prf, ZTD) 1,4
    type(gpsprofile)             :: prf
    type(gpsdiff)   , intent(out):: ZTD(:)

    real(dp)        , parameter  :: Avog  = 6.02214e26_dp
    real(dp)        , parameter  :: Boltz = 1.38065e-23_dp
    real(dp)        , parameter  :: mwDAir= 28.966_dp
    real(dp)        , parameter  :: Rd=Avog*Boltz/mwDAir
    
    integer(i4)                  :: i, ngpslev
    type(gpsdiff)                :: p, t, q, tv
    real(dp)                     :: g
    
    call gpspre(prf)
    call gpstem(prf)
    call gpsQ  (prf)
    ngpslev=prf%ngpslev
    do i = 1, ngpslev
       tv = (1._dp+0.6069_dp*prf%qst(i))*prf%tst(i)
       g  = gpsgravityabove(prf%rLat, prf%gst(i)%Var)
       ZTD(i)=prf%rst(i)*Rd*tv/(g*prf%pst(i))
    enddo
  end subroutine gpsZTD


  subroutine gpssat(prf, sat),3
    type(gpsprofile)             :: prf
    type(gpsdiff)   , intent(out):: sat(:)

    real(dp)        , parameter  :: lnPw0 = 1.809926_dp;
    real(dp)        , parameter  :: Aliq  = -5.350334e3_dp;
    real(dp)        , parameter  :: Asol  = -6.149077e3_dp;
    real(dp)        , parameter  :: InvT0 = 1._dp / 273.15_dp;
    
    integer(i4)                  :: i, ngpslev
    type(gpsdiff)                :: p, t, q, pw
    type(gpsdiff)                :: ExcInvT, lnPw, wvp
    
    call gpspre(prf)
    call gpstem(prf)
    call gpsQ  (prf)
    ngpslev=prf%ngpslev

    do i = 1, ngpslev
       p  = prf%pst(i)
       t  = prf%tst(i)
       q  = prf%qst(i)
       pw = 1.6069_dp*q*p/(1._dp+0.6069_dp*q)

       ExcInvT = 1._dp / t - InvT0;
       if ( t%Var .GE. 273.15_dp ) then
          lnPw = lnPw0 + Aliq * ExcInvT
       else
          lnPw = lnPw0 + Asol * ExcInvT
       endif
       wvp = 1.e2_dp * exp( lnPw )

       sat(i) = pw / wvp
    enddo
  end subroutine gpssat
  
end module modgps05refstruct