MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
MathUtil.f90
Go to the documentation of this file.
2  use kindmodule, only: dp, i4b, lgp
3  use errorutilmodule, only: pstop
7 
8  implicit none
9  private
10  public :: f1d, is_close, mod_offset, zero_ch, zero_br, &
12 
13  interface mod_offset
14  module procedure :: mod_offset_int, mod_offset_dbl
15  end interface mod_offset
16 
17  interface
18  function f1d(x) result(fx)
19  import dp
20  real(dp), intent(in) :: x
21  real(dp) :: fx
22  end function
23  end interface
24 
25 contains
26 
27  !> @brief Check if a real value is approximately equal to another.
28  !!
29  !! By default the determination is symmetric in a and b, as in
30  !! Python's math.isclose, with relative difference scaled by a
31  !! factor of the larger absolute value of a and b. The formula
32  !! is: abs(a - b) <= max(rtol * max(abs(a), abs(b)), atol).
33  !!
34  !! If symmetric is set to false the test is asymmetric in a and
35  !! b, with b taken to be the reference value, and the alternate
36  !! formula (abs(a - b) <= (atol + rtol * abs(b))) is used. This
37  !! is the approach taken by numpy.allclose.
38  !!
39  !! Defaults for rtol and atol are DSAME and DZERO, respectively.
40  !! If a or b are near 0 (especially if either is 0), an absolute
41  !! tolerance suitable for the particular case should be provided.
42  !! For a justification of a zero absolute tolerance default see:
43  !! https://peps.python.org/pep-0485/#absolute-tolerance-default
44  !<
45  pure logical function is_close(a, b, rtol, atol, symmetric)
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
86  end function is_close
87 
88  !> @brief Modulo with offset for integer values.
89  pure function mod_offset_int(a, n, d) result(mo)
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)
104  end function mod_offset_int
105 
106  !> @brief Modulo with offset for double precision values.
107  pure function mod_offset_dbl(a, n, d) result(mo)
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)
122  end function mod_offset_dbl
123 
124  !> @brief Compute zeros on an interval using Chadrupatla's method
125  !!
126  !! A zero of the function f{x} is computed in the interval (x0, x1)
127  !! given tolerance epsa using Chandrupatla's method. FORTRAN code based
128  !! generally on pseudocode in Scherer, POJ (2013) "Computational Physics:
129  !! Simulation of Classical and Quantum Systems," 2nd ed., Springer, New York.
130  !!
131  !<
132  function zero_ch(x0, x1, f, epsa) result(z)
133  ! -- dummy
134  real(dp) :: x0, x1
135  procedure(f1d), pointer, intent(in) :: f
136  real(dp) :: epsa
137  real(dp) :: z
138  ! -- local
139  real(dp) :: epsm
140  real(dp) :: a, b, c, t
141  real(dp) :: aminusb, cminusb
142  real(dp) :: fa, fb, fc, fm, ft
143  real(dp) :: faminusfb, fcminusfb
144  real(dp) :: phi, philo, phihi
145  real(dp) :: racb, rcab, rbca
146  real(dp) :: tol, tl, tlc
147  real(dp) :: xi, xm, xt
148 
149  epsm = epsilon(x0)
150  b = x0
151  a = x1
152  c = x1
153  aminusb = a - b
154  fb = f(b)
155  fa = f(a)
156  fc = f(c)
157  t = 5d-1
158 
159  ! check whether a or b is the solution
160  if (fa == dzero) then
161  z = a
162  return
163  else if (fb == dzero) then
164  z = b
165  return
166  end if
167 
168  do while (.true.)
169  xt = a - t * aminusb
170  ft = f(xt)
171  if (sign(ft, fa) == ft) then
172  c = a
173  fc = fa
174  a = xt
175  fa = ft
176  else
177  c = b
178  b = a
179  a = xt
180  fc = fb
181  fb = fa
182  fa = ft
183  end if
184  aminusb = a - b
185  cminusb = c - b
186  faminusfb = fa - fb
187  fcminusfb = fc - fb
188  xm = a
189  fm = fa
190  if (dabs(fb) < dabs(fa)) then
191  xm = b
192  fm = fb
193  end if
194  tol = dtwo * epsm * dabs(xm) + epsa
195  tl = tol / dabs(cminusb)
196  if ((tl > 5d-1) .or. (fm == dzero)) then
197  z = xm
198  return
199  end if
200  xi = aminusb / cminusb
201  phi = faminusfb / fcminusfb
202  philo = done - dsqrt(done - xi)
203  phihi = dsqrt(xi)
204  if ((phi > philo) .and. (phi < phihi)) then
205  racb = fa / fcminusfb
206  rcab = fc / faminusfb
207  rbca = fb / (fc - fa)
208  t = racb * (rcab - rbca * (c - a) / aminusb)
209  if (t < tl) then
210  t = tl
211  else
212  tlc = done - tl
213  if (t > tlc) then
214  t = tlc
215  end if
216  end if
217  else
218  t = 5d-1
219  end if
220  end do
221  end function
222 
223  !> @brief Compute a zero of the function f(x) in the interval (x0, x1).
224  !!
225  !! A zero of the function f(x) is computed in the interval ax,bx.
226  !!
227  !! Input:
228  !!
229  !! ax left endpoint of initial interval
230  !! bx right endpoint of initial interval
231  !! f function subprogram which evaluates f(x) for any x in
232  !! the interval ax,bx
233  !! tol desired length of the interval of uncertainty of the
234  !! final result (.ge.0.)
235  !!
236  !! Output:
237  !!
238  !! zero_br abscissa approximating a zero of f in the interval ax,bx
239  !!
240  !! it is assumed that f(ax) and f(bx) have opposite signs
241  !! this is checked, and an error message is printed if this is not
242  !! satisfied. zero_br returns a zero x in the given interval
243  !! ax,bx to within a tolerance 4*macheps*abs(x)+tol, where macheps is
244  !! the relative machine precision defined as the smallest representable
245  !! number such that 1.+macheps .gt. 1.
246  !! this function subprogram is a slightly modified translation of
247  !! the algol 60 procedure zero given in richard brent, algorithms for
248  !! minimization without derivatives, prentice-hall, inc. (1973).
249  !<
250  function zero_br(ax, bx, f, tol) result(z)
251  ! -- dummy
252  real(dp) :: ax, bx
253  procedure(f1d), pointer, intent(in) :: f
254  real(dp) :: tol
255  real(dp) :: z
256  ! -- local
257  real(dp) :: eps
258  real(dp) :: a, b, c, d, e, s, p, q
259  real(dp) :: fa, fb, fc, r, tol1, xm
260  logical(LGP) :: rs
261 
262  eps = epsilon(ax)
263  tol1 = eps + done
264 
265  a = ax
266  b = bx
267  fa = f(a)
268  fb = f(b)
269 
270  ! check if a or b is the solution
271  if (fa == dzero) then
272  z = a
273  return
274  else if (fb == dzero) then
275  z = b
276  return
277  end if
278 
279  ! check that f(ax) and f(bx) have opposite sign
280  if (fa * (fb / dabs(fb)) .ge. dzero) &
281  call pstop(1, 'f(ax) and f(bx) do not have different signs,')
282 
283  rs = .true. ! var reset
284  do while (.true.)
285  if (rs) then
286  c = a
287  fc = fa
288  d = b - a
289  e = d
290  end if
291 
292  if (.not. (dabs(fc) .ge. dabs(fb))) then
293  a = b
294  b = c
295  c = a
296  fa = fb
297  fb = fc
298  fc = fa
299  end if
300 
301  tol1 = dtwo * eps * dabs(b) + dhalf * tol
302  xm = dhalf * (c - b)
303  if ((dabs(xm) .le. tol1) .or. (fb .eq. dzero)) then
304  z = b
305  return
306  end if
307 
308  ! see if a bisection is forced
309  if ((dabs(e) .ge. tol1) .and. (dabs(fa) .gt. dabs(fb))) then
310  s = fb / fa
311  if (a .ne. c) then
312  ! inverse quadratic interpolation
313  q = fa / fc
314  r = fb / fc
315  p = s * (dtwo * xm * q * (q - r) - (b - a) * (r - done))
316  q = (q - done) * (r - done) * (s - done)
317  else
318  ! linear interpolation
319  p = dtwo * xm * s
320  q = done - s
321  end if
322 
323  if (p .le. dzero) then
324  p = -p
325  else
326  q = -q
327  end if
328  s = e
329  e = d
330  if (((dtwo * p) .ge. (dthree * xm * q - dabs(tol1 * q))) .or. &
331  (p .ge. dabs(dhalf * s * q))) then
332  d = xm
333  e = d
334  else
335  d = p / q
336  end if
337  else
338  d = xm
339  e = d
340  end if
341 
342  a = b
343  fa = fb
344 
345  if (dabs(d) .le. tol1) then
346  if (xm .le. dzero) then
347  b = b - tol1
348  else
349  b = b + tol1
350  end if
351  else
352  b = b + d
353  end if
354 
355  fb = f(b)
356  rs = (fb * (fc / dabs(fc))) .gt. dzero
357  end do
358  end function zero_br
359 
360  !> @brief Calculate a numerical perturbation given the value of x
361  !!
362  !! Calculate a perturbation value to use for a numerical derivative
363  !! calculation. Taken from the book "Solving Nonlinear Equations with
364  !! Newton's Method" 2003, by C.T. Kelley. Method also used in the
365  !! SWR Process for MODFLOW-2005.
366  !<
367  function get_perturbation(x) result(res)
368  use constantsmodule, only: dprecsqrt, done
369  real(dp), intent(in) :: x !< value that will be perturbed by the result
370  real(dp) :: res !< calculated perturbation value
371  res = dprecsqrt * max(abs(x), done) * sign(done, x)
372  end function get_perturbation
373 
374 end module mathutilmodule
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dsame
real constant for values that are considered the same based on machine precision
Definition: Constants.f90:121
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:44
integer(i4b), parameter lenhugeline
maximum length of a huge line
Definition: Constants.f90:16
real(dp), parameter dprecsqrt
Definition: Constants.f90:120
real(dp), parameter dhalf
real constant 1/2
Definition: Constants.f90:67
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:64
real(dp), parameter dprec
real constant machine precision
Definition: Constants.f90:119
integer(i4b), parameter maxcharlen
maximum length of char string
Definition: Constants.f90:46
@ vsummary
write summary output
Definition: Constants.f90:187
real(dp), parameter dtwo
real constant 2
Definition: Constants.f90:78
real(dp), parameter dthree
real constant 3
Definition: Constants.f90:79
real(dp), parameter done
real constant 1
Definition: Constants.f90:75
subroutine pstop(status, message)
Stop the program, optionally specifying an error status code.
Definition: ErrorUtil.f90:24
This module defines variable data types.
Definition: kind.f90:8
real(dp) function, public zero_ch(x0, x1, f, epsa)
Compute zeros on an interval using Chadrupatla's method.
Definition: MathUtil.f90:133
pure integer(i4b) function mod_offset_int(a, n, d)
Modulo with offset for integer values.
Definition: MathUtil.f90:90
pure real(dp) function mod_offset_dbl(a, n, d)
Modulo with offset for double precision values.
Definition: MathUtil.f90:108
real(dp) function, public get_perturbation(x)
Calculate a numerical perturbation given the value of x.
Definition: MathUtil.f90:368
real(dp) function, public zero_br(ax, bx, f, tol)
Compute a zero of the function f(x) in the interval (x0, x1).
Definition: MathUtil.f90:251
pure logical function, public is_close(a, b, rtol, atol, symmetric)
Check if a real value is approximately equal to another.
Definition: MathUtil.f90:46