36 integer(I4B),
intent(in) :: npts
37 real(dp),
dimension(npts),
intent(in) :: stations
43 w = stations(npts) - stations(1)
58 integer(I4B),
intent(in) :: npts
59 real(dp),
dimension(npts),
intent(in) :: stations
60 real(dp),
dimension(npts),
intent(in) :: heights
61 real(dp),
intent(in) :: d
65 real(dp),
dimension(npts - 1) :: widths
87 integer(I4B),
intent(in) :: n
88 integer(I4B),
intent(in) :: npts
89 real(dp),
dimension(npts),
intent(in) :: heights
90 real(dp),
intent(in) :: d
91 logical,
intent(in) :: leftface
97 if (heights(n - 1) > d)
then
99 else if (heights(n - 1) > heights(n))
then
100 vwf = heights(n - 1) - heights(n)
103 if (heights(n + 2) > d)
then
104 vwf = d - heights(n + 1)
105 else if (heights(n + 2) > heights(n + 1))
then
106 vwf = heights(n + 2) - heights(n + 1)
120 integer(I4B),
intent(in) :: npts
121 real(dp),
dimension(npts),
intent(in) :: stations
122 real(dp),
dimension(npts),
intent(in) :: heights
123 real(dp),
intent(in) :: d
127 real(dp),
dimension(npts - 1) :: perimeters
137 p = p + perimeters(n)
150 integer(I4B),
intent(in) :: npts
151 real(dp),
dimension(npts),
intent(in) :: stations
152 real(dp),
dimension(npts),
intent(in) :: heights
153 real(dp),
intent(in) :: d
157 real(dp),
dimension(npts - 1) :: areas
180 integer(I4B),
intent(in) :: npts
181 real(dp),
dimension(npts),
intent(in) :: stations
182 real(dp),
dimension(npts),
intent(in) :: heights
183 real(dp),
intent(in) :: d
189 real(dp),
dimension(npts - 1) :: areas
190 real(dp),
dimension(npts - 1) :: perimeters
202 p = p + perimeters(n)
230 roughness, conv_fact, slope, d)
result(q)
232 integer(I4B),
intent(in) :: npts
233 real(dp),
dimension(npts),
intent(in) :: stations
234 real(dp),
dimension(npts),
intent(in) :: heights
235 real(dp),
dimension(npts),
intent(in) :: roughfracs
236 real(dp),
intent(in) :: roughness
237 real(dp),
intent(in) :: conv_fact
238 real(dp),
intent(in) :: slope
239 real(dp),
intent(in) :: d
242 real(dp) :: conveyance
248 real(dp),
dimension(npts - 1) :: areas
249 real(dp),
dimension(npts - 1) :: perimeters
264 p = p + perimeters(n)
276 r = roughness * roughfracs(n)
277 if (p * r >
dzero)
then
280 conveyance = conveyance + a * rh**
dtwothirds / r
283 q = conv_fact * conveyance * sqrt(slope)
297 integer(I4B),
intent(in) :: npts
298 real(DP),
dimension(npts),
intent(in) :: stations
299 real(DP),
dimension(npts),
intent(in) :: heights
300 logical,
dimension(npts - 1),
intent(inout) :: leftv
301 logical,
dimension(npts - 1),
intent(inout) :: rightv
313 if (stations(n - 1) == stations(n) .and. heights(n - 1) > heights(n))
then
318 if (n < npts - 1)
then
319 if (stations(n + 2) == stations(n + 1) .and. &
320 heights(n + 2) > heights(n + 1))
then
336 integer(I4B),
intent(in) :: npts
337 real(DP),
dimension(npts),
intent(in) :: stations
338 real(DP),
dimension(npts),
intent(in) :: heights
339 real(DP),
intent(in) :: d
340 real(DP),
dimension(npts - 1),
intent(inout) :: p
351 logical,
dimension(npts - 1) :: leftv, rightv
375 if (xlen >
dzero)
then
383 dlen = min(d, dmax) - dmin
388 p(n) = sqrt(xlen**
dtwo + dlen**
dtwo)
400 if (n < npts - 1)
then
417 integer(I4B),
intent(in) :: npts
418 real(DP),
dimension(npts),
intent(in) :: stations
419 real(DP),
dimension(npts),
intent(in) :: heights
420 real(DP),
intent(in) :: d
421 real(DP),
dimension(npts - 1),
intent(inout) :: a
449 if (xlen >
dzero)
then
453 a(n) = xlen * (d - dmax)
457 if (dmax /= dmin .and. d > dmin)
then
459 a(n) = a(n) +
dhalf * (d - dmin) * xlen
461 a(n) = a(n) +
dhalf * (dmax - dmin) * xlen
477 integer(I4B),
intent(in) :: npts
478 real(DP),
dimension(npts),
intent(in) :: stations
479 real(DP),
dimension(npts),
intent(in) :: heights
480 real(DP),
intent(in) :: d
481 real(DP),
dimension(npts - 1),
intent(inout) :: w
521 real(dp),
intent(inout) :: x0
522 real(dp),
intent(inout) :: x1
523 real(dp),
intent(in) :: d0
524 real(dp),
intent(in) :: d1
525 real(dp),
intent(inout) :: dmax
526 real(dp),
intent(inout) :: dmin
527 real(dp),
intent(in) :: d
547 else if (d < dmax)
then
550 if (abs(dlen) >
dzero)
then
556 dx = (d - d1) * slope
561 dx = (d - d0) * slope
This module contains simulation constants.
real(dp), parameter dtwothirds
real constant 2/3
real(dp), parameter dhalf
real constant 1/2
real(dp), parameter dzero
real constant zero
real(dp), parameter dtwo
real constant 2
real(dp), parameter done
real constant 1
This module contains stateless sfr subroutines and functions.
real(dp) function, public get_wetted_topwidth(npts, stations, heights, d)
Calculate the wetted top width for a reach.
real(dp) function get_wet_vert_face(n, npts, heights, d, leftface)
Calculate wetted vertical height.
subroutine get_wetted_topwidths(npts, stations, heights, d, w)
Calculate the wetted top widths for each line segment.
subroutine get_cross_section_areas(npts, stations, heights, d, a)
Calculate the cross-sectional areas for each line segment.
subroutine get_wetted_perimeters(npts, stations, heights, d, p)
Calculate the wetted perimeters for each line segment.
pure subroutine get_wetted_station(x0, x1, d0, d1, dmax, dmin, d)
Calculate the station values for the wetted portion of the cross-section.
real(dp) function, public get_wetted_perimeter(npts, stations, heights, d)
Calculate the wetted perimeter for a reach.
real(dp) function, public get_cross_section_area(npts, stations, heights, d)
Calculate the cross-sectional area for a reach.
real(dp) function, public get_saturated_topwidth(npts, stations)
Calculate the saturated top width for a reach.
subroutine determine_vert_neighbors(npts, stations, heights, leftv, rightv)
Determine vertical segments.
real(dp) function, public get_mannings_section(npts, stations, heights, roughfracs, roughness, conv_fact, slope, d)
Calculate the manning's discharge for a reach.
real(dp) function, public get_hydraulic_radius(npts, stations, heights, d)
Calculate the hydraulic radius for a reach.
This module defines variable data types.