25 use mo_kind,
only : i4, dp
26 use mo_message,
only : message, error_message
63 use mo_append,
only : append
71 integer(i4),
intent(in) :: idomain
73 integer(i4) :: nrows1, ncols1
75 integer(i4) :: nrows11, ncols11
79 integer(i4) :: icc, jcc
81 integer(i4) :: iu, id, jl, jr
84 integer(i4),
dimension(:, :),
allocatable :: l11id_on_l1
87 integer(i4),
dimension(:, :),
allocatable :: l1id_on_l11
90 integer(i4),
dimension(:, :),
allocatable :: dummy_2d_id
92 real(dp) :: cellfactorrbyh
94 integer(i4) :: cellfactorrbyh_inv
97 nrows1 =
level1(idomain)%nrows
98 nrows11 =
level11(idomain)%nrows
99 ncols1 =
level1(idomain)%ncols
100 ncols11 =
level11(idomain)%ncols
103 allocate (l11id_on_l1(nrows1, ncols1))
104 allocate (l1id_on_l11(nrows11, ncols11))
109 cellfactorrbyh =
level11(idomain)%cellsize /
level1(idomain)%cellsize
113 if (cellfactorrbyh .lt. 1._dp)
then
114 allocate (dummy_2d_id(nrows1, ncols1))
116 cellfactorrbyh_inv = int(1. / cellfactorrbyh, i4)
120 if(.not.
level1(idomain)%mask(icc, jcc)) cycle
123 iu = (icc - 1) * cellfactorrbyh_inv + 1
124 id = min(icc * cellfactorrbyh_inv, nrows11)
125 jl = (jcc - 1) * cellfactorrbyh_inv + 1
126 jr = min(jcc * cellfactorrbyh_inv, ncols11)
128 l1id_on_l11(iu : id, jl : jr) = merge(dummy_2d_id(icc, jcc),
nodata_i4,
level11(idomain)%mask(iu : id, jl : jr))
132 allocate (dummy_2d_id(nrows11, ncols11))
138 if(.not.
level11(idomain)%mask(icc, jcc)) cycle
142 iu = (icc - 1) * nint(cellfactorrbyh, i4) + 1
143 id = icc * nint(cellfactorrbyh, i4)
144 jl = (jcc - 1) * nint(cellfactorrbyh, i4) + 1
145 jr = jcc * nint(cellfactorrbyh, i4)
149 if(id > nrows1) id = nrows1
151 if(jr > ncols1) jr = ncols1
154 l11id_on_l1(iu : id, jl : jr) = merge(dummy_2d_id(icc, jcc),
nodata_i4,
level1(idomain)%mask(iu : id, jl : jr))
164 deallocate(l11id_on_l1, l1id_on_l11, dummy_2d_id)
229 use mo_append,
only : append
235 use mo_string_utils,
only : num2str
240 integer(i4),
intent(in) :: idomain
242 integer(i4) :: ncells0
244 integer(i4) :: nrows0, ncols0
246 integer(i4) :: s0, e0
248 integer(i4) :: nrows11, ncols11
251 integer(i4) :: nnodes
253 integer(i4) :: ii, jj, kk, ic, jc
255 integer(i4) :: iu, id
257 integer(i4) :: jl, jr
259 integer(i4) :: irow, jcol
261 integer(i4),
dimension(:, :),
allocatable :: id0
263 integer(i4),
dimension(:, :),
allocatable :: fdir0
265 integer(i4),
dimension(:, :),
allocatable :: facc0
267 integer(i4),
dimension(:, :),
allocatable :: fdir11
270 integer(i4),
dimension(:),
allocatable :: rowout
273 integer(i4),
dimension(:),
allocatable :: colout
275 integer(i4),
dimension(:, :),
allocatable :: drasc0
278 integer(i4),
dimension(:, :),
allocatable :: oloc
282 integer(i4) :: faccmax, idmax
285 integer(i4) :: noutlet
288 integer(i4) :: old_noutlet
293 type(
grid),
pointer :: level0_idomain => null()
316 nrows0 = level0_idomain%nrows
317 ncols0 = level0_idomain%ncols
318 ncells0 = level0_idomain%ncells
319 nrows11 =
level11(idomain)%nrows
320 ncols11 =
level11(idomain)%ncols
321 nnodes =
level11(idomain)%ncells
322 s0 = level0_idomain%iStart
323 e0 = level0_idomain%iEnd
325 allocate (id0(nrows0, ncols0))
326 allocate (facc0(nrows0, ncols0))
327 allocate (fdir0(nrows0, ncols0))
328 allocate (drasc0(nrows0, ncols0))
329 allocate (fdir11(nrows11, ncols11))
330 allocate (rowout(nnodes))
331 allocate (colout(nnodes))
332 allocate (oloc(1, 2))
346 id0(:, :) = unpack(level0_idomain%Id, level0_idomain%mask,
nodata_i4)
351 IF(ncells0 .EQ. nnodes)
THEN
352 oloc(1, :) = maxloc(facc0, level0_idomain%mask)
353 kk =
l0_l11_remap(idomain)%lowres_id_on_highres(oloc(1, 1), oloc(1, 2))
355 if(ncells0 .EQ. 1)
then
356 fdir11(1, 1) = fdir0(oloc(1, 1), oloc(1, 2))
358 fdir11(:, :) = fdir0(:, :)
360 fdir11(
level11(idomain)%CellCoor(kk, 1),
level11(idomain)%CellCoor(kk, 2)) = 0
363 ii =
level11(idomain)%CellCoor(kk, 1)
364 jj =
level11(idomain)%CellCoor(kk, 2)
369 ii = level0_idomain%CellCoor(kk, 1)
370 jj = level0_idomain%CellCoor(kk, 2)
380 irow = level0_idomain%CellCoor(ii, 1)
381 jcol = level0_idomain%CellCoor(ii, 2)
385 if ((irow .le. 0_i4) .or. (irow .gt. nrows0) .or. &
386 (jcol .le. 0_i4) .or. (jcol .gt. ncols0))
then
389 if (fdir0(irow, jcol) .le. 0) is_outlet = .true.
393 noutlet = noutlet + 1_i4
395 if (noutlet .eq. 1)
then
396 oloc(1, :) = level0_idomain%CellCoor(ii, :)
398 call append(oloc, level0_idomain%CellCoor(ii : ii, :))
401 kk =
l0_l11_remap(idomain)%lowres_id_on_highres(oloc(noutlet, 1), oloc(noutlet, 2))
402 drasc0(oloc(noutlet, 1), oloc(noutlet, 2)) = kk
409 if (maxval(facc0(iu : id, jl : jr)) .eq. facc0(oloc(noutlet, 1), oloc(noutlet, 2)))
then
411 rowout(kk) = oloc(noutlet, 1)
412 colout(kk) = oloc(noutlet, 2)
413 fdir11(
level11(idomain)%CellCoor(kk, 1),
level11(idomain)%CellCoor(kk, 2)) = 0
423 if (rowout(kk) > 0) cycle
425 ic =
level11(idomain)%CellCoor(kk, 1)
426 jc =
level11(idomain)%CellCoor(kk, 2)
439 if ((facc0(iu, jj) > faccmax) .and. &
440 (fdir0(iu, jj) == 32 .or. &
441 fdir0(iu, jj) == 64 .or. &
442 fdir0(iu, jj) == 128))
then
443 faccmax = facc0(iu, jj)
451 if ((facc0(ii, jr) > faccmax) .and. &
452 (fdir0(ii, jr) == 1 .or. &
453 fdir0(ii, jr) == 2 .or. &
454 fdir0(ii, jr) == 128))
then
455 faccmax = facc0(ii, jr)
463 if ((facc0(id, jj) > faccmax) .and. &
464 (fdir0(id, jj) == 2 .or. &
465 fdir0(id, jj) == 4 .or. &
466 fdir0(id, jj) == 8))
then
467 faccmax = facc0(id, jj)
475 if ((facc0(ii, jl) > faccmax) .and. &
476 (fdir0(ii, jl) == 8 .or. &
477 fdir0(ii, jl) == 16 .or. &
478 fdir0(ii, jl) == 32))
then
479 faccmax = facc0(ii, jl)
486 ii = level0_idomain%CellCoor(idmax, 1)
487 jj = level0_idomain%CellCoor(idmax, 2)
493 if (ii == iu .and. jj == jl)
then
494 select case (fdir0(ii, jj))
502 elseif (ii == iu .and. jj == jr)
then
503 select case (fdir0(ii, jj))
511 elseif (ii == id .and. jj == jl)
then
512 select case (fdir0(ii, jj))
520 elseif (ii == id .and. jj == jr)
then
521 select case (fdir0(ii, jj))
541 call error_message(
'Error L11_flow_direction: side = -1')
562 if (idomain .eq. 1)
then
563 call append(
l0_drasc, pack(drasc0(:, :), level0_idomain%mask))
564 else if (
domainmeta%L0DataFrom(idomain) == idomain)
then
565 call append(
l0_drasc, pack(drasc0(:, :), level0_idomain%mask))
570 old_noutlet =
size(
domain_mrm(idomain)%L0_rowOutlet, dim = 1)
571 if (noutlet .le. old_noutlet)
then
572 domain_mrm(idomain)%L0_rowOutlet(: noutlet) = oloc(:, 1)
573 domain_mrm(idomain)%L0_colOutlet(: noutlet) = oloc(:, 2)
576 domain_mrm(idomain)%L0_rowOutlet(: old_noutlet) = oloc(: old_noutlet, 1)
577 domain_mrm(idomain)%L0_colOutlet(: old_noutlet) = oloc(: old_noutlet, 2)
580 call append(
domain_mrm(idomain)%L0_rowOutlet, oloc(old_noutlet + 1 :, 1))
581 call append(
domain_mrm(idomain)%L0_colOutlet, oloc(old_noutlet + 1 :, 2))
592 call message(
' Domain: ' // num2str(idomain,
'(i3)'))
593 call message(
' Number of outlets found at Level 0:.. ' // num2str(noutlet,
'(i7)'))
594 call message(
' Number of outlets found at Level 11:. ' // num2str(count(fdir11 .eq. 0_i4),
'(i7)'))
597 deallocate(fdir0, facc0, fdir11, rowout, colout, drasc0)
632 use mo_append,
only : append
639 integer(i4),
intent(in) :: idomain
641 integer(i4),
dimension(:, :),
allocatable :: fdir11
643 integer(i4),
dimension(:, :),
allocatable :: dummy_2d_id
645 integer(i4) :: jj, kk, ic, jc
647 integer(i4) :: fn, tn
649 integer(i4),
dimension(:),
allocatable :: nlinkfromn, nlinkton
655 allocate (nlinkfromn(
level11(idomain)%nCells))
656 allocate (nlinkton(
level11(idomain)%nCells))
658 allocate (dummy_2d_id(
level11(idomain)%nrows,
level11(idomain)%ncols))
675 do kk = 1,
level11(idomain)%nCells
676 ic =
level11(idomain)%CellCoor(kk, 1)
677 jc =
level11(idomain)%CellCoor(kk, 2)
680 tn = dummy_2d_id(ic, jc)
693 call append(
l11_ton, nlinkton(:))
696 deallocate (fdir11, nlinkfromn, nlinkton)
730 use mo_append,
only : append
738 integer(i4),
intent(in) :: idomain
740 integer(i4) :: nlinks
743 integer(i4),
dimension(:),
allocatable :: nlinkfromn
746 integer(i4),
dimension(:),
allocatable :: nlinkton
749 integer(i4),
dimension(:),
allocatable :: nlinkrorder
752 integer(i4),
dimension(:),
allocatable :: nlinklabel
755 logical,
dimension(:),
allocatable :: nlinksink
758 integer(i4),
dimension(:),
allocatable :: netperm
760 integer(i4) :: ii, jj, kk
770 allocate (nlinkfromn(
level11(idomain)%nCells))
771 allocate (nlinkton(
level11(idomain)%nCells))
772 allocate (nlinkrorder(
level11(idomain)%nCells))
773 allocate (nlinklabel(
level11(idomain)%nCells))
774 allocate (nlinksink(
level11(idomain)%nCells))
775 allocate (netperm(
level11(idomain)%nCells))
779 nlinkrorder(1 : nlinks) = 1
781 nlinklabel(1 : nlinks) = 0
783 nlinksink(:) = .false.
787 if(
level11(idomain)%nCells .GT. 1)
then
792 loop1 :
do ii = 1, nlinks
793 loop2 :
do jj = 1, nlinks
794 if (jj == ii) cycle loop2
795 if (nlinkfromn(ii) == nlinkton(jj))
then
798 if (nlinkrorder(ii) == -9) cycle loop1
804 if (nlinkrorder(ii) == 1)
then
811 do while (minval(nlinkrorder(1 : nlinks)) < 0)
813 loop3 :
do ii = 1, nlinks
814 if (.NOT. nlinkrorder(ii) == -9) cycle loop3
816 loop4 :
do jj = 1, nlinks
817 if (jj == ii .OR. nlinkfromn(ii) /= nlinkton(jj))
then
819 else if (.NOT. (nlinkfromn(ii) == nlinkton(jj) .AND. nlinkrorder(jj) > 0))
then
835 if (
l11_fdir(
level11(idomain)%iStart + nlinkton(ii) - 1_i4) .eq. 0_i4)
then
836 nlinksink(ii) = .true.
841 where(nlinksink) nlinklabel = 2
845 netperm(nlinkrorder(ii)) = ii
861 deallocate (nlinkfromn, nlinkton, nlinkrorder, nlinklabel, nlinksink, netperm)
893 use mo_append,
only : append
899 use mo_string_utils,
only : num2str
904 integer(i4),
intent(in) :: idomain
906 integer(i4) :: nlinks
909 integer(i4),
dimension(:),
allocatable :: rowout
912 integer(i4),
dimension(:),
allocatable :: colout
914 integer(i4),
dimension(:),
allocatable :: nlinkfromn
916 integer(i4),
dimension(:),
allocatable :: netperm
918 integer(i4),
dimension(:),
allocatable :: nlinkfromrow
920 integer(i4),
dimension(:),
allocatable :: nlinkfromcol
922 integer(i4),
dimension(:),
allocatable :: nlinktorow
924 integer(i4),
dimension(:),
allocatable :: nlinktocol
926 integer(i4),
dimension(:, :),
allocatable :: fdir0
928 integer(i4),
dimension(:, :),
allocatable :: drasc0
930 integer(i4) :: ii, rr, kk, s0, e0
932 integer(i4) :: inode, irow, jcol, prevrow, prevcol
935 integer(i4),
dimension(:, :),
allocatable :: oloc
938 integer(i4) :: noutlets
943 type(
grid),
pointer :: level0_idomain => null()
947 s0 = level0_idomain%iStart
948 e0 = level0_idomain%iEnd
952 nlinks =
level11(idomain)%nCells - noutlets
957 allocate (rowout(
level11(idomain)%nCells))
958 allocate (colout(
level11(idomain)%nCells))
959 allocate (nlinkfromn(
level11(idomain)%nCells))
960 allocate (netperm(
level11(idomain)%nCells))
961 allocate (nlinkfromrow(
level11(idomain)%nCells))
962 allocate (nlinkfromcol(
level11(idomain)%nCells))
963 allocate (nlinktorow(
level11(idomain)%nCells))
964 allocate (nlinktocol(
level11(idomain)%nCells))
965 allocate (fdir0(level0_idomain%nrows, level0_idomain%ncols))
966 allocate (drasc0(level0_idomain%nrows, level0_idomain%ncols))
981 if(
level11(idomain)%nCells .GT. 1)
then
993 allocate(oloc(noutlets, 2))
994 oloc(:, 1) =
domain_mrm(idomain)%L0_rowOutlet(: noutlets)
995 oloc(:, 2) =
domain_mrm(idomain)%L0_colOutlet(: noutlets)
1001 inode = nlinkfromn(ii)
1002 irow = rowout(inode)
1003 jcol = colout(inode)
1006 nlinkfromrow(ii) = irow
1007 nlinkfromcol(ii) = jcol
1012 if (irow .eq. oloc(kk, 1) .and. jcol .eq. oloc(kk, 2)) is_outlet = .true.
1017 nlinktorow(ii) = irow
1018 nlinktocol(ii) = jcol
1022 do while (.not. (drasc0(irow, jcol) > 0))
1028 if (irow .eq. oloc(kk, 1) .and. jcol .eq. oloc(kk, 2))
exit
1030 if (prevrow .eq. irow .and. prevcol .eq. jcol)
then
1031 call error_message(
'Something went wrong during L11_link_location, ', &
1032 'movedownonecell got stuck in infinite loop at cell (', num2str(irow),
' ', &
1037 nlinktorow(ii) = irow
1038 nlinktocol(ii) = jcol
1050 call append(
l11_frow, nlinkfromrow(:))
1051 call append(
l11_fcol, nlinkfromcol(:))
1052 call append(
l11_trow, nlinktorow(:))
1053 call append(
l11_tcol, nlinktocol(:))
1056 deallocate (rowout, colout, nlinkfromn, netperm, nlinkfromrow, &
1057 nlinkfromcol, nlinktorow, nlinktocol, fdir0, drasc0)
1094 use mo_append,
only : append
1104 integer(i4),
intent(in) :: idomain
1106 integer(i4),
dimension(:, :),
allocatable :: drasc0
1108 integer(i4),
dimension(:, :),
allocatable :: fdir0
1110 integer(i4),
dimension(:, :),
allocatable :: gaugeloc0
1112 integer(i4),
dimension(:, :),
allocatable :: inflowgaugeloc0
1114 integer(i4),
dimension(:, :),
allocatable :: dracell0
1116 integer(i4) :: ii, jj, kk, ll, s0, e0
1120 integer(i4) :: irow, jcol
1122 type(
grid),
pointer :: level0_idomain => null()
1126 s0 = level0_idomain%iStart
1127 e0 = level0_idomain%iEnd
1131 allocate (drasc0(level0_idomain%nrows, level0_idomain%ncols))
1132 allocate (fdir0(level0_idomain%nrows, level0_idomain%ncols))
1133 allocate (gaugeloc0(level0_idomain%nrows, level0_idomain%ncols))
1134 allocate (inflowgaugeloc0(level0_idomain%nrows, level0_idomain%ncols))
1135 allocate (dracell0(level0_idomain%nrows, level0_idomain%ncols))
1144 drasc0(:, :) = unpack(
l0_drasc(s0 : e0), &
1146 fdir0(:, :) = unpack(
l0_fdir(s0 : e0), &
1153 do kk = 1, level0_idomain%nCells
1154 ii = level0_idomain%CellCoor(kk, 1)
1155 jj = level0_idomain%CellCoor(kk, 2)
1156 isc = drasc0(ii, jj)
1160 do while (.NOT. isc > 0)
1163 isc = drasc0(irow, jcol)
1165 dracell0(ii, jj) = isc
1169 if (gaugeloc0(ii, jj) .NE.
nodata_i4)
then
1173 if (
domain_mrm(idomain)%gaugeIdList(ll) .EQ. gaugeloc0(ii, jj))
then
1179 if (inflowgaugeloc0(ii, jj) .NE.
nodata_i4)
then
1183 if (
domain_mrm(idomain)%InflowGaugeIdList(ll) .EQ. inflowgaugeloc0(ii, jj)) &
1195 if (idomain .eq. 1)
then
1196 call append(
l0_dracell, pack(dracell0(:, :), level0_idomain%mask))
1197 else if (
domainmeta%L0DataFrom(idomain) == idomain)
then
1198 call append(
l0_dracell, pack(dracell0(:, :), level0_idomain%mask))
1202 deallocate (drasc0, fdir0, gaugeloc0, dracell0)
1239 use mo_append,
only : append
1246 use mo_percentile,
only: percentile
1251 integer(i4),
intent(in) :: idomain
1253 integer(i4) :: nlinks
1255 integer(i4),
dimension(:, :),
allocatable :: id0
1257 integer(i4),
dimension(:, :),
allocatable :: fdir0
1259 real(dp),
dimension(:, :),
allocatable :: elev0
1261 real(dp),
dimension(:, :),
allocatable :: cellarea0
1263 integer(i4),
dimension(:, :),
allocatable :: streamnet0
1265 integer(i4),
dimension(:, :),
allocatable :: floodplain0
1268 integer(i4),
dimension(:),
allocatable :: netperm
1270 integer(i4),
dimension(:),
allocatable :: nlinkfromrow
1272 integer(i4),
dimension(:),
allocatable :: nlinkfromcol
1274 integer(i4),
dimension(:),
allocatable :: nlinktorow
1276 integer(i4),
dimension(:),
allocatable :: nlinktocol
1278 real(dp),
dimension(:),
allocatable :: nlinklength
1280 real(dp),
dimension(:),
allocatable :: nlinkafloodplain
1282 real(dp),
dimension(:),
allocatable :: nlinkslope
1284 integer(i4) :: ii, rr, ns, s0, e0
1286 integer(i4) :: frow, fcol
1288 integer(i4) :: fid, tid
1290 integer(i4),
dimension(:, :),
allocatable :: stack, append_chunk
1292 integer(i4),
dimension(:),
allocatable :: dummy_1d
1296 integer(i4),
dimension(:, :),
allocatable :: nodata_i4_tmp
1298 real(dp),
dimension(:, :),
allocatable :: nodata_dp_tmp
1300 type(
grid),
pointer :: level0_idomain => null()
1304 s0 = level0_idomain%iStart
1305 e0 = level0_idomain%iEnd
1310 allocate (id0(level0_idomain%nrows, level0_idomain%ncols))
1311 allocate (elev0(level0_idomain%nrows, level0_idomain%ncols))
1312 allocate (fdir0(level0_idomain%nrows, level0_idomain%ncols))
1313 allocate (cellarea0(level0_idomain%nrows, level0_idomain%ncols))
1314 allocate (streamnet0(level0_idomain%nrows, level0_idomain%ncols))
1315 allocate (floodplain0(level0_idomain%nrows, level0_idomain%ncols))
1319 allocate (stack(
level11(idomain)%nCells, 2))
1320 allocate (dummy_1d(2))
1321 allocate (append_chunk(8, 2))
1322 allocate (netperm(
level11(idomain)%nCells))
1323 allocate (nlinkfromrow(
level11(idomain)%nCells))
1324 allocate (nlinkfromcol(
level11(idomain)%nCells))
1325 allocate (nlinktorow(
level11(idomain)%nCells))
1326 allocate (nlinktocol(
level11(idomain)%nCells))
1327 allocate (nlinklength(
level11(idomain)%nCells))
1328 allocate (nlinkafloodplain(
level11(idomain)%nCells))
1329 allocate (nlinkslope(
level11(idomain)%nCells))
1331 allocate (nodata_i4_tmp(level0_idomain%nrows, level0_idomain%ncols))
1332 allocate (nodata_dp_tmp(level0_idomain%nrows, level0_idomain%ncols))
1357 if(
level11(idomain)%nCells .GT. 1)
then
1359 id0(:, :) = unpack(level0_idomain%Id, level0_idomain%mask, nodata_i4_tmp)
1360 elev0(:, :) = unpack(
l0_elev(s0 : e0), &
1361 level0_idomain%mask, nodata_dp_tmp)
1362 fdir0(:, :) = unpack(
l0_fdir(s0 : e0), &
1363 level0_idomain%mask, nodata_i4_tmp)
1364 cellarea0(:, :) = unpack(level0_idomain%CellArea, level0_idomain%mask, nodata_dp_tmp)
1380 frow = nlinkfromrow(ii)
1381 fcol = nlinkfromcol(ii)
1384 streamnet0(frow, fcol) = ii
1385 floodplain0(frow, fcol) = ii
1393 nlinkslope(ii) = elev0(frow, fcol)
1395 fid = id0(frow, fcol)
1396 tid = id0(nlinktorow(ii), nlinktocol(ii))
1398 do while (.NOT. (fid == tid))
1405 if (ns + 8 .gt.
size(stack, 1))
then
1406 call append(stack, append_chunk)
1408 call moveup(elev0, fdir0, frow, fcol, stack, ns)
1413 dummy_1d = stack(1, :)
1414 stack(:
size(stack, dim = 1) - 1, :) = stack(2 :, :)
1415 stack(
size(stack, dim = 1), :) = dummy_1d
1417 if (stack(1, 1) > 0 .and. stack(1, 2) > 0) floodplain0(stack(1, 1), stack(1, 2)) = ii
1418 ns = count(stack > 0) / 2
1424 streamnet0(frow, fcol) = ii
1425 floodplain0(frow, fcol) = ii
1426 fid = id0(frow, fcol)
1432 nlinklength(ii) = nlinklength(ii) + length
1437 nlinkslope(ii) = (nlinkslope(ii) - elev0(frow, fcol)) / nlinklength(ii)
1439 if (nlinkslope(ii) < 0.0001_dp) nlinkslope(ii) = 0.0001_dp
1442 nlinkafloodplain(ii) = sum(cellarea0(:, :), mask = (floodplain0(:, :) == ii))
1452 if (count(nlinklength(:) .ge. 0._dp) .gt. 2)
then
1453 length = percentile(pack(nlinklength(:), nlinklength(:) .ge. 0._dp), 40._dp)
1454 nlinklength(:) = merge(nlinklength(:), length, (nlinklength(:) .gt. length))
1465 if (idomain .eq. 1)
then
1466 call append(
l0_streamnet, pack(streamnet0(:, :), level0_idomain%mask))
1467 call append(
l0_floodplain, pack(floodplain0(:, :), level0_idomain%mask))
1468 else if (
domainmeta%L0DataFrom(idomain) == idomain)
then
1469 call append(
l0_streamnet, pack(streamnet0(:, :), level0_idomain%mask))
1470 call append(
l0_floodplain, pack(floodplain0(:, :), level0_idomain%mask))
1481 id0, elev0, fdir0, streamnet0, floodplain0, &
1482 cellarea0, stack, netperm, nlinkfromrow, nlinkfromcol, nlinktorow, nlinktocol, &
1483 nlinklength, nlinkafloodplain, nlinkslope, dummy_1d)
1484 deallocate(nodata_i4_tmp, nodata_dp_tmp)
1521 use mo_append,
only : append
1531 integer(i4),
intent(in) :: lcclassimp
1533 logical,
intent(in) :: do_init
1535 integer(i4) :: nlinks
1537 real(dp),
dimension(:),
pointer :: nlinkafloodplain => null()
1539 real(dp),
dimension(:,:),
allocatable :: temp_array
1541 integer(i4) :: ii, idomain, iilc, s0, e0
1543 type(
grid),
pointer :: level0_idomain => null()
1553 s0 = level0_idomain%iStart
1554 e0 = level0_idomain%iEnd
1560 if(nlinks .GT. 0)
then
1562 temp_array(ii, iilc) = sum(level0_idomain%CellArea(:), &
1564 / nlinkafloodplain(ii)
1570 deallocate(temp_array)
1604 subroutine moveup(elev0, fDir0, fi, fj, ss, nn)
1607 use mo_utils,
only : ge, le
1611 real(dp),
dimension(:, :),
allocatable,
intent(IN) :: elev0
1613 integer(i4),
dimension(:, :),
allocatable,
intent(IN) :: fDir0
1616 integer(i4),
intent(IN) :: fi, fj
1618 integer(i4),
dimension(:, :),
intent(INOUT) :: ss
1620 integer(i4),
intent(INOUT) :: nn
1622 integer(i4) :: ii, jj, ip, im, jp, jm
1624 integer(i4) :: nrows, ncols
1634 nrows =
size(fdir0, 1)
1635 ncols =
size(fdir0, 2)
1638 if (jp <= ncols)
then
1639 if ((fdir0(ii, jp) == 16) .and. &
1640 (le((elev0(ii, jp) - elev0(fi, fj)),
deltah)) .and. &
1641 (ge((elev0(ii, jp) - elev0(fi, fj)), 0.0_dp)) &
1651 if ((ip <= nrows) .and. &
1653 if ((fdir0(ip, jp) == 32) .and. &
1654 (le((elev0(ip, jp) - elev0(fi, fj)),
deltah)) .and. &
1655 (ge((elev0(ii, jp) - elev0(fi, fj)), 0.0_dp)) &
1665 if ((ip <= nrows) .and. &
1667 if ((fdir0(ip, jj) == 64) .and. &
1668 (le((elev0(ip, jj) - elev0(fi, fj)),
deltah)) .and. &
1669 (ge((elev0(ii, jp) - elev0(fi, fj)), 0.0_dp)) &
1679 if ((ip <= nrows) .and. &
1680 (jp <= ncols) .and. &
1682 if ((fdir0(ip, jm) == 128) .and. &
1683 (le((elev0(ip, jm) - elev0(fi, fj)),
deltah)) .and. &
1684 (ge((elev0(ii, jp) - elev0(fi, fj)), 0.0_dp)) &
1694 if ((jm >= 1) .and. &
1696 if ((fdir0(ii, jm) == 1) .and. &
1697 (le((elev0(ii, jm) - elev0(fi, fj)),
deltah)) .and. &
1698 (ge((elev0(ii, jp) - elev0(fi, fj)), 0.0_dp)) &
1708 if ((im >= 1) .and. &
1709 (jp <= ncols) .and. &
1711 if ((fdir0(im, jm) == 2) .and. &
1712 (le((elev0(im, jm) - elev0(fi, fj)),
deltah)) .and. &
1713 (ge((elev0(ii, jp) - elev0(fi, fj)), 0.0_dp)) &
1723 if ((im >= 1) .and. &
1725 if ((fdir0(im, jj) == 4) .and. &
1726 (le((elev0(im, jj) - elev0(fi, fj)),
deltah)) .and. &
1727 (ge((elev0(ii, jp) - elev0(fi, fj)), 0.0_dp)) &
1737 if ((im >= 1) .and. &
1739 if ((fdir0(im, jp) == 8) .and. &
1740 (le((elev0(im, jp) - elev0(fi, fj)),
deltah)) .and. &
1741 (ge((elev0(ii, jp) - elev0(fi, fj)), 0.0_dp)) &
1780 integer(i4),
intent(IN) :: fDir
1782 integer(i4),
intent(INOUT) :: iRow, jCol
1840 subroutine celllength(iDomain, fDir, iRow, jCol, iCoorSystem, length)
1844 use mo_constants,
only : sqrt2_dp
1848 integer(i4),
intent(IN) :: iDomain
1850 integer(i4),
intent(IN) :: fDir
1852 integer(i4),
intent(IN) :: iRow
1854 integer(i4),
intent(IN) :: jCol
1856 integer(i4),
intent(IN) :: iCoorSystem
1858 real(dp),
intent(OUT) :: length
1860 integer(i4) :: iRow_to, jCol_to
1862 real(dp) :: lat_1, long_1, lat_2, long_2
1864 type(
grid),
pointer :: level0_iDomain => null()
1870 IF(icoorsystem .EQ. 0)
THEN
1878 length = length * level0_idomain%cellsize
1881 ELSE IF(icoorsystem .EQ. 1)
THEN
1889 lat_1 = level0_idomain%yllcorner + real((level0_idomain%ncols - jcol), dp) * level0_idomain%cellsize + &
1890 0.5_dp * level0_idomain%cellsize
1891 long_1 = level0_idomain%xllcorner + real((irow - 1), dp) * level0_idomain%cellsize + &
1892 0.5_dp * level0_idomain%cellsize
1894 lat_2 = level0_idomain%yllcorner + real((level0_idomain%ncols - jcol_to), dp) * level0_idomain%cellsize + &
1895 0.5_dp * level0_idomain%cellsize
1896 long_2 = level0_idomain%xllcorner + real((irow_to - 1), dp) * level0_idomain%cellsize + &
1897 0.5_dp * level0_idomain%cellsize
1937 use mo_constants,
only : radiusearth_dp, twopi_dp
1942 real(dp),
intent(in) :: lat1, long1, lat2, long2
1945 real(dp),
intent(out) :: distance_out
1966 dtor = twopi_dp / 360.0_dp
1967 theta1 = dtor * long1
1969 theta2 = dtor * long2
1972 term1 = cos(phi1) * cos(theta1) * cos(phi2) * cos(theta2)
1973 term2 = cos(phi1) * sin(theta1) * cos(phi2) * sin(theta2)
1974 term3 = sin(phi1) * sin(phi2)
1975 temp = term1 + term2 + term3
1976 if(temp .GT. 1.0_dp) temp = 1.0_dp
1978 distance_out = radiusearth_dp * acos(temp);
2038 use mo_append,
only : append
2043 integer(i4),
intent(in) :: idomain
2044 real(dp),
dimension(:,:),
allocatable :: facc11
2045 integer(i4) :: ii, jj
2046 integer(i4) :: s11, e11
2047 integer(i4) :: nrows11, ncols11
2048 integer(i4),
dimension(:,:),
allocatable :: fdir11
2049 logical,
dimension(:,:),
allocatable :: mask11
2052 nrows11 =
level11(idomain)%nrows
2053 ncols11 =
level11(idomain)%ncols
2056 mask11 =
level11(idomain)%mask
2059 allocate(fdir11(nrows11, ncols11))
2060 allocate(facc11(nrows11, ncols11))
2074 if (fdir11(ii,jj) .eq. 0)
then
2086 call append(
l11_facc, pack(facc11(:,:),mask11))
2089 deallocate(fdir11, facc11, mask11)
2097 integer(i4),
intent(in) :: fdir(:,:)
2098 real(dp),
intent (inout) :: facc(:,:)
2099 integer(i4),
intent(in) :: ii, jj
2100 integer(i4),
intent(in) :: nrow, ncol
2112 if (jj+1 .le. ncol)
then
2113 if (fdir(ii,jj+1) .eq. 16_i4)
then
2115 facc(ii,jj) = facc(ii,jj) + facc(ii,jj+1)
2119 if ((ii+1 .le. nrow) .and. (jj+1 .le. ncol))
then
2120 if (fdir(ii+1,jj+1) .eq. 32_i4)
then
2122 facc(ii,jj) = facc(ii,jj) + facc(ii+1,jj+1)
2126 if (ii+1 .le. nrow)
then
2127 if (fdir(ii+1,jj) .eq. 64_i4)
then
2129 facc(ii,jj) = facc(ii,jj) + facc(ii+1,jj)
2133 if ((ii+1 .le. nrow) .and. (jj-1 .ge. 1))
then
2134 if (fdir(ii+1,jj-1) .eq. 128_i4)
then
2136 facc(ii,jj) = facc(ii,jj) + facc(ii+1,jj-1)
2140 if (jj-1 .ge. 1)
then
2141 if (fdir(ii,jj-1) .eq. 1_i4)
then
2143 facc(ii,jj) = facc(ii,jj) + facc(ii,jj-1)
2147 if ((ii-1 .ge. 1) .and. (jj-1 .ge. 1))
then
2148 if (fdir(ii-1,jj-1) .eq. 2_i4)
then
2150 facc(ii,jj) = facc(ii,jj) + facc(ii-1,jj-1)
2154 if (ii-1 .ge. 1)
then
2155 if (fdir(ii-1,jj) .eq. 4_i4)
then
2157 facc(ii,jj) = facc(ii,jj) + facc(ii-1,jj)
2161 if ((ii-1 .ge. 1) .and. (jj+1 .le. ncol))
then
2162 if (fdir11(ii-1,jj+1) .eq. 8_i4)
then
2164 facc(ii,jj) = facc(ii,jj) + facc(ii-1,jj+1)
2223 use mo_mad,
only: mad
2224 use mo_append,
only: append
2247 integer(i4),
intent(in) :: idomain
2248 real(dp),
dimension(:),
intent(in) :: param
2251 integer(i4) :: ncells0
2252 integer(i4) :: nnodes
2253 integer(i4) :: nlinks
2254 integer(i4) :: nrows0, ncols0
2255 integer(i4) :: istart0, iend0
2256 integer(i4) :: nrows11, ncols11
2257 integer(i4) :: istart11, iend11
2258 logical,
dimension(:,:),
allocatable :: mask0
2259 integer(i4),
dimension(:,:),
allocatable :: id0
2260 integer(i4),
dimension(:,:),
allocatable :: fdir0
2261 integer(i4),
dimension(:,:),
allocatable :: facc0
2262 real(dp),
dimension(:,:),
allocatable :: slope0
2263 real(dp),
dimension(:),
allocatable :: slope_tmp
2264 real(dp),
dimension(:,:),
allocatable :: cellarea0
2265 integer(i4),
dimension(:),
allocatable :: netperm
2266 integer(i4),
dimension(:),
allocatable :: nlinkfromrow
2267 integer(i4),
dimension(:),
allocatable :: nlinkfromcol
2268 integer(i4),
dimension(:),
allocatable :: nlinktorow
2269 integer(i4),
dimension(:),
allocatable :: nlinktocol
2270 integer(i4) :: ii, rr, ns
2271 integer(i4) :: frow, fcol
2272 integer(i4) :: fid, tid
2273 real(dp),
dimension(:),
allocatable :: stack, append_chunk
2274 integer(i4),
dimension(:),
allocatable :: dummy_1d
2276 real(dp) :: l0_link_slope
2277 real(dp),
dimension(:),
allocatable :: celerity11
2278 real(dp),
dimension(:,:),
allocatable :: celerity0
2280 integer(i4),
dimension(:,:),
allocatable :: nodata_i4_tmp
2281 real(dp),
dimension(:,:),
allocatable :: nodata_dp_tmp
2282 logical,
dimension(:),
allocatable :: slopemask0
2284 type(
grid),
pointer :: level0_idomain
2288 nrows0 = level0_idomain%nrows
2289 ncols0 = level0_idomain%ncols
2290 ncells0 = level0_idomain%ncells
2291 istart0 = level0_idomain%iStart
2292 iend0 = level0_idomain%iEnd
2293 mask0 = level0_idomain%mask
2296 istart11 =
level11(idomain)%iStart
2297 iend11 =
level11(idomain)%iEnd
2298 nrows11 =
level11(idomain)%nrows
2299 ncols11 =
level11(idomain)%ncols
2300 nnodes =
level11(idomain)%ncells
2305 allocate ( id0( nrows0, ncols0 ) )
2306 allocate ( slope0( nrows0, ncols0 ) )
2307 allocate ( fdir0( nrows0, ncols0 ) )
2308 allocate ( facc0( nrows0, ncols0 ) )
2309 allocate ( cellarea0( nrows0, ncols0 ) )
2310 allocate ( celerity0( nrows0, ncols0 ) )
2311 allocate ( slopemask0( ncells0 ) )
2315 allocate ( stack( 1 ) )
2316 allocate ( append_chunk( 1 ) )
2317 allocate ( dummy_1d( 2 ))
2318 allocate ( netperm( nnodes ) )
2319 allocate ( nlinkfromrow( nnodes ) )
2320 allocate ( nlinkfromcol( nnodes ) )
2321 allocate ( nlinktorow( nnodes ) )
2322 allocate ( nlinktocol( nnodes ) )
2323 allocate ( celerity11( nnodes ) )
2324 allocate ( slope_tmp( nnodes ) )
2326 allocate (nodata_i4_tmp( nrows0, ncols0 ) )
2327 allocate (nodata_dp_tmp( nrows0, ncols0 ) )
2345 slopemask0(:) = .false.
2352 if(nnodes .GT. 1)
then
2354 id0(:,:) = unpack(level0_idomain%Id(1:ncells0), mask0, nodata_i4_tmp)
2355 fdir0(:,:) = unpack(
l0_fdir(istart0:iend0), mask0, nodata_i4_tmp)
2356 facc0(:,:) = unpack(
l0_facc(istart0:iend0), mask0, nodata_i4_tmp)
2357 cellarea0(:,:) = unpack(level0_idomain%cellarea(1:ncells0), mask0, nodata_dp_tmp)
2360 slope_tmp =
l0_slope(istart0:iend0)
2361 where ( slope_tmp .lt. 0.1_dp ) slope_tmp = 0.1_dp
2366 if( count(slopemask0) .GT. 1)
then
2367 slope_tmp = mad(arr = slope_tmp, z = 2.25_dp, mask = slopemask0, tout=
'u', mval=0.1_dp)
2369 slope0(:,:) = unpack(slope_tmp, mask0, nodata_dp_tmp )
2373 nlinkfromrow(:) =
l11_frow( istart11 : iend11 )
2374 nlinkfromcol(:) =
l11_fcol( istart11 : iend11 )
2375 nlinktorow(:) =
l11_trow( istart11 : iend11 )
2376 nlinktocol(:) =
l11_tcol( istart11 : iend11 )
2381 frow = nlinkfromrow(ii)
2382 fcol = nlinkfromcol(ii)
2388 fid = id0( frow, fcol )
2389 tid = id0( nlinktorow(ii) , nlinktocol(ii) )
2391 l0_link_slope = slope0(frow, fcol) / 100._dp
2394 stack(ns) = param(1) * sqrt(l0_link_slope)
2396 celerity0(frow, fcol) = stack(ns)
2398 fid = id0(frow, fcol)
2399 if( .NOT. (fid == tid))
then
2400 call append(stack, append_chunk)
2408 celerity11(ii) =
size(stack) / sum(1/stack(:))
2418 celerity11(:) = 1._dp
2424 l0_celerity(istart0:iend0) = pack(celerity0(:,:), mask0)
2428 mask0, id0, slope_tmp, slopemask0, &
2429 slope0, fdir0, cellarea0, &
2430 stack, netperm, nlinkfromrow, nlinkfromcol, nlinktorow, nlinktocol)
recursive subroutine calculate_l11_flow_accumulation(fdir, facc, ii, jj, nrow, ncol)
Provides constants commonly used by mHM, mRM and MPR.
real(dp), parameter, public nodata_dp
integer(i4), parameter, public nodata_i4
Provides common types needed by mHM, mRM and/or mpr.
Provides structures needed by mHM, mRM and/or mpr.
type(domain_meta), public domainmeta
integer(i4), public nlcoverscene
real(dp), dimension(:), allocatable, public l0_elev
integer(i4), public iflag_cordinate_sys
integer(i4), dimension(:, :), allocatable, public l0_lcover
integer(i4), dimension(nprocesses, 3), public processmatrix
type(grid), dimension(:), allocatable, target, public level1
type(grid), dimension(:), allocatable, target, public level0
Global variables for mpr only.
real(dp), dimension(:), allocatable, public l0_slope
Provides mRM specific constants.
real(dp), parameter, public deltah
Global variables for mRM only.
type(gridremapper), dimension(:), allocatable, public l0_l11_remap
integer(i4), dimension(:), allocatable, public l0_floodplain
real(dp), dimension(:), allocatable, public l11_facc
integer(i4), dimension(:), allocatable, public l11_netperm
real(dp), dimension(:, :), allocatable, public l11_nlinkfracfpimp
integer(i4), dimension(:), allocatable, public l11_l1_id
integer(i4), dimension(:), allocatable, public l1_l11_id
real(dp), dimension(:), allocatable, public l0_celerity
real(dp), dimension(:), allocatable, target, public l11_afloodplain
integer(i4), dimension(:), allocatable, public l0_fdir
integer(i4), dimension(:), allocatable, public l11_fcol
type(sink_cells_t), dimension(:), allocatable sink_cells
sink cell ids for each domain
integer(i4), dimension(:), allocatable, public l0_dracell
integer(i4), dimension(:), allocatable, public l11_label
integer(i4), dimension(:), allocatable, public l0_streamnet
integer(i4), dimension(:), allocatable, public l11_fromn
real(dp), dimension(:), allocatable, public l11_length
integer(i4), dimension(:), allocatable, public l11_rowout
integer(i4), dimension(:), allocatable, public l11_ton
type(domaininfo_mrm), dimension(:), allocatable, target, public domain_mrm
integer(i4), dimension(:), allocatable, public l11_trow
type(grid), dimension(:), allocatable, target, public level11
integer(i4), dimension(:), allocatable, public l11_fdir
integer(i4), dimension(:), allocatable, public l0_facc
integer(i4), dimension(:), allocatable, public l11_tcol
integer(i4), dimension(:), allocatable, public l11_frow
real(dp), dimension(:), allocatable, public l11_slope
logical, dimension(:), allocatable, public l11_sink
integer(i4), dimension(:), allocatable, public l0_inflowgaugeloc
real(dp), dimension(:), allocatable, public l11_celerity
integer(i4), dimension(:), allocatable, public l11_noutlets
integer(i4), dimension(:), allocatable, public l0_drasc
integer(i4), dimension(:), allocatable, public l11_rorder
integer(i4), dimension(:), allocatable, public l0_gaugeloc
integer(i4), dimension(:), allocatable, public l11_colout
Startup drainage network for mHM.
subroutine, public l11_stream_features(idomain)
Stream features (stream network and floodplain)
subroutine, public l11_flow_direction(idomain)
Determine the flow direction of the upscaled river network at level L11.
subroutine celllength(idomain, fdir, irow, jcol, icoorsystem, length)
TODO: add description.
subroutine moveup(elev0, fdir0, fi, fj, ss, nn)
TODO: add description.
subroutine, public get_distance_two_lat_lon_points(lat1, long1, lat2, long2, distance_out)
estimate distance in [m] between two points in a lat-lon
subroutine, public l11_fraction_sealed_floodplain(lcclassimp, do_init)
Fraction of the flood plain with impervious cover.
subroutine, public l11_set_drain_outlet_gauges(idomain)
Draining cell identification and Set gauging node.
subroutine, public l11_set_network_topology(idomain)
Set network topology.
subroutine, public l11_routing_order(idomain)
Find routing order, headwater cells and sink.
subroutine movedownonecell(fdir, irow, jcol)
TODO: add description.
subroutine, public l11_link_location(idomain)
Estimate the LO (row,col) location for each routing link at level L11.
subroutine, public l11_calc_celerity(idomain, param)
L11 celerity based on L0 elevation and L0 fAcc.
subroutine, public l11_l1_mapping(idomain)
TODO: add description.
subroutine, public l11_flow_accumulation(idomain)
Calculates L11 flow accumulation per grid cell.