MODFLOW 6  version 6.7.0.dev0
USGS Modular Hydrologic Model
SwfCxsUtils.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, lgp
14  use constantsmodule, only: dzero, dhalf, dtwothirds, done, dtwo, &
15  dem5, dp999, dthree
16  use smoothingmodule, only: schsmooth
17 
18  implicit none
19  private
20  public :: calc_depth_from_q
21  public :: calc_qman
22  public :: get_saturated_topwidth
23  public :: get_wetted_topwidth
24  public :: get_wetted_perimeter
25  public :: get_cross_section_area
26  public :: get_hydraulic_radius
27  public :: get_hydraulic_radius_xf
28  public :: get_mannings_section
29  public :: calc_composite_roughness
30  public :: get_conveyance
31  public :: get_composite_conveyance
32 
33 contains
34 
35  !> @brief Calculate the depth at the midpoint of a irregular cross-section
36  !<
37  function calc_depth_from_q(qrch, width, rough, slope, &
38  cxs_xf, cxs_h, cxs_rf, unitconv, &
39  icalcmeth) result(depth)
40  ! -- dummy variables
41  real(dp), intent(in) :: qrch !< streamflow
42  real(dp), intent(in) :: width !< reach width
43  real(dp), intent(in) :: rough !< mannings reach roughness coefficient
44  real(dp), intent(in) :: slope !< reach bottom slope
45  real(dp), dimension(:), intent(in) :: cxs_xf ! xfraction distances for this cross section
46  real(dp), dimension(:), intent(in) :: cxs_h ! heights for this cross section
47  real(dp), dimension(:), intent(in) :: cxs_rf ! mannings fractions for this cross section
48  real(dp), intent(in) :: unitconv !< conversion unit numerator in mannings equation
49  integer(I4B), intent(in) :: icalcmeth !< manning calculation method
50 
51  ! -- return
52  real(dp) :: depth !< stream depth at midpoint of reach
53 
54  if (icalcmeth == 1) then
55  ! slower but robust bisection method
56  depth = calc_depth_from_q_bisect(qrch, width, rough, slope, &
57  cxs_xf, cxs_h, cxs_rf, unitconv, &
58  icalcmeth)
59  else
60  ! faster but less forgiving newton-raphson method
61  depth = calc_depth_from_q_nr(qrch, width, rough, slope, &
62  cxs_xf, cxs_h, cxs_rf, unitconv, &
63  icalcmeth)
64  end if
65  end function calc_depth_from_q
66 
67  !> @brief Calculate the depth at the midpoint of a irregular cross-section
68  !!
69  !! Method to calculate the depth at the midpoint of a reach with a
70  !! irregular cross-section using bisection.
71  !!
72  !<
73  function calc_depth_from_q_bisect(qrch, width, rough, slope, &
74  cxs_xf, cxs_h, cxs_rf, unitconv, &
75  icalcmeth) result(depth)
76  ! -- dummy variables
77  !class(SfrType) :: this !< SfrType object
78  !integer(I4B), intent(in) :: n !< reach number
79  real(dp), intent(in) :: qrch !< streamflow
80  real(dp), intent(in) :: width !< reach width
81  real(dp), intent(in) :: rough !< mannings reach roughness coefficient
82  real(dp), intent(in) :: slope !< reach bottom slope
83  real(dp), dimension(:), intent(in) :: cxs_xf ! xfraction distances for this cross section
84  real(dp), dimension(:), intent(in) :: cxs_h ! heights for this cross section
85  real(dp), dimension(:), intent(in) :: cxs_rf ! mannings fractions for this cross section
86  real(dp), intent(in) :: unitconv !< conversion unit numerator in mannings equation
87  integer(I4B), intent(in) :: icalcmeth !< manning calculation method
88 
89  ! -- return
90  real(dp) :: depth !< stream depth at midpoint of reach
91 
92  ! -- local variables
93  integer(I4B) :: maxiter = 100
94  integer(I4B) :: iter
95  real(dp) :: dmaxchg = dem5
96  real(dp) :: a
97  real(dp) :: b
98  real(dp) :: c
99  real(dp) :: f_a
100  real(dp) :: f_c
101 
102  ! constrain the bisection range
103  a = dzero
104  b = maxval(cxs_h) * dtwo
105  !
106  ! -- bisection iteration
107  bisectiter: do iter = 1, maxiter
108 
109  c = (a + b) / dtwo
110  f_c = calc_qman(c, width, rough, slope, &
111  cxs_xf, cxs_h, cxs_rf, unitconv, icalcmeth) - qrch
112  if (f_c == dzero .or. (b - a) / dtwo < dmaxchg) then
113  depth = c
114  return
115  end if
116 
117  f_a = calc_qman(a, width, rough, slope, &
118  cxs_xf, cxs_h, cxs_rf, unitconv, icalcmeth) - qrch
119  if (sign(done, f_c) == sign(done, f_a)) then
120  a = c
121  else
122  b = c
123  end if
124 
125  end do bisectiter
126  !
127  ! TODO: raise an error
128  print *, "bisection routine failed"
129  end function calc_depth_from_q_bisect
130 
131  !> @brief Calculate the depth at the midpoint of a irregular cross-section
132  !!
133  !! Method to calculate the depth at the midpoint of a reach with a
134  !! irregular cross-section using Newton-Raphson.
135  !!
136  !<
137  function calc_depth_from_q_nr(qrch, width, rough, slope, &
138  cxs_xf, cxs_h, cxs_rf, unitconv, &
139  icalcmeth) result(depth)
140  ! -- dummy variables
141  !class(SfrType) :: this !< SfrType object
142  !integer(I4B), intent(in) :: n !< reach number
143  real(dp), intent(in) :: qrch !< streamflow
144  real(dp), intent(in) :: width !< reach width
145  real(dp), intent(in) :: rough !< mannings reach roughness coefficient
146  real(dp), intent(in) :: slope !< reach bottom slope
147  real(dp), dimension(:), intent(in) :: cxs_xf ! xfraction distances for this cross section
148  real(dp), dimension(:), intent(in) :: cxs_h ! heights for this cross section
149  real(dp), dimension(:), intent(in) :: cxs_rf ! mannings fractions for this cross section
150  real(dp), intent(in) :: unitconv !< conversion unit numerator in mannings equation
151  integer(I4B), intent(in) :: icalcmeth !< manning calculation method
152 
153  ! -- return
154  real(dp) :: depth !< stream depth at midpoint of reach
155 
156  ! -- local variables
157  integer(I4B) :: maxiter = 100
158  integer(I4B) :: iter
159  real(dp) :: dmaxchg = dem5
160  real(dp) :: deps = dem5 * dp999
161  real(dp) :: perturbation
162  real(dp) :: q0
163  real(dp) :: q1
164  real(dp) :: dq
165  real(dp) :: derv
166  real(dp) :: d
167  real(dp) :: dd
168  real(dp) :: residual
169  !
170  ! -- initialize variables
171  perturbation = deps * dtwo
172  d = dzero
173  q0 = dzero
174  residual = q0 - qrch
175  !
176  ! -- Newton-Raphson iteration
177  nriter: do iter = 1, maxiter
178  !call this%sfr_calc_qman(n, d + perturbation, q1)
179  q1 = calc_qman(d + perturbation, width, rough, slope, &
180  cxs_xf, cxs_h, cxs_rf, unitconv, icalcmeth)
181  dq = (q1 - q0)
182  if (dq /= dzero) then
183  derv = perturbation / (q1 - q0)
184  else
185  derv = dzero
186  end if
187  dd = derv * residual
188  d = d - dd
189  !call this%sfr_calc_qman(n, d, q0)
190  q0 = calc_qman(d, width, rough, slope, &
191  cxs_xf, cxs_h, cxs_rf, unitconv, icalcmeth)
192  residual = q0 - qrch
193  !
194  ! -- check for convergence
195  if (abs(dd) < dmaxchg) then
196  exit nriter
197  end if
198 
199  end do nriter
200  depth = d
201  end function calc_depth_from_q_nr
202 
203  !> @brief Calculate streamflow using Manning's equation
204  !!
205  !! Method to calculate the streamflow using Manning's equation for a
206  !! single reach defined by a cross section.
207  !!
208  !<
209  function calc_qman(depth, width, rough, slope, &
210  cxs_xf, cxs_h, cxs_rf, unitconv, icalcmeth) result(qman)
211 
212  ! -- dummy variables
213  real(dp), intent(in) :: depth !< reach depth
214  real(dp), intent(in) :: width !< reach width
215  real(dp), intent(in) :: rough !< mannings reach roughness coefficient
216  real(dp), intent(in) :: slope !< reach bottom slope
217  real(dp), dimension(:), intent(in) :: cxs_xf ! xfraction distances for this cross section
218  real(dp), dimension(:), intent(in) :: cxs_h ! heights for this cross section
219  real(dp), dimension(:), intent(in) :: cxs_rf ! mannings fractions for this cross section
220  real(dp), intent(in) :: unitconv !< conversion unit numerator in mannings equation
221  integer(I4B), intent(in) :: icalcmeth ! 0 is composite linear, 1 is by section, 2 is composite nonlinear
222 
223  ! return value
224  real(dp) :: qman !< calculated mannings flow
225 
226  ! -- local variables
227  integer(I4B) :: linmeth
228 
229  select case (icalcmeth)
230  case (0) ! composite linear
231  linmeth = 0
232  qman = calc_qman_composite(depth, width, rough, slope, &
233  cxs_xf, cxs_h, cxs_rf, unitconv, &
234  linmeth)
235  case (1) ! by section
236  qman = calc_qman_by_section(depth, width, rough, slope, &
237  cxs_xf, cxs_h, cxs_rf, unitconv)
238  case (2) ! composite nonlinear
239  linmeth = 1
240  qman = calc_qman_composite(depth, width, rough, slope, &
241  cxs_xf, cxs_h, cxs_rf, unitconv, &
242  linmeth)
243  end select
244  end function calc_qman
245 
246  function calc_qman_composite(depth, width, rough, slope, &
247  cxs_xf, cxs_h, cxs_rf, unitconv, &
248  linmeth) result(qman)
249 
250  ! -- dummy variables
251  real(dp), intent(in) :: depth !< reach depth
252  real(dp), intent(in) :: width !< reach width
253  real(dp), intent(in) :: rough !< mannings reach roughness coefficient
254  real(dp), intent(in) :: slope !< reach bottom slope
255  real(dp), dimension(:), intent(in) :: cxs_xf ! xfraction distances for this cross section
256  real(dp), dimension(:), intent(in) :: cxs_h ! heights for this cross section
257  real(dp), dimension(:), intent(in) :: cxs_rf ! mannings fractions for this cross section
258  real(dp), intent(in) :: unitconv !< conversion unit numerator in mannings equation
259  integer(I4B), intent(in) :: linmeth !< method for composite (0) linear (1) nonlinear
260 
261  ! return value
262  real(dp) :: qman !< calculated mannings flow
263 
264  ! -- local variables
265  integer(I4B) :: npts
266  real(dp) :: sat
267  real(dp) :: derv
268  real(dp) :: aw
269  real(dp) :: wp
270  real(dp) :: rh
271  real(dp) :: rough_composite
272  !
273  ! -- initialize variables
274  qman = dzero
275  !
276  ! -- calculate Manning's discharge for non-zero depths
277  if (depth > dzero) then
278 
279  npts = size(cxs_xf)
280  rough_composite = calc_composite_roughness(npts, depth, width, &
281  rough, slope, &
282  cxs_xf, cxs_h, cxs_rf, &
283  linmeth)
284  wp = get_wetted_perimeter(npts, cxs_xf, cxs_h, width, depth)
285  aw = get_cross_section_area(npts, cxs_xf, cxs_h, width, depth)
286  if (wp > dzero) then
287  rh = aw / wp
288  else
289  rh = dzero
290  end if
291  qman = unitconv * aw * (rh**dtwothirds) * sqrt(slope) / rough
292  call schsmooth(depth, sat, derv)
293  qman = sat * qman
294  end if
295  end function calc_qman_composite
296 
297  function calc_composite_roughness(npts, depth, width, rough, slope, &
298  cxs_xf, cxs_h, cxs_rf, &
299  linmeth) result(rc)
300  ! -- dummy variables
301  integer(I4B), intent(in) :: npts
302  real(dp), intent(in) :: depth !< reach depth
303  real(dp), intent(in) :: width !< reach width
304  real(dp), intent(in) :: rough !< mannings reach roughness coefficient
305  real(dp), intent(in) :: slope !< reach bottom slope
306  real(dp), dimension(:), intent(in) :: cxs_xf ! xfraction distances for this cross section
307  real(dp), dimension(:), intent(in) :: cxs_h ! heights for this cross section
308  real(dp), dimension(:), intent(in) :: cxs_rf ! mannings fractions for this cross section
309  integer(I4B), intent(in) :: linmeth !< method for composite calculation; linear 0 or nonlinear 1
310 
311  ! return value
312  real(dp) :: rc !< composite roughness
313 
314  ! -- local variables
315  integer(I4B) :: n
316  real(dp) :: wp
317  real(dp) :: exp1
318  real(dp) :: exp2
319  real(dp), dimension(npts) :: stations
320  real(dp), dimension(npts - 1) :: perimeters
321  !
322  ! stations
323  do n = 1, npts
324  stations(n) = cxs_xf(n) * width
325  end do
326  !
327  ! -- select method
328  select case (linmeth)
329  case (0) ! linear
330  exp1 = done
331  exp2 = done
332  case (1) ! nonlinear
333  exp1 = 1.5d0
334  exp2 = dtwo / dthree
335  end select
336  !
337  ! -- initialize variables
338  rc = rough
339  !
340  ! -- calculate composite roughness
341  if (depth > dzero) then
342 
343  if (npts > 1) then
344  call get_wetted_perimeters(npts, stations, cxs_h, depth, perimeters)
345  wp = dzero
346  rc = dzero
347  do n = 1, npts - 1
348  rc = rc + (rough * cxs_rf(n))**exp1 * perimeters(n)
349  wp = wp + perimeters(n)
350  end do
351  if (wp > dzero) then
352  rc = (rc / wp)**exp2
353  else
354  rc = rough
355  end if
356  end if
357 
358  end if
359  end function calc_composite_roughness
360 
361  function calc_qman_by_section(depth, width, rough, slope, &
362  cxs_xf, cxs_h, cxs_rf, unitconv) result(qman)
363 
364  ! -- dummy variables
365  real(dp), intent(in) :: depth !< reach depth
366  real(dp), intent(in) :: width !< reach width
367  real(dp), intent(in) :: rough !< mannings reach roughness coefficient
368  real(dp), intent(in) :: slope !< reach bottom slope
369  real(dp), dimension(:), intent(in) :: cxs_xf ! xfraction distances for this cross section
370  real(dp), dimension(:), intent(in) :: cxs_h ! heights for this cross section
371  real(dp), dimension(:), intent(in) :: cxs_rf ! mannings fractions for this cross section
372  real(dp), intent(in) :: unitconv !< conversion unit numerator in mannings equation
373 
374  ! return value
375  real(dp) :: qman !< calculated mannings flow
376 
377  ! -- local variables
378  integer(I4B) :: npts
379  real(dp) :: sat
380  real(dp) :: derv
381  real(dp) :: aw
382  real(dp) :: wp
383  real(dp) :: rh
384  !
385  ! -- initialize variables
386  qman = dzero
387  !
388  ! -- calculate Manning's discharge for non-zero depths
389  if (depth > dzero) then
390  npts = size(cxs_xf)
391  !
392  ! -- set constant terms for Manning's equation
393  call schsmooth(depth, sat, derv)
394  !
395  ! -- calculate the mannings coefficient that is a
396  ! function of depth
397  if (npts > 1) then
398  !
399  ! -- call function that calculates flow for an
400  ! n-point cross section
401  qman = get_mannings_section(npts, &
402  cxs_xf, &
403  cxs_h, &
404  cxs_rf, &
405  rough, &
406  unitconv, &
407  width, &
408  slope, &
409  depth)
410  else
411 
412  ! hydraulically wide channel (only 1 point is defined)
413  aw = width * depth
414  wp = width
415  if (wp > dzero) then
416  rh = aw / wp
417  else
418  rh = dzero
419  end if
420  qman = unitconv * aw * (rh**dtwothirds) * sqrt(slope) / rough
421  end if
422  !
423  ! -- calculate stream flow
424  qman = sat * qman
425  end if
426  end function calc_qman_by_section
427 
428  !> @brief Calculate the saturated top width for a reach
429  !!
430  !! Function to calculate the maximum top width for a reach using the
431  !! cross-section station data.
432  !!
433  !! @return w saturated top width
434  !<
435  function get_saturated_topwidth(npts, stations) result(w)
436  ! -- dummy variables
437  integer(I4B), intent(in) :: npts !< number of station-height data for a reach
438  real(dp), dimension(npts), intent(in) :: stations !< cross-section station distances (x-distance)
439  ! -- local variables
440  real(dp) :: w
441  !
442  ! -- calculate the saturated top width
443  if (npts > 1) then
444  w = stations(npts) - stations(1)
445  else
446  w = stations(1)
447  end if
448  end function get_saturated_topwidth
449 
450  !> @brief Calculate the wetted top width for a reach
451  !!
452  !! Function to calculate the wetted top width for a reach using the
453  !! cross-section station-height data given a passed depth.
454  !!
455  !! @return w wetted top width
456  !<
457  function get_wetted_topwidth(npts, xfraction, heights, width, d) result(w)
458  ! -- dummy variables
459  integer(I4B), intent(in) :: npts !< number of station-height data for a reach
460  real(dp), dimension(npts), intent(in) :: xfraction !< cross-section station distances (x-distance)
461  real(dp), dimension(npts), intent(in) :: heights !< cross-section height data
462  real(dp), intent(in) :: width !< reach width
463  real(dp), intent(in) :: d !< depth to evaluate cross-section
464  ! -- local variables
465  integer(I4B) :: n
466  real(dp) :: w
467  real(dp), dimension(npts) :: stations
468  real(dp), dimension(npts - 1) :: widths
469  !
470  ! -- calculate station from xfractions and width
471  do n = 1, npts
472  stations(n) = xfraction(n) * width
473  end do
474  !
475  ! -- initialize the wetted perimeter for the reach
476  w = dzero
477  !
478  ! -- calculate the wetted top width for each line segment
479  call get_wetted_topwidths(npts, stations, heights, d, widths)
480  !
481  ! -- calculate the wetted top widths
482  do n = 1, npts - 1
483  w = w + widths(n)
484  end do
485  end function get_wetted_topwidth
486 
487  !> @brief Calculate the wetted perimeter for a reach
488  !!
489  !! Function to calculate the wetted perimeter for a reach using the
490  !! cross-section station-height data given a passed depth.
491  !!
492  !! @return p wetted perimeter
493  !<
494  function get_wetted_perimeter(npts, xfraction, heights, width, d) result(p)
495  ! -- dummy variables
496  integer(I4B), intent(in) :: npts !< number of station-height data for a reach
497  real(dp), dimension(npts), intent(in) :: xfraction !< cross-section x fractions
498  real(dp), dimension(npts), intent(in) :: heights !< cross-section height data
499  real(dp), intent(in) :: width !< cross section width
500  real(dp), intent(in) :: d !< depth to evaluate cross-section
501  ! -- local variables
502  integer(I4B) :: n
503  real(dp) :: p
504  real(dp), dimension(npts) :: stations
505  real(dp), dimension(npts - 1) :: perimeters
506  !
507  ! -- initialize stations
508  do n = 1, npts
509  stations(n) = xfraction(n) * width
510  end do
511  !
512  ! -- initialize the wetted perimeter for the reach
513  p = dzero
514  !
515  ! -- calculate the wetted perimeter for each line segment
516  call get_wetted_perimeters(npts, stations, heights, d, perimeters)
517  !
518  ! -- calculate the wetted perimenter
519  do n = 1, npts - 1
520  p = p + perimeters(n)
521  end do
522  end function get_wetted_perimeter
523 
524  !> @brief Calculate the cross-sectional area for a reach
525  !!
526  !! Function to calculate the cross-sectional area for a reach using
527  !! the cross-section station-height data given a passed depth.
528  !!
529  !! @return a cross-sectional area
530  !<
531  function get_cross_section_area(npts, xfraction, heights, width, d) result(a)
532  ! -- dummy variables
533  integer(I4B), intent(in) :: npts !< number of station-height data for a reach
534  real(dp), dimension(npts), intent(in) :: xfraction !< cross-section station distances (x-distance)
535  real(dp), dimension(npts), intent(in) :: heights !< cross-section height data
536  real(dp), intent(in) :: width !< reach width
537  real(dp), intent(in) :: d !< depth to evaluate cross-section
538  ! -- local variables
539  integer(I4B) :: n
540  real(dp) :: a
541  real(dp), dimension(npts) :: stations
542  real(dp), dimension(npts - 1) :: areas
543  !
544  ! -- calculate station from xfractions and width
545  do n = 1, npts
546  stations(n) = xfraction(n) * width
547  end do
548  !
549  ! -- initialize the area
550  a = dzero
551  !
552  ! -- calculate the cross-sectional area for each line segment
553  call get_cross_section_areas(npts, stations, heights, d, areas)
554  !
555  ! -- calculate the cross-sectional area
556  do n = 1, npts - 1
557  a = a + areas(n)
558  end do
559  end function get_cross_section_area
560 
561  !> @brief Calculate conveyance
562  !!
563  !! Return a calculated conveyance. Calculation
564  !! depends on whether the channel is a rectangular channel
565  !< with vertical sides or something else.
566  function get_conveyance(npts, xfraction, heights, cxs_rf, &
567  width, rough, d) result(c)
568  ! -- dummy variables
569  integer(I4B), intent(in) :: npts !< number of station-height data for a reach
570  real(dp), dimension(npts), intent(in) :: xfraction !< cross-section station distances (x-distance)
571  real(dp), dimension(npts), intent(in) :: heights !< cross-section height data
572  real(dp), dimension(:), intent(in) :: cxs_rf ! mannings fractions for this cross section
573  real(dp), intent(in) :: width !< reach width
574  real(dp), intent(in) :: rough !< reach roughness
575  real(dp), intent(in) :: d !< depth to evaluate cross-section
576  ! result
577  real(dp) :: c
578  if (is_rectangular(xfraction)) then
579  c = get_rectangular_conveyance(npts, xfraction, heights, cxs_rf, &
580  width, rough, d)
581  else
582  c = get_composite_conveyance(npts, xfraction, heights, cxs_rf, &
583  width, rough, d)
584  end if
585  end function get_conveyance
586 
587  !> @brief Calculate conveyance for rectangular channel
588  !<
589  function get_rectangular_conveyance(npts, xfraction, heights, cxs_rf, &
590  width, rough, d) result(c)
591  ! -- dummy variables
592  integer(I4B), intent(in) :: npts !< number of station-height data for a reach
593  real(dp), dimension(npts), intent(in) :: xfraction !< cross-section station distances (x-distance)
594  real(dp), dimension(npts), intent(in) :: heights !< cross-section height data
595  real(dp), dimension(:), intent(in) :: cxs_rf ! mannings fractions for this cross section
596  real(dp), intent(in) :: width !< reach width
597  real(dp), intent(in) :: rough !< reach roughness
598  real(dp), intent(in) :: d !< depth to evaluate cross-section
599  ! -- local variables
600  real(dp) :: c
601  real(dp) :: a
602  real(dp) :: p
603  real(dp) :: ravg
604 
605  ! find cross sectional area and wetted perimeter
606  a = get_cross_section_area(npts, xfraction, heights, width, d)
607  p = get_wetted_perimeter(npts, xfraction, heights, width, d)
608 
609  ! calculate roughness or determine an average
610  if (has_uniform_resistance(cxs_rf)) then
611  ravg = rough * cxs_rf(1)
612  else
613  ravg = calc_composite_roughness(npts, d, width, rough, dzero, &
614  xfraction, heights, cxs_rf, 0)
615  end if
616 
617  ! make conveyance calculation
618  c = dzero
619  if (p > dzero) then
620  c = a / ravg * (a / p)**dtwothirds
621  end if
622  end function get_rectangular_conveyance
623 
624  !> @brief Determine if cross section is rectangular
625  !<
626  function is_rectangular(xfraction)
627  ! dummy
628  real(dp), dimension(:), intent(in) :: xfraction !< cross-section station distances (x-distance)
629  ! result
630  logical(LGP) :: is_rectangular
631  is_rectangular = .false.
632  if (size(xfraction) == 4) then
633  if (xfraction(1) == xfraction(2) .and. &
634  xfraction(3) == xfraction(4)) then
635  is_rectangular = .true.
636  end if
637  end if
638  end function is_rectangular
639 
640  !> @brief Determine if roughness is uniform for the section
641  !<
642  function has_uniform_resistance(cxs_rf)
643  ! dummy
644  real(dp), dimension(:), intent(in) :: cxs_rf !< cross-section station distances (x-distance)
645  ! result
646  logical(LGP) :: has_uniform_resistance
647  ! local
648  real(dp) :: rmin
649  real(dp) :: rmax
650  rmin = minval(cxs_rf(1:3))
651  rmax = maxval(cxs_rf(1:3))
652  has_uniform_resistance = (rmin == rmax)
653  end function has_uniform_resistance
654 
655  !> @brief Calculate composite conveyance
656  !!
657  !! Function to calculate the composite conveyance.
658  !! This is based on the approach in HEC-RAS, where
659  !! a conveyance is calculated for each line segment
660  !! in the cross section, and then summed to produce
661  !! a total conveyance.
662  !!
663  !<
664  function get_composite_conveyance(npts, xfraction, heights, cxs_rf, &
665  width, rough, d) result(c)
666  ! -- dummy variables
667  integer(I4B), intent(in) :: npts !< number of station-height data for a reach
668  real(dp), dimension(npts), intent(in) :: xfraction !< cross-section station distances (x-distance)
669  real(dp), dimension(npts), intent(in) :: heights !< cross-section height data
670  real(dp), dimension(:), intent(in) :: cxs_rf ! mannings fractions for this cross section
671  real(dp), intent(in) :: width !< reach width
672  real(dp), intent(in) :: rough !< reach roughness
673  real(dp), intent(in) :: d !< depth to evaluate cross-section
674  ! -- local variables
675  integer(I4B) :: n
676  real(dp) :: c
677  real(dp) :: p
678  real(dp) :: rc
679  real(dp) :: rh
680  real(dp) :: cn
681  real(dp), dimension(npts) :: stations
682  real(dp), dimension(npts - 1) :: areas
683  real(dp), dimension(npts - 1) :: perimeters
684  !
685  ! -- calculate station from xfractions and width
686  do n = 1, npts
687  stations(n) = xfraction(n) * width
688  end do
689  !
690  ! -- calculate the cross-sectional area for each line segment
691  call get_cross_section_areas(npts, stations, heights, d, areas)
692  !
693  ! -- calculate the wetted perimeter for each line segment
694  call get_wetted_perimeters(npts, stations, heights, d, perimeters)
695  !
696  ! -- calculate the composite conveyance
697  c = dzero
698  do n = 1, npts - 1
699  rc = cxs_rf(n) * rough
700  p = perimeters(n)
701  cn = dzero
702  if (p > dzero) then
703  rh = areas(n) / p
704  cn = areas(n) / rc * rh**dtwothirds
705  end if
706  c = c + cn
707  end do
708  end function get_composite_conveyance
709 
710  !> @brief Calculate the hydraulic radius for a reach
711  !!
712  !! Function to calculate the hydraulic radius for a reach using
713  !! the cross-section station-height data given a passed depth.
714  !!
715  !! @return r hydraulic radius
716  !<
717  function get_hydraulic_radius(npts, stations, heights, d) result(r)
718  ! -- dummy variables
719  integer(I4B), intent(in) :: npts !< number of station-height data for a reach
720  real(dp), dimension(npts), intent(in) :: stations !< cross-section station distances (x-distance)
721  real(dp), dimension(npts), intent(in) :: heights !< cross-section height data
722  real(dp), intent(in) :: d !< depth to evaluate cross-section
723  ! -- local variables
724  integer(I4B) :: n
725  real(dp) :: r
726  real(dp) :: p
727  real(dp) :: a
728  real(dp), dimension(npts - 1) :: areas
729  real(dp), dimension(npts - 1) :: perimeters
730  !
731  ! -- initialize the hydraulic radius, perimeter, and area
732  r = dzero
733  p = dzero
734  a = dzero
735  !
736  ! -- calculate the wetted perimeter for each line segment
737  call get_wetted_perimeters(npts, stations, heights, d, perimeters)
738  !
739  ! -- calculate the wetted perimenter
740  do n = 1, npts - 1
741  p = p + perimeters(n)
742  end do
743  !
744  ! -- calculate the hydraulic radius only if the perimeter is non-zero
745  if (p > dzero) then
746  !
747  ! -- calculate the cross-sectional area for each line segment
748  call get_cross_section_areas(npts, stations, heights, d, areas)
749  !
750  ! -- calculate the cross-sectional area
751  do n = 1, npts - 1
752  a = a + areas(n)
753  end do
754  !
755  ! -- calculate the hydraulic radius
756  r = a / p
757  end if
758  end function get_hydraulic_radius
759 
760  !> @brief Calculate the hydraulic radius for a reach
761  !!
762  !! Function to calculate the hydraulic radius for a reach using
763  !! the cross-section xfraction-height data given a passed depth.
764  !! This is different from get_hydraulic_radius as it requires
765  !! xfraction and width instead of station.
766  !!
767  !! @return r hydraulic radius
768  !<
769  function get_hydraulic_radius_xf(npts, xfraction, heights, width, d) result(r)
770  ! -- dummy variables
771  integer(I4B), intent(in) :: npts !< number of station-height data for a reach
772  real(dp), dimension(npts), intent(in) :: xfraction !< cross-section station distances (x-distance)
773  real(dp), dimension(npts), intent(in) :: heights !< cross-section height data
774  real(dp), intent(in) :: width !< reach width
775  real(dp), intent(in) :: d !< depth to evaluate cross-section
776  ! -- local variables
777  integer(I4B) :: n
778  real(dp) :: r
779  real(dp), dimension(npts) :: stations
780  !
781  ! -- calculate station from xfractions and width
782  do n = 1, npts
783  stations(n) = xfraction(n) * width
784  end do
785  !
786  ! -- calculate hydraulic radius
787  r = get_hydraulic_radius(npts, stations, heights, d)
788  end function get_hydraulic_radius_xf
789 
790  !> @brief Calculate the manning's discharge for a reach
791  !!
792  !! Function to calculate the mannings discharge for a reach
793  !! by calculating the discharge for each section, which can
794  !! have a unique Manning's coefficient given a passed depth.
795  !!
796  !! @return q reach discharge
797  !<
798  function get_mannings_section(npts, xfraction, heights, roughfracs, &
799  roughness, conv_fact, width, slope, d) result(q)
800  ! -- dummy variables
801  integer(I4B), intent(in) :: npts !< number of station-height data for a reach
802  real(dp), dimension(npts), intent(in) :: xfraction !< cross-section station distances (x-distance)
803  real(dp), dimension(npts), intent(in) :: heights !< cross-section height data
804  real(dp), dimension(npts), intent(in) :: roughfracs !< cross-section Mannings roughness fraction data
805  real(dp), intent(in) :: roughness !< base reach roughness
806  real(dp), intent(in) :: conv_fact !< unit conversion factor
807  real(dp), intent(in) :: width !< reach width
808  real(dp), intent(in) :: slope !< reach slope
809  real(dp), intent(in) :: d !< depth to evaluate cross-section
810  ! -- local variables
811  integer(I4B) :: n
812  real(dp) :: q
813  real(dp) :: rh
814  real(dp) :: r
815  real(dp) :: p
816  real(dp) :: a
817  real(dp), dimension(npts) :: stations
818  real(dp), dimension(npts - 1) :: areas
819  real(dp), dimension(npts - 1) :: perimeters
820  !
821  ! -- initialize the hydraulic radius, perimeter, and area
822  q = dzero
823  rh = dzero
824  r = dzero
825  p = dzero
826  a = dzero
827  !
828  ! -- calculate station from xfractions and width
829  do n = 1, npts
830  stations(n) = xfraction(n) * width
831  end do
832  !
833  ! -- calculate the wetted perimeter for each line segment
834  call get_wetted_perimeters(npts, stations, heights, d, perimeters)
835  !
836  ! -- calculate the wetted perimenter
837  do n = 1, npts - 1
838  p = p + perimeters(n)
839  end do
840  !
841  ! -- calculate the hydraulic radius only if the perimeter is non-zero
842  if (p > dzero) then
843  !
844  ! -- calculate the cross-sectional area for each line segment
845  call get_cross_section_areas(npts, stations, heights, d, areas)
846  !
847  ! -- calculate the cross-sectional area
848  do n = 1, npts - 1
849  p = perimeters(n)
850  r = roughness * roughfracs(n)
851  if (p * r > dzero) then
852  a = areas(n)
853  rh = a / p
854  q = q + conv_fact * a * rh**dtwothirds * sqrt(slope) / r
855  end if
856  end do
857  end if
858  end function get_mannings_section
859 
860  ! -- private functions and subroutines
861 
862  !> @brief Calculate the wetted perimeters for each line segment
863  !!
864  !! Subroutine to calculate the wetted perimeter for each line segment
865  !! that defines the reach using the cross-section station-height
866  !! data given a passed depth.
867  !!
868  !<
869  subroutine get_wetted_perimeters(npts, stations, heights, d, p)
870  ! -- dummy variables
871  integer(I4B), intent(in) :: npts !< number of station-height data for a reach
872  real(DP), dimension(npts), intent(in) :: stations !< cross-section station distances (x-distance)
873  real(DP), dimension(npts), intent(in) :: heights !< cross-section height data
874  real(DP), intent(in) :: d !< depth to evaluate cross-section
875  real(DP), dimension(npts - 1), intent(inout) :: p !< wetted perimeter for each line segment
876  ! -- local variables
877  integer(I4B) :: n
878  real(DP) :: x0
879  real(DP) :: x1
880  real(DP) :: d0
881  real(DP) :: d1
882  real(DP) :: dmax
883  real(DP) :: dmin
884  real(DP) :: xlen
885  real(DP) :: dlen
886  !
887  ! -- iterate over the station-height data
888  do n = 1, npts - 1
889  !
890  ! -- initialize the wetted perimeter
891  p(n) = dzero
892  !
893  ! -- initialize station-height data for segment
894  x0 = stations(n)
895  x1 = stations(n + 1)
896  d0 = heights(n)
897  d1 = heights(n + 1)
898  !
899  ! -- get the start and end station position of the wetted segment
900  call get_wetted_station(x0, x1, d0, d1, dmax, dmin, d)
901  !
902  ! -- calculate the wetted perimeter for the segment
903  xlen = x1 - x0
904  dlen = dzero
905  if (xlen > dzero) then
906  if (d > dmax) then
907  dlen = dmax - dmin
908  else
909  dlen = d - dmin
910  end if
911  else
912  if (d > dmin) then
913  dlen = min(d, dmax) - dmin
914  else
915  dlen = dzero
916  end if
917  end if
918  p(n) = sqrt(xlen**dtwo + dlen**dtwo)
919  end do
920  end subroutine get_wetted_perimeters
921 
922  !> @brief Calculate the cross-sectional areas for each line segment
923  !!
924  !! Subroutine to calculate the cross-sectional area for each line segment
925  !! that defines the reach using the cross-section station-height
926  !! data given a passed depth.
927  !!
928  !<
929  subroutine get_cross_section_areas(npts, stations, heights, d, a)
930  ! -- dummy variables
931  integer(I4B), intent(in) :: npts !< number of station-height data for a reach
932  real(DP), dimension(npts), intent(in) :: stations !< cross-section station distances (x-distance)
933  real(DP), dimension(npts), intent(in) :: heights !< cross-section height data
934  real(DP), intent(in) :: d !< depth to evaluate cross-section
935  real(DP), dimension(npts - 1), intent(inout) :: a !< cross-sectional area for each line segment
936  ! -- local variables
937  integer(I4B) :: n
938  real(DP) :: x0
939  real(DP) :: x1
940  real(DP) :: d0
941  real(DP) :: d1
942  real(DP) :: dmax
943  real(DP) :: dmin
944  real(DP) :: xlen
945  !
946  ! -- iterate over the station-height data
947  do n = 1, npts - 1
948  !
949  ! -- initialize the cross-sectional area
950  a(n) = dzero
951  !
952  ! -- initialize station-height data for segment
953  x0 = stations(n)
954  x1 = stations(n + 1)
955  d0 = heights(n)
956  d1 = heights(n + 1)
957  !
958  ! -- get the start and end station position of the wetted segment
959  call get_wetted_station(x0, x1, d0, d1, dmax, dmin, d)
960  !
961  ! -- calculate the cross-sectional area for the segment
962  xlen = x1 - x0
963  if (xlen > dzero) then
964  !
965  ! -- add the area above dmax
966  if (d > dmax) then
967  a(n) = xlen * (d - dmax)
968  end if
969  !
970  ! -- add the area below dmax
971  if (dmax /= dmin .and. d > dmin) then
972  if (d < dmax) then
973  a(n) = a(n) + dhalf * (d - dmin) * xlen
974  else
975  a(n) = a(n) + dhalf * (dmax - dmin) * xlen
976  end if
977  end if
978  end if
979  end do
980  end subroutine get_cross_section_areas
981 
982  !> @brief Calculate the wetted top widths for each line segment
983  !!
984  !! Subroutine to calculate the wetted top width for each line segment
985  !! that defines the reach using the cross-section station-height
986  !! data given a passed depth.
987  !!
988  !<
989  subroutine get_wetted_topwidths(npts, stations, heights, d, w)
990  ! -- dummy variables
991  integer(I4B), intent(in) :: npts !< number of station-height data for a reach
992  real(DP), dimension(npts), intent(in) :: stations !< cross-section station distances (x-distance)
993  real(DP), dimension(npts), intent(in) :: heights !< cross-section height data
994  real(DP), intent(in) :: d !< depth to evaluate cross-section
995  real(DP), dimension(npts - 1), intent(inout) :: w !< wetted top widths for each line segment
996  ! -- local variables
997  integer(I4B) :: n
998  real(DP) :: x0
999  real(DP) :: x1
1000  real(DP) :: d0
1001  real(DP) :: d1
1002  real(DP) :: dmax
1003  real(DP) :: dmin
1004  !
1005  ! -- iterate over the station-height data
1006  do n = 1, npts - 1
1007  !
1008  ! -- initialize station-height data for segment
1009  x0 = stations(n)
1010  x1 = stations(n + 1)
1011  d0 = heights(n)
1012  d1 = heights(n + 1)
1013  !
1014  ! -- get the start and end station position of the wetted segment
1015  call get_wetted_station(x0, x1, d0, d1, dmax, dmin, d)
1016  !
1017  ! -- calculate the wetted top width for the segment
1018  w(n) = x1 - x0
1019  end do
1020  end subroutine get_wetted_topwidths
1021 
1022  !> @brief Calculate the station values for the wetted portion of the cross-section
1023  !!
1024  !! Subroutine to calculate the station values that define the extent of the
1025  !! wetted portion of the cross section for a line segment. The left (x0) and
1026  !! right (x1) station positions are altered if the passed depth is less
1027  !! than the maximum line segment depth. If the line segment is dry the left
1028  !! and right station are equal. Otherwise the wetted station values are equal
1029  !! to the full line segment or smaller if the passed depth is less than
1030  !! the maximum line segment depth.
1031  !!
1032  !<
1033  pure subroutine get_wetted_station(x0, x1, d0, d1, dmax, dmin, d)
1034  ! -- dummy variables
1035  real(dp), intent(inout) :: x0 !< left station position
1036  real(dp), intent(inout) :: x1 !< right station position
1037  real(dp), intent(in) :: d0 !< depth at the left station
1038  real(dp), intent(in) :: d1 !< depth at the right station
1039  real(dp), intent(inout) :: dmax !< maximum depth
1040  real(dp), intent(inout) :: dmin !< minimum depth
1041  real(dp), intent(in) :: d !< depth to evaluate cross-section
1042  ! -- local variables
1043  real(dp) :: xlen
1044  real(dp) :: dlen
1045  real(dp) :: slope
1046  real(dp) :: dx
1047  real(dp) :: xt
1048  real(dp) :: xt0
1049  real(dp) :: xt1
1050  !
1051  ! -- calculate the minimum and maximum depth
1052  dmin = min(d0, d1)
1053  dmax = max(d0, d1)
1054  !
1055  ! -- if d is less than or equal to the minimum value the
1056  ! station length (xlen) is zero
1057  if (d <= dmin) then
1058  x1 = x0
1059  ! -- if d is between dmin and dmax station length is less
1060  ! than d1 - d0
1061  else if (d < dmax) then
1062  xlen = x1 - x0
1063  dlen = d1 - d0
1064  if (abs(dlen) > dzero) then
1065  slope = xlen / dlen
1066  else
1067  slope = dzero
1068  end if
1069  if (d0 > d1) then
1070  dx = (d - d1) * slope
1071  xt = x1 + dx
1072  xt0 = xt
1073  xt1 = x1
1074  else
1075  dx = (d - d0) * slope
1076  xt = x0 + dx
1077  xt0 = x0
1078  xt1 = xt
1079  end if
1080  x0 = xt0
1081  x1 = xt1
1082  end if
1083  end subroutine get_wetted_station
1084 
1085 end module swfcxsutilsmodule
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dtwothirds
real constant 2/3
Definition: Constants.f90:70
real(dp), parameter dp999
real constant 999/1000
Definition: Constants.f90:74
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 dem5
real constant 1e-5
Definition: Constants.f90:108
real(dp), parameter dtwo
real constant 2
Definition: Constants.f90:79
real(dp), parameter dthree
real constant 3
Definition: Constants.f90:80
real(dp), parameter done
real constant 1
Definition: Constants.f90:76
This module defines variable data types.
Definition: kind.f90:8
subroutine schsmooth(d, smooth, dwdh)
@ brief sChSmooth
This module contains stateless sfr subroutines and functions.
Definition: SwfCxsUtils.f90:11
real(dp) function, public calc_qman(depth, width, rough, slope, cxs_xf, cxs_h, cxs_rf, unitconv, icalcmeth)
Calculate streamflow using Manning's equation.
real(dp) function, public get_composite_conveyance(npts, xfraction, heights, cxs_rf, width, rough, d)
Calculate composite conveyance.
real(dp) function, public get_mannings_section(npts, xfraction, heights, roughfracs, roughness, conv_fact, width, slope, d)
Calculate the manning's discharge for a reach.
real(dp) function, public get_hydraulic_radius_xf(npts, xfraction, heights, width, d)
Calculate the hydraulic radius for a reach.
subroutine get_cross_section_areas(npts, stations, heights, d, a)
Calculate the cross-sectional areas for each line segment.
logical(lgp) function has_uniform_resistance(cxs_rf)
Determine if roughness is uniform for the section.
logical(lgp) function is_rectangular(xfraction)
Determine if cross section is rectangular.
real(dp) function, public calc_composite_roughness(npts, depth, width, rough, slope, cxs_xf, cxs_h, cxs_rf, linmeth)
real(dp) function, public calc_depth_from_q(qrch, width, rough, slope, cxs_xf, cxs_h, cxs_rf, unitconv, icalcmeth)
Calculate the depth at the midpoint of a irregular cross-section.
Definition: SwfCxsUtils.f90:40
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_hydraulic_radius(npts, stations, heights, d)
Calculate the hydraulic radius for a reach.
real(dp) function calc_qman_by_section(depth, width, rough, slope, cxs_xf, cxs_h, cxs_rf, unitconv)
subroutine get_wetted_topwidths(npts, stations, heights, d, w)
Calculate the wetted top widths for each line segment.
real(dp) function, public get_cross_section_area(npts, xfraction, heights, width, d)
Calculate the cross-sectional area for a reach.
real(dp) function, public get_wetted_topwidth(npts, xfraction, heights, width, d)
Calculate the wetted top width for a reach.
real(dp) function calc_qman_composite(depth, width, rough, slope, cxs_xf, cxs_h, cxs_rf, unitconv, linmeth)
real(dp) function calc_depth_from_q_nr(qrch, width, rough, slope, cxs_xf, cxs_h, cxs_rf, unitconv, icalcmeth)
Calculate the depth at the midpoint of a irregular cross-section.
real(dp) function, public get_conveyance(npts, xfraction, heights, cxs_rf, width, rough, d)
Calculate conveyance.
subroutine get_wetted_perimeters(npts, stations, heights, d, p)
Calculate the wetted perimeters for each line segment.
real(dp) function, public get_wetted_perimeter(npts, xfraction, heights, width, d)
Calculate the wetted perimeter for a reach.
real(dp) function get_rectangular_conveyance(npts, xfraction, heights, cxs_rf, width, rough, d)
Calculate conveyance for rectangular channel.
real(dp) function, public get_saturated_topwidth(npts, stations)
Calculate the saturated top width for a reach.
real(dp) function calc_depth_from_q_bisect(qrch, width, rough, slope, cxs_xf, cxs_h, cxs_rf, unitconv, icalcmeth)
Calculate the depth at the midpoint of a irregular cross-section.
Definition: SwfCxsUtils.f90:76