Compute the gradually-varied flow profile of a prismatic channel.

compute_profile(So, n, Q, y0, Cm, g, B, SS, z0 = 0, x0 = 0, stepdist,
  totaldist)

Arguments

So

Channel slope [\(L L^{-1}\)].

n

Manning's roughness coefficient.

Q

Flow rate [\(L^3 T^{-1}\)].

y0

The water depth at the control section [\(L\)].

Cm

Unit conversion coefficient for Manning's equation. For SI units, Cm = 1.

g

Gravitational acceleration [\(L T^{-2}\)].

B

Channel bottom width [\(L\)].

SS

Channel sideslope [\(L L^{-1}\)].

z0

Elevation reference datum at control section [\(L\)]. Default is 0.

x0

Distance reference at control section [\(L\)]. Default is 0.

stepdist

The spatial interval used in the Standard step method [\(L\)].

totaldist

The total distance upstream (or downstream) to compute the profile [\(L\)].

Value

data.frame with columns:

x

Along-channel distance.

z

Elevation.

y

Flow depth.

v

Flow velocity.

A

Flow area.

Sf

Friction slope.

E

Total energy.

Fr

Froude Number.

Details

Computes the longitudinal water surface profile of a prismatic channel using the standard step method by solving the non-linear ODE $$\frac{dy}{dx} = \frac{S_0 - S_f}{1 - Fr^2}$$ The standard-step method operates by stepping along the channel by a constant distance interval, starting from a cross-section where the flow depth is known (the control section). The flow depth is computed at the adjacent cross-section (target section). The computed value at the target is then used as the basis for computing flow depth at the next cross-section, i.e. the previous target section becomes the new control section for each step. A Newton-Raphson scheme is used each step to compute the flow depth and friction slope. Technically, the average friction slope of the control and target section is used to compute the flow depth at the target section.

Examples

# example M1 profile compute_profile(0.001, 0.045, 250, 2.7, 1.486, 32.2, 100, 0, stepdist = 10, totaldist = 3000)
#> M1 profile specified. Computing upstream profile
#> Simulation type: #> Gradually-varied flow #> #> Call: #> compute_profile(So = 0.001, n = 0.045, Q = 250, y0 = 2.7, Cm = 1.486, #> g = 32.2, B = 100, SS = 0, stepdist = 10, totaldist = 3000) #> #> #> Channel geometry: #> So:0.001 #> n :0.045 #> B :100 #> SS:0 #> #> Model specification: #> delta.x :10 #> channel.length:3000 #> normal.depth :1.71 #> critical.depth:0.579 #> #> Data: #> x : num [1:301] 0 -10 -20 -30 -40 -50 -60 -70 -80 -90 ... #> z : num [1:301] 0 0.01 0.02 0.03 0.04 0.05 0.06 0.07 0.08 0.09 ... #> y : num [1:301] 2.7 2.69 2.68 2.68 2.67 ... #> v : num [1:301] 0.926 0.929 0.931 0.934 0.937 ... #> A : num [1:301] 270 269 268 268 267 ... #> Sf: num [1:301] 0.000224 0.000226 0.000229 0.000231 0.000233 ... #> E : num [1:301] 2.71 2.72 2.72 2.72 2.72 ... #> Fr: num [1:301] 0.0993 0.0997 0.1002 0.1006 0.101 ...
# example M2 profile compute_profile(0.001, 0.045, 250, 0.64, 1.486, 32.2, 100, 0, stepdist = 10, totaldist = 3000)
#> M2 profile specified. Computing upstream profile
#> Simulation type: #> Gradually-varied flow #> #> Call: #> compute_profile(So = 0.001, n = 0.045, Q = 250, y0 = 0.64, Cm = 1.486, #> g = 32.2, B = 100, SS = 0, stepdist = 10, totaldist = 3000) #> #> #> Channel geometry: #> So:0.001 #> n :0.045 #> B :100 #> SS:0 #> #> Model specification: #> delta.x :10 #> channel.length:3000 #> normal.depth :1.71 #> critical.depth:0.579 #> #> Data: #> x : num [1:301] 0 -10 -20 -30 -40 -50 -60 -70 -80 -90 ... #> z : num [1:301] 0 0.01 0.02 0.03 0.04 0.05 0.06 0.07 0.08 0.09 ... #> y : num [1:301] 0.64 0.92 0.995 1.05 1.094 ... #> v : num [1:301] 3.91 2.72 2.51 2.38 2.29 ... #> A : num [1:301] 64 92 99.5 105 109.4 ... #> Sf: num [1:301] 0.0258 0.00775 0.00598 0.005 0.00437 ... #> E : num [1:301] 0.877 1.045 1.113 1.168 1.215 ... #> Fr: num [1:301] 0.86 0.499 0.444 0.409 0.385 ...
# example S2 profile compute_profile(0.005, 0.01, 250, 2.65, 1.486, 32.2, 10, 0, stepdist = 10, totaldist = 2000)
#> S2 profile specified. Computing downstream profile
#> Simulation type: #> Gradually-varied flow #> #> Call: #> compute_profile(So = 0.005, n = 0.01, Q = 250, y0 = 2.65, Cm = 1.486, #> g = 32.2, B = 10, SS = 0, stepdist = 10, totaldist = 2000) #> #> #> Channel geometry: #> So:0.005 #> n :0.01 #> B :10 #> SS:0 #> #> Model specification: #> delta.x :10 #> channel.length:2000 #> normal.depth :1.92 #> critical.depth:2.69 #> #> Data: #> x : num [1:201] 0 10 20 30 40 50 60 70 80 90 ... #> z : num [1:201] 0 -0.05 -0.1 -0.15 -0.2 -0.25 -0.3 -0.35 -0.4 -0.45 ... #> y : num [1:201] 2.65 2.47 2.4 2.35 2.31 ... #> v : num [1:201] 9.43 10.12 10.42 10.65 10.84 ... #> A : num [1:201] 26.5 24.7 24 23.5 23.1 ... #> Sf: num [1:201] 0.00194 0.00237 0.00259 0.00275 0.00289 ... #> E : num [1:201] 4.03 4.01 3.99 3.96 3.93 ... #> Fr: num [1:201] 1.02 1.13 1.19 1.23 1.26 ...
# example S3 profile compute_profile(0.005, 0.01, 250, 0.5, 1.486, 32.2, 10, 0, stepdist = 10, totaldist = 2000)
#> S3 profile specified. Computing downstream profile
#> Simulation type: #> Gradually-varied flow #> #> Call: #> compute_profile(So = 0.005, n = 0.01, Q = 250, y0 = 0.5, Cm = 1.486, #> g = 32.2, B = 10, SS = 0, stepdist = 10, totaldist = 2000) #> #> #> Channel geometry: #> So:0.005 #> n :0.01 #> B :10 #> SS:0 #> #> Model specification: #> delta.x :10 #> channel.length:2000 #> normal.depth :1.92 #> critical.depth:2.69 #> #> Data: #> x : num [1:201] 0 10 20 30 40 50 60 70 80 90 ... #> z : num [1:201] 0 -0.05 -0.1 -0.15 -0.2 -0.25 -0.3 -0.35 -0.4 -0.45 ... #> y : num [1:201] 0.5 0.521 0.541 0.561 0.581 ... #> v : num [1:201] 50 48 46.2 44.5 43 ... #> A : num [1:201] 5 5.21 5.41 5.61 5.81 ... #> Sf: num [1:201] 0.324 0.285 0.252 0.224 0.2 ... #> E : num [1:201] 39.3 36.3 33.6 31.2 29.1 ... #> Fr: num [1:201] 12.46 11.73 11.07 10.48 9.94 ...