55 allocate (method%name)
59 method%delegates = .true.
66 deallocate (this%name)
70 subroutine load_disv(this, particle, next_level, submethod)
78 integer(I4B),
intent(in) :: next_level
79 class(
methodtype),
pointer,
intent(inout) :: submethod
86 select type (cell => this%cell)
88 ic = particle%idomain(next_level)
89 call this%load_cell_defn(ic, cell%defn)
90 if (this%fmi%ibdgwfsat0(ic) == 0)
then
94 trackctl=this%trackctl, &
95 tracktimes=this%tracktimes)
97 else if (particle%ifrctrn > 0)
then
101 trackctl=this%trackctl, &
102 tracktimes=this%tracktimes)
104 else if (cell%defn%can_be_rect)
then
110 trackctl=this%trackctl, &
111 tracktimes=this%tracktimes)
113 else if (cell%defn%can_be_quad)
then
119 trackctl=this%trackctl, &
120 tracktimes=this%tracktimes)
126 trackctl=this%trackctl, &
127 tracktimes=this%tracktimes)
141 integer(I4B) :: inface
146 integer(I4B) :: idiag
151 select type (dis => this%fmi%dis)
153 inface = particle%iboundary(2)
154 idiag = dis%con%ia(cell%defn%icell)
155 inbr = cell%defn%facenbr(inface)
157 ic = dis%con%ja(ipos)
158 icu = dis%get_nodeuser(ic)
159 call get_jk(icu, dis%ncpl, dis%nlay, icpl, ilay)
166 if (ic == particle%icp .and. inface == 7 .and. ilay < particle%ilay)
then
167 particle%advancing = .false.
168 particle%idomain(2) = particle%icp
170 particle%izone = particle%izp
171 call this%save(particle, reason=3)
174 particle%icp = particle%idomain(2)
175 particle%izp = particle%izone
178 particle%idomain(2) = ic
183 call this%map_neighbor(cell%defn, inface, z)
185 particle%iboundary(2) = inface
186 particle%idomain(3:) = 0
187 particle%iboundary(3:) = 0
200 integer(I4B) :: idiag
203 idiag = this%fmi%dis%con%ia(cell%defn%icell)
204 inbr = cell%defn%facenbr(particle%iboundary(2))
208 this%flowja(ipos) = this%flowja(ipos) -
done
211 this%flowja(this%fmi%dis%con%isym(ipos)) &
212 = this%flowja(this%fmi%dis%con%isym(ipos)) +
done
224 select type (c => this%cell)
230 if (cell%defn%facenbr(particle%iboundary(2)) .eq. 0)
then
232 particle%advancing = .false.
233 call this%save(particle, reason=3)
236 call this%load_particle(cell, particle)
237 if (.not. particle%advancing)
return
240 call this%update_flowja(cell, particle)
250 integer(I4B),
intent(inout) :: inface
251 double precision,
intent(inout) :: z
254 integer(I4B) :: npolyvertsin
256 integer(I4B) :: npolyverts
258 integer(I4B) :: inbrnbr
269 inbr = defn%facenbr(inface)
270 if (inbr .eq. 0)
then
277 j = this%fmi%dis%con%ia(icin)
278 ic = this%fmi%dis%con%ja(j + inbr)
279 call this%load_cell_defn(ic, this%neighbor)
281 npolyvertsin = defn%npolyverts
282 npolyverts = this%neighbor%npolyverts
283 if (inface .eq. npolyvertsin + 2)
then
285 inface = npolyverts + 3
286 else if (inface .eq. npolyvertsin + 3)
then
288 inface = npolyverts + 2
291 j = this%fmi%dis%con%ia(ic)
292 do m = 1, npolyverts + 3
293 inbrnbr = this%neighbor%facenbr(m)
294 if (this%fmi%dis%con%ja(j + inbrnbr) .eq. icin)
then
302 zrel = (z - botfrom) / (topfrom - botfrom)
303 top = this%fmi%dis%top(ic)
304 bot = this%fmi%dis%bot(ic)
305 sat = this%fmi%gwfsat(ic)
306 z = bot + zrel * sat * (top - bot)
315 real(DP),
intent(in) :: tmax
317 call this%track(particle, 1, tmax)
324 integer(I4B),
intent(in) :: ic
327 call this%load_properties(ic, defn)
328 call this%load_polygon(defn)
329 call this%load_neighbors(defn)
330 call this%load_indicators(defn)
331 call this%load_flows(defn)
338 integer(I4B),
intent(in) :: ic
344 integer(I4B) :: icu, icpl, ilay
347 defn%iatop =
get_iatop(this%fmi%dis%get_ncpl(), &
348 this%fmi%dis%get_nodeuser(ic))
349 top = this%fmi%dis%top(ic)
350 bot = this%fmi%dis%bot(ic)
351 sat = this%fmi%gwfsat(ic)
352 top = bot + sat * (top - bot)
356 defn%porosity = this%porosity(ic)
357 defn%retfactor = this%retfactor(ic)
358 defn%izone = this%izone(ic)
359 select type (dis => this%fmi%dis)
361 icu = dis%get_nodeuser(ic)
362 call get_jk(icu, dis%ncpl, dis%nlay, icpl, ilay)
372 call this%fmi%dis%get_polyverts( &
376 defn%npolyverts =
size(defn%polyvert, dim=2) - 1
396 integer(I4B) :: istart1
397 integer(I4B) :: istart2
398 integer(I4B) :: istop1
399 integer(I4B) :: istop2
400 integer(I4B) :: isharedface
402 integer(I4B) :: nfaces
403 integer(I4B) :: nslots
406 nfaces = defn%npolyverts + 3
407 nslots =
size(defn%facenbr)
408 if (nslots < nfaces)
call expandarray(defn%facenbr, nfaces - nslots)
410 select type (dis => this%fmi%dis)
415 icu1 = dis%get_nodeuser(ic1)
416 ncpl = dis%get_ncpl()
417 call get_jk(icu1, ncpl, dis%nlay, j1, k1)
418 istart1 = dis%iavert(j1)
419 istop1 = dis%iavert(j1 + 1) - 1
420 do iloc = 1, dis%con%ia(ic1 + 1) - dis%con%ia(ic1) - 1
421 ipos = dis%con%ia(ic1) + iloc
423 if (dis%con%mask(ipos) == 0) cycle
424 ic2 = dis%con%ja(ipos)
425 icu2 = dis%get_nodeuser(ic2)
426 call get_jk(icu2, ncpl, dis%nlay, j2, k2)
427 istart2 = dis%iavert(j2)
428 istop2 = dis%iavert(j2 + 1) - 1
430 dis%javert(istart2:istop2), &
432 if (isharedface /= 0)
then
434 defn%facenbr(isharedface) = int(iloc, 1)
438 defn%facenbr(defn%npolyverts + 2) = int(iloc, 1)
439 else if (k2 < k1)
then
441 defn%facenbr(defn%npolyverts + 3) = int(iloc, 1)
443 call pstop(1,
"k2 should be <> k1, since no shared edge face")
449 defn%facenbr(defn%npolyverts + 1) = defn%facenbr(1)
460 integer(I4B) :: nfaces
461 integer(I4B) :: nslots
464 nfaces = defn%npolyverts + 3
465 nslots =
size(defn%faceflow)
466 if (nslots < nfaces)
call expandarray(defn%faceflow, nfaces - nslots)
472 defn%faceflow =
dzero
474 call this%load_boundary_flows_to_defn_poly(defn)
475 call this%load_face_flows_to_defn_poly(defn)
478 defn%distflow = this%fmi%SourceFlows(defn%icell) + &
479 this%fmi%SinkFlows(defn%icell) + &
480 this%fmi%StorageFlows(defn%icell)
483 if (this%fmi%SinkFlows(defn%icell) .ne.
dzero)
then
495 integer(I4B) :: m, n, nfaces
498 nfaces = defn%npolyverts + 3
502 q = this%fmi%gwfflowja(this%fmi%dis%con%ia(defn%icell) + n)
503 defn%faceflow(m) = defn%faceflow(m) + q
505 if (defn%faceflow(m) <
dzero) defn%inoexitface = 0
516 integer(I4B) :: ioffset
520 ioffset = (defn%icell - 1) * 10
521 defn%faceflow(1) = defn%faceflow(1) + &
522 this%fmi%BoundaryFlows(ioffset + 4)
523 defn%faceflow(2) = defn%faceflow(2) + &
524 this%fmi%BoundaryFlows(ioffset + 2)
525 defn%faceflow(3) = defn%faceflow(3) + &
526 this%fmi%BoundaryFlows(ioffset + 3)
527 defn%faceflow(4) = defn%faceflow(4) + &
528 this%fmi%BoundaryFlows(ioffset + 1)
529 defn%faceflow(5) = defn%faceflow(1)
530 defn%faceflow(6) = defn%faceflow(6) + &
531 this%fmi%BoundaryFlows(ioffset + 9)
532 defn%faceflow(7) = defn%faceflow(7) + &
533 this%fmi%BoundaryFlows(ioffset + 10)
546 integer(I4B) :: ioffset
550 integer(I4B) :: mdiff
552 integer(I4B) :: irectvert(5)
554 ioffset = (defn%icell - 1) * 10
560 else if (n .eq. 4)
then
565 qbf = this%fmi%BoundaryFlows(ioffset + nbf)
567 do m = 1, defn%npolyverts
568 if (.not. defn%ispv180(m))
then
573 irectvert(5) = irectvert(1)
575 m2 = irectvert(n + 1)
576 if (m2 .lt. m1) m2 = m2 + defn%npolyverts
578 if (mdiff .eq. 1)
then
580 defn%faceflow(m1) = defn%faceflow(m1) + qbf
584 defn%faceflow(m1) = defn%faceflow(m1) + qbf
585 defn%faceflow(m1 + 1) = defn%faceflow(m1 + 1) + qbf
589 m = defn%npolyverts + 1
590 defn%faceflow(m) = defn%faceflow(1)
593 defn%faceflow(m) = defn%faceflow(m) + &
594 this%fmi%BoundaryFlows(ioffset + 9)
597 defn%faceflow(m) = defn%faceflow(m) + &
598 this%fmi%BoundaryFlows(ioffset + 10)
610 integer(I4B) :: npolyverts
611 integer(I4B) :: ioffset
615 npolyverts = defn%npolyverts
618 do iv = 1, npolyverts
619 defn%faceflow(iv) = &
620 defn%faceflow(iv) + &
621 this%fmi%BoundaryFlows(ioffset + iv)
623 defn%faceflow(npolyverts + 1) = defn%faceflow(1)
624 defn%faceflow(npolyverts + 2) = &
625 defn%faceflow(npolyverts + 2) + &
627 defn%faceflow(npolyverts + 3) = &
628 defn%faceflow(npolyverts + 3) + &
642 integer(I4B) :: npolyverts
648 integer(I4B) :: num90
649 integer(I4B) :: num180
660 s0mag, s2x, s2y, s2mag, sinang, cosang, dotprod
664 npolyverts = defn%npolyverts
666 if (
size(defn%ispv180) < npolyverts + 3) &
669 defn%ispv180(1:npolyverts + 1) = .false.
670 defn%can_be_rect = .false.
671 defn%can_be_quad = .false.
683 else if (m1 .eq. npolyverts)
then
690 x0 = defn%polyvert(1, m0)
691 y0 = defn%polyvert(2, m0)
692 x1 = defn%polyvert(1, m1)
693 y1 = defn%polyvert(2, m1)
694 x2 = defn%polyvert(1, m2)
695 y2 = defn%polyvert(2, m2)
698 s0mag = dsqrt(s0x * s0x + s0y * s0y)
701 s2mag = dsqrt(s2x * s2x + s2y * s2y)
702 sinang = (s0x * s2y - s0y * s2x) / (s0mag * s2mag)
703 cosang = dsqrt(dabs(
done - (sinang * sinang)))
704 if (dabs(sinang) .lt. epsang)
then
705 dotprod = s0x * s2x + s0y * s2y
706 if (dotprod .lt.
dzero)
then
709 defn%ispv180(m) = .true.
712 if (dabs(cosang) .lt. epsang) num90 = num90 + 1
718 defn%ispv180(npolyverts + 1) = defn%ispv180(1)
721 if (num90 .eq. 4)
then
722 if (num180 .eq. 0)
then
723 defn%can_be_rect = .true.
725 defn%can_be_quad = .true.
subroutine, public create_defn(cellDefn)
Create a new cell definition object.
integer(i4b) function, public get_iatop(ncpl, icu)
Get the index corresponding to top elevation of a cell in the grid. This is -1 if the cell is in the ...
integer(i4b), parameter, public max_poly_cells
subroutine, public create_cell_poly(cell)
Create a new polygonal cell.
subroutine, public cell_poly_to_rect(poly, rect)
Convert CellPoly representation to CellRect. Assumes the conversion is possible.
subroutine, public cell_poly_to_quad(poly, quad)
Convert CellPoly representation to CellRectQuad. Assumes the conversion is possible.
This module contains simulation constants.
real(dp), parameter dzero
real constant zero
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 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.
Cell-level tracking methods.
type(methodcellpollocktype), pointer, public method_cell_plck
type(methodcellpasstobottype), pointer, public method_cell_ptb
type(methodcellpollockquadtype), pointer, public method_cell_quad
type(methodcellternarytype), pointer, public method_cell_tern
subroutine, public create_method_disv(method)
Create a new vertex grid (DISV) tracking method.
subroutine load_boundary_flows_to_defn_rect(this, defn)
Load boundary flows from the grid into a rectangular cell. Assumes cell index and number of vertices ...
subroutine load_neighbors(this, defn)
Loads face neighbors to cell definition from the grid Assumes cell index and number of vertices are a...
subroutine map_neighbor(this, defn, inface, z)
Map location on cell face to shared cell face of neighbor.
subroutine load_indicators(this, defn)
Load 180-degree vertex indicator array and set flags indicating how cell can be represented....
subroutine load_boundary_flows_to_defn_rect_quad(this, defn)
Load boundary flows from the grid into rectangular quadcell. Assumes cell index and number of vertice...
subroutine update_flowja(this, cell, particle)
subroutine load_disv(this, particle, next_level, submethod)
Load the cell and the tracking method.
subroutine load_face_flows_to_defn_poly(this, defn)
subroutine load_polygon(this, defn)
subroutine load_particle(this, cell, particle)
subroutine load_properties(this, ic, defn)
Loads cell properties to cell definition from the grid.
subroutine load_flows(this, defn)
Load flows into the cell definition. These include face, boundary and net distributed flows....
subroutine load_boundary_flows_to_defn_poly(this, defn)
Load boundary flows from the grid into a polygonal cell. Assumes cell index and number of vertices ar...
subroutine load_cell_defn(this, ic, defn)
Load cell definition from the grid.
subroutine apply_disv(this, particle, tmax)
Apply the DISV tracking method to a particle.
subroutine pass_disv(this, particle)
Pass a particle to the next cell, if there is one.
Particle tracking strategies.
Base grid cell definition.
Base type for grid cells of a concrete type. Contains a cell-definition which is information shared b...
Vertex grid discretization.
Base type for particle tracking methods.
Particle tracked by the PRT model.