Skip to contents

heads() computes the hydraulic head at the given x and y coordinates for an aem object.

omega() computes the complex potential for an aem or element object at the given x and y coordinates.

potential() computes the discharge potential for an aem or element object at the given x and y coordinates.

streamfunction() computes the stream function for an aem or element object at the given x and y coordinates.

Usage

heads(aem, x, y, as.grid = FALSE, na.below = TRUE, ...)

omega(...)

potential(...)

streamfunction(...)

# S3 method for class 'aem'
omega(aem, x, y, as.grid = FALSE, ...)

# S3 method for class 'aem'
potential(aem, x, y, as.grid = FALSE, ...)

# S3 method for class 'aem'
streamfunction(aem, x, y, as.grid = FALSE, ...)

# S3 method for class 'element'
omega(element, x, y, ...)

# S3 method for class 'element'
potential(element, x, y, ...)

# S3 method for class 'element'
streamfunction(element, x, y, ...)

Arguments

aem

aem object.

x

numeric x coordinates to evaluate the variable at.

y

numeric y coordinates to evaluate the variable at.

as.grid

logical, should a matrix be returned? Defaults to FALSE. See details.

na.below

logical indicating if calculated head values below the aquifer base should be set to NA. Defaults to TRUE. See potential_to_head().

...

ignored

element

analytic element of class element.

Value

For heads(), a vector of length(x) (equal to length(y)) with the hydraulic head values at x and y. If as.grid = TRUE, a matrix of dimensions c(length(y), length(x)) described by marginal vectors x and y containing the hydraulic head values at the grid points. The heads are computed from potential() and the aquifer parameters using potential_to_head().

For omega(), the same as for heads() but containing the complex potential values evaluated at x and y.

For potential(), the same as for heads() but containing the discharge potential values evaluated at x and y, which are the real components of omega().

For streamfunction(), the same as for heads() but containing the stream function values evaluated at x and y, which are the imaginary components of omega().

Details

heads() should not to be confused with utils::head(), which returns the first part of an object.

Examples

w <- well(xw = 55, yw = 0, Q = 200)
uf <- uniformflow(gradient = 0.002, angle = -45, TR = 100)
rf <- constant(xc = -1000, yc = 1000, hc = 10)
ml <- aem(k = 10, top = 10, base = -15, n = 0.2, w, uf, rf)

xg <- seq(-100, 100, length = 5)
yg <- seq(-75, 75, length = 3)

# Hydraulic heads
heads(ml, c(50, 0), c(25, -25))
#> [1] 8.280544 8.398209
heads(ml, xg, yg, as.grid = TRUE)
#>          [,1]     [,2]     [,3]     [,4]     [,5]
#> [1,] 8.660110 8.591307 8.517038 8.458073 8.448488
#> [2,] 8.601056 8.518423 8.400544 8.041392 8.312647
#> [3,] 8.570281 8.501215 8.426661 8.367468 8.357846

# do not confuse heads() with utils::head, which will give an error
try(
head(ml, c(50, 0), c(25, -25))
)
#> Error in head.default(ml, c(50, 0), c(25, -25)) : 
#>   invalid 'n' - must have length one when dim(x) is NULL, got 2

# Complex potential
omega(ml, c(50, 0), c(25, -25))
#> [1] 2709.919+45.67669i 2737.381-82.88449i

# Discharge potential
potential(ml, c(50, 0), c(25, -25))
#> [1] 2709.919 2737.381

# Stream function
streamfunction(ml, c(50, 0), c(25, -25))
#> [1]  45.67669 -82.88449

# For elements
omega(w, c(50, 0), c(-25, 25))
#> [1] 103.0842-56.28330i 130.5466+86.42003i

potential(w, c(50, 0), c(-25, 25))
#> [1] 103.0842 130.5466

streamfunction(w, c(50, 0), c(-25, 25))
#> [1] -56.28330  86.42003