erosion_enritchment_ratio Subroutine

public subroutine erosion_enritchment_ratio(j, da, da9, flu, pr, precip, rp, xqd, yd)

Arguments

Type IntentOptional AttributesName
integer :: j
real(kind=dp), intent(in) :: da
real(kind=dp), intent(in) :: da9
real(kind=dp), intent(in), dimension(:):: flu
real(kind=dp), intent(in) :: pr
real(kind=dp), intent(in) :: precip
real(kind=dp), intent(in) :: rp
real(kind=dp), intent(in) :: xqd
real(kind=dp), intent(in) :: yd

Called by

proc~~erosion_enritchment_ratio~~CalledByGraph proc~erosion_enritchment_ratio erosion_enritchment_ratio proc~runsubbasin runsubbasin proc~runsubbasin->proc~erosion_enritchment_ratio 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 erosion_enritchment_ratio(j, da, da9, flu, pr, precip, rp, xqd, yd)
    !**** PURPOSE: THIS SUBROUTINE CALCULATES ENRICHMENT RATIO FOR SUBBASIN
    !**** CALLED IN: SUBBASIN
    !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
    !     PARAMETERS & VARIABLES
    !
    !      >>>>> COMMON PARAMETERS & VARIABLES
    !      da = area of the basin, km2
    !      da9 = 100. * da = basin area in ha, from readbas
    !      er = enrichment ratio for subbasin
    !      flu(j) = fraction of subbasin area in the basin
    !      parsz(5, j) = particle size distribution, calc in subbasin
    !      pct(5, j) = delivery ratios (part. size distr. for sediments)
    !      pr = peak runoff rate, m3/sec., in peakq
    !      precip = precipitation, mm
    !      rp = alpha for rainfall, the fraction of total rainfall
    !                   occuring during 0.5h
    !      xqd = daily surface runoff, mm, calc in subbasin
    !      yd = daily soil loss, t, calc in ysed
    !      >>>>>

    !      >>>>> STATIC PARAMETERS
    !      bet = local par
    !      cy = local par
    !      dia = local par
    !      dr = local par
    !      durf = local par
    !      jj = local par
    !      nsz = local par
    !      rep = local par
    !      rinf = local par
    !      rp1 = local par
    !      tot = local par
    !      x1 = local par
    !      x2 = local par
    !      >>>>>
    !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~

    !**** Include common parameters

    real(dp), intent(in) :: da
    real(dp), intent(in) :: da9
    real(dp), dimension(:), intent(in) :: flu
    real(dp), intent(in) :: pr
    real(dp), intent(in) :: precip
    real(dp), intent(in) :: rp
    real(dp), intent(in) :: xqd
    real(dp), intent(in) :: yd
    integer j, jj
    real(dp) bet, cy, dr, durf, rep, rinf, rp1, tot, x1, x2

    rp1 = - 2. * precip *log(1. - rp)
    durf = 4.605 * precip / rp1
    rinf = (precip - xqd) / durf
    rep = (rp1 - rinf) * da * flu(j) / 3.6
    dr = (pr / (rep + 1.e-10)) ** .56
    bet =log(dr) / 4.47
    tot = 0.

    do jj = 1, nsz
      pct(jj, j) = parsz(jj, j) * exp(bet * dia(jj))
      tot = tot + pct(jj, j)
    end do

    do jj = 1, nsz
      pct(jj, j) = pct(jj, j) / tot
    end do

    cy = 100000. * yd / (da9 * xqd + 1.e-6)
    x2 = - log10(dr) / 2.699
    x1 = 1. / .25 ** x2

    !**** CALC ENRICHMENT RATIO FOR SUBBASIN
    if (cy > 1.e-6) then
      er = .78 * cy ** (- .2468)
    else
      er = 0.
    endif
    if (er .lt. 1.) er = 1.
    if (er .gt. 3.5) er = 3.5

    return
  end subroutine erosion_enritchment_ratio