64 allocate (method%name)
68 method%delegates = .true.
81 deallocate (this%name)
84 call this%method_cell_plck%deallocate()
85 deallocate (this%method_cell_plck)
86 call this%method_cell_quad%deallocate()
87 deallocate (this%method_cell_quad)
88 call this%method_cell_tern%deallocate()
89 deallocate (this%method_cell_tern)
90 call this%method_cell_ptb%deallocate()
91 deallocate (this%method_cell_ptb)
95 subroutine load_disv(this, particle, next_level, submethod)
103 integer(I4B),
intent(in) :: next_level
104 class(
methodtype),
pointer,
intent(inout) :: submethod
111 select type (cell => this%cell)
113 ic = particle%itrdomain(next_level)
114 call this%load_cell_defn(ic, cell%defn)
115 if (this%fmi%ibdgwfsat0(ic) == 0)
then
118 call this%method_cell_ptb%init( &
121 events=this%events, &
122 tracktimes=this%tracktimes)
123 submethod => this%method_cell_ptb
124 else if (particle%frctrn)
then
126 call this%method_cell_tern%init( &
129 events=this%events, &
130 tracktimes=this%tracktimes)
131 submethod => this%method_cell_tern
132 else if (cell%defn%can_be_rect)
then
137 call this%method_cell_plck%init( &
140 events=this%events, &
141 tracktimes=this%tracktimes)
142 submethod => this%method_cell_plck
143 else if (cell%defn%can_be_quad)
then
148 call this%method_cell_quad%init( &
151 events=this%events, &
152 tracktimes=this%tracktimes)
153 submethod => this%method_cell_quad
156 call this%method_cell_tern%init( &
159 events=this%events, &
160 tracktimes=this%tracktimes)
161 submethod => this%method_cell_tern
176 integer(I4B) :: inface
181 integer(I4B) :: idiag
186 select type (dis => this%fmi%dis)
189 idiag = dis%con%ia(cell%defn%icell)
190 inbr = cell%defn%facenbr(inface)
192 ic = dis%con%ja(ipos)
193 icu = dis%get_nodeuser(ic)
194 call get_jk(icu, dis%ncpl, dis%nlay, icpl, ilay)
201 call this%map_neighbor(cell%defn, inface, z)
218 integer(I4B) :: idiag
221 idiag = this%fmi%dis%con%ia(cell%defn%icell)
226 this%flowja(ipos) = this%flowja(ipos) -
done
229 this%flowja(this%fmi%dis%con%isym(ipos)) &
230 = this%flowja(this%fmi%dis%con%isym(ipos)) +
done
243 integer(I4B) :: iface
244 logical(LGP) :: no_neighbors
246 select type (c => this%cell)
250 no_neighbors = cell%defn%facenbr(iface) == 0
253 if (no_neighbors)
then
258 call this%load_particle(cell, particle)
259 call this%update_flowja(cell, particle)
268 integer(I4B),
intent(inout) :: inface
269 double precision,
intent(inout) :: z
272 integer(I4B) :: npolyvertsin
274 integer(I4B) :: npolyverts
276 integer(I4B) :: inbrnbr
287 inbr = defn%facenbr(inface)
288 if (inbr .eq. 0)
then
295 j = this%fmi%dis%con%ia(icin)
296 ic = this%fmi%dis%con%ja(j + inbr)
297 call this%load_cell_defn(ic, this%neighbor)
299 npolyvertsin = defn%npolyverts
300 npolyverts = this%neighbor%npolyverts
301 if (inface .eq. npolyvertsin + 2)
then
303 inface = npolyverts + 3
304 else if (inface .eq. npolyvertsin + 3)
then
306 inface = npolyverts + 2
309 j = this%fmi%dis%con%ia(ic)
310 do m = 1, npolyverts + 3
311 inbrnbr = this%neighbor%facenbr(m)
312 if (this%fmi%dis%con%ja(j + inbrnbr) .eq. icin)
then
320 zrel = (z - botfrom) / (topfrom - botfrom)
321 top = this%fmi%dis%top(ic)
322 bot = this%fmi%dis%bot(ic)
323 sat = this%fmi%gwfsat(ic)
324 z = bot + zrel * sat * (top - bot)
333 real(DP),
intent(in) :: tmax
335 call this%track(particle, 1, tmax)
342 integer(I4B),
intent(in) :: ic
345 call this%load_cell_properties(ic, defn)
346 call this%load_cell_polygon(defn)
347 call this%load_cell_neighbors(defn)
348 call this%load_cell_saturation_status(defn)
349 call this%load_cell_flags(defn)
350 call this%load_cell_flows(defn)
357 integer(I4B),
intent(in) :: ic
363 integer(I4B) :: icu, icpl, ilay
366 defn%iatop =
get_iatop(this%fmi%dis%get_ncpl(), &
367 this%fmi%dis%get_nodeuser(ic))
368 top = this%fmi%dis%top(ic)
369 bot = this%fmi%dis%bot(ic)
370 sat = this%fmi%gwfsat(ic)
371 top = bot + sat * (top - bot)
375 defn%porosity = this%porosity(ic)
376 defn%retfactor = this%retfactor(ic)
377 defn%izone = this%izone(ic)
378 select type (dis => this%fmi%dis)
380 icu = dis%get_nodeuser(ic)
381 call get_jk(icu, dis%ncpl, dis%nlay, icpl, ilay)
392 call this%fmi%dis%get_polyverts( &
396 defn%npolyverts =
size(defn%polyvert, dim=2) - 1
416 integer(I4B) :: istart1
417 integer(I4B) :: istart2
418 integer(I4B) :: istop1
419 integer(I4B) :: istop2
420 integer(I4B) :: isharedface
422 integer(I4B) :: nfaces
423 integer(I4B) :: nslots
426 nfaces = defn%npolyverts + 3
427 nslots =
size(defn%facenbr)
428 if (nslots < nfaces)
call expandarray(defn%facenbr, nfaces - nslots)
430 select type (dis => this%fmi%dis)
435 icu1 = dis%get_nodeuser(ic1)
436 ncpl = dis%get_ncpl()
437 call get_jk(icu1, ncpl, dis%nlay, j1, k1)
438 istart1 = dis%iavert(j1)
439 istop1 = dis%iavert(j1 + 1) - 1
440 do iloc = 1, dis%con%ia(ic1 + 1) - dis%con%ia(ic1) - 1
441 ipos = dis%con%ia(ic1) + iloc
443 if (dis%con%mask(ipos) == 0) cycle
444 ic2 = dis%con%ja(ipos)
445 icu2 = dis%get_nodeuser(ic2)
446 call get_jk(icu2, ncpl, dis%nlay, j2, k2)
447 istart2 = dis%iavert(j2)
448 istop2 = dis%iavert(j2 + 1) - 1
450 dis%javert(istart2:istop2), &
452 if (isharedface /= 0)
then
454 defn%facenbr(isharedface) = int(iloc, 1)
458 defn%facenbr(defn%npolyverts + 2) = int(iloc, 1)
459 else if (k2 < k1)
then
461 defn%facenbr(defn%npolyverts + 3) = int(iloc, 1)
463 call pstop(1,
"k2 should be <> k1, since no shared edge face")
469 defn%facenbr(defn%npolyverts + 1) = defn%facenbr(1)
480 integer(I4B) :: nfaces, nslots
483 nfaces = defn%npolyverts + 3
484 nslots =
size(defn%faceflow)
485 if (nslots < nfaces)
call expandarray(defn%faceflow, nfaces - nslots)
491 defn%faceflow =
dzero
492 call this%load_cell_boundary_flows(defn)
493 call this%load_cell_face_flows(defn)
494 call this%cap_cell_wt_flow(defn)
495 call this%load_cell_no_exit_face(defn)
498 defn%distflow = this%fmi%SourceFlows(defn%icell) + &
499 this%fmi%SinkFlows(defn%icell) + &
500 this%fmi%StorageFlows(defn%icell)
503 if (this%fmi%SinkFlows(defn%icell) .ne.
dzero)
then
515 integer(I4B) :: m, n, nfaces
518 nfaces = defn%npolyverts + 3
522 q = this%fmi%gwfflowja(this%fmi%dis%con%ia(defn%icell) + n)
523 defn%faceflow(m) = defn%faceflow(m) + q
536 integer(I4B) :: ic, iv, npolyverts
539 npolyverts = defn%npolyverts
540 do iv = 1, npolyverts
541 defn%faceflow(iv) = &
542 defn%faceflow(iv) + &
543 this%fmi%BoundaryFlows(ic, iv)
545 defn%faceflow(npolyverts + 1) = defn%faceflow(1)
546 defn%faceflow(npolyverts + 2) = &
547 defn%faceflow(npolyverts + 2) + &
548 this%fmi%BoundaryFlows(ic, this%fmi%max_faces - 1)
549 defn%faceflow(npolyverts + 3) = &
550 defn%faceflow(npolyverts + 3) + &
551 this%fmi%BoundaryFlows(ic, this%fmi%max_faces)
564 integer(I4B) :: npolyverts
570 integer(I4B) :: num90
571 integer(I4B) :: num180
582 s0mag, s2x, s2y, s2mag, sinang, cosang, dotprod
586 npolyverts = defn%npolyverts
588 if (
size(defn%ispv180) < npolyverts + 3) &
591 defn%ispv180(1:npolyverts + 1) = .false.
592 defn%can_be_rect = .false.
593 defn%can_be_quad = .false.
605 else if (m1 .eq. npolyverts)
then
612 x0 = defn%polyvert(1, m0)
613 y0 = defn%polyvert(2, m0)
614 x1 = defn%polyvert(1, m1)
615 y1 = defn%polyvert(2, m1)
616 x2 = defn%polyvert(1, m2)
617 y2 = defn%polyvert(2, m2)
620 s0mag = dsqrt(s0x * s0x + s0y * s0y)
623 s2mag = dsqrt(s2x * s2x + s2y * s2y)
624 sinang = (s0x * s2y - s0y * s2x) / (s0mag * s2mag)
625 cosang = dsqrt(dabs(
done - (sinang * sinang)))
626 if (dabs(sinang) .lt. epsang)
then
627 dotprod = s0x * s2x + s0y * s2y
628 if (dotprod .lt.
dzero)
then
631 defn%ispv180(m) = .true.
634 if (dabs(cosang) .lt. epsang) num90 = num90 + 1
640 defn%ispv180(npolyverts + 1) = defn%ispv180(1)
643 if (num90 .eq. 4)
then
644 if (num180 .eq. 0)
then
645 defn%can_be_rect = .true.
647 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 ...
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.
subroutine, public create_method_cell_ptb(method)
Create a new pass-to-bottom tracking method.
procedure subroutine, public create_method_cell_pollock(method)
Create a tracking method.
procedure subroutine, public create_method_cell_quad(method)
Create a new Pollock quad-refined cell method.
subroutine, public create_method_cell_ternary(method)
Create a tracking method.
subroutine, public create_method_disv(method)
Create a new vertex grid (DISV) tracking method.
subroutine load_cell_boundary_flows(this, defn)
Load boundary flows from the grid into a polygonal cell. Assumes cell index and number of vertices ar...
subroutine load_cell_polygon(this, defn)
subroutine load_cell_flows(this, defn)
Load flows into the cell definition. These include face, boundary and net distributed flows....
subroutine load_cell_flags(this, defn)
Load 180-degree vertex indicator array and set flags indicating how cell can be represented....
subroutine map_neighbor(this, defn, inface, z)
Map location on cell face to shared cell face of neighbor.
subroutine load_cell_properties(this, ic, defn)
Loads cell properties to cell definition from the grid.
subroutine update_flowja(this, cell, particle)
subroutine load_cell_face_flows(this, defn)
subroutine load_cell_neighbors(this, defn)
Loads face neighbors to cell definition from the grid Assumes cell index and number of vertices are a...
subroutine load_disv(this, particle, next_level, submethod)
Load the cell and the tracking method.
subroutine load_particle(this, cell, particle)
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.
@, public level_subfeature
@, public terminate
particle terminated
@ term_boundary
terminated at a boundary face
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.