subbasin_read_climate Subroutine

public subroutine subbasin_read_climate(humi, mb, ra, subp, tmn, tmx, tx)

Uses

  • proc~~subbasin_read_climate~~UsesGraph proc~subbasin_read_climate subbasin_read_climate module~input input proc~subbasin_read_climate->module~input module~utilities utilities module~input->module~utilities

Arguments

Type IntentOptional AttributesName
real(kind=dp), intent(inout), dimension(:):: humi
integer, intent(inout) :: mb
real(kind=dp), intent(inout), dimension(:):: ra
real(kind=dp), intent(inout), dimension(:):: subp
real(kind=dp), intent(inout), dimension(:):: tmn
real(kind=dp), intent(inout), dimension(:):: tmx
real(kind=dp), intent(inout), dimension(:):: tx

Calls

proc~~subbasin_read_climate~~CallsGraph proc~subbasin_read_climate subbasin_read_climate proc~read_real_column read_real_column proc~subbasin_read_climate->proc~read_real_column proc~input_type_conversion_error input_type_conversion_error proc~read_real_column->proc~input_type_conversion_error proc~move_lines move_lines proc~read_real_column->proc~move_lines proc~log_error log_error proc~read_real_column->proc~log_error proc~check_range check_range proc~read_real_column->proc~check_range proc~read_csv_item read_csv_item proc~read_real_column->proc~read_csv_item proc~header_column_index header_column_index proc~read_real_column->proc~header_column_index proc~input_error_column_not_found input_error_column_not_found proc~read_real_column->proc~input_error_column_not_found proc~input_type_conversion_error->proc~log_error proc~log_message log_message proc~log_error->proc~log_message proc~check_range->proc~log_error proc~log_warn log_warn proc~check_range->proc~log_warn proc~out_of_range_error out_of_range_error proc~check_range->proc~out_of_range_error proc~header_column_index->proc~move_lines proc~header_column_index->proc~input_error_column_not_found proc~input_error_column_not_found->proc~log_error proc~log_warn->proc~log_message proc~out_of_range_error->proc~log_error proc~to_string to_string proc~out_of_range_error->proc~to_string 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~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~~subbasin_read_climate~~CalledByGraph proc~subbasin_read_climate subbasin_read_climate proc~time_process_day time_process_day proc~time_process_day->proc~subbasin_read_climate 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


Source Code

  subroutine subbasin_read_climate(humi, mb, ra, subp, tmn, tmx, tx)

    !**** PURPOSE: this subroutine reads climate data
    !              & hydrological data from runoff.dat
    !**** CALLED IN: MAIN
    !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
    !      PARAMETERS & VARIABLES
    !
    !      >>>>> COMMON PARAMETERS & VARIABLES
    !      daylen(j) = day length in subbasin
    !      humi(j) = air humidity in the subbasin, %
    !      ida = current day in the year (Julian day)
    !      mb = number of subbasins
    !block pit = parameter for estimation of the day length
    !      precip = precipitation in the current subbasin, mm
    !      ra(j) = radiation, J/cm^2
    !      subp(j) = daily precipitation in the subbasin, mm
    !      tmn(j) = daily min temp in the subbasin, degree C
    !      tmx(j) = daily max temp in the subbasin, degree C
    !      tx(j) = daily aver temp in the subbasin, degree C
    !      ylc(j) = cos(lat()/clt), lat() - latitude, clt=57.296,
    !                  used to calc rmx in evap
    !      yls(j) = sin(lat()/clt), lat() - latitude, clt=57.296,
    !                  used to calc rmx in evap
    !      >>>>>

    !      >>>>> STATIC PARAMETERS
    !      ch = parameter for estim. day length
    !      h = parameter for estim. day length
    !      id = day
    !      iiyr = year
    !      imn = month
    !      k = count parameter
    !      sd = parameter for estim. day length
    !      xi = day (real(dp) number)
    !      >>>>>
    !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~

    use input, only : read_real_column

    real(dp), dimension(:), intent(inout) :: humi
    integer, intent(inout) :: mb
    real(dp), dimension(:), intent(inout) :: ra
    real(dp), dimension(:), intent(inout) :: subp
    real(dp), dimension(:), intent(inout) :: tmn
    real(dp), dimension(:), intent(inout) :: tmx
    real(dp), dimension(:), intent(inout) :: tx
    integer k, fid
    real(dp), dimension(mb) :: mx, mn

    !**** READ CLIMATE DATA & MEASURED DISCHARGE (current day)
    fid = subbasin_climate_file_id
    call read_real_column(fid, array=tx, index=tmean_column_ix)  ! no skip
    call read_real_column(fid, array=tmn, index=tmin_column_ix, skip=-mb)
    call read_real_column(fid, array=tmx, index=tmax_column_ix, skip=-mb)
    call read_real_column(fid, array=subp, index=precipitation_column_ix, skip=-mb)
    if (radiation_column_ix > 0) then
      call read_real_column(fid, array=ra, index=radiation_column_ix, skip=-mb)
    else
      ra = 1.
    end if
    if (humidity_column_ix > 0) then
      call read_real_column(fid, array=humi, index=humidity_column_ix, skip=-mb)
    else
      humi = - 999.9
    end if

    ! In order to avoid errors if Tmax < Tmin the following conversion is done
    do k = 1, mb
      mx(k) = max(tmx(k), tmn(k))
      mn(k) = min(tmx(k), tmn(k))
      tmx(k) = mx(k)
      tmn(k) = mn(k)
    end do
  end subroutine subbasin_read_climate