MODFLOW 6  version 6.7.0.dev0
USGS Modular Hydrologic Model
SfrCrossSectionUtils.f90
Go to the documentation of this file.
1 !> @brief This module contains stateless sfr subroutines and functions
2 !!
3 !! This module contains the functions to calculate the wetted perimeter
4 !! and cross-sectional area for a reach cross-section that are used in
5 !! the streamflow routing (SFR) package. It also contains subroutines to
6 !! calculate the wetted perimeter and cross-sectional area for each
7 !! line segment in the cross-section. This module does not depend on the
8 !! SFR package.
9 !!
10 !<
12 
13  use kindmodule, only: dp, i4b
15 
16  implicit none
17  private
18  public :: get_saturated_topwidth
19  public :: get_wetted_topwidth
20  public :: get_wetted_perimeter
21  public :: get_cross_section_area
22  public :: get_hydraulic_radius
23  public :: get_mannings_section
24 
25 contains
26 
27  !> @brief Calculate the saturated top width for a reach
28  !!
29  !! Function to calculate the maximum top width for a reach using the
30  !! cross-section station data.
31  !!
32  !! @return w saturated top width
33  !<
34  function get_saturated_topwidth(npts, stations) result(w)
35  ! -- dummy variables
36  integer(I4B), intent(in) :: npts !< number of station-height data for a reach
37  real(dp), dimension(npts), intent(in) :: stations !< cross-section station distances (x-distance)
38  ! -- local variables
39  real(dp) :: w
40  !
41  ! -- calculate the saturated top width
42  if (npts > 1) then
43  w = stations(npts) - stations(1)
44  else
45  w = stations(1)
46  end if
47  end function get_saturated_topwidth
48 
49  !> @brief Calculate the wetted top width for a reach
50  !!
51  !! Function to calculate the wetted top width for a reach using the
52  !! cross-section station-height data given a passed depth.
53  !!
54  !! @return w wetted top width
55  !<
56  function get_wetted_topwidth(npts, stations, heights, d) result(w)
57  ! -- dummy variables
58  integer(I4B), intent(in) :: npts !< number of station-height data for a reach
59  real(dp), dimension(npts), intent(in) :: stations !< cross-section station distances (x-distance)
60  real(dp), dimension(npts), intent(in) :: heights !< cross-section height data
61  real(dp), intent(in) :: d !< depth to evaluate cross-section
62  ! -- local variables
63  integer(I4B) :: n
64  real(dp) :: w
65  real(dp), dimension(npts - 1) :: widths
66  !
67  ! -- initialize the wetted perimeter for the reach
68  w = dzero
69  !
70  ! -- calculate the wetted top width for each line segment
71  call get_wetted_topwidths(npts, stations, heights, d, widths)
72  !
73  ! -- calculate the wetted top widths
74  do n = 1, npts - 1
75  w = w + widths(n)
76  end do
77  end function get_wetted_topwidth
78 
79  !> @brief Calculate wetted vertical height
80  !!
81  !! For segments flanked by vertically-oriented neighboring segments,
82  !! return the length of the submerged, vertically-oriented, neighboring face
83  !!
84  !<
85  function get_wet_vert_face(n, npts, heights, d, leftface) result(vwf)
86  ! -- dummy
87  integer(I4B), intent(in) :: n !< index to be evaluated
88  integer(I4B), intent(in) :: npts !< length of heights vector
89  real(dp), dimension(npts), intent(in) :: heights !< cross-section height data
90  real(dp), intent(in) :: d
91  logical, intent(in) :: leftface
92  ! -- local
93  real(dp) :: vwf !< vertically wetted face length
94  !
95  ! -- calculate the vertically-oriented wetted face length
96  if (leftface) then
97  if (heights(n - 1) > d) then
98  vwf = d - heights(n)
99  else if (heights(n - 1) > heights(n)) then
100  vwf = heights(n - 1) - heights(n)
101  end if
102  else
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)
107  end if
108  end if
109  end function get_wet_vert_face
110 
111  !> @brief Calculate the wetted perimeter for a reach
112  !!
113  !! Function to calculate the wetted perimeter for a reach using the
114  !! cross-section station-height data given a passed depth.
115  !!
116  !! @return p wetted perimeter
117  !<
118  function get_wetted_perimeter(npts, stations, heights, d) result(p)
119  ! -- dummy variables
120  integer(I4B), intent(in) :: npts !< number of station-height data for a reach
121  real(dp), dimension(npts), intent(in) :: stations !< cross-section station distances (x-distance)
122  real(dp), dimension(npts), intent(in) :: heights !< cross-section height data
123  real(dp), intent(in) :: d !< depth to evaluate cross-section
124  ! -- local variables
125  integer(I4B) :: n
126  real(dp) :: p
127  real(dp), dimension(npts - 1) :: perimeters
128  !
129  ! -- initialize the wetted perimeter for the reach
130  p = dzero
131  !
132  ! -- calculate the wetted perimeter for each line segment
133  call get_wetted_perimeters(npts, stations, heights, d, perimeters)
134  !
135  ! -- calculate the wetted perimenter
136  do n = 1, npts - 1
137  p = p + perimeters(n)
138  end do
139  end function get_wetted_perimeter
140 
141  !> @brief Calculate the cross-sectional area for a reach
142  !!
143  !! Function to calculate the cross-sectional area for a reach using
144  !! the cross-section station-height data given a passed depth.
145  !!
146  !! @return a cross-sectional area
147  !<
148  function get_cross_section_area(npts, stations, heights, d) result(a)
149  ! -- dummy variables
150  integer(I4B), intent(in) :: npts !< number of station-height data for a reach
151  real(dp), dimension(npts), intent(in) :: stations !< cross-section station distances (x-distance)
152  real(dp), dimension(npts), intent(in) :: heights !< cross-section height data
153  real(dp), intent(in) :: d !< depth to evaluate cross-section
154  ! -- local variables
155  integer(I4B) :: n
156  real(dp) :: a
157  real(dp), dimension(npts - 1) :: areas
158  !
159  ! -- initialize the area
160  a = dzero
161  !
162  ! -- calculate the cross-sectional area for each line segment
163  call get_cross_section_areas(npts, stations, heights, d, areas)
164  !
165  ! -- calculate the cross-sectional area
166  do n = 1, npts - 1
167  a = a + areas(n)
168  end do
169  end function get_cross_section_area
170 
171  !> @brief Calculate the hydraulic radius for a reach
172  !!
173  !! Function to calculate the hydraulic radius for a reach using
174  !! the cross-section station-height data given a passed depth.
175  !!
176  !! @return r hydraulic radius
177  !<
178  function get_hydraulic_radius(npts, stations, heights, d) result(r)
179  ! -- dummy variables
180  integer(I4B), intent(in) :: npts !< number of station-height data for a reach
181  real(dp), dimension(npts), intent(in) :: stations !< cross-section station distances (x-distance)
182  real(dp), dimension(npts), intent(in) :: heights !< cross-section height data
183  real(dp), intent(in) :: d !< depth to evaluate cross-section
184  ! -- local variables
185  integer(I4B) :: n
186  real(dp) :: r
187  real(dp) :: p
188  real(dp) :: a
189  real(dp), dimension(npts - 1) :: areas
190  real(dp), dimension(npts - 1) :: perimeters
191  !
192  ! -- initialize the hydraulic radius, perimeter, and area
193  r = dzero
194  p = dzero
195  a = dzero
196  !
197  ! -- calculate the wetted perimeter for each line segment
198  call get_wetted_perimeters(npts, stations, heights, d, perimeters)
199  !
200  ! -- calculate the wetted perimenter
201  do n = 1, npts - 1
202  p = p + perimeters(n)
203  end do
204  !
205  ! -- calculate the hydraulic radius only if the perimeter is non-zero
206  if (p > dzero) then
207  !
208  ! -- calculate the cross-sectional area for each line segment
209  call get_cross_section_areas(npts, stations, heights, d, areas)
210  !
211  ! -- calculate the cross-sectional area
212  do n = 1, npts - 1
213  a = a + areas(n)
214  end do
215  !
216  ! -- calculate the hydraulic radius
217  r = a / p
218  end if
219  end function get_hydraulic_radius
220 
221  !> @brief Calculate the manning's discharge for a reach
222  !!
223  !! Function to calculate the mannings discharge for a reach
224  !! by calculating the discharge for each section, which can
225  !! have a unique Manning's coefficient given a passed depth.
226  !!
227  !! @return q reach discharge
228  !<
229  function get_mannings_section(npts, stations, heights, roughfracs, &
230  roughness, conv_fact, slope, d) result(q)
231  ! -- dummy variables
232  integer(I4B), intent(in) :: npts !< number of station-height data for a reach
233  real(dp), dimension(npts), intent(in) :: stations !< cross-section station distances (x-distance)
234  real(dp), dimension(npts), intent(in) :: heights !< cross-section height data
235  real(dp), dimension(npts), intent(in) :: roughfracs !< cross-section Mannings roughness fraction data
236  real(dp), intent(in) :: roughness !< base reach roughness
237  real(dp), intent(in) :: conv_fact !< unit conversion factor
238  real(dp), intent(in) :: slope !< reach slope
239  real(dp), intent(in) :: d !< depth to evaluate cross-section
240  ! -- local variables
241  integer(I4B) :: n
242  real(dp) :: conveyance
243  real(dp) :: q
244  real(dp) :: rh
245  real(dp) :: r
246  real(dp) :: p
247  real(dp) :: a
248  real(dp), dimension(npts - 1) :: areas
249  real(dp), dimension(npts - 1) :: perimeters
250  !
251  ! -- initialize the conveyance, hydraulic radius, perimeter, and area
252  conveyance = dzero
253  q = dzero
254  rh = dzero
255  r = dzero
256  p = dzero
257  a = dzero
258  !
259  ! -- calculate the wetted perimeter for each line segment
260  call get_wetted_perimeters(npts, stations, heights, d, perimeters)
261  !
262  ! -- calculate the wetted perimenter
263  do n = 1, npts - 1
264  p = p + perimeters(n)
265  end do
266  !
267  ! -- calculate the hydraulic radius only if the perimeter is non-zero
268  if (p > dzero) then
269  !
270  ! -- calculate the cross-sectional area for each line segment
271  call get_cross_section_areas(npts, stations, heights, d, areas)
272  !
273  ! -- calculate the cross-sectional area
274  do n = 1, npts - 1
275  p = perimeters(n)
276  r = roughness * roughfracs(n)
277  if (p * r > dzero) then
278  a = areas(n)
279  rh = a / p
280  conveyance = conveyance + a * rh**dtwothirds / r
281  end if
282  end do
283  q = conv_fact * conveyance * sqrt(slope)
284  end if
285  end function get_mannings_section
286 
287  ! -- private functions and subroutines
288 
289  !> @brief Determine vertical segments
290  !!
291  !! Subroutine to cycle through each segment (npts - 1) and determine
292  !! whether neighboring segments are vertically-oriented.
293  !!
294  !<
295  subroutine determine_vert_neighbors(npts, stations, heights, leftv, rightv)
296  ! -- dummy
297  integer(I4B), intent(in) :: npts !< number of station-height data for a reach
298  real(DP), dimension(npts), intent(in) :: stations !< cross-section station distances (x-distance)
299  real(DP), dimension(npts), intent(in) :: heights !< cross-section height data
300  logical, dimension(npts - 1), intent(inout) :: leftv
301  logical, dimension(npts - 1), intent(inout) :: rightv
302  ! -- local
303  integer(I4B) :: n
304  !
305  ! -- default neighboring segments to false unless determined otherwise
306  ! o 2 pt x-section has 1 segment (no neighbors to eval)
307  ! o 3+ pt x-section has at the very least one neighbor to eval
308  do n = 1, npts - 1
309  leftv(n) = .false.
310  rightv(n) = .false.
311  ! -- left neighbor
312  if (n > 1) then
313  if (stations(n - 1) == stations(n) .and. heights(n - 1) > heights(n)) then
314  leftv(n) = .true.
315  end if
316  end if
317  ! -- right neighbor
318  if (n < npts - 1) then
319  if (stations(n + 2) == stations(n + 1) .and. &
320  heights(n + 2) > heights(n + 1)) then
321  rightv(n) = .true.
322  end if
323  end if
324  end do
325  end subroutine determine_vert_neighbors
326 
327  !> @brief Calculate the wetted perimeters for each line segment
328  !!
329  !! Subroutine to calculate the wetted perimeter for each line segment
330  !! that defines the reach using the cross-section station-height
331  !! data given a passed depth.
332  !!
333  !<
334  subroutine get_wetted_perimeters(npts, stations, heights, d, p)
335  ! -- dummy variables
336  integer(I4B), intent(in) :: npts !< number of station-height data for a reach
337  real(DP), dimension(npts), intent(in) :: stations !< cross-section station distances (x-distance)
338  real(DP), dimension(npts), intent(in) :: heights !< cross-section height data
339  real(DP), intent(in) :: d !< depth to evaluate cross-section
340  real(DP), dimension(npts - 1), intent(inout) :: p !< wetted perimeter for each line segment
341  ! -- local variables
342  integer(I4B) :: n
343  real(DP) :: x0
344  real(DP) :: x1
345  real(DP) :: d0
346  real(DP) :: d1
347  real(DP) :: dmax
348  real(DP) :: dmin
349  real(DP) :: xlen
350  real(DP) :: dlen
351  logical, dimension(npts - 1) :: leftv, rightv
352  !
353  ! -- set neighbor status
354  call determine_vert_neighbors(npts, stations, heights, leftv, rightv)
355  !
356  ! -- iterate over the station-height data
357  do n = 1, npts - 1
358  !
359  ! -- initialize the wetted perimeter
360  p(n) = dzero
361  !
362  ! -- initialize station-height data for segment
363  x0 = stations(n)
364  x1 = stations(n + 1)
365  d0 = heights(n)
366  d1 = heights(n + 1)
367  !
368  ! -- get the start and end station position of the wetted segment
369  call get_wetted_station(x0, x1, d0, d1, dmax, dmin, d)
370  !
371  ! -- calculate the wetted perimeter for the segment
372  ! - bottom wetted length
373  xlen = x1 - x0
374  dlen = dzero
375  if (xlen > dzero) then
376  if (d > dmax) then
377  dlen = dmax - dmin
378  else
379  dlen = d - dmin
380  end if
381  else
382  if (d > dmin) then
383  dlen = min(d, dmax) - dmin
384  else
385  dlen = dzero
386  end if
387  end if
388  p(n) = sqrt(xlen**dtwo + dlen**dtwo)
389  !
390  ! -- if neighboring segments are vertical, account for their
391  ! contribution to wetted perimeter
392  !
393  ! left neighbor (if applicable)
394  if (n > 1) then
395  if (leftv(n)) then
396  p(n) = p(n) + get_wet_vert_face(n, npts, heights, d, .true.)
397  end if
398  end if
399  ! right neighbor (if applicable)
400  if (n < npts - 1) then
401  if (rightv(n)) then
402  p(n) = p(n) + get_wet_vert_face(n, npts, heights, d, .false.)
403  end if
404  end if
405  end do
406  end subroutine get_wetted_perimeters
407 
408  !> @brief Calculate the cross-sectional areas for each line segment
409  !!
410  !! Subroutine to calculate the cross-sectional area for each line segment
411  !! that defines the reach using the cross-section station-height
412  !! data given a passed depth.
413  !!
414  !<
415  subroutine get_cross_section_areas(npts, stations, heights, d, a)
416  ! -- dummy variables
417  integer(I4B), intent(in) :: npts !< number of station-height data for a reach
418  real(DP), dimension(npts), intent(in) :: stations !< cross-section station distances (x-distance)
419  real(DP), dimension(npts), intent(in) :: heights !< cross-section height data
420  real(DP), intent(in) :: d !< depth to evaluate cross-section
421  real(DP), dimension(npts - 1), intent(inout) :: a !< cross-sectional area for each line segment
422  ! -- local variables
423  integer(I4B) :: n
424  real(DP) :: x0
425  real(DP) :: x1
426  real(DP) :: d0
427  real(DP) :: d1
428  real(DP) :: dmax
429  real(DP) :: dmin
430  real(DP) :: xlen
431  !
432  ! -- iterate over the station-height data
433  do n = 1, npts - 1
434  !
435  ! -- initialize the cross-sectional area
436  a(n) = dzero
437  !
438  ! -- initialize station-height data for segment
439  x0 = stations(n)
440  x1 = stations(n + 1)
441  d0 = heights(n)
442  d1 = heights(n + 1)
443  !
444  ! -- get the start and end station position of the wetted segment
445  call get_wetted_station(x0, x1, d0, d1, dmax, dmin, d)
446  !
447  ! -- calculate the cross-sectional area for the segment
448  xlen = x1 - x0
449  if (xlen > dzero) then
450  !
451  ! -- add the area above dmax
452  if (d > dmax) then
453  a(n) = xlen * (d - dmax)
454  end if
455  !
456  ! -- add the area below dmax
457  if (dmax /= dmin .and. d > dmin) then
458  if (d < dmax) then
459  a(n) = a(n) + dhalf * (d - dmin) * xlen
460  else
461  a(n) = a(n) + dhalf * (dmax - dmin) * xlen
462  end if
463  end if
464  end if
465  end do
466  end subroutine get_cross_section_areas
467 
468  !> @brief Calculate the wetted top widths for each line segment
469  !!
470  !! Subroutine to calculate the wetted top width for each line segment
471  !! that defines the reach using the cross-section station-height
472  !! data given a passed depth.
473  !!
474  !<
475  subroutine get_wetted_topwidths(npts, stations, heights, d, w)
476  ! -- dummy variables
477  integer(I4B), intent(in) :: npts !< number of station-height data for a reach
478  real(DP), dimension(npts), intent(in) :: stations !< cross-section station distances (x-distance)
479  real(DP), dimension(npts), intent(in) :: heights !< cross-section height data
480  real(DP), intent(in) :: d !< depth to evaluate cross-section
481  real(DP), dimension(npts - 1), intent(inout) :: w !< wetted top widths for each line segment
482  ! -- local variables
483  integer(I4B) :: n
484  real(DP) :: x0
485  real(DP) :: x1
486  real(DP) :: d0
487  real(DP) :: d1
488  real(DP) :: dmax
489  real(DP) :: dmin
490  !
491  ! -- iterate over the station-height data
492  do n = 1, npts - 1
493  !
494  ! -- initialize station-height data for segment
495  x0 = stations(n)
496  x1 = stations(n + 1)
497  d0 = heights(n)
498  d1 = heights(n + 1)
499  !
500  ! -- get the start and end station position of the wetted segment
501  call get_wetted_station(x0, x1, d0, d1, dmax, dmin, d)
502  !
503  ! -- calculate the wetted top width for the segment
504  w(n) = x1 - x0
505  end do
506  end subroutine get_wetted_topwidths
507 
508  !> @brief Calculate the station values for the wetted portion of the cross-section
509  !!
510  !! Subroutine to calculate the station values that define the extent of the
511  !! wetted portion of the cross section for a line segment. The left (x0) and
512  !! right (x1) station positions are altered if the passed depth is less
513  !! than the maximum line segment depth. If the line segment is dry the left
514  !! and right station are equal. Otherwise the wetted station values are equal
515  !! to the full line segment or smaller if the passed depth is less than
516  !! the maximum line segment depth.
517  !!
518  !<
519  pure subroutine get_wetted_station(x0, x1, d0, d1, dmax, dmin, d)
520  ! -- dummy variables
521  real(dp), intent(inout) :: x0 !< left station position
522  real(dp), intent(inout) :: x1 !< right station position
523  real(dp), intent(in) :: d0 !< depth at the left station
524  real(dp), intent(in) :: d1 !< depth at the right station
525  real(dp), intent(inout) :: dmax !< maximum depth
526  real(dp), intent(inout) :: dmin !< minimum depth
527  real(dp), intent(in) :: d !< depth to evaluate cross-section
528  ! -- local variables
529  real(dp) :: xlen
530  real(dp) :: dlen
531  real(dp) :: slope
532  real(dp) :: dx
533  real(dp) :: xt
534  real(dp) :: xt0
535  real(dp) :: xt1
536  !
537  ! -- calculate the minimum and maximum depth
538  dmin = min(d0, d1)
539  dmax = max(d0, d1)
540  !
541  ! -- if d is less than or equal to the minimum value the
542  ! station length (xlen) is zero
543  if (d <= dmin) then
544  x1 = x0
545  ! -- if d is between dmin and dmax station length is less
546  ! than d1 - d0
547  else if (d < dmax) then
548  xlen = x1 - x0
549  dlen = d1 - d0
550  if (abs(dlen) > dzero) then
551  slope = xlen / dlen
552  else
553  slope = dzero
554  end if
555  if (d0 > d1) then
556  dx = (d - d1) * slope
557  xt = x1 + dx
558  xt0 = xt
559  xt1 = x1
560  else
561  dx = (d - d0) * slope
562  xt = x0 + dx
563  xt0 = x0
564  xt1 = xt
565  end if
566  x0 = xt0
567  x1 = xt1
568  end if
569  end subroutine get_wetted_station
570 
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dtwothirds
real constant 2/3
Definition: Constants.f90:70
real(dp), parameter dhalf
real constant 1/2
Definition: Constants.f90:68
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
real(dp), parameter dtwo
real constant 2
Definition: Constants.f90:79
real(dp), parameter done
real constant 1
Definition: Constants.f90:76
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.
Definition: kind.f90:8