17 real(dp),
intent(in) :: x, a, b
18 between = ((x >= a .and. x <= b) .or. (x <= a .and. x >= b))
28 real(dp),
intent(in) :: x
29 real(dp),
intent(in) :: y
30 real(dp),
allocatable,
intent(in) :: poly(:, :)
32 integer(I4B) :: i, ii, num_verts
33 real(dp) :: xa, xb, ya, yb, c = 0.0_dp
36 num_verts =
size(poly, 2)
37 xa = poly(1, num_verts)
38 ya = poly(2, num_verts)
40 do i = 0, num_verts - 1
41 ii = mod(i, num_verts) + 1
45 if ((x == xa .and. y == ya) .or. &
46 (x == xb .and. y == yb))
then
50 else if (ya == yb .and. &
56 else if (
between(y, ya, yb))
then
57 if ((y == ya .and. yb >= ya) .or. &
58 (y == yb .and. ya >= yb))
then
64 c = (xa - x) * (yb - y) - (xb - x) * (ya - y)
69 else if ((ya < yb) .eqv. (c > 0))
then
82 function get_node(ilay, irow, icol, nlay, nrow, ncol)
83 integer(I4B),
intent(in) :: ilay, irow, icol, nlay, nrow, ncol
86 if (nlay > 0 .and. nrow > 0 .and. ncol > 0 .and. &
87 ilay > 0 .and. ilay <= nlay .and. &
88 irow > 0 .and. irow <= nrow .and. &
89 icol > 0 .and. icol <= ncol)
then
91 icol + ncol * (irow - 1) + (ilay - 1) * nrow * ncol
99 subroutine get_ijk(nodenumber, nrow, ncol, nlay, irow, icol, ilay)
101 integer(I4B),
intent(in) :: nodenumber
102 integer(I4B),
intent(in) :: nrow
103 integer(I4B),
intent(in) :: ncol
104 integer(I4B),
intent(in) :: nlay
105 integer(I4B),
intent(out) :: irow
106 integer(I4B),
intent(out) :: icol
107 integer(I4B),
intent(out) :: ilay
109 integer(I4B) :: nodes
112 nodes = nlay * nrow * ncol
113 if (nodenumber < 1 .or. nodenumber > nodes)
then
118 ilay = (nodenumber - 1) / (ncol * nrow) + 1
119 ij = nodenumber - (ilay - 1) * ncol * nrow
120 irow = (ij - 1) / ncol + 1
121 icol = ij - (irow - 1) * ncol
127 subroutine get_jk(nodenumber, ncpl, nlay, icpl, ilay)
129 integer(I4B),
intent(in) :: nodenumber
130 integer(I4B),
intent(in) :: ncpl
131 integer(I4B),
intent(in) :: nlay
132 integer(I4B),
intent(out) :: icpl
133 integer(I4B),
intent(out) :: ilay
135 integer(I4B) :: nodes
138 if (nodenumber < 1 .or. nodenumber > nodes)
then
142 ilay = (nodenumber - 1) / ncpl + 1
143 icpl = nodenumber - (ilay - 1) * ncpl
148 pure function skew(v, s, invert)
result(res)
150 real(dp),
intent(in) :: v(2)
151 real(dp),
intent(in) :: s(3)
152 logical(LGP),
intent(in),
optional :: invert
155 logical(LGP) :: linvert
156 real(dp) :: sxx, sxy, syy
159 if (
present(invert))
then
168 if (.not. linvert)
then
169 res(1) = sxx * v(1) + sxy * v(2)
173 res(1) = (v(1) - sxy * res(2)) / sxx
180 xorigin, yorigin, zorigin, &
184 real(dp) :: xin, yin, zin
185 real(dp) :: xout, yout, zout
186 real(dp),
optional :: xorigin, yorigin, zorigin
187 real(dp),
optional :: sinrot, cosrot
188 logical(LGP),
optional :: invert
190 logical(LGP) :: ltranslate, lrotate, linvert
192 real(dp) :: lxorigin, lyorigin, lzorigin
193 real(dp) :: lsinrot, lcosrot
196 call defaults(lxorigin, lyorigin, lzorigin, &
197 lsinrot, lcosrot, linvert, &
198 ltranslate, lrotate, &
199 xorigin, yorigin, zorigin, &
200 sinrot, cosrot, invert)
203 if (.not. linvert)
then
206 xout = xin - lxorigin
207 yout = yin - lyorigin
208 zout = zin - lzorigin
217 xout = x * lcosrot + y * lsinrot
218 yout = -x * lsinrot + y * lcosrot
223 x = xin * lcosrot - yin * lsinrot
224 y = xin * lsinrot + yin * lcosrot
232 zout = zin + lzorigin
238 subroutine compose(xorigin, yorigin, zorigin, &
240 xorigin_new, yorigin_new, zorigin_new, &
241 sinrot_new, cosrot_new, &
244 real(dp) :: xorigin, yorigin, zorigin
245 real(dp) :: sinrot, cosrot
246 real(dp),
optional :: xorigin_new, yorigin_new, zorigin_new
247 real(dp),
optional :: sinrot_new, cosrot_new
248 logical(LGP),
optional :: invert
250 logical(LGP) :: ltranslate, lrotate, linvert
251 real(dp) :: xorigin_add, yorigin_add, zorigin_add
252 real(dp) :: sinrot_add, cosrot_add
253 real(dp) :: x0, y0, z0, s0, c0
256 call defaults(xorigin_add, yorigin_add, zorigin_add, &
257 sinrot_add, cosrot_add, linvert, &
258 ltranslate, lrotate, &
259 xorigin_new, yorigin_new, zorigin_new, &
260 sinrot_new, cosrot_new, invert)
270 if (.not. linvert)
then
277 call transform(xorigin_add, yorigin_add, zorigin_add, &
278 xorigin, yorigin, zorigin, &
279 x0, y0, z0, s0, c0, .true.)
285 sinrot = cosrot_add * s0 + sinrot_add * c0
286 cosrot = cosrot_add * c0 - sinrot_add * s0
296 call transform(-xorigin_add, -yorigin_add, zorigin_add, &
297 x0, y0, z0, xorigin, yorigin, zorigin, &
298 -sinrot_add, cosrot_add, .true.)
300 xorigin = c0 * x0 - s0 * y0
301 yorigin = s0 * x0 + c0 * y0
308 sinrot = cosrot_add * s0 - sinrot_add * c0
309 cosrot = cosrot_add * c0 + sinrot_add * s0
317 invert, translate, rotate, &
318 xorigin_opt, yorigin_opt, zorigin_opt, &
319 sinrot_opt, cosrot_opt, invert_opt)
321 real(DP) :: xorigin, yorigin, zorigin
322 real(DP) :: sinrot, cosrot
323 logical(LGP) :: invert, translate, rotate
324 real(DP),
optional :: xorigin_opt, yorigin_opt, zorigin_opt
325 real(DP),
optional :: sinrot_opt, cosrot_opt
326 logical(LGP),
optional :: invert_opt
330 if (
present(xorigin_opt))
then
331 xorigin = xorigin_opt
335 if (
present(yorigin_opt))
then
336 yorigin = yorigin_opt
340 if (
present(zorigin_opt))
then
341 zorigin = zorigin_opt
347 if (
present(sinrot_opt))
then
349 if (
present(cosrot_opt))
then
354 cosrot = dsqrt(
done - sinrot * sinrot)
357 else if (
present(cosrot_opt))
then
361 sinrot = dsqrt(
done - cosrot * cosrot)
365 if (
present(invert_opt)) invert = invert_opt
369 function area(xv, yv, cw)
result(a)
371 real(dp),
dimension(:),
intent(in) :: xv
372 real(dp),
dimension(:),
intent(in) :: yv
373 logical(LGP),
intent(in),
optional :: cw
378 if (
present(cw))
then
388 a = -
dhalf * sum(xv(:) * cshift(yv(:), s) - cshift(xv(:), s) * yv(:))
403 integer(I4B),
dimension(:) :: iverts1
404 integer(I4B),
dimension(:) :: iverts2
405 integer(I4B),
intent(out) :: iface
408 integer(I4B) :: il1, iil1
409 integer(I4B) :: il2, iil2
410 logical(LGP) :: found
411 logical(LGP) :: wrapped
417 wrapped = iverts1(1) == iverts1(nv1)
424 outerloop:
do il1 = 1, nv1 - 1
426 if (iverts1(il1) == iverts2(il2))
then
431 if (wrapped) iil2 = iil2 - 1
435 if (iverts1(iil1) == iverts2(iil2))
then
444 if (wrapped) iil1 = iil1 - 1
448 if (iverts1(iil1) == iverts2(iil2))
then
467 real(dp),
intent(inout) :: alpha
468 real(dp),
intent(inout) :: beta
469 real(dp),
intent(out) :: gamma
470 real(dp),
intent(in),
optional :: pad
477 if (
present(pad))
then
480 call pstop(1,
"pad must be between 0 and 1/3, inclusive")
485 gamma =
done - alpha - beta
489 if (alpha < lolimit)
then
493 gamma =
done - alpha - beta
497 if (beta < lolimit)
then
505 else if (beta > hilimit)
then
517 if (beta < lolimit)
then
521 gamma =
done - alpha - beta
525 if (alpha < lolimit)
then
533 else if (alpha > hilimit)
then
545 if (gamma < lolimit)
then
548 delta =
dhalf * (lolimit - gamma)
550 alpha = alpha - delta
This module contains simulation constants.
real(dp), parameter dsame
real constant for values that are considered the same based on machine precision
real(dp), parameter dep3
real constant 1000
real(dp), parameter donethird
real constant 1/3
real(dp), parameter dhalf
real constant 1/2
real(dp), parameter dzero
real constant zero
real(dp), parameter dtwo
real constant 2
real(dp), parameter done
real constant 1
subroutine pstop(status, message)
Stop the program, optionally specifying an error status code.
subroutine, public shared_face(iverts1, iverts2, iface)
Find the lateral face shared by two cells.
subroutine, public transform(xin, yin, zin, xout, yout, zout, xorigin, yorigin, zorigin, sinrot, cosrot, invert)
Apply a 3D translation and optional 2D rotation to coordinates.
integer(i4b) function, public get_node(ilay, irow, icol, nlay, nrow, ncol)
Get node number, given layer, row, and column indices for a structured grid. If any argument is inval...
subroutine defaults(xorigin, yorigin, zorigin, sinrot, cosrot, invert, translate, rotate, xorigin_opt, yorigin_opt, zorigin_opt, sinrot_opt, cosrot_opt, invert_opt)
Process arguments and set defaults. Internal use only.
logical function, public between(x, a, b)
Check if a value is between two other values (inclusive).
subroutine, public clamp_bary(alpha, beta, gamma, pad)
Clamp barycentric coordinates to the interior of a triangle, with optional padding some minimum dista...
real(dp) function, public area(xv, yv, cw)
Calculate polygon area, with vertices given in CW or CCW order.
subroutine, public compose(xorigin, yorigin, zorigin, sinrot, cosrot, xorigin_new, yorigin_new, zorigin_new, sinrot_new, cosrot_new, invert)
Apply a 3D translation and 2D rotation to an existing transformation.
logical function, public point_in_polygon(x, y, poly)
Check if a point is within a polygon.
subroutine, public get_ijk(nodenumber, nrow, ncol, nlay, irow, icol, ilay)
Get row, column and layer indices from node number and grid dimensions. If nodenumber is invalid,...
pure real(dp) function, dimension(2), public skew(v, s, invert)
Skew a 2D vector along the x-axis.
subroutine, public get_jk(nodenumber, ncpl, nlay, icpl, ilay)
Get layer index and within-layer node index from node number and grid dimensions. If nodenumber is in...
This module defines variable data types.