MODFLOW 6  version 6.7.0.dev1
USGS Modular Hydrologic Model
swfcxsutilsmodule Module Reference

This module contains stateless sfr subroutines and functions. More...

Functions/Subroutines

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. More...
 
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. More...
 
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. More...
 
real(dp) function, public calc_qman (depth, width, rough, slope, cxs_xf, cxs_h, cxs_rf, unitconv, icalcmeth)
 Calculate streamflow using Manning's equation. More...
 
real(dp) function calc_qman_composite (depth, width, rough, slope, cxs_xf, cxs_h, cxs_rf, unitconv, linmeth)
 
real(dp) function, public calc_composite_roughness (npts, depth, width, rough, slope, cxs_xf, cxs_h, cxs_rf, linmeth)
 
real(dp) function calc_qman_by_section (depth, width, rough, slope, cxs_xf, cxs_h, cxs_rf, unitconv)
 
real(dp) function, public get_saturated_topwidth (npts, xfraction, width)
 Calculate the saturated top width for a reach. More...
 
real(dp) function, public get_wetted_topwidth (npts, xfraction, heights, width, d)
 Calculate the wetted top width for a reach. More...
 
real(dp) function, public get_wetted_perimeter (npts, xfraction, heights, width, d)
 Calculate the wetted perimeter for a reach. More...
 
real(dp) function, public get_cross_section_area (npts, xfraction, heights, width, d)
 Calculate the cross-sectional area for a reach. More...
 
real(dp) function, public get_conveyance (npts, xfraction, heights, cxs_rf, width, rough, d)
 Calculate conveyance. More...
 
real(dp) function get_rectangular_conveyance (npts, xfraction, heights, cxs_rf, width, rough, d)
 Calculate conveyance for rectangular channel. More...
 
logical(lgp) function is_rectangular (xfraction)
 Determine if cross section is rectangular. More...
 
logical(lgp) function has_uniform_resistance (cxs_rf)
 Determine if roughness is uniform for the section. More...
 
real(dp) function, public get_composite_conveyance (npts, xfraction, heights, cxs_rf, width, rough, d)
 Calculate composite conveyance. More...
 
real(dp) function, public get_hydraulic_radius (npts, stations, heights, d)
 Calculate the hydraulic radius for a reach. More...
 
real(dp) function, public get_hydraulic_radius_xf (npts, xfraction, heights, width, d)
 Calculate the hydraulic radius for a reach. More...
 
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. More...
 
subroutine get_wetted_perimeters (npts, stations, heights, d, p)
 Calculate the wetted perimeters for each line segment. More...
 
subroutine get_cross_section_areas (npts, stations, heights, d, a)
 Calculate the cross-sectional areas for each line segment. More...
 
subroutine get_wetted_topwidths (npts, stations, heights, d, w)
 Calculate the wetted top widths for each line segment. More...
 
pure subroutine get_wetted_station (x0, x1, d0, d1, dmax, dmin, d)
 Calculate the station values for the wetted portion of the cross-section. More...
 

Detailed Description

This module contains the functions to calculate the wetted perimeter and cross-sectional area for a reach cross-section that are used in the streamflow routing (SFR) package. It also contains subroutines to calculate the wetted perimeter and cross-sectional area for each line segment in the cross-section. This module does not depend on the SFR package.

Function/Subroutine Documentation

◆ calc_composite_roughness()

real(dp) function, public swfcxsutilsmodule::calc_composite_roughness ( integer(i4b), intent(in)  npts,
real(dp), intent(in)  depth,
real(dp), intent(in)  width,
real(dp), intent(in)  rough,
real(dp), intent(in)  slope,
real(dp), dimension(:), intent(in)  cxs_xf,
real(dp), dimension(:), intent(in)  cxs_h,
real(dp), dimension(:), intent(in)  cxs_rf,
integer(i4b), intent(in)  linmeth 
)
Parameters
[in]depthreach depth
[in]widthreach width
[in]roughmannings reach roughness coefficient
[in]slopereach bottom slope
[in]linmethmethod for composite calculation; linear 0 or nonlinear 1
Returns
composite roughness

Definition at line 297 of file SwfCxsUtils.f90.

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
Here is the call graph for this function:
Here is the caller graph for this function:

◆ calc_depth_from_q()

real(dp) function, public swfcxsutilsmodule::calc_depth_from_q ( real(dp), intent(in)  qrch,
real(dp), intent(in)  width,
real(dp), intent(in)  rough,
real(dp), intent(in)  slope,
real(dp), dimension(:), intent(in)  cxs_xf,
real(dp), dimension(:), intent(in)  cxs_h,
real(dp), dimension(:), intent(in)  cxs_rf,
real(dp), intent(in)  unitconv,
integer(i4b), intent(in)  icalcmeth 
)
Parameters
[in]qrchstreamflow
[in]widthreach width
[in]roughmannings reach roughness coefficient
[in]slopereach bottom slope
[in]unitconvconversion unit numerator in mannings equation
[in]icalcmethmanning calculation method
Returns
stream depth at midpoint of reach

Definition at line 37 of file SwfCxsUtils.f90.

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
Here is the call graph for this function:

◆ calc_depth_from_q_bisect()

real(dp) function swfcxsutilsmodule::calc_depth_from_q_bisect ( real(dp), intent(in)  qrch,
real(dp), intent(in)  width,
real(dp), intent(in)  rough,
real(dp), intent(in)  slope,
real(dp), dimension(:), intent(in)  cxs_xf,
real(dp), dimension(:), intent(in)  cxs_h,
real(dp), dimension(:), intent(in)  cxs_rf,
real(dp), intent(in)  unitconv,
integer(i4b), intent(in)  icalcmeth 
)
private

Method to calculate the depth at the midpoint of a reach with a irregular cross-section using bisection.

Parameters
[in]qrchstreamflow
[in]widthreach width
[in]roughmannings reach roughness coefficient
[in]slopereach bottom slope
[in]unitconvconversion unit numerator in mannings equation
[in]icalcmethmanning calculation method
Returns
stream depth at midpoint of reach

Definition at line 73 of file SwfCxsUtils.f90.

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"
Here is the call graph for this function:
Here is the caller graph for this function:

◆ calc_depth_from_q_nr()

real(dp) function swfcxsutilsmodule::calc_depth_from_q_nr ( real(dp), intent(in)  qrch,
real(dp), intent(in)  width,
real(dp), intent(in)  rough,
real(dp), intent(in)  slope,
real(dp), dimension(:), intent(in)  cxs_xf,
real(dp), dimension(:), intent(in)  cxs_h,
real(dp), dimension(:), intent(in)  cxs_rf,
real(dp), intent(in)  unitconv,
integer(i4b), intent(in)  icalcmeth 
)
private

Method to calculate the depth at the midpoint of a reach with a irregular cross-section using Newton-Raphson.

Parameters
[in]qrchstreamflow
[in]widthreach width
[in]roughmannings reach roughness coefficient
[in]slopereach bottom slope
[in]unitconvconversion unit numerator in mannings equation
[in]icalcmethmanning calculation method
Returns
stream depth at midpoint of reach

Definition at line 137 of file SwfCxsUtils.f90.

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
Here is the call graph for this function:
Here is the caller graph for this function:

◆ calc_qman()

real(dp) function, public swfcxsutilsmodule::calc_qman ( real(dp), intent(in)  depth,
real(dp), intent(in)  width,
real(dp), intent(in)  rough,
real(dp), intent(in)  slope,
real(dp), dimension(:), intent(in)  cxs_xf,
real(dp), dimension(:), intent(in)  cxs_h,
real(dp), dimension(:), intent(in)  cxs_rf,
real(dp), intent(in)  unitconv,
integer(i4b), intent(in)  icalcmeth 
)

Method to calculate the streamflow using Manning's equation for a single reach defined by a cross section.

Parameters
[in]depthreach depth
[in]widthreach width
[in]roughmannings reach roughness coefficient
[in]slopereach bottom slope
[in]unitconvconversion unit numerator in mannings equation
Returns
calculated mannings flow

Definition at line 209 of file SwfCxsUtils.f90.

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
Here is the call graph for this function:
Here is the caller graph for this function:

◆ calc_qman_by_section()

real(dp) function swfcxsutilsmodule::calc_qman_by_section ( real(dp), intent(in)  depth,
real(dp), intent(in)  width,
real(dp), intent(in)  rough,
real(dp), intent(in)  slope,
real(dp), dimension(:), intent(in)  cxs_xf,
real(dp), dimension(:), intent(in)  cxs_h,
real(dp), dimension(:), intent(in)  cxs_rf,
real(dp), intent(in)  unitconv 
)
private
Parameters
[in]depthreach depth
[in]widthreach width
[in]roughmannings reach roughness coefficient
[in]slopereach bottom slope
[in]unitconvconversion unit numerator in mannings equation
Returns
calculated mannings flow

Definition at line 361 of file SwfCxsUtils.f90.

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
Here is the call graph for this function:
Here is the caller graph for this function:

◆ calc_qman_composite()

real(dp) function swfcxsutilsmodule::calc_qman_composite ( real(dp), intent(in)  depth,
real(dp), intent(in)  width,
real(dp), intent(in)  rough,
real(dp), intent(in)  slope,
real(dp), dimension(:), intent(in)  cxs_xf,
real(dp), dimension(:), intent(in)  cxs_h,
real(dp), dimension(:), intent(in)  cxs_rf,
real(dp), intent(in)  unitconv,
integer(i4b), intent(in)  linmeth 
)
private
Parameters
[in]depthreach depth
[in]widthreach width
[in]roughmannings reach roughness coefficient
[in]slopereach bottom slope
[in]unitconvconversion unit numerator in mannings equation
[in]linmethmethod for composite (0) linear (1) nonlinear
Returns
calculated mannings flow

Definition at line 246 of file SwfCxsUtils.f90.

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
Here is the call graph for this function:
Here is the caller graph for this function:

◆ get_composite_conveyance()

real(dp) function, public swfcxsutilsmodule::get_composite_conveyance ( integer(i4b), intent(in)  npts,
real(dp), dimension(npts), intent(in)  xfraction,
real(dp), dimension(npts), intent(in)  heights,
real(dp), dimension(:), intent(in)  cxs_rf,
real(dp), intent(in)  width,
real(dp), intent(in)  rough,
real(dp), intent(in)  d 
)

Function to calculate the composite conveyance. This is based on the approach in HEC-RAS, where a conveyance is calculated for each line segment in the cross section, and then summed to produce a total conveyance.

Parameters
[in]nptsnumber of station-height data for a reach
[in]xfractioncross-section station distances (x-distance)
[in]heightscross-section height data
[in]widthreach width
[in]roughreach roughness
[in]ddepth to evaluate cross-section

Definition at line 672 of file SwfCxsUtils.f90.

674  ! -- dummy variables
675  integer(I4B), intent(in) :: npts !< number of station-height data for a reach
676  real(DP), dimension(npts), intent(in) :: xfraction !< cross-section station distances (x-distance)
677  real(DP), dimension(npts), intent(in) :: heights !< cross-section height data
678  real(DP), dimension(:), intent(in) :: cxs_rf ! mannings fractions for this cross section
679  real(DP), intent(in) :: width !< reach width
680  real(DP), intent(in) :: rough !< reach roughness
681  real(DP), intent(in) :: d !< depth to evaluate cross-section
682  ! -- local variables
683  integer(I4B) :: n
684  real(DP) :: c
685  real(DP) :: p
686  real(DP) :: rc
687  real(DP) :: rh
688  real(DP) :: cn
689  real(DP), dimension(npts) :: stations
690  real(DP), dimension(npts - 1) :: areas
691  real(DP), dimension(npts - 1) :: perimeters
692  !
693  ! -- calculate station from xfractions and width
694  do n = 1, npts
695  stations(n) = xfraction(n) * width
696  end do
697  !
698  ! -- calculate the cross-sectional area for each line segment
699  call get_cross_section_areas(npts, stations, heights, d, areas)
700  !
701  ! -- calculate the wetted perimeter for each line segment
702  call get_wetted_perimeters(npts, stations, heights, d, perimeters)
703  !
704  ! -- calculate the composite conveyance
705  c = dzero
706  do n = 1, npts - 1
707  rc = cxs_rf(n) * rough
708  p = perimeters(n)
709  cn = dzero
710  if (p > dzero) then
711  rh = areas(n) / p
712  cn = areas(n) / rc * rh**dtwothirds
713  end if
714  c = c + cn
715  end do
Here is the call graph for this function:
Here is the caller graph for this function:

◆ get_conveyance()

real(dp) function, public swfcxsutilsmodule::get_conveyance ( integer(i4b), intent(in)  npts,
real(dp), dimension(npts), intent(in)  xfraction,
real(dp), dimension(npts), intent(in)  heights,
real(dp), dimension(:), intent(in)  cxs_rf,
real(dp), intent(in)  width,
real(dp), intent(in)  rough,
real(dp), intent(in)  d 
)

Return a calculated conveyance. Calculation depends on whether the channel is a rectangular channel

Parameters
[in]nptsnumber of station-height data for a reach
[in]xfractioncross-section station distances (x-distance)
[in]heightscross-section height data
[in]widthreach width
[in]roughreach roughness
[in]ddepth to evaluate cross-section

Definition at line 574 of file SwfCxsUtils.f90.

576  ! -- dummy variables
577  integer(I4B), intent(in) :: npts !< number of station-height data for a reach
578  real(DP), dimension(npts), intent(in) :: xfraction !< cross-section station distances (x-distance)
579  real(DP), dimension(npts), intent(in) :: heights !< cross-section height data
580  real(DP), dimension(:), intent(in) :: cxs_rf ! mannings fractions for this cross section
581  real(DP), intent(in) :: width !< reach width
582  real(DP), intent(in) :: rough !< reach roughness
583  real(DP), intent(in) :: d !< depth to evaluate cross-section
584  ! result
585  real(DP) :: c
586  if (is_rectangular(xfraction)) then
587  c = get_rectangular_conveyance(npts, xfraction, heights, cxs_rf, &
588  width, rough, d)
589  else
590  c = get_composite_conveyance(npts, xfraction, heights, cxs_rf, &
591  width, rough, d)
592  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ get_cross_section_area()

real(dp) function, public swfcxsutilsmodule::get_cross_section_area ( integer(i4b), intent(in)  npts,
real(dp), dimension(npts), intent(in)  xfraction,
real(dp), dimension(npts), intent(in)  heights,
real(dp), intent(in)  width,
real(dp), intent(in)  d 
)

Function to calculate the cross-sectional area for a reach using the cross-section station-height data given a passed depth.

Returns
a cross-sectional area
Parameters
[in]nptsnumber of station-height data for a reach
[in]xfractioncross-section station distances (x-distance)
[in]heightscross-section height data
[in]widthreach width
[in]ddepth to evaluate cross-section

Definition at line 539 of file SwfCxsUtils.f90.

540  ! -- dummy variables
541  integer(I4B), intent(in) :: npts !< number of station-height data for a reach
542  real(DP), dimension(npts), intent(in) :: xfraction !< cross-section station distances (x-distance)
543  real(DP), dimension(npts), intent(in) :: heights !< cross-section height data
544  real(DP), intent(in) :: width !< reach width
545  real(DP), intent(in) :: d !< depth to evaluate cross-section
546  ! -- local variables
547  integer(I4B) :: n
548  real(DP) :: a
549  real(DP), dimension(npts) :: stations
550  real(DP), dimension(npts - 1) :: areas
551  !
552  ! -- calculate station from xfractions and width
553  do n = 1, npts
554  stations(n) = xfraction(n) * width
555  end do
556  !
557  ! -- initialize the area
558  a = dzero
559  !
560  ! -- calculate the cross-sectional area for each line segment
561  call get_cross_section_areas(npts, stations, heights, d, areas)
562  !
563  ! -- calculate the cross-sectional area
564  do n = 1, npts - 1
565  a = a + areas(n)
566  end do
Here is the call graph for this function:
Here is the caller graph for this function:

◆ get_cross_section_areas()

subroutine swfcxsutilsmodule::get_cross_section_areas ( integer(i4b), intent(in)  npts,
real(dp), dimension(npts), intent(in)  stations,
real(dp), dimension(npts), intent(in)  heights,
real(dp), intent(in)  d,
real(dp), dimension(npts - 1), intent(inout)  a 
)
private

Subroutine to calculate the cross-sectional area for each line segment that defines the reach using the cross-section station-height data given a passed depth.

Parameters
[in]nptsnumber of station-height data for a reach
[in]stationscross-section station distances (x-distance)
[in]heightscross-section height data
[in]ddepth to evaluate cross-section
[in,out]across-sectional area for each line segment

Definition at line 937 of file SwfCxsUtils.f90.

938  ! -- dummy variables
939  integer(I4B), intent(in) :: npts !< number of station-height data for a reach
940  real(DP), dimension(npts), intent(in) :: stations !< cross-section station distances (x-distance)
941  real(DP), dimension(npts), intent(in) :: heights !< cross-section height data
942  real(DP), intent(in) :: d !< depth to evaluate cross-section
943  real(DP), dimension(npts - 1), intent(inout) :: a !< cross-sectional area for each line segment
944  ! -- local variables
945  integer(I4B) :: n
946  real(DP) :: x0
947  real(DP) :: x1
948  real(DP) :: d0
949  real(DP) :: d1
950  real(DP) :: dmax
951  real(DP) :: dmin
952  real(DP) :: xlen
953  !
954  ! -- iterate over the station-height data
955  do n = 1, npts - 1
956  !
957  ! -- initialize the cross-sectional area
958  a(n) = dzero
959  !
960  ! -- initialize station-height data for segment
961  x0 = stations(n)
962  x1 = stations(n + 1)
963  d0 = heights(n)
964  d1 = heights(n + 1)
965  !
966  ! -- get the start and end station position of the wetted segment
967  call get_wetted_station(x0, x1, d0, d1, dmax, dmin, d)
968  !
969  ! -- calculate the cross-sectional area for the segment
970  xlen = x1 - x0
971  if (xlen > dzero) then
972  !
973  ! -- add the area above dmax
974  if (d > dmax) then
975  a(n) = xlen * (d - dmax)
976  end if
977  !
978  ! -- add the area below dmax
979  if (dmax /= dmin .and. d > dmin) then
980  if (d < dmax) then
981  a(n) = a(n) + dhalf * (d - dmin) * xlen
982  else
983  a(n) = a(n) + dhalf * (dmax - dmin) * xlen
984  end if
985  end if
986  end if
987  end do
Here is the call graph for this function:
Here is the caller graph for this function:

◆ get_hydraulic_radius()

real(dp) function, public swfcxsutilsmodule::get_hydraulic_radius ( integer(i4b), intent(in)  npts,
real(dp), dimension(npts), intent(in)  stations,
real(dp), dimension(npts), intent(in)  heights,
real(dp), intent(in)  d 
)

Function to calculate the hydraulic radius for a reach using the cross-section station-height data given a passed depth.

Returns
r hydraulic radius
Parameters
[in]nptsnumber of station-height data for a reach
[in]stationscross-section station distances (x-distance)
[in]heightscross-section height data
[in]ddepth to evaluate cross-section

Definition at line 725 of file SwfCxsUtils.f90.

726  ! -- dummy variables
727  integer(I4B), intent(in) :: npts !< number of station-height data for a reach
728  real(DP), dimension(npts), intent(in) :: stations !< cross-section station distances (x-distance)
729  real(DP), dimension(npts), intent(in) :: heights !< cross-section height data
730  real(DP), intent(in) :: d !< depth to evaluate cross-section
731  ! -- local variables
732  integer(I4B) :: n
733  real(DP) :: r
734  real(DP) :: p
735  real(DP) :: a
736  real(DP), dimension(npts - 1) :: areas
737  real(DP), dimension(npts - 1) :: perimeters
738  !
739  ! -- initialize the hydraulic radius, perimeter, and area
740  r = dzero
741  p = dzero
742  a = dzero
743  !
744  ! -- calculate the wetted perimeter for each line segment
745  call get_wetted_perimeters(npts, stations, heights, d, perimeters)
746  !
747  ! -- calculate the wetted perimenter
748  do n = 1, npts - 1
749  p = p + perimeters(n)
750  end do
751  !
752  ! -- calculate the hydraulic radius only if the perimeter is non-zero
753  if (p > dzero) then
754  !
755  ! -- calculate the cross-sectional area for each line segment
756  call get_cross_section_areas(npts, stations, heights, d, areas)
757  !
758  ! -- calculate the cross-sectional area
759  do n = 1, npts - 1
760  a = a + areas(n)
761  end do
762  !
763  ! -- calculate the hydraulic radius
764  r = a / p
765  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ get_hydraulic_radius_xf()

real(dp) function, public swfcxsutilsmodule::get_hydraulic_radius_xf ( integer(i4b), intent(in)  npts,
real(dp), dimension(npts), intent(in)  xfraction,
real(dp), dimension(npts), intent(in)  heights,
real(dp), intent(in)  width,
real(dp), intent(in)  d 
)

Function to calculate the hydraulic radius for a reach using the cross-section xfraction-height data given a passed depth. This is different from get_hydraulic_radius as it requires xfraction and width instead of station.

Returns
r hydraulic radius
Parameters
[in]nptsnumber of station-height data for a reach
[in]xfractioncross-section station distances (x-distance)
[in]heightscross-section height data
[in]widthreach width
[in]ddepth to evaluate cross-section

Definition at line 777 of file SwfCxsUtils.f90.

778  ! -- dummy variables
779  integer(I4B), intent(in) :: npts !< number of station-height data for a reach
780  real(DP), dimension(npts), intent(in) :: xfraction !< cross-section station distances (x-distance)
781  real(DP), dimension(npts), intent(in) :: heights !< cross-section height data
782  real(DP), intent(in) :: width !< reach width
783  real(DP), intent(in) :: d !< depth to evaluate cross-section
784  ! -- local variables
785  integer(I4B) :: n
786  real(DP) :: r
787  real(DP), dimension(npts) :: stations
788  !
789  ! -- calculate station from xfractions and width
790  do n = 1, npts
791  stations(n) = xfraction(n) * width
792  end do
793  !
794  ! -- calculate hydraulic radius
795  r = get_hydraulic_radius(npts, stations, heights, d)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ get_mannings_section()

real(dp) function, public swfcxsutilsmodule::get_mannings_section ( integer(i4b), intent(in)  npts,
real(dp), dimension(npts), intent(in)  xfraction,
real(dp), dimension(npts), intent(in)  heights,
real(dp), dimension(npts), intent(in)  roughfracs,
real(dp), intent(in)  roughness,
real(dp), intent(in)  conv_fact,
real(dp), intent(in)  width,
real(dp), intent(in)  slope,
real(dp), intent(in)  d 
)

Function to calculate the mannings discharge for a reach by calculating the discharge for each section, which can have a unique Manning's coefficient given a passed depth.

Returns
q reach discharge
Parameters
[in]nptsnumber of station-height data for a reach
[in]xfractioncross-section station distances (x-distance)
[in]heightscross-section height data
[in]roughfracscross-section Mannings roughness fraction data
[in]roughnessbase reach roughness
[in]conv_factunit conversion factor
[in]widthreach width
[in]slopereach slope
[in]ddepth to evaluate cross-section

Definition at line 806 of file SwfCxsUtils.f90.

808  ! -- dummy variables
809  integer(I4B), intent(in) :: npts !< number of station-height data for a reach
810  real(DP), dimension(npts), intent(in) :: xfraction !< cross-section station distances (x-distance)
811  real(DP), dimension(npts), intent(in) :: heights !< cross-section height data
812  real(DP), dimension(npts), intent(in) :: roughfracs !< cross-section Mannings roughness fraction data
813  real(DP), intent(in) :: roughness !< base reach roughness
814  real(DP), intent(in) :: conv_fact !< unit conversion factor
815  real(DP), intent(in) :: width !< reach width
816  real(DP), intent(in) :: slope !< reach slope
817  real(DP), intent(in) :: d !< depth to evaluate cross-section
818  ! -- local variables
819  integer(I4B) :: n
820  real(DP) :: q
821  real(DP) :: rh
822  real(DP) :: r
823  real(DP) :: p
824  real(DP) :: a
825  real(DP), dimension(npts) :: stations
826  real(DP), dimension(npts - 1) :: areas
827  real(DP), dimension(npts - 1) :: perimeters
828  !
829  ! -- initialize the hydraulic radius, perimeter, and area
830  q = dzero
831  rh = dzero
832  r = dzero
833  p = dzero
834  a = dzero
835  !
836  ! -- calculate station from xfractions and width
837  do n = 1, npts
838  stations(n) = xfraction(n) * width
839  end do
840  !
841  ! -- calculate the wetted perimeter for each line segment
842  call get_wetted_perimeters(npts, stations, heights, d, perimeters)
843  !
844  ! -- calculate the wetted perimenter
845  do n = 1, npts - 1
846  p = p + perimeters(n)
847  end do
848  !
849  ! -- calculate the hydraulic radius only if the perimeter is non-zero
850  if (p > dzero) then
851  !
852  ! -- calculate the cross-sectional area for each line segment
853  call get_cross_section_areas(npts, stations, heights, d, areas)
854  !
855  ! -- calculate the cross-sectional area
856  do n = 1, npts - 1
857  p = perimeters(n)
858  r = roughness * roughfracs(n)
859  if (p * r > dzero) then
860  a = areas(n)
861  rh = a / p
862  q = q + conv_fact * a * rh**dtwothirds * sqrt(slope) / r
863  end if
864  end do
865  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ get_rectangular_conveyance()

real(dp) function swfcxsutilsmodule::get_rectangular_conveyance ( integer(i4b), intent(in)  npts,
real(dp), dimension(npts), intent(in)  xfraction,
real(dp), dimension(npts), intent(in)  heights,
real(dp), dimension(:), intent(in)  cxs_rf,
real(dp), intent(in)  width,
real(dp), intent(in)  rough,
real(dp), intent(in)  d 
)
private
Parameters
[in]nptsnumber of station-height data for a reach
[in]xfractioncross-section station distances (x-distance)
[in]heightscross-section height data
[in]widthreach width
[in]roughreach roughness
[in]ddepth to evaluate cross-section

Definition at line 597 of file SwfCxsUtils.f90.

599  ! -- dummy variables
600  integer(I4B), intent(in) :: npts !< number of station-height data for a reach
601  real(DP), dimension(npts), intent(in) :: xfraction !< cross-section station distances (x-distance)
602  real(DP), dimension(npts), intent(in) :: heights !< cross-section height data
603  real(DP), dimension(:), intent(in) :: cxs_rf ! mannings fractions for this cross section
604  real(DP), intent(in) :: width !< reach width
605  real(DP), intent(in) :: rough !< reach roughness
606  real(DP), intent(in) :: d !< depth to evaluate cross-section
607  ! -- local variables
608  real(DP) :: c
609  real(DP) :: a
610  real(DP) :: p
611  real(DP) :: ravg
612 
613  ! find cross sectional area and wetted perimeter
614  a = get_cross_section_area(npts, xfraction, heights, width, d)
615  p = get_wetted_perimeter(npts, xfraction, heights, width, d)
616 
617  ! calculate roughness or determine an average
618  if (has_uniform_resistance(cxs_rf)) then
619  ravg = rough * cxs_rf(1)
620  else
621  ravg = calc_composite_roughness(npts, d, width, rough, dzero, &
622  xfraction, heights, cxs_rf, 0)
623  end if
624 
625  ! make conveyance calculation
626  c = dzero
627  if (p > dzero) then
628  c = a / ravg * (a / p)**dtwothirds
629  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ get_saturated_topwidth()

real(dp) function, public swfcxsutilsmodule::get_saturated_topwidth ( integer(i4b), intent(in)  npts,
real(dp), dimension(npts), intent(in)  xfraction,
real(dp), intent(in)  width 
)

Function to calculate the maximum top width for a reach using the cross-section station data.

Returns
w saturated top width
Parameters
[in]nptsnumber of station-height data for a reach
[in]xfractioncross-section station fractional distances (x-distance)

Definition at line 435 of file SwfCxsUtils.f90.

436  ! -- dummy variables
437  integer(I4B), intent(in) :: npts !< number of station-height data for a reach
438  real(DP), dimension(npts), intent(in) :: xfraction !< cross-section station fractional distances (x-distance)
439  real(DP), intent(in) :: width
440  ! -- local variables
441  integer(I4B) :: n
442  real(DP) :: w
443  real(DP), dimension(npts) :: stations
444  !
445  ! -- calculate station from xfractions and width
446  do n = 1, npts
447  stations(n) = xfraction(n) * width
448  end do
449  !
450  ! -- calculate the saturated top width
451  if (npts > 1) then
452  w = stations(npts) - stations(1)
453  else
454  w = stations(1)
455  end if
Here is the caller graph for this function:

◆ get_wetted_perimeter()

real(dp) function, public swfcxsutilsmodule::get_wetted_perimeter ( integer(i4b), intent(in)  npts,
real(dp), dimension(npts), intent(in)  xfraction,
real(dp), dimension(npts), intent(in)  heights,
real(dp), intent(in)  width,
real(dp), intent(in)  d 
)

Function to calculate the wetted perimeter for a reach using the cross-section station-height data given a passed depth.

Returns
p wetted perimeter
Parameters
[in]nptsnumber of station-height data for a reach
[in]xfractioncross-section x fractions
[in]heightscross-section height data
[in]widthcross section width
[in]ddepth to evaluate cross-section

Definition at line 502 of file SwfCxsUtils.f90.

503  ! -- dummy variables
504  integer(I4B), intent(in) :: npts !< number of station-height data for a reach
505  real(DP), dimension(npts), intent(in) :: xfraction !< cross-section x fractions
506  real(DP), dimension(npts), intent(in) :: heights !< cross-section height data
507  real(DP), intent(in) :: width !< cross section width
508  real(DP), intent(in) :: d !< depth to evaluate cross-section
509  ! -- local variables
510  integer(I4B) :: n
511  real(DP) :: p
512  real(DP), dimension(npts) :: stations
513  real(DP), dimension(npts - 1) :: perimeters
514  !
515  ! -- initialize stations
516  do n = 1, npts
517  stations(n) = xfraction(n) * width
518  end do
519  !
520  ! -- initialize the wetted perimeter for the reach
521  p = dzero
522  !
523  ! -- calculate the wetted perimeter for each line segment
524  call get_wetted_perimeters(npts, stations, heights, d, perimeters)
525  !
526  ! -- calculate the wetted perimenter
527  do n = 1, npts - 1
528  p = p + perimeters(n)
529  end do
Here is the call graph for this function:
Here is the caller graph for this function:

◆ get_wetted_perimeters()

subroutine swfcxsutilsmodule::get_wetted_perimeters ( integer(i4b), intent(in)  npts,
real(dp), dimension(npts), intent(in)  stations,
real(dp), dimension(npts), intent(in)  heights,
real(dp), intent(in)  d,
real(dp), dimension(npts - 1), intent(inout)  p 
)
private

Subroutine to calculate the wetted perimeter for each line segment that defines the reach using the cross-section station-height data given a passed depth.

Parameters
[in]nptsnumber of station-height data for a reach
[in]stationscross-section station distances (x-distance)
[in]heightscross-section height data
[in]ddepth to evaluate cross-section
[in,out]pwetted perimeter for each line segment

Definition at line 877 of file SwfCxsUtils.f90.

878  ! -- dummy variables
879  integer(I4B), intent(in) :: npts !< number of station-height data for a reach
880  real(DP), dimension(npts), intent(in) :: stations !< cross-section station distances (x-distance)
881  real(DP), dimension(npts), intent(in) :: heights !< cross-section height data
882  real(DP), intent(in) :: d !< depth to evaluate cross-section
883  real(DP), dimension(npts - 1), intent(inout) :: p !< wetted perimeter for each line segment
884  ! -- local variables
885  integer(I4B) :: n
886  real(DP) :: x0
887  real(DP) :: x1
888  real(DP) :: d0
889  real(DP) :: d1
890  real(DP) :: dmax
891  real(DP) :: dmin
892  real(DP) :: xlen
893  real(DP) :: dlen
894  !
895  ! -- iterate over the station-height data
896  do n = 1, npts - 1
897  !
898  ! -- initialize the wetted perimeter
899  p(n) = dzero
900  !
901  ! -- initialize station-height data for segment
902  x0 = stations(n)
903  x1 = stations(n + 1)
904  d0 = heights(n)
905  d1 = heights(n + 1)
906  !
907  ! -- get the start and end station position of the wetted segment
908  call get_wetted_station(x0, x1, d0, d1, dmax, dmin, d)
909  !
910  ! -- calculate the wetted perimeter for the segment
911  xlen = x1 - x0
912  dlen = dzero
913  if (xlen > dzero) then
914  if (d > dmax) then
915  dlen = dmax - dmin
916  else
917  dlen = d - dmin
918  end if
919  else
920  if (d > dmin) then
921  dlen = min(d, dmax) - dmin
922  else
923  dlen = dzero
924  end if
925  end if
926  p(n) = sqrt(xlen**dtwo + dlen**dtwo)
927  end do
Here is the call graph for this function:
Here is the caller graph for this function:

◆ get_wetted_station()

pure subroutine swfcxsutilsmodule::get_wetted_station ( real(dp), intent(inout)  x0,
real(dp), intent(inout)  x1,
real(dp), intent(in)  d0,
real(dp), intent(in)  d1,
real(dp), intent(inout)  dmax,
real(dp), intent(inout)  dmin,
real(dp), intent(in)  d 
)
private

Subroutine to calculate the station values that define the extent of the wetted portion of the cross section for a line segment. The left (x0) and right (x1) station positions are altered if the passed depth is less than the maximum line segment depth. If the line segment is dry the left and right station are equal. Otherwise the wetted station values are equal to the full line segment or smaller if the passed depth is less than the maximum line segment depth.

Parameters
[in,out]x0left station position
[in,out]x1right station position
[in]d0depth at the left station
[in]d1depth at the right station
[in,out]dmaxmaximum depth
[in,out]dminminimum depth
[in]ddepth to evaluate cross-section

Definition at line 1041 of file SwfCxsUtils.f90.

1042  ! -- dummy variables
1043  real(DP), intent(inout) :: x0 !< left station position
1044  real(DP), intent(inout) :: x1 !< right station position
1045  real(DP), intent(in) :: d0 !< depth at the left station
1046  real(DP), intent(in) :: d1 !< depth at the right station
1047  real(DP), intent(inout) :: dmax !< maximum depth
1048  real(DP), intent(inout) :: dmin !< minimum depth
1049  real(DP), intent(in) :: d !< depth to evaluate cross-section
1050  ! -- local variables
1051  real(DP) :: xlen
1052  real(DP) :: dlen
1053  real(DP) :: slope
1054  real(DP) :: xt
1055  !
1056  ! -- calculate the minimum and maximum depth
1057  dmin = min(d0, d1)
1058  dmax = max(d0, d1)
1059  !
1060  ! -- calculate x0 and x1
1061  if (d <= dmin) then
1062  ! -- if d is less than or equal to the minimum value the
1063  ! station length (xlen) is zero
1064  x1 = x0
1065  else if (d >= dmax) then
1066  ! x0 and x1 unchanged (full wetted width)
1067  continue
1068  else
1069  ! -- if d is between dmin and dmax, station length is less
1070  ! than d1 - d0
1071  xlen = x1 - x0
1072  dlen = d1 - d0
1073  ! because of preceding checks
1074  ! we know dmin<d<dmax, dmax>dmin, dlen > 0
1075  ! no need to check for dlen == 0
1076  slope = xlen / dlen
1077  xt = x0 + slope * (d - d0)
1078  if (d0 > d1) then
1079  ! x1 unchanged
1080  x0 = xt
1081  else
1082  !x0 unchanged
1083  x1 = xt
1084  end if
1085  end if
Here is the caller graph for this function:

◆ get_wetted_topwidth()

real(dp) function, public swfcxsutilsmodule::get_wetted_topwidth ( integer(i4b), intent(in)  npts,
real(dp), dimension(npts), intent(in)  xfraction,
real(dp), dimension(npts), intent(in)  heights,
real(dp), intent(in)  width,
real(dp), intent(in)  d 
)

Function to calculate the wetted top width for a reach using the cross-section station-height data given a passed depth.

Returns
w wetted top width
Parameters
[in]nptsnumber of station-height data for a reach
[in]xfractioncross-section station distances (x-distance)
[in]heightscross-section height data
[in]widthreach width
[in]ddepth to evaluate cross-section

Definition at line 465 of file SwfCxsUtils.f90.

466  ! -- dummy variables
467  integer(I4B), intent(in) :: npts !< number of station-height data for a reach
468  real(DP), dimension(npts), intent(in) :: xfraction !< cross-section station distances (x-distance)
469  real(DP), dimension(npts), intent(in) :: heights !< cross-section height data
470  real(DP), intent(in) :: width !< reach width
471  real(DP), intent(in) :: d !< depth to evaluate cross-section
472  ! -- local variables
473  integer(I4B) :: n
474  real(DP) :: w
475  real(DP), dimension(npts) :: stations
476  real(DP), dimension(npts - 1) :: widths
477  !
478  ! -- calculate station from xfractions and width
479  do n = 1, npts
480  stations(n) = xfraction(n) * width
481  end do
482  !
483  ! -- initialize the wetted perimeter for the reach
484  w = dzero
485  !
486  ! -- calculate the wetted top width for each line segment
487  call get_wetted_topwidths(npts, stations, heights, d, widths)
488  !
489  ! -- calculate the wetted top widths
490  do n = 1, npts - 1
491  w = w + widths(n)
492  end do
Here is the call graph for this function:
Here is the caller graph for this function:

◆ get_wetted_topwidths()

subroutine swfcxsutilsmodule::get_wetted_topwidths ( integer(i4b), intent(in)  npts,
real(dp), dimension(npts), intent(in)  stations,
real(dp), dimension(npts), intent(in)  heights,
real(dp), intent(in)  d,
real(dp), dimension(npts - 1), intent(inout)  w 
)
private

Subroutine to calculate the wetted top width for each line segment that defines the reach using the cross-section station-height data given a passed depth.

Parameters
[in]nptsnumber of station-height data for a reach
[in]stationscross-section station distances (x-distance)
[in]heightscross-section height data
[in]ddepth to evaluate cross-section
[in,out]wwetted top widths for each line segment

Definition at line 997 of file SwfCxsUtils.f90.

998  ! -- dummy variables
999  integer(I4B), intent(in) :: npts !< number of station-height data for a reach
1000  real(DP), dimension(npts), intent(in) :: stations !< cross-section station distances (x-distance)
1001  real(DP), dimension(npts), intent(in) :: heights !< cross-section height data
1002  real(DP), intent(in) :: d !< depth to evaluate cross-section
1003  real(DP), dimension(npts - 1), intent(inout) :: w !< wetted top widths for each line segment
1004  ! -- local variables
1005  integer(I4B) :: n
1006  real(DP) :: x0
1007  real(DP) :: x1
1008  real(DP) :: d0
1009  real(DP) :: d1
1010  real(DP) :: dmax
1011  real(DP) :: dmin
1012  !
1013  ! -- iterate over the station-height data
1014  do n = 1, npts - 1
1015  !
1016  ! -- initialize station-height data for segment
1017  x0 = stations(n)
1018  x1 = stations(n + 1)
1019  d0 = heights(n)
1020  d1 = heights(n + 1)
1021  !
1022  ! -- get the start and end station position of the wetted segment
1023  call get_wetted_station(x0, x1, d0, d1, dmax, dmin, d)
1024  !
1025  ! -- calculate the wetted top width for the segment
1026  w(n) = x1 - x0
1027  end do
Here is the call graph for this function:
Here is the caller graph for this function:

◆ has_uniform_resistance()

logical(lgp) function swfcxsutilsmodule::has_uniform_resistance ( real(dp), dimension(:), intent(in)  cxs_rf)
private
Parameters
[in]cxs_rfcross-section station distances (x-distance)

Definition at line 650 of file SwfCxsUtils.f90.

651  ! dummy
652  real(DP), dimension(:), intent(in) :: cxs_rf !< cross-section station distances (x-distance)
653  ! result
654  logical(LGP) :: has_uniform_resistance
655  ! local
656  real(DP) :: rmin
657  real(DP) :: rmax
658  rmin = minval(cxs_rf(1:3))
659  rmax = maxval(cxs_rf(1:3))
660  has_uniform_resistance = (rmin == rmax)
Here is the caller graph for this function:

◆ is_rectangular()

logical(lgp) function swfcxsutilsmodule::is_rectangular ( real(dp), dimension(:), intent(in)  xfraction)
private
Parameters
[in]xfractioncross-section station distances (x-distance)

Definition at line 634 of file SwfCxsUtils.f90.

635  ! dummy
636  real(DP), dimension(:), intent(in) :: xfraction !< cross-section station distances (x-distance)
637  ! result
638  logical(LGP) :: is_rectangular
639  is_rectangular = .false.
640  if (size(xfraction) == 4) then
641  if (xfraction(1) == xfraction(2) .and. &
642  xfraction(3) == xfraction(4)) then
643  is_rectangular = .true.
644  end if
645  end if
Here is the caller graph for this function: