MODFLOW 6  version 6.7.0.dev0
USGS Modular Hydrologic Model
mathutilmodule Module Reference

Data Types

interface  mod_offset
 
interface  f1d
 

Functions/Subroutines

pure logical function, public is_close (a, b, rtol, atol, symmetric)
 Check if a real value is approximately equal to another. More...
 
pure integer(i4b) function mod_offset_int (a, n, d)
 Modulo with offset for integer values. More...
 
pure real(dp) function mod_offset_dbl (a, n, d)
 Modulo with offset for double precision values. More...
 
real(dp) function, public zero_ch (x0, x1, f, epsa)
 Compute zeros on an interval using Chadrupatla's method. More...
 
real(dp) function, public zero_br (ax, bx, f, tol)
 Compute a zero of the function f(x) in the interval (x0, x1). More...
 
real(dp) function, public get_perturbation (x)
 Calculate a numerical perturbation given the value of x. More...
 
pure real(dp) function, dimension(:), allocatable, public arange (start, stop, step)
 Return reals separated by the given step over the given interval. More...
 
pure real(dp) function, dimension(:), allocatable, public linspace (start, stop, num)
 Return evenly spaced reals over the given interval. More...
 

Function/Subroutine Documentation

◆ arange()

pure real(dp) function, dimension(:), allocatable, public mathutilmodule::arange ( real(dp), intent(in)  start,
real(dp), intent(in)  stop,
real(dp), intent(in), optional  step 
)

This function is designed to behave like NumPy's arange. Adapted from https://stackoverflow.com/a/65839162/6514033.

Definition at line 383 of file MathUtil.f90.

384  ! dummy
385  real(DP), intent(in) :: start, stop
386  real(DP), intent(in), optional :: step
387  real(DP), allocatable :: array(:)
388  ! local
389  real(DP) :: lstep
390  integer(I4B) :: i, n
391 
392  if (present(step)) then
393  lstep = step
394  else
395  lstep = done
396  end if
397 
398  if (step <= 0) then
399  allocate (array(0))
400  else
401  n = ceiling((stop - start) / step)
402  allocate (array(n))
403  do i = 0, n - 1
404  array(i + 1) = start + step * i
405  end do
406  end if
Here is the caller graph for this function:

◆ get_perturbation()

real(dp) function, public mathutilmodule::get_perturbation ( real(dp), intent(in)  x)

Calculate a perturbation value to use for a numerical derivative calculation. Taken from the book "Solving Nonlinear Equations with Newton's Method" 2003, by C.T. Kelley. Method also used in the SWR Process for MODFLOW-2005.

Parameters
[in]xvalue that will be perturbed by the result
Returns
calculated perturbation value

Definition at line 371 of file MathUtil.f90.

372  use constantsmodule, only: dprecsqrt, done
373  real(DP), intent(in) :: x !< value that will be perturbed by the result
374  real(DP) :: res !< calculated perturbation value
375  res = dprecsqrt * max(abs(x), done) * sign(done, x)
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dprecsqrt
Definition: Constants.f90:121
real(dp), parameter done
real constant 1
Definition: Constants.f90:76
Here is the caller graph for this function:

◆ is_close()

pure logical function, public mathutilmodule::is_close ( real(dp), intent(in)  a,
real(dp), intent(in)  b,
real(dp), intent(in), optional  rtol,
real(dp), intent(in), optional  atol,
logical(lgp), intent(in), optional  symmetric 
)

By default the determination is symmetric in a and b, as in Python's math.isclose, with relative difference scaled by a factor of the larger absolute value of a and b. The formula is: abs(a - b) <= max(rtol * max(abs(a), abs(b)), atol).

If symmetric is set to false the test is asymmetric in a and b, with b taken to be the reference value, and the alternate formula (abs(a - b) <= (atol + rtol * abs(b))) is used. This is the approach taken by numpy.allclose.

Defaults for rtol and atol are DSAME and DZERO, respectively. If a or b are near 0 (especially if either is 0), an absolute tolerance suitable for the particular case should be provided. For a justification of a zero absolute tolerance default see: https://peps.python.org/pep-0485/#absolute-tolerance-default

Parameters
[in]afirst real
[in]bsecond real (reference value if asymmetric)
[in]rtolrelative tolerance (default=DSAME)
[in]atolabsolute tolerance (default=DZERO)
[in]symmetrictoggle (a)symmetric comparison

Definition at line 45 of file MathUtil.f90.

46  ! dummy
47  real(DP), intent(in) :: a !< first real
48  real(DP), intent(in) :: b !< second real (reference value if asymmetric)
49  real(DP), intent(in), optional :: rtol !< relative tolerance (default=DSAME)
50  real(DP), intent(in), optional :: atol !< absolute tolerance (default=DZERO)
51  logical(LGP), intent(in), optional :: symmetric !< toggle (a)symmetric comparison
52  ! local
53  real(DP) :: lrtol, latol
54  logical(LGP) :: lsymmetric
55 
56  ! check for exact equality
57  if (a == b) then
58  is_close = .true.
59  return
60  end if
61 
62  ! process optional arguments
63  if (.not. present(rtol)) then
64  lrtol = dsame
65  else
66  lrtol = rtol
67  end if
68  if (.not. present(atol)) then
69  latol = dzero
70  else
71  latol = atol
72  end if
73  if (.not. present(symmetric)) then
74  lsymmetric = .true.
75  else
76  lsymmetric = symmetric
77  end if
78 
79  if (lsymmetric) then
80  ! "weak" symmetric test, https://peps.python.org/pep-0485/#which-symmetric-test
81  is_close = abs(a - b) <= max(lrtol * max(abs(a), abs(b)), latol)
82  else
83  ! asymmetric, https://numpy.org/doc/stable/reference/generated/numpy.isclose.html
84  is_close = (abs(a - b) <= (latol + lrtol * abs(b)))
85  end if
Here is the caller graph for this function:

◆ linspace()

pure real(dp) function, dimension(:), allocatable, public mathutilmodule::linspace ( real(dp), intent(in)  start,
real(dp), intent(in)  stop,
integer(i4b), intent(in)  num 
)

This function is designed to behave like NumPy's linspace. Adapted from https://stackoverflow.com/a/57211848/6514033.

Definition at line 414 of file MathUtil.f90.

415  ! dummy
416  real(DP), intent(in) :: start, stop
417  integer(I4B), intent(in) :: num
418  real(DP), allocatable :: array(:)
419  ! local
420  real(DP) :: range
421  integer(I4B) :: i
422 
423  allocate (array(num))
424  range = stop - start
425 
426  if (num == 0) return
427  if (num == 1) then
428  array(1) = start
429  return
430  end if
431  do i = 1, num
432  array(i) = start + range * (i - 1) / (num - 1)
433  end do

◆ mod_offset_dbl()

pure real(dp) function mathutilmodule::mod_offset_dbl ( real(dp), intent(in)  a,
real(dp), intent(in)  n,
real(dp), intent(in), optional  d 
)
private
Parameters
[in]adividend
[in]ndivisor
[in]doffset

Definition at line 107 of file MathUtil.f90.

108  ! -- dummy
109  real(DP), intent(in) :: a !< dividend
110  real(DP), intent(in) :: n !< divisor
111  real(DP), intent(in), optional :: d !< offset
112  real(DP) :: mo
113  ! -- local
114  real(DP) :: ld
115 
116  if (present(d)) then
117  ld = d
118  else
119  ld = 0
120  end if
121  mo = a - n * floor((a - ld) / n)

◆ mod_offset_int()

pure integer(i4b) function mathutilmodule::mod_offset_int ( integer(i4b), intent(in)  a,
integer(i4b), intent(in)  n,
integer(i4b), intent(in), optional  d 
)
private
Parameters
[in]adividend
[in]ndivisor
[in]doffset

Definition at line 89 of file MathUtil.f90.

90  ! -- dummy
91  integer(I4B), intent(in) :: a !< dividend
92  integer(I4B), intent(in) :: n !< divisor
93  integer(I4B), intent(in), optional :: d !< offset
94  integer(I4B) :: mo
95  ! -- local
96  integer(I4B) :: ld
97 
98  if (present(d)) then
99  ld = d
100  else
101  ld = 0
102  end if
103  mo = a - n * floor(real(a - ld) / n)

◆ zero_br()

real(dp) function, public mathutilmodule::zero_br ( real(dp)  ax,
real(dp)  bx,
procedure(f1d), intent(in), pointer  f,
real(dp)  tol 
)

A zero of the function f(x) is computed in the interval ax,bx.

Input:

ax left endpoint of initial interval bx right endpoint of initial interval f function subprogram which evaluates f(x) for any x in the interval ax,bx tol desired length of the interval of uncertainty of the final result (.ge.0.)

Output:

zero_br abscissa approximating a zero of f in the interval ax,bx

it is assumed  that   f(ax)   and   f(bx)   have  opposite  signs

this is checked, and an error message is printed if this is not satisfied. zero_br returns a zero x in the given interval ax,bx to within a tolerance 4*macheps*abs(x)+tol, where macheps is the relative machine precision defined as the smallest representable number such that 1.+macheps .gt. 1. this function subprogram is a slightly modified translation of the algol 60 procedure zero given in richard brent, algorithms for minimization without derivatives, prentice-hall, inc. (1973).

This subroutine was obtained by the authors of MODFLOW 6 from netlib.org with the understanding that it is freely available. It has subsequently undergone minor modification to suit the current application.

Definition at line 254 of file MathUtil.f90.

255  ! -- dummy
256  real(DP) :: ax, bx
257  procedure(f1d), pointer, intent(in) :: f
258  real(DP) :: tol
259  real(DP) :: z
260  ! -- local
261  real(DP) :: eps
262  real(DP) :: a, b, c, d, e, s, p, q
263  real(DP) :: fa, fb, fc, r, tol1, xm
264  logical(LGP) :: rs
265 
266  eps = epsilon(ax)
267  tol1 = eps + done
268 
269  a = ax
270  b = bx
271  fa = f(a)
272  fb = f(b)
273 
274  ! check if a or b is the solution
275  if (fa == dzero) then
276  z = a
277  return
278  else if (fb == dzero) then
279  z = b
280  return
281  end if
282 
283  ! check that f(ax) and f(bx) have opposite sign
284  if (fa * (fb / dabs(fb)) .ge. dzero) &
285  call pstop(1, 'f(ax) and f(bx) do not have different signs,')
286 
287  rs = .true. ! var reset
288  do while (.true.)
289  if (rs) then
290  c = a
291  fc = fa
292  d = b - a
293  e = d
294  end if
295 
296  if (.not. (dabs(fc) .ge. dabs(fb))) then
297  a = b
298  b = c
299  c = a
300  fa = fb
301  fb = fc
302  fc = fa
303  end if
304 
305  tol1 = dtwo * eps * dabs(b) + dhalf * tol
306  xm = dhalf * (c - b)
307  if ((dabs(xm) .le. tol1) .or. (fb .eq. dzero)) then
308  z = b
309  return
310  end if
311 
312  ! see if a bisection is forced
313  if ((dabs(e) .ge. tol1) .and. (dabs(fa) .gt. dabs(fb))) then
314  s = fb / fa
315  if (a .ne. c) then
316  ! inverse quadratic interpolation
317  q = fa / fc
318  r = fb / fc
319  p = s * (dtwo * xm * q * (q - r) - (b - a) * (r - done))
320  q = (q - done) * (r - done) * (s - done)
321  else
322  ! linear interpolation
323  p = dtwo * xm * s
324  q = done - s
325  end if
326 
327  if (p .le. dzero) then
328  p = -p
329  else
330  q = -q
331  end if
332  s = e
333  e = d
334  if (((dtwo * p) .ge. (dthree * xm * q - dabs(tol1 * q))) .or. &
335  (p .ge. dabs(dhalf * s * q))) then
336  d = xm
337  e = d
338  else
339  d = p / q
340  end if
341  else
342  d = xm
343  e = d
344  end if
345 
346  a = b
347  fa = fb
348 
349  if (dabs(d) .le. tol1) then
350  if (xm .le. dzero) then
351  b = b - tol1
352  else
353  b = b + tol1
354  end if
355  else
356  b = b + d
357  end if
358 
359  fb = f(b)
360  rs = (fb * (fc / dabs(fc))) .gt. dzero
361  end do
Here is the call graph for this function:
Here is the caller graph for this function:

◆ zero_ch()

real(dp) function, public mathutilmodule::zero_ch ( real(dp)  x0,
real(dp)  x1,
procedure(f1d), intent(in), pointer  f,
real(dp)  epsa 
)

A zero of the function f{x} is computed in the interval (x0, x1) given tolerance epsa using Chandrupatla's method. FORTRAN code based generally on pseudocode in Scherer, POJ (2013) "Computational Physics: Simulation of Classical and Quantum Systems," 2nd ed., Springer, New York.

Definition at line 131 of file MathUtil.f90.

132  ! -- dummy
133  real(DP) :: x0, x1
134  procedure(f1d), pointer, intent(in) :: f
135  real(DP) :: epsa
136  real(DP) :: z
137  ! -- local
138  real(DP) :: epsm
139  real(DP) :: a, b, c, t
140  real(DP) :: aminusb, cminusb
141  real(DP) :: fa, fb, fc, fm, ft
142  real(DP) :: faminusfb, fcminusfb
143  real(DP) :: phi, philo, phihi
144  real(DP) :: racb, rcab, rbca
145  real(DP) :: tol, tl, tlc
146  real(DP) :: xi, xm, xt
147 
148  epsm = epsilon(x0)
149  b = x0
150  a = x1
151  c = x1
152  aminusb = a - b
153  fb = f(b)
154  fa = f(a)
155  fc = f(c)
156  t = 5d-1
157 
158  ! check whether a or b is the solution
159  if (fa == dzero) then
160  z = a
161  return
162  else if (fb == dzero) then
163  z = b
164  return
165  end if
166 
167  do while (.true.)
168  xt = a - t * aminusb
169  ft = f(xt)
170  if (sign(ft, fa) == ft) then
171  c = a
172  fc = fa
173  a = xt
174  fa = ft
175  else
176  c = b
177  b = a
178  a = xt
179  fc = fb
180  fb = fa
181  fa = ft
182  end if
183  aminusb = a - b
184  cminusb = c - b
185  faminusfb = fa - fb
186  fcminusfb = fc - fb
187  xm = a
188  fm = fa
189  if (dabs(fb) < dabs(fa)) then
190  xm = b
191  fm = fb
192  end if
193  tol = dtwo * epsm * dabs(xm) + epsa
194  tl = tol / dabs(cminusb)
195  if ((tl > 5d-1) .or. (fm == dzero)) then
196  z = xm
197  return
198  end if
199  xi = aminusb / cminusb
200  phi = faminusfb / fcminusfb
201  philo = done - dsqrt(done - xi)
202  phihi = dsqrt(xi)
203  if ((phi > philo) .and. (phi < phihi)) then
204  racb = fa / fcminusfb
205  rcab = fc / faminusfb
206  rbca = fb / (fc - fa)
207  t = racb * (rcab - rbca * (c - a) / aminusb)
208  if (t < tl) then
209  t = tl
210  else
211  tlc = done - tl
212  if (t > tlc) then
213  t = tlc
214  end if
215  end if
216  else
217  t = 5d-1
218  end if
219  end do
Here is the caller graph for this function: