MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
DisvGeom.f90
Go to the documentation of this file.
1 module disvgeom
2 
3  use kindmodule, only: dp, i4b
4  use geomutilmodule, only: get_node, get_jk
5  implicit none
6  private
7  public :: disvgeomtype
8  public :: line_unit_vector
9 
11 
12  integer(I4B) :: k
13  integer(I4B) :: j
14  integer(I4B) :: nodeusr
15  integer(I4B) :: nodered
16  integer(I4B) :: nlay
17  integer(I4B) :: ncpl
18  logical :: reduced
19  integer(I4B) :: nodes ! number of reduced nodes; nodes = nlay *ncpl when grid is NOT reduced
20  real(dp) :: top
21  real(dp) :: bot
22  real(dp), pointer, dimension(:) :: top_grid => null()
23  real(dp), pointer, dimension(:) :: bot_grid => null()
24  integer(I4B), pointer, dimension(:) :: iavert => null()
25  integer(I4B), pointer, dimension(:) :: javert => null()
26  real(dp), pointer, dimension(:, :) :: vertex_grid => null()
27  real(dp), pointer, dimension(:, :) :: cellxy_grid => null()
28  integer(I4B), pointer, dimension(:, :) :: nodereduced => null() ! nodered = nodereduced(nodeusr)
29  integer(I4B), pointer, dimension(:) :: nodeuser => null() ! nodeusr = nodesuser(nodered)
30 
31  contains
32 
33  procedure :: init
34  generic :: set => set_kj, set_nodered
35  procedure :: set_kj
36  procedure :: set_nodered
37  procedure :: cell_setup
38  procedure :: cprops
39  procedure :: edge_normal
40  procedure :: connection_vector
41  procedure :: shares_edge
42  procedure :: get_area
43 
44  end type disvgeomtype
45 
46 contains
47 
48  !> @brief Initialize
49  !<
50  subroutine init(this, nlay, ncpl, nodes, top_grid, bot_grid, iavert, &
51  javert, vertex_grid, cellxy_grid, nodereduced, nodeuser)
52  ! -- dummy
53  class(disvgeomtype) :: this
54  integer(I4B), intent(in) :: nlay
55  integer(I4B), intent(in) :: ncpl
56  integer(I4B), intent(in) :: nodes
57  real(DP), dimension(nodes), target :: top_grid
58  real(DP), dimension(nodes), target :: bot_grid
59  integer(I4B), dimension(:), target :: iavert
60  integer(I4B), dimension(:), target :: javert
61  real(DP), dimension(:, :), target :: vertex_grid
62  real(DP), dimension(:, :), target :: cellxy_grid
63  integer(I4B), dimension(ncpl, nlay), target :: nodereduced
64  integer(I4B), dimension(nodes), target :: nodeuser
65  ! -- local
66  integer(I4B) :: nodesuser
67  !
68  this%nlay = nlay
69  this%ncpl = ncpl
70  this%nodes = nodes
71  this%top_grid => top_grid
72  this%bot_grid => bot_grid
73  this%iavert => iavert
74  this%javert => javert
75  this%vertex_grid => vertex_grid
76  this%cellxy_grid => cellxy_grid
77  this%nodereduced => nodereduced
78  this%nodeuser => nodeuser
79  nodesuser = ncpl * nlay
80  !
81  if (nodes < nodesuser) then
82  this%reduced = .true.
83  else
84  this%reduced = .false.
85  end if
86  end subroutine init
87 
88  !> @brief Set node IDs
89  !<
90  subroutine set_kj(this, k, j)
91  ! -- dummy
92  class(disvgeomtype) :: this
93  integer(I4B), intent(in) :: k
94  integer(I4B), intent(in) :: j
95  !
96  this%k = k
97  this%j = j
98  this%nodeusr = get_node(k, 1, j, this%nlay, 1, this%ncpl)
99  if (this%reduced) then
100  this%nodered = this%nodereduced(k, j)
101  else
102  this%nodered = this%nodeusr
103  end if
104  call this%cell_setup()
105  !
106  ! -- Return
107  return
108  end subroutine set_kj
109 
110  !> @brief Set reduced node number
111  !<
112  subroutine set_nodered(this, nodered)
113  ! -- dummy
114  class(disvgeomtype) :: this
115  integer(I4B), intent(in) :: nodered
116  !
117  this%nodered = nodered
118  !
119  if (this%reduced) then
120  this%nodeusr = this%nodeuser(nodered)
121  else
122  this%nodeusr = nodered
123  end if
124  !
125  call get_jk(this%nodeusr, this%ncpl, this%nlay, this%j, this%k)
126  call this%cell_setup()
127  !
128  ! -- Return
129  return
130  end subroutine set_nodered
131 
132  !> @brief Set top and bottom elevations of grid cell
133  !<
134  subroutine cell_setup(this)
135  ! -- dummy
136  class(disvgeomtype) :: this
137  !
138  this%top = this%top_grid(this%nodered)
139  this%bot = this%bot_grid(this%nodered)
140  end subroutine cell_setup
141 
142  subroutine cprops(this, cell2, hwva, cl1, cl2, ax, ihc)
143  ! -- module
144  use constantsmodule, only: dzero, dhalf, done
145  ! -- dummy
146  class(disvgeomtype) :: this
147  type(disvgeomtype) :: cell2
148  real(DP), intent(out) :: hwva
149  real(DP), intent(out) :: cl1
150  real(DP), intent(out) :: cl2
151  real(DP), intent(out) :: ax
152  integer(I4B), intent(out) :: ihc
153  ! -- local
154  integer(I4B) :: ivert1, ivert2
155  integer(I4B) :: istart1, istart2, istop1, istop2
156  real(DP) :: x0, y0, x1, y1, x2, y2
157  !
158  if (this%j == cell2%j) then
159  !
160  ! -- Cells share same j index, so must be a vertical connection
161  ihc = 0
162  hwva = this%get_area()
163  cl1 = dhalf * (this%top - this%bot)
164  cl2 = dhalf * (cell2%top - cell2%bot)
165  ax = dzero
166  else
167  !
168  ! -- Must be horizontal connection
169  ihc = 1
170  istart1 = this%iavert(this%j)
171  istop1 = this%iavert(this%j + 1) - 1
172  istart2 = cell2%iavert(cell2%j)
173  istop2 = this%iavert(cell2%j + 1) - 1
174  call shared_edge(this%javert(istart1:istop1), &
175  this%javert(istart2:istop2), &
176  ivert1, ivert2)
177  if (ivert1 == 0 .or. ivert2 == 0) then
178  !
179  ! -- Cells do not share an edge
180  hwva = dzero
181  cl1 = done
182  cl2 = done
183  else
184  x1 = this%vertex_grid(1, ivert1)
185  y1 = this%vertex_grid(2, ivert1)
186  x2 = this%vertex_grid(1, ivert2)
187  y2 = this%vertex_grid(2, ivert2)
188  hwva = distance(x1, y1, x2, y2)
189  !
190  ! -- cl1
191  x0 = this%cellxy_grid(1, this%j)
192  y0 = this%cellxy_grid(2, this%j)
193  cl1 = distance_normal(x0, y0, x1, y1, x2, y2)
194  !
195  ! -- cl2
196  x0 = this%cellxy_grid(1, cell2%j)
197  y0 = this%cellxy_grid(2, cell2%j)
198  cl2 = distance_normal(x0, y0, x1, y1, x2, y2)
199  !
200  ! -- anglex
201  x1 = this%vertex_grid(1, ivert1)
202  y1 = this%vertex_grid(2, ivert1)
203  x2 = this%vertex_grid(1, ivert2)
204  y2 = this%vertex_grid(2, ivert2)
205  ax = anglex(x1, y1, x2, y2)
206  end if
207  end if
208  !
209  ! -- Return
210  return
211  end subroutine cprops
212 
213  !> @brief Return the x and y components of an outward normal facing vector
214  !<
215  subroutine edge_normal(this, cell2, xcomp, ycomp)
216  ! -- module
217  use constantsmodule, only: dzero, dhalf, done
218  ! -- dummy
219  class(disvgeomtype) :: this
220  type(disvgeomtype) :: cell2
221  real(DP), intent(out) :: xcomp
222  real(DP), intent(out) :: ycomp
223  ! -- local
224  integer(I4B) :: ivert1, ivert2
225  integer(I4B) :: istart1, istart2, istop1, istop2
226  real(DP) :: x1, y1, x2, y2
227  !
228  istart1 = this%iavert(this%j)
229  istop1 = this%iavert(this%j + 1) - 1
230  istart2 = cell2%iavert(cell2%j)
231  istop2 = this%iavert(cell2%j + 1) - 1
232  call shared_edge(this%javert(istart1:istop1), &
233  this%javert(istart2:istop2), &
234  ivert1, ivert2)
235  x1 = this%vertex_grid(1, ivert1)
236  y1 = this%vertex_grid(2, ivert1)
237  x2 = this%vertex_grid(1, ivert2)
238  y2 = this%vertex_grid(2, ivert2)
239  !
240  call line_unit_normal(x1, y1, x2, y2, xcomp, ycomp)
241  !
242  ! -- Return
243  return
244  end subroutine edge_normal
245 
246  !> @brief Return the x y and z components of a unit vector that points from
247  !! from the center of this to the center of cell2, and the straight-line
248  !! connection length
249  !<
250  subroutine connection_vector(this, cell2, nozee, satn, satm, xcomp, &
251  ycomp, zcomp, conlen)
252  ! -- module
253  use constantsmodule, only: dzero, dhalf, done
254  ! -- dummy
255  class(disvgeomtype) :: this
256  type(disvgeomtype) :: cell2
257  logical, intent(in) :: nozee
258  real(DP), intent(in) :: satn
259  real(DP), intent(in) :: satm
260  real(DP), intent(out) :: xcomp
261  real(DP), intent(out) :: ycomp
262  real(DP), intent(out) :: zcomp
263  real(DP), intent(out) :: conlen
264  ! -- local
265  real(DP) :: x1, y1, z1, x2, y2, z2
266  !
267  x1 = this%cellxy_grid(1, this%j)
268  y1 = this%cellxy_grid(2, this%j)
269  x2 = this%cellxy_grid(1, cell2%j)
270  y2 = this%cellxy_grid(2, cell2%j)
271  if (nozee) then
272  z1 = dzero
273  z2 = dzero
274  else
275  z1 = this%bot + dhalf * satn * (this%top - this%bot)
276  z2 = cell2%bot + dhalf * satm * (cell2%top - cell2%bot)
277  end if
278  !
279  call line_unit_vector(x1, y1, z1, x2, y2, z2, xcomp, ycomp, zcomp, &
280  conlen)
281  !
282  ! -- Return
283  return
284  end subroutine connection_vector
285 
286  !> @brief Return true if this shares a horizontal edge with cell2
287  !<
288  function shares_edge(this, cell2) result(l)
289  ! -- dummy
290  class(disvgeomtype) :: this
291  type(disvgeomtype) :: cell2
292  ! -- return
293  logical l
294  ! -- local
295  integer(I4B) :: istart1, istop1, istart2, istop2
296  integer(I4B) :: ivert1, ivert2
297  !
298  istart1 = this%iavert(this%j)
299  istop1 = this%iavert(this%j + 1) - 1
300  istart2 = cell2%iavert(cell2%j)
301  istop2 = this%iavert(cell2%j + 1) - 1
302  call shared_edge(this%javert(istart1:istop1), &
303  this%javert(istart2:istop2), &
304  ivert1, ivert2)
305  l = .true.
306  if (ivert1 == 0 .or. ivert2 == 0) then
307  l = .false.
308  end if
309  !
310  ! -- Return
311  return
312  end function shares_edge
313 
314  !> @brief Find two common vertices shared by cell1 and cell2.
315  !!
316  !! Return 0 if there are no shared edges. Proceed forward through ivlist1
317  !! and backward through ivlist2 as a clockwise face in cell1 must correspond
318  !! to a counter clockwise face in cell2.
319  !<
320  subroutine shared_edge(ivlist1, ivlist2, ivert1, ivert2)
321  ! -- dummy
322  integer(I4B), dimension(:) :: ivlist1
323  integer(I4B), dimension(:) :: ivlist2
324  integer(I4B), intent(out) :: ivert1
325  integer(I4B), intent(out) :: ivert2
326  ! -- local
327  integer(I4B) :: nv1
328  integer(I4B) :: nv2
329  integer(I4B) :: il1
330  integer(I4B) :: il2
331  logical :: found
332  !
333  found = .false.
334  nv1 = size(ivlist1)
335  nv2 = size(ivlist2)
336  ivert1 = 0
337  ivert2 = 0
338  outerloop: do il1 = 1, nv1 - 1
339  do il2 = nv2, 2, -1
340  if (ivlist1(il1) == ivlist2(il2) .and. &
341  ivlist1(il1 + 1) == ivlist2(il2 - 1)) then
342  found = .true.
343  ivert1 = ivlist1(il1)
344  ivert2 = ivlist1(il1 + 1)
345  exit outerloop
346  end if
347  end do
348  if (found) exit
349  end do outerloop
350  end subroutine shared_edge
351 
352  !> @brief Calculate and return the area of the cell
353  !!
354  !! a = 1/2 *[(x1*y2 + x2*y3 + x3*y4 + ... + xn*y1) -
355  !! (x2*y1 + x3*y2 + x4*y3 + ... + x1*yn)]
356  !<
357  function get_area(this) result(area)
358  ! -- module
359  use constantsmodule, only: dzero, dhalf
360  ! -- dummy
361  class(disvgeomtype) :: this
362  ! -- return
363  real(dp) :: area
364  ! -- local
365  integer(I4B) :: ivert
366  integer(I4B) :: nvert
367  integer(I4B) :: icount
368  integer(I4B) :: iv1
369  real(dp) :: x
370  real(dp) :: y
371  real(dp) :: x1
372  real(dp) :: y1
373  !
374  area = dzero
375  nvert = this%iavert(this%j + 1) - this%iavert(this%j)
376  icount = 1
377  iv1 = this%javert(this%iavert(this%j))
378  x1 = this%vertex_grid(1, iv1)
379  y1 = this%vertex_grid(2, iv1)
380  do ivert = this%iavert(this%j), this%iavert(this%j + 1) - 1
381  x = this%vertex_grid(1, this%javert(ivert))
382  if (icount < nvert) then
383  y = this%vertex_grid(2, this%javert(ivert + 1))
384  else
385  y = this%vertex_grid(2, this%javert(this%iavert(this%j)))
386  end if
387  area = area + (x - x1) * (y - y1)
388  icount = icount + 1
389  end do
390  !
391  icount = 1
392  do ivert = this%iavert(this%j), this%iavert(this%j + 1) - 1
393  y = this%vertex_grid(2, this%javert(ivert))
394  if (icount < nvert) then
395  x = this%vertex_grid(1, this%javert(ivert + 1))
396  else
397  x = this%vertex_grid(1, this%javert(this%iavert(this%j)))
398  end if
399  area = area - (x - x1) * (y - y1)
400  icount = icount + 1
401  end do
402  !
403  area = abs(area) * dhalf
404  !
405  ! -- Return
406  return
407  end function get_area
408 
409  !> @brief Calculate the angle that the x-axis makes with a line that is
410  !! normal to the two points.
411  !!
412  !! This assumes that vertices are numbered clockwise so that the angle is for
413  !! the normal outward of cell n.
414  !<
415  function anglex(x1, y1, x2, y2) result(ax)
416  ! -- modules
417  use constantsmodule, only: dzero, dtwo, dpi
418  ! -- dummy
419  real(dp), intent(in) :: x1
420  real(dp), intent(in) :: x2
421  real(dp), intent(in) :: y1
422  real(dp), intent(in) :: y2
423  ! -- return
424  real(dp) :: ax
425  ! -- local
426  real(dp) :: dx
427  real(dp) :: dy
428  !
429  dx = x2 - x1
430  dy = y2 - y1
431  ax = atan2(dx, -dy)
432  if (ax < dzero) ax = dtwo * dpi + ax
433  !
434  ! -- Return
435  return
436  end function anglex
437 
438  !> @brief Calculate distance between two points
439  !<
440  function distance(x1, y1, x2, y2) result(d)
441  ! -- dummy
442  real(dp), intent(in) :: x1
443  real(dp), intent(in) :: x2
444  real(dp), intent(in) :: y1
445  real(dp), intent(in) :: y2
446  ! -- return
447  real(dp) :: d
448  !
449  d = (x1 - x2)**2 + (y1 - y2)**2
450  d = sqrt(d)
451  !
452  ! -- Return
453  return
454  end function distance
455 
456  !> @brief Calculate normal distance from point (x0, y0) to line defined by
457  !! two points, (x1, y1), (x2, y2).
458  !<
459  function distance_normal(x0, y0, x1, y1, x2, y2) result(d)
460  ! -- dummy
461  real(dp), intent(in) :: x0
462  real(dp), intent(in) :: y0
463  real(dp), intent(in) :: x1
464  real(dp), intent(in) :: y1
465  real(dp), intent(in) :: x2
466  real(dp), intent(in) :: y2
467  ! -- return
468  real(dp) :: d
469  !
470  d = abs((x2 - x1) * (y1 - y0) - (x1 - x0) * (y2 - y1))
471  d = d / distance(x1, y1, x2, y2)
472  !
473  ! -- Return
474  return
475  end function distance_normal
476 
477  !> @brief Calculate the normal vector components (xcomp and ycomp) for a line
478  !! defined by two points, (x0, y0), (x1, y1).
479  !<
480  subroutine line_unit_normal(x0, y0, x1, y1, xcomp, ycomp)
481  ! -- dummy
482  real(DP), intent(in) :: x0
483  real(DP), intent(in) :: y0
484  real(DP), intent(in) :: x1
485  real(DP), intent(in) :: y1
486  real(DP), intent(out) :: xcomp
487  real(DP), intent(out) :: ycomp
488  ! -- local
489  real(DP) :: dx, dy, vmag
490  !
491  dx = x1 - x0
492  dy = y1 - y0
493  vmag = sqrt(dx**2 + dy**2)
494  xcomp = -dy / vmag
495  ycomp = dx / vmag
496  !
497  ! -- Return
498  return
499  end subroutine line_unit_normal
500 
501  !> @brief Calculate the vector components (xcomp, ycomp, and zcomp) for a
502  !! line defined by two points, (x0, y0, z0), (x1, y1, z1).
503  !!
504  !! Also return the magnitude of the original vector, vmag.
505  !<
506  subroutine line_unit_vector(x0, y0, z0, x1, y1, z1, &
507  xcomp, ycomp, zcomp, vmag)
508  ! -- dummy
509  real(dp), intent(in) :: x0
510  real(dp), intent(in) :: y0
511  real(dp), intent(in) :: z0
512  real(dp), intent(in) :: x1
513  real(dp), intent(in) :: y1
514  real(dp), intent(in) :: z1
515  real(dp), intent(out) :: xcomp
516  real(dp), intent(out) :: ycomp
517  real(dp), intent(out) :: zcomp
518  real(dp) :: vmag
519  ! -- local
520  real(dp) :: dx, dy, dz
521  !
522  dx = x1 - x0
523  dy = y1 - y0
524  dz = z1 - z0
525  vmag = sqrt(dx**2 + dy**2 + dz**2)
526  xcomp = dx / vmag
527  ycomp = dy / vmag
528  zcomp = dz / vmag
529  !
530  ! -- Return
531  return
532  end subroutine line_unit_vector
533 
534 end module disvgeom
subroutine init()
Definition: GridSorting.f90:24
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dhalf
real constant 1/2
Definition: Constants.f90:67
real(dp), parameter dpi
real constant
Definition: Constants.f90:127
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:64
real(dp), parameter dtwo
real constant 2
Definition: Constants.f90:78
real(dp), parameter done
real constant 1
Definition: Constants.f90:75
subroutine cprops(this, cell2, hwva, cl1, cl2, ax, ihc)
Definition: DisvGeom.f90:143
real(dp) function anglex(x1, y1, x2, y2)
Calculate the angle that the x-axis makes with a line that is normal to the two points.
Definition: DisvGeom.f90:416
subroutine line_unit_normal(x0, y0, x1, y1, xcomp, ycomp)
Calculate the normal vector components (xcomp and ycomp) for a line defined by two points,...
Definition: DisvGeom.f90:481
subroutine cell_setup(this)
Set top and bottom elevations of grid cell.
Definition: DisvGeom.f90:135
subroutine shared_edge(ivlist1, ivlist2, ivert1, ivert2)
Find two common vertices shared by cell1 and cell2.
Definition: DisvGeom.f90:321
real(dp) function distance_normal(x0, y0, x1, y1, x2, y2)
Calculate normal distance from point (x0, y0) to line defined by two points, (x1, y1),...
Definition: DisvGeom.f90:460
subroutine set_nodered(this, nodered)
Set reduced node number.
Definition: DisvGeom.f90:113
logical function shares_edge(this, cell2)
Return true if this shares a horizontal edge with cell2.
Definition: DisvGeom.f90:289
subroutine edge_normal(this, cell2, xcomp, ycomp)
Return the x and y components of an outward normal facing vector.
Definition: DisvGeom.f90:216
real(dp) function get_area(this)
Calculate and return the area of the cell.
Definition: DisvGeom.f90:358
subroutine connection_vector(this, cell2, nozee, satn, satm, xcomp, ycomp, zcomp, conlen)
Return the x y and z components of a unit vector that points from from the center of this to the cent...
Definition: DisvGeom.f90:252
subroutine, public line_unit_vector(x0, y0, z0, x1, y1, z1, xcomp, ycomp, zcomp, vmag)
Calculate the vector components (xcomp, ycomp, and zcomp) for a line defined by two points,...
Definition: DisvGeom.f90:508
real(dp) function distance(x1, y1, x2, y2)
Calculate distance between two points.
Definition: DisvGeom.f90:441
subroutine set_kj(this, k, j)
Set node IDs.
Definition: DisvGeom.f90:91
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...
Definition: GeomUtil.f90:81
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...
Definition: GeomUtil.f90:126
This module defines variable data types.
Definition: kind.f90:8