evapotranspiration_radiation Subroutine

public subroutine evapotranspiration_radiation(ida, mb, tmx)

Arguments

Type IntentOptional AttributesName
integer, intent(in) :: ida
integer, intent(in) :: mb
real(kind=dp), intent(in), dimension(:):: tmx

Calls

proc~~evapotranspiration_radiation~~CallsGraph proc~evapotranspiration_radiation evapotranspiration_radiation proc~log_warn log_warn proc~evapotranspiration_radiation->proc~log_warn proc~log_message log_message proc~log_warn->proc~log_message proc~log_write log_write proc~log_message->proc~log_write proc~log_format_message log_format_message proc~log_message->proc~log_format_message proc~to_string to_string proc~log_write->proc~to_string proc~date_time_str date_time_str proc~log_format_message->proc~date_time_str proc~colourise colourise proc~log_format_message->proc~colourise proc~string_index string_index proc~colourise->proc~string_index

Called by

proc~~evapotranspiration_radiation~~CalledByGraph proc~evapotranspiration_radiation evapotranspiration_radiation proc~time_process_day time_process_day proc~time_process_day->proc~evapotranspiration_radiation 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 evapotranspiration_radiation(ida, mb, tmx)

    integer, intent(in) :: ida
    integer, intent(in) :: mb
    real(dp), dimension(:), intent(in) :: tmx
    integer :: k
    real(dp) :: gb, dr, ds, ws2, ram_1, radi, ek2

    gb = 0.
    ! ek1 = empirical constant
    !ek1=0.16 ! used in West Africa (Niger Basin)
    !ek1 = 0.115 ! used in East Africa (Blue Nile Basin), old value
    !ec1 = 0.145 ! used in East Africa (Blue Nile Basin), for new SWIM trunk version
    ek2 = 0. !empirische Konstante 1

    dr = 1 + (.033 * cos(((2. * 3.141592) / 365.) * ida)) !Inverse relative Distanz Erde - Sonne
    ds = .409 * sin(((2. * 3.141592) / 365.) * ida - 1.39) ! Deklination der Sonne

    do k = 1, mb
      ! Note, if latitude is larger 67 degrees North, this function does not work!
      if (lat(k) > 66.5 ) then
        lat(k) = 66.5
        call log_warn("evapotranspiration_radiation", "Latitude is > 66.5, "// &
          "in subbasin (N) set to 66.5", i1=k)
      end if
      if (lat(k) < - 66.5 ) then
        lat(k) = - 66.5
        call log_warn("evapotranspiration_radiation", "Latitude is < -66.5, "// &
          "in subbasin (N) set to -66.5", i1=k)
      end if
      gb = (3.1415 / 180.) * lat(k) !Breite in Bogenwinkel
      ws2 = acos((- tan(gb)) * tan(ds)) !Stundenwinkel beim Sonnenuntergang
      ram_1 = (ws2 * sin(gb) * sin(ds)) + (cos(gb) * cos(ds) * sin(ws2))
      radi = ((24. * 60.) / 3.141592) * .0820 * dr * ram_1
      radi = (radi * 100.) ! Formel in MJ / m2 / d
      ra(k) = ec1 * radi * ((tmx(k) - tmn(k)) ** .5) + ek2 ! Formel von Hargreaves(1985)
    end do

  end subroutine evapotranspiration_radiation