14 integer(I4B) :: nodeusr
15 integer(I4B) :: nodered
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()
29 integer(I4B),
pointer,
dimension(:) :: nodeuser => null()
50 subroutine init(this, nlay, ncpl, nodes, top_grid, bot_grid, iavert, &
51 javert, vertex_grid, cellxy_grid, nodereduced, nodeuser)
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
66 integer(I4B) :: nodesuser
71 this%top_grid => top_grid
72 this%bot_grid => bot_grid
75 this%vertex_grid => vertex_grid
76 this%cellxy_grid => cellxy_grid
77 this%nodereduced => nodereduced
78 this%nodeuser => nodeuser
79 nodesuser = ncpl * nlay
81 if (nodes < nodesuser)
then
84 this%reduced = .false.
93 integer(I4B),
intent(in) :: k
94 integer(I4B),
intent(in) :: 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)
102 this%nodered = this%nodeusr
104 call this%cell_setup()
112 integer(I4B),
intent(in) :: nodered
114 this%nodered = nodered
116 if (this%reduced)
then
117 this%nodeusr = this%nodeuser(nodered)
119 this%nodeusr = nodered
122 call get_jk(this%nodeusr, this%ncpl, this%nlay, this%j, this%k)
123 call this%cell_setup()
132 this%top = this%top_grid(this%nodered)
133 this%bot = this%bot_grid(this%nodered)
136 subroutine cprops(this, cell2, hwva, cl1, cl2, ax, ihc)
142 real(DP),
intent(out) :: hwva
143 real(DP),
intent(out) :: cl1
144 real(DP),
intent(out) :: cl2
145 real(DP),
intent(out) :: ax
146 integer(I4B),
intent(out) :: ihc
148 integer(I4B) :: ivert1, ivert2
149 integer(I4B) :: istart1, istart2, istop1, istop2
150 real(DP) :: x0, y0, x1, y1, x2, y2
152 if (this%j == cell2%j)
then
156 hwva = this%get_area()
157 cl1 =
dhalf * (this%top - this%bot)
158 cl2 =
dhalf * (cell2%top - cell2%bot)
164 istart1 = this%iavert(this%j)
165 istop1 = this%iavert(this%j + 1) - 1
166 istart2 = cell2%iavert(cell2%j)
167 istop2 = this%iavert(cell2%j + 1) - 1
169 this%javert(istart2:istop2), &
171 if (ivert1 == 0 .or. ivert2 == 0)
then
178 x1 = this%vertex_grid(1, ivert1)
179 y1 = this%vertex_grid(2, ivert1)
180 x2 = this%vertex_grid(1, ivert2)
181 y2 = this%vertex_grid(2, ivert2)
185 x0 = this%cellxy_grid(1, this%j)
186 y0 = this%cellxy_grid(2, this%j)
190 x0 = this%cellxy_grid(1, cell2%j)
191 y0 = this%cellxy_grid(2, cell2%j)
195 x1 = this%vertex_grid(1, ivert1)
196 y1 = this%vertex_grid(2, ivert1)
197 x2 = this%vertex_grid(1, ivert2)
198 y2 = this%vertex_grid(2, ivert2)
199 ax =
anglex(x1, y1, x2, y2)
212 real(DP),
intent(out) :: xcomp
213 real(DP),
intent(out) :: ycomp
215 integer(I4B) :: ivert1, ivert2
216 integer(I4B) :: istart1, istart2, istop1, istop2
217 real(DP) :: x1, y1, x2, y2
219 istart1 = this%iavert(this%j)
220 istop1 = this%iavert(this%j + 1) - 1
221 istart2 = cell2%iavert(cell2%j)
222 istop2 = this%iavert(cell2%j + 1) - 1
224 this%javert(istart2:istop2), &
226 x1 = this%vertex_grid(1, ivert1)
227 y1 = this%vertex_grid(2, ivert1)
228 x2 = this%vertex_grid(1, ivert2)
229 y2 = this%vertex_grid(2, ivert2)
239 ycomp, zcomp, conlen)
245 logical,
intent(in) :: nozee
246 real(DP),
intent(in) :: satn
247 real(DP),
intent(in) :: satm
248 real(DP),
intent(out) :: xcomp
249 real(DP),
intent(out) :: ycomp
250 real(DP),
intent(out) :: zcomp
251 real(DP),
intent(out) :: conlen
253 real(DP) :: x1, y1, z1, x2, y2, z2
255 x1 = this%cellxy_grid(1, this%j)
256 y1 = this%cellxy_grid(2, this%j)
257 x2 = this%cellxy_grid(1, cell2%j)
258 y2 = this%cellxy_grid(2, cell2%j)
263 z1 = this%bot +
dhalf * satn * (this%top - this%bot)
264 z2 = cell2%bot +
dhalf * satm * (cell2%top - cell2%bot)
280 integer(I4B) :: istart1, istop1, istart2, istop2
281 integer(I4B) :: ivert1, ivert2
283 istart1 = this%iavert(this%j)
284 istop1 = this%iavert(this%j + 1) - 1
285 istart2 = cell2%iavert(cell2%j)
286 istop2 = this%iavert(cell2%j + 1) - 1
288 this%javert(istart2:istop2), &
291 if (ivert1 == 0 .or. ivert2 == 0)
then
304 integer(I4B),
dimension(:) :: ivlist1
305 integer(I4B),
dimension(:) :: ivlist2
306 integer(I4B),
intent(out) :: ivert1
307 integer(I4B),
intent(out) :: ivert2
320 outerloop:
do il1 = 1, nv1 - 1
322 if (ivlist1(il1) == ivlist2(il2) .and. &
323 ivlist1(il1 + 1) == ivlist2(il2 - 1))
then
325 ivert1 = ivlist1(il1)
326 ivert2 = ivlist1(il1 + 1)
347 integer(I4B) :: ivert
348 integer(I4B) :: nvert
349 integer(I4B) :: icount
357 nvert = this%iavert(this%j + 1) - this%iavert(this%j)
359 iv1 = this%javert(this%iavert(this%j))
360 x1 = this%vertex_grid(1, iv1)
361 y1 = this%vertex_grid(2, iv1)
362 do ivert = this%iavert(this%j), this%iavert(this%j + 1) - 1
363 x = this%vertex_grid(1, this%javert(ivert))
364 if (icount < nvert)
then
365 y = this%vertex_grid(2, this%javert(ivert + 1))
367 y = this%vertex_grid(2, this%javert(this%iavert(this%j)))
369 area = area + (x - x1) * (y - y1)
374 do ivert = this%iavert(this%j), this%iavert(this%j + 1) - 1
375 y = this%vertex_grid(2, this%javert(ivert))
376 if (icount < nvert)
then
377 x = this%vertex_grid(1, this%javert(ivert + 1))
379 x = this%vertex_grid(1, this%javert(this%iavert(this%j)))
381 area = area - (x - x1) * (y - y1)
385 area = abs(area) *
dhalf
394 function anglex(x1, y1, x2, y2)
result(ax)
398 real(dp),
intent(in) :: x1
399 real(dp),
intent(in) :: x2
400 real(dp),
intent(in) :: y1
401 real(dp),
intent(in) :: y2
418 real(dp),
intent(in) :: x1
419 real(dp),
intent(in) :: x2
420 real(dp),
intent(in) :: y1
421 real(dp),
intent(in) :: y2
425 d = (x1 - x2)**2 + (y1 - y2)**2
434 real(dp),
intent(in) :: x0
435 real(dp),
intent(in) :: y0
436 real(dp),
intent(in) :: x1
437 real(dp),
intent(in) :: y1
438 real(dp),
intent(in) :: x2
439 real(dp),
intent(in) :: y2
443 d = abs((x2 - x1) * (y1 - y0) - (x1 - x0) * (y2 - y1))
452 real(DP),
intent(in) :: x0
453 real(DP),
intent(in) :: y0
454 real(DP),
intent(in) :: x1
455 real(DP),
intent(in) :: y1
456 real(DP),
intent(out) :: xcomp
457 real(DP),
intent(out) :: ycomp
459 real(DP) :: dx, dy, vmag
463 vmag = sqrt(dx**2 + dy**2)
474 xcomp, ycomp, zcomp, vmag)
476 real(dp),
intent(in) :: x0
477 real(dp),
intent(in) :: y0
478 real(dp),
intent(in) :: z0
479 real(dp),
intent(in) :: x1
480 real(dp),
intent(in) :: y1
481 real(dp),
intent(in) :: z1
482 real(dp),
intent(out) :: xcomp
483 real(dp),
intent(out) :: ycomp
484 real(dp),
intent(out) :: zcomp
487 real(dp) :: dx, dy, dz
492 vmag = sqrt(dx**2 + dy**2 + dz**2)
This module contains simulation constants.
real(dp), parameter dhalf
real constant 1/2
real(dp), parameter dpi
real constant
real(dp), parameter dzero
real constant zero
real(dp), parameter dtwo
real constant 2
real(dp), parameter done
real constant 1
subroutine cprops(this, cell2, hwva, cl1, cl2, ax, ihc)
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.
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,...
subroutine cell_setup(this)
Set top and bottom elevations of grid cell.
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),...
subroutine set_nodered(this, nodered)
Set reduced node number.
logical function shares_edge(this, cell2)
Return true if this shares a horizontal edge with cell2.
subroutine, public shared_edge(ivlist1, ivlist2, ivert1, ivert2)
Find two common vertices shared by cell1 and cell2.
subroutine edge_normal(this, cell2, xcomp, ycomp)
Return the x and y components of an outward normal facing vector.
real(dp) function get_area(this)
Calculate and return the area of the cell.
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...
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,...
real(dp) function distance(x1, y1, x2, y2)
Calculate distance between two points.
subroutine set_kj(this, k, j)
Set node IDs.
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, 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.