56 allocate (method%type)
60 method%delegates = .true.
67 deallocate (this%type)
71 subroutine load_disv(this, particle, next_level, submethod)
79 integer(I4B),
intent(in) :: next_level
80 class(
methodtype),
pointer,
intent(inout) :: submethod
87 select type (cell => this%cell)
90 ic = particle%idomain(next_level)
91 call this%load_cell_defn(ic, cell%defn)
92 if (this%fmi%ibdgwfsat0(ic) == 0)
then
97 trackfilectl=this%trackfilectl, &
98 tracktimes=this%tracktimes)
102 if (particle%ifrctrn > 0)
then
105 trackfilectl=this%trackfilectl, &
106 tracktimes=this%tracktimes)
108 else if (cell%defn%can_be_rect)
then
113 trackfilectl=this%trackfilectl, &
114 tracktimes=this%tracktimes)
116 else if (cell%defn%can_be_quad)
then
121 trackfilectl=this%trackfilectl, &
122 tracktimes=this%tracktimes)
127 trackfilectl=this%trackfilectl, &
128 tracktimes=this%tracktimes)
143 integer(I4B) :: inface
148 integer(I4B) :: idiag
153 select type (dis => this%fmi%dis)
156 inface = particle%iboundary(2)
157 idiag = dis%con%ia(cell%defn%icell)
158 inbr = cell%defn%facenbr(inface)
160 ic = dis%con%ja(ipos)
161 icu = dis%get_nodeuser(ic)
162 call get_jk(icu, dis%ncpl, dis%nlay, icpl, ilay)
163 particle%idomain(2) = ic
169 call this%map_neighbor(cell%defn, inface, z)
170 particle%iboundary(2) = inface
171 particle%idomain(3:) = 0
172 particle%iboundary(3:) = 0
185 integer(I4B) :: idiag
188 idiag = this%fmi%dis%con%ia(cell%defn%icell)
189 inbr = cell%defn%facenbr(particle%iboundary(2))
193 this%flowja(ipos) = this%flowja(ipos) -
done
196 this%flowja(this%fmi%dis%con%isym(ipos)) &
197 = this%flowja(this%fmi%dis%con%isym(ipos)) +
done
209 select type (c => this%cell)
215 if (cell%defn%facenbr(particle%iboundary(2)) .eq. 0)
then
217 particle%advancing = .false.
218 call this%save(particle, reason=3)
222 call this%load_particle(cell, particle)
223 call this%update_flowja(cell, particle)
233 integer(I4B),
intent(inout) :: inface
234 double precision,
intent(inout) :: z
237 integer(I4B) :: npolyvertsin
239 integer(I4B) :: npolyverts
241 integer(I4B) :: inbrnbr
252 inbr = defn%facenbr(inface)
253 if (inbr .eq. 0)
then
260 j = this%fmi%dis%con%ia(icin)
261 ic = this%fmi%dis%con%ja(j + inbr)
262 call this%load_cell_defn(ic, this%neighbor)
264 npolyvertsin = defn%npolyverts
265 npolyverts = this%neighbor%npolyverts
266 if (inface .eq. npolyvertsin + 2)
then
268 inface = npolyverts + 3
269 else if (inface .eq. npolyvertsin + 3)
then
271 inface = npolyverts + 2
274 j = this%fmi%dis%con%ia(ic)
275 do m = 1, npolyverts + 3
276 inbrnbr = this%neighbor%facenbr(m)
277 if (this%fmi%dis%con%ja(j + inbrnbr) .eq. icin)
then
285 zrel = (z - botfrom) / (topfrom - botfrom)
286 top = this%fmi%dis%top(ic)
287 bot = this%fmi%dis%bot(ic)
288 sat = this%fmi%gwfsat(ic)
289 z = bot + zrel * sat * (top - bot)
298 real(DP),
intent(in) :: tmax
299 call this%track(particle, 1, tmax)
306 integer(I4B),
intent(in) :: ic
310 call this%load_properties(ic, defn)
313 call this%load_polygon(defn)
316 call this%load_neighbors(defn)
319 call this%load_indicators(defn)
322 call this%load_flows(defn)
329 integer(I4B),
intent(in) :: ic
337 defn%iatop =
get_iatop(this%fmi%dis%get_ncpl(), &
338 this%fmi%dis%get_nodeuser(ic))
339 top = this%fmi%dis%top(ic)
340 bot = this%fmi%dis%bot(ic)
341 sat = this%fmi%gwfsat(ic)
342 top = bot + sat * (top - bot)
346 defn%porosity = this%porosity(ic)
347 defn%retfactor = this%retfactor(ic)
348 defn%izone = this%izone(ic)
357 call this%fmi%dis%get_polyverts( &
361 defn%npolyverts =
size(defn%polyvert, dim=2) - 1
382 integer(I4B) :: istart1
383 integer(I4B) :: istart2
384 integer(I4B) :: istop1
385 integer(I4B) :: istop2
386 integer(I4B) :: isharedface
388 integer(I4B) :: nfaces
389 integer(I4B) :: nslots
392 nfaces = defn%npolyverts + 3
393 nslots =
size(defn%facenbr)
394 if (nslots < nfaces)
call expandarray(defn%facenbr, nfaces - nslots)
396 select type (dis => this%fmi%dis)
401 icu1 = dis%get_nodeuser(ic1)
402 ncpl = dis%get_ncpl()
403 call get_jk(icu1, ncpl, dis%nlay, j1, k1)
404 istart1 = dis%iavert(j1)
405 istop1 = dis%iavert(j1 + 1) - 1
406 do iloc = 1, dis%con%ia(ic1 + 1) - dis%con%ia(ic1) - 1
407 ipos = dis%con%ia(ic1) + iloc
409 if (dis%con%mask(ipos) == 0) cycle
410 ic2 = dis%con%ja(ipos)
411 icu2 = dis%get_nodeuser(ic2)
412 call get_jk(icu2, ncpl, dis%nlay, j2, k2)
413 istart2 = dis%iavert(j2)
414 istop2 = dis%iavert(j2 + 1) - 1
416 dis%javert(istart2:istop2), &
418 if (isharedface /= 0)
then
420 defn%facenbr(isharedface) = int(iloc, 1)
424 defn%facenbr(defn%npolyverts + 2) = int(iloc, 1)
425 else if (k2 < k1)
then
427 defn%facenbr(defn%npolyverts + 3) = int(iloc, 1)
429 call pstop(1,
"k2 should be <> k1, since no shared edge face")
435 defn%facenbr(defn%npolyverts + 1) = defn%facenbr(1)
447 integer(I4B) :: nfaces
448 integer(I4B) :: nslots
451 nfaces = defn%npolyverts + 3
452 nslots =
size(defn%faceflow)
453 if (nslots < nfaces)
call expandarray(defn%faceflow, nfaces - nslots)
459 defn%faceflow =
dzero
461 call this%load_boundary_flows_to_defn_poly(defn)
462 call this%load_face_flows_to_defn_poly(defn)
465 defn%distflow = this%fmi%SourceFlows(defn%icell) + &
466 this%fmi%SinkFlows(defn%icell) + &
467 this%fmi%StorageFlows(defn%icell)
470 if (this%fmi%SinkFlows(defn%icell) .ne.
dzero)
then
483 integer(I4B) :: m, n, nfaces
486 nfaces = defn%npolyverts + 3
490 q = this%fmi%gwfflowja(this%fmi%dis%con%ia(defn%icell) + n)
491 defn%faceflow(m) = defn%faceflow(m) + q
493 if (defn%faceflow(m) <
dzero) defn%inoexitface = 0
504 integer(I4B) :: ioffset
508 ioffset = (defn%icell - 1) * 10
509 defn%faceflow(1) = defn%faceflow(1) + &
510 this%fmi%BoundaryFlows(ioffset + 4)
511 defn%faceflow(2) = defn%faceflow(2) + &
512 this%fmi%BoundaryFlows(ioffset + 2)
513 defn%faceflow(3) = defn%faceflow(3) + &
514 this%fmi%BoundaryFlows(ioffset + 3)
515 defn%faceflow(4) = defn%faceflow(4) + &
516 this%fmi%BoundaryFlows(ioffset + 1)
517 defn%faceflow(5) = defn%faceflow(1)
518 defn%faceflow(6) = defn%faceflow(6) + &
519 this%fmi%BoundaryFlows(ioffset + 9)
520 defn%faceflow(7) = defn%faceflow(7) + &
521 this%fmi%BoundaryFlows(ioffset + 10)
535 integer(I4B) :: ioffset
539 integer(I4B) :: mdiff
541 integer(I4B) :: irectvert(5)
543 ioffset = (defn%icell - 1) * 10
549 else if (n .eq. 4)
then
554 qbf = this%fmi%BoundaryFlows(ioffset + nbf)
556 do m = 1, defn%npolyverts
557 if (.not. defn%ispv180(m))
then
562 irectvert(5) = irectvert(1)
564 m2 = irectvert(n + 1)
565 if (m2 .lt. m1) m2 = m2 + defn%npolyverts
567 if (mdiff .eq. 1)
then
569 defn%faceflow(m1) = defn%faceflow(m1) + qbf
573 defn%faceflow(m1) = defn%faceflow(m1) + qbf
574 defn%faceflow(m1 + 1) = defn%faceflow(m1 + 1) + qbf
578 m = defn%npolyverts + 1
579 defn%faceflow(m) = defn%faceflow(1)
582 defn%faceflow(m) = defn%faceflow(m) + &
583 this%fmi%BoundaryFlows(ioffset + 9)
586 defn%faceflow(m) = defn%faceflow(m) + &
587 this%fmi%BoundaryFlows(ioffset + 10)
599 integer(I4B) :: npolyverts
600 integer(I4B) :: ioffset
604 npolyverts = defn%npolyverts
607 do iv = 1, npolyverts
608 defn%faceflow(iv) = &
609 defn%faceflow(iv) + &
610 this%fmi%BoundaryFlows(ioffset + iv)
612 defn%faceflow(npolyverts + 1) = defn%faceflow(1)
613 defn%faceflow(npolyverts + 2) = &
614 defn%faceflow(npolyverts + 2) + &
616 defn%faceflow(npolyverts + 3) = &
617 defn%faceflow(npolyverts + 3) + &
630 integer(I4B) :: npolyverts
636 integer(I4B) :: num90
637 integer(I4B) :: num180
648 s0mag, s2x, s2y, s2mag, sinang, cosang, dotprod
652 npolyverts = defn%npolyverts
655 if (
size(defn%ispv180) < npolyverts + 3) &
660 defn%ispv180(1:npolyverts + 1) = .false.
661 defn%can_be_rect = .false.
662 defn%can_be_quad = .false.
673 else if (m1 .eq. npolyverts)
then
680 x0 = defn%polyvert(1, m0)
681 y0 = defn%polyvert(2, m0)
682 x1 = defn%polyvert(1, m1)
683 y1 = defn%polyvert(2, m1)
684 x2 = defn%polyvert(1, m2)
685 y2 = defn%polyvert(2, m2)
688 s0mag = dsqrt(s0x * s0x + s0y * s0y)
691 s2mag = dsqrt(s2x * s2x + s2y * s2y)
692 sinang = (s0x * s2y - s0y * s2x) / (s0mag * s2mag)
693 cosang = dsqrt(dabs(
done - (sinang * sinang)))
694 if (dabs(sinang) .lt. epsang)
then
695 dotprod = s0x * s2x + s0y * s2y
696 if (dotprod .lt.
dzero)
then
699 defn%ispv180(m) = .true.
702 if (dabs(cosang) .lt. epsang) num90 = num90 + 1
708 defn%ispv180(npolyverts + 1) = defn%ispv180(1)
710 if (num90 .eq. 4)
then
711 if (num180 .eq. 0)
then
712 defn%can_be_rect = .true.
714 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-grid method.
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.
Manages particle track (i.e. pathline) files.