24 USE mo_kind,
ONLY : i4, dp
105 read_states, processCase, global_routing_param, L1_total_runoff, L1_areaCell, L1_L11_Id, &
106 L11_areaCell, L11_L1_Id, L11_netPerm, L11_fromN, L11_toN, L11_nOutlets, timestep, tsRoutFactor, &
107 nNodes, nInflowGauges, InflowGaugeIndexList, InflowGaugeHeadwater, InflowGaugeNodeList, &
108 InflowDischarge, nGauges, gaugeIndexList, gaugeNodeList, map_flag, L11_length, L11_slope, &
109 L11_FracFPimp, L11_C1, L11_C2, L11_qOut, L11_qTIN, L11_qTR, L11_qMod, GaugeDischarge &
112 use mo_constants,
only : t0_dp
121 logical,
intent(in) :: read_states
123 integer(i4),
intent(in) :: processcase
125 real(dp),
dimension(:),
intent(in) :: global_routing_param
127 real(dp),
dimension(:),
intent(in) :: l1_total_runoff
129 real(dp),
dimension(:),
intent(in) :: l1_areacell
131 integer(i4),
dimension(:),
intent(in) ::
l1_l11_id
135 integer(i4),
dimension(:),
intent(in) ::
l11_l1_id
137 integer(i4),
dimension(:),
intent(in) ::
l11_netperm
139 integer(i4),
dimension(:),
intent(in) ::
l11_fromn
141 integer(i4),
dimension(:),
intent(in) ::
l11_ton
145 integer(i4),
intent(in) :: timestep
147 real(dp),
intent(in) :: tsroutfactor
149 integer(i4),
intent(in) :: nnodes
151 integer(i4),
intent(in) :: ninflowgauges
153 integer(i4),
dimension(:),
intent(in) :: inflowgaugeindexlist
155 logical,
dimension(:),
intent(in) :: inflowgaugeheadwater
157 integer(i4),
dimension(:),
intent(in) :: inflowgaugenodelist
159 real(dp),
dimension(:),
intent(in) :: inflowdischarge
161 integer(i4),
intent(in) :: ngauges
163 integer(i4),
dimension(:),
intent(in) :: gaugeindexlist
165 integer(i4),
dimension(:),
intent(in) :: gaugenodelist
167 logical,
intent(in) :: map_flag
169 real(dp),
dimension(:),
intent(in) ::
l11_length
171 real(dp),
dimension(:),
intent(in) ::
l11_slope
173 real(dp),
dimension(:),
intent(in) :: l11_fracfpimp
175 real(dp),
dimension(:),
intent(inout) ::
l11_c1
177 real(dp),
dimension(:),
intent(inout) ::
l11_c2
179 real(dp),
dimension(:),
intent(inout) ::
l11_qout
181 real(dp),
dimension(:, :),
intent(inout) ::
l11_qtin
183 real(dp),
dimension(:, :),
intent(inout) ::
l11_qtr
185 real(dp),
dimension(:),
intent(inout) ::
l11_qmod
187 real(dp),
dimension(:),
intent(inout) :: gaugedischarge
189 integer(i4) :: s11, e11
193 integer(i4) :: rout_loop
195 real(dp),
dimension(size(L11_qMod, dim = 1)) :: l11_qacc
196 real(dp),
dimension(:),
allocatable ::
l11_e_acc
211 if (processcase .eq. 1_i4 .AND. (.not. read_states))
then
213 if (nnodes .GT. 1)
then
214 call reg_rout(global_routing_param, &
224 rout_loop = max(1_i4, nint(1._dp / tsroutfactor))
228 l1_total_runoff, l1_areacell,
l1_l11_id, &
236 inflowgaugeindexlist, &
237 inflowgaugeheadwater, &
238 inflowgaugenodelist, &
243 if(nnodes .GT. 1)
then
258 inflowgaugeheadwater, &
259 inflowgaugenodelist, &
276 inflowgaugeheadwater, &
277 inflowgaugenodelist, &
284 l11_qmod = l11_qacc / real(rout_loop, dp)
300 gaugedischarge(gaugeindexlist(gg)) =
l11_qmod(gaugenodelist(gg))
380 subroutine l11_routing(nNodes, nLinks, netPerm, netLink_fromN, netLink_toN, netLink_C1, netLink_C2, netNode_qOUT, &
381 nInflowGauges, InflowHeadwater, InflowNodeList, netNode_qTIN, netNode_qTR, netNode_Qmod)
385 integer(i4),
intent(in) :: nNodes
387 integer(i4),
intent(in) :: nLinks
389 integer(i4),
dimension(:),
intent(in) :: netPerm
391 integer(i4),
dimension(:),
intent(in) :: netLink_fromN
393 integer(i4),
dimension(:),
intent(in) :: netLink_toN
395 real(dp),
dimension(:),
intent(in) :: netLink_C1
397 real(dp),
dimension(:),
intent(in) :: netLink_C2
399 real(dp),
dimension(:),
intent(in) :: netNode_qOUT
401 integer(i4),
intent(in) :: nInflowGauges
403 logical,
dimension(:),
intent(in) :: InflowHeadwater
405 integer(i4),
dimension(:),
intent(in) :: InflowNodeList
407 real(dp),
dimension(:, :),
intent(inout) :: netNode_qTIN
409 real(dp),
dimension(:, :),
intent(inout) :: netNode_qTR
411 real(dp),
dimension(nNodes),
intent(out) :: netNode_Qmod
413 integer(i4) :: g, i, k, iNode, tNode
415 integer(i4),
parameter :: IT = 2
417 integer(i4),
parameter :: IT1 = 1
428 netnode_qtin(:, it) = 0.0_dp
438 inode = netlink_fromn(i)
439 tnode = netlink_ton(i)
442 netnode_qtin(inode, it) = netnode_qtin(inode, it) + netnode_qout(inode)
445 netnode_qtr(inode, it) = netnode_qtr(inode, it1) &
446 + netlink_c1(i) * (netnode_qtin(inode, it1) - netnode_qtr(inode, it1)) &
447 + netlink_c2(i) * (netnode_qtin(inode, it) - netnode_qtin(inode, it1))
450 if (ninflowgauges .GT. 0)
then
451 do g = 1, ninflowgauges
453 if ((tnode == inflownodelist(g)) .AND. (.NOT. inflowheadwater(g))) netnode_qtr(inode, it) = 0.0_dp
458 netnode_qtin(tnode, it) = netnode_qtin(tnode, it) + netnode_qtr(inode, it)
466 tnode = netlink_ton(netperm(nlinks))
467 netnode_qtin(tnode, it) = netnode_qtin(tnode, it) + netnode_qout(tnode)
475 netnode_qmod(1 : nnodes) = netnode_qtin(1 : nnodes, it)
477 netnode_qtr(1 : nnodes, it1) = netnode_qtr(1 : nnodes, it)
478 netnode_qtin(1 : nnodes, it1) = netnode_qtin(1 : nnodes, it)
Global variables for mRM only.
real(dp), dimension(:, :), allocatable, public l11_qtin
integer(i4), dimension(:), allocatable, public l11_netperm
real(dp), dimension(:), allocatable, public l11_qout
integer(i4), dimension(:), allocatable, public l11_l1_id
integer(i4), dimension(:), allocatable, public l1_l11_id
type(riv_temp_type), public riv_temp_pcs
This is a container for the river temperature routing process (pcs)
integer(i4), dimension(:), allocatable, public l11_fromn
real(dp), dimension(:), allocatable, public l11_length
real(dp), dimension(:), allocatable, public l11_qmod
integer(i4), dimension(:), allocatable, public l11_ton
real(dp), dimension(:), allocatable, public l11_c1
real(dp), dimension(:), allocatable, public l11_slope
integer(i4), dimension(:), allocatable, public l11_noutlets
real(dp), dimension(:, :), allocatable, public l11_qtr
real(dp), dimension(:), allocatable, public l11_areacell
real(dp), dimension(:), allocatable, public l11_c2
Perform Multiscale Parameter Regionalization on Routing Parameters.
subroutine, public reg_rout(param, length, slope, ffpimp, ts, c1, c2)
Regionalized routing.
Performs pre-processing for routing for mHM at level L11.
subroutine l11_e_acc(qall, efecarea, l1_l11_id, l11_areacell, l11_l1_id, ts, map_flag, qacc)
temperature energy accumulation at L11.
subroutine, public add_inflow(ninflowgauges, inflowindexlist, inflowheadwater, inflownodelist, qinflow, qout)
Adds inflow discharge to the runoff produced at the cell where the inflow is occurring.
subroutine, public l11_runoff_acc(qall, efecarea, l1_l11_id, l11_areacell, l11_l1_id, ts, map_flag, qacc)
total runoff accumulation at L11.
Performs runoff routing for mHM at level L11.
subroutine, public mrm_routing(read_states, processcase, global_routing_param, l1_total_runoff, l1_areacell, l1_l11_id, l11_areacell, l11_l1_id, l11_netperm, l11_fromn, l11_ton, l11_noutlets, timestep, tsroutfactor, nnodes, ninflowgauges, inflowgaugeindexlist, inflowgaugeheadwater, inflowgaugenodelist, inflowdischarge, ngauges, gaugeindexlist, gaugenodelist, map_flag, l11_length, l11_slope, l11_fracfpimp, l11_c1, l11_c2, l11_qout, l11_qtin, l11_qtr, l11_qmod, gaugedischarge)
route water given runoff
subroutine l11_routing(nnodes, nlinks, netperm, netlink_fromn, netlink_ton, netlink_c1, netlink_c2, netnode_qout, ninflowgauges, inflowheadwater, inflownodelist, netnode_qtin, netnode_qtr, netnode_qmod)
Performs runoff routing for mHM at L11 upscaled network (Routing Network).