soil_curve_transmission_losses Subroutine

public subroutine soil_curve_transmission_losses(j, chk, chl, chw)

Arguments

Type IntentOptional AttributesName
integer :: j
real(kind=dp), intent(in), dimension(:, :):: chk
real(kind=dp), intent(in), dimension(:, :):: chl
real(kind=dp), intent(in), dimension(:, :):: chw

Called by

proc~~soil_curve_transmission_losses~~CalledByGraph proc~soil_curve_transmission_losses soil_curve_transmission_losses proc~runsubbasin runsubbasin proc~runsubbasin->proc~soil_curve_transmission_losses proc~time_process_day time_process_day proc~time_process_day->proc~runsubbasin proc~time_process_month time_process_month proc~time_process_month->proc~time_process_day proc~time_process_years time_process_years proc~time_process_years->proc~time_process_month program~swim swim program~swim->proc~time_process_years

Contents


Source Code

  subroutine soil_curve_transmission_losses(j, chk, chl, chw)
    !**** PURPOSE: THIS SUBROUTINE COMPUTES CHANNEL TRANSMISSION LOSSES
    !              and recalculates surface runoff xqd and peak runoff rate pr
    !**** CALLED IN: SUBBASIN
    !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
    !     PARAMETERS & VARIABLES
    !
    !      >>>>> COMMON PARAMETERS & VARIABLES
    !      chk(2, j) = effective hydraulic conductivity of main channel, mm/h
    !      chl(2, j) = main channel length, km
    !      chw(2, j) = average width of main channel, m
    !      dur = flow duration, h
    !      ida = current day
    !      itran = switch code to print from tran()
    !      iv = 1 in subbasin (reach 1 as outlet)
    !      pr = peak runoff rate, m3/sec.
    !      q1 = remember old xqd, mm, to compare with new one in subbasin
    !      vo = runoff volume = xqd * da * flu(j) * 1000, m3,
    !                 calc in subbasin
    !      xqd = surface runoff, mm, recalculation
    !      >>>>>

    !      >>>>> STATIC PARAMETERS
    !      a = local par
    !      ak = local par
    !      axw = local par
    !      b = local par
    !      bxw = local par
    !      pr1 = local par
    !      pxw = local par
    !      vol = local par
    !      xx = local par
    !      yy = local par
    !      zz = local par
    !      >>>>>
    !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
    !**** Include common parameters

    real(dp), dimension(:, :), intent(in) :: chk
    real(dp), dimension(:, :), intent(in) :: chl
    real(dp), dimension(:, :), intent(in) :: chw
    integer j
    real(dp) a, ak, axw, b, bxw, pr1, pxw, vol, xx, yy, zz

    !      q1 = xqd
    yy = vo / (xqd * 1000.)
    a = - .0001831 * chk(iv, j) * dur
    xx = .2649 * chk(iv, j) * dur / vo

    !**** RECALC xqd, pr
    !     iv = 1 (from subbasin)
    if (xx .lt. 1.) then
      ak = - 1.09 *log(1. - xx)
      b = exp(- ak)
      if (.not. ((1.-b) .le. 1.e-20)) then
        zz = - 2.04 * ak * chw(iv, j) * chl(iv, j)
        if (zz .ge. - 30.) then
          bxw = exp(zz)
          axw = (a / (1. - b)) * (1. - bxw)
          pxw = (- axw / bxw)
          vol = vo / 1233.5
          xqd = 0.
          if (vol .gt. pxw) xqd = axw + bxw * vol
          xqd = xqd * 1.234 / yy
          pr1 = 35.3 * pr
          pr = 12.1 * axw / dur - (1. - bxw) * vol + bxw * pr1
          pr = pr / 35.3
          if (pr .lt. 0.) pr = 0.
        else
          xqd = 0.
        end if
      end if
    end if

    return
  end subroutine soil_curve_transmission_losses