hydrotope_subbasin_division Subroutine

public subroutine hydrotope_subbasin_division(mb, neap, sbar)

Uses

  • proc~~hydrotope_subbasin_division~~UsesGraph proc~hydrotope_subbasin_division hydrotope_subbasin_division module~utilities utilities proc~hydrotope_subbasin_division->module~utilities

Arguments

Type IntentOptional AttributesName
integer, intent(in) :: mb
integer, intent(inout) :: neap(:)
real(kind=dp), intent(inout) :: sbar(:)

Calls

proc~~hydrotope_subbasin_division~~CallsGraph proc~hydrotope_subbasin_division hydrotope_subbasin_division proc~log_info log_info proc~hydrotope_subbasin_division->proc~log_info proc~log_error log_error proc~hydrotope_subbasin_division->proc~log_error proc~log_message log_message proc~log_info->proc~log_message proc~log_error->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~~hydrotope_subbasin_division~~CalledByGraph proc~hydrotope_subbasin_division hydrotope_subbasin_division proc~hydrotope_initialise hydrotope_initialise proc~hydrotope_initialise->proc~hydrotope_subbasin_division proc~initialise initialise proc~initialise->proc~hydrotope_initialise program~swim swim program~swim->proc~initialise

Contents


Source Code

  subroutine hydrotope_subbasin_division(mb, neap, sbar)
    use utilities, only : log_error
    integer, intent(in) :: mb
    integer, intent(inout) :: neap(:)
    real(dp), intent(inout) :: sbar(:)
    integer i, j, jea

    ! Count max number of hydrotopes per subbasin
    meap = 0
    do i = 1, hydrotope_csv_size
      j = hydrotope_subbasin_id(i)
      if (j > 0) then
        neap(j) = neap(j) + 1 ! number of HRUs per subbasin
        sbar(j) = sbar(j) + hydrotope_area(i) ! total subbasin area in m^2
        if (j.gt.mb) call log_error("hydrotope_subbasin_division", &
          "subbasin_id in hydrotope input out of range.")
        if (neap(j) .gt. meap) meap = neap(j)
      else
        exit
      end if ! (j > 0)
    end do ! i = 1, mbeap

    ! Fill frar and mstruct sparse arrays
    allocate(frar(mb, meap))
    frar = 0.
    allocate(mstruc(mb, meap, 7))
    mstruc = 0

    ! fraction of hydrotop area of subbasin area
    i = 1
    do j = 1, mb
      frar(j, :neap(j)) = hydrotope_area(i:i+neap(j)-1) ! area in m^2 !!! recalculated below as fraction of subbasin
      mstruc(j, :neap(j), 1) = landuse_ids(i:i+neap(j)-1) ! land use ID
      mstruc(j, :neap(j), 2) = soil_ids(i:i+neap(j)-1) ! soil type ID
      mstruc(j, :neap(j), 3) = wetland(i:i+neap(j)-1) ! wetland
      mstruc(j, :neap(j), 4) = crop_management_id(i:i+neap(j)-1) ! LU management (this parameter has disappeared from hydrotope.csv)
      mstruc(j, :neap(j), 5) = elevations(i:i+neap(j)-1) ! HRU mean elevation
      mstruc(j, :neap(j), 6) = glaciers(i:i+neap(j)-1) ! initial glacier depth
      mstruc(j, :neap(j), 7) = irrigations(i:i+neap(j)-1) ! initial glacier depth
      i = i + neap(j)
      do jea = 1, neap(j)
        frar(j, jea) = frar(j, jea) / sbar(j)
      end do
    end do

    call log_info("hydrotope_subbasin_division", &
      "Max. number HRUs in subbasin:", int=meap)
  end subroutine hydrotope_subbasin_division