MODFLOW 6  version 6.7.0.dev0
USGS Modular Hydrologic Model
methoddisvmodule Module Reference

Data Types

type  methoddisvtype
 

Functions/Subroutines

subroutine, public create_method_disv (method)
 Create a new vertex grid (DISV) tracking method. More...
 
subroutine deallocate (this)
 Destroy the tracking method. More...
 
subroutine load_disv (this, particle, next_level, submethod)
 Load the cell and the tracking method. More...
 
subroutine load_particle (this, cell, particle)
 
subroutine update_flowja (this, cell, particle)
 
subroutine pass_disv (this, particle)
 Pass a particle to the next cell, if there is one. More...
 
subroutine map_neighbor (this, defn, inface, z)
 Map location on cell face to shared cell face of neighbor. More...
 
subroutine apply_disv (this, particle, tmax)
 Apply the DISV tracking method to a particle. More...
 
subroutine load_cell_defn (this, ic, defn)
 Load cell definition from the grid. More...
 
subroutine load_properties (this, ic, defn)
 Loads cell properties to cell definition from the grid. More...
 
subroutine load_polygon (this, defn)
 
subroutine load_neighbors (this, defn)
 Loads face neighbors to cell definition from the grid Assumes cell index and number of vertices are already loaded. More...
 
subroutine load_flows (this, defn)
 Load flows into the cell definition. These include face, boundary and net distributed flows. Assumes cell index and number of vertices are already loaded. More...
 
subroutine load_face_flows_to_defn_poly (this, defn)
 
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 are already loaded. More...
 
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 vertices are already loaded. More...
 
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 are already loaded. More...
 
subroutine load_indicators (this, defn)
 Load 180-degree vertex indicator array and set flags indicating how cell can be represented. Assumes cell index and number of vertices are already loaded. More...
 

Function/Subroutine Documentation

◆ apply_disv()

subroutine methoddisvmodule::apply_disv ( class(methoddisvtype), intent(inout)  this,
type(particletype), intent(inout), pointer  particle,
real(dp), intent(in)  tmax 
)
private

Definition at line 312 of file MethodDisv.f90.

313  class(MethodDisvType), intent(inout) :: this
314  type(ParticleType), pointer, intent(inout) :: particle
315  real(DP), intent(in) :: tmax
316 
317  call this%track(particle, 1, tmax)

◆ create_method_disv()

subroutine, public methoddisvmodule::create_method_disv ( type(methoddisvtype), pointer  method)

Definition at line 48 of file MethodDisv.f90.

49  ! dummy
50  type(MethodDisvType), pointer :: method
51  ! local
52  type(CellPolyType), pointer :: cell
53 
54  allocate (method)
55  allocate (method%name)
56  call create_cell_poly(cell)
57  method%cell => cell
58  method%name = "disv"
59  method%delegates = .true.
60  call create_defn(method%neighbor)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ deallocate()

subroutine methoddisvmodule::deallocate ( class(methoddisvtype), intent(inout)  this)
private

Definition at line 64 of file MethodDisv.f90.

65  class(MethodDisvType), intent(inout) :: this
66  deallocate (this%name)

◆ load_boundary_flows_to_defn_poly()

subroutine methoddisvmodule::load_boundary_flows_to_defn_poly ( class(methoddisvtype), intent(inout)  this,
type(celldefntype), intent(inout)  defn 
)
private

Definition at line 604 of file MethodDisv.f90.

605  ! dummy
606  class(MethodDisvType), intent(inout) :: this
607  type(CellDefnType), intent(inout) :: defn
608  ! local
609  integer(I4B) :: ic
610  integer(I4B) :: npolyverts
611  integer(I4B) :: ioffset
612  integer(I4B) :: iv
613 
614  ic = defn%icell
615  npolyverts = defn%npolyverts
616 
617  ioffset = (ic - 1) * max_poly_cells
618  do iv = 1, npolyverts
619  defn%faceflow(iv) = &
620  defn%faceflow(iv) + &
621  this%fmi%BoundaryFlows(ioffset + iv)
622  end do
623  defn%faceflow(npolyverts + 1) = defn%faceflow(1)
624  defn%faceflow(npolyverts + 2) = &
625  defn%faceflow(npolyverts + 2) + &
626  this%fmi%BoundaryFlows(ioffset + max_poly_cells - 1)
627  defn%faceflow(npolyverts + 3) = &
628  defn%faceflow(npolyverts + 3) + &
629  this%fmi%BoundaryFlows(ioffset + max_poly_cells)
630 

◆ load_boundary_flows_to_defn_rect()

subroutine methoddisvmodule::load_boundary_flows_to_defn_rect ( class(methoddisvtype), intent(inout)  this,
type(celldefntype), intent(inout)  defn 
)
private

Definition at line 511 of file MethodDisv.f90.

512  ! dummy
513  class(MethodDisvType), intent(inout) :: this
514  type(CellDefnType), intent(inout) :: defn
515  ! local
516  integer(I4B) :: ioffset
517 
518  ! assignment of BoundaryFlows to faceflow below assumes clockwise
519  ! ordering of faces, with face 1 being the "western" face
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)

◆ load_boundary_flows_to_defn_rect_quad()

subroutine methoddisvmodule::load_boundary_flows_to_defn_rect_quad ( class(methoddisvtype), intent(inout)  this,
type(celldefntype), intent(inout)  defn 
)
private

Definition at line 538 of file MethodDisv.f90.

539  ! dummy
540  class(MethodDisvType), intent(inout) :: this
541  type(CellDefnType), intent(inout) :: defn
542  ! local
543  integer(I4B) :: m
544  integer(I4B) :: n
545  integer(I4B) :: nn
546  integer(I4B) :: ioffset
547  integer(I4B) :: nbf
548  integer(I4B) :: m1
549  integer(I4B) :: m2
550  integer(I4B) :: mdiff
551  real(DP) :: qbf
552  integer(I4B) :: irectvert(5)
553 
554  ioffset = (defn%icell - 1) * 10
555 
556  ! Polygon faces in positions 1 through npolyverts
557  do n = 1, 4
558  if (n .eq. 2) then
559  nbf = 4
560  else if (n .eq. 4) then
561  nbf = 1
562  else
563  nbf = n
564  end if
565  qbf = this%fmi%BoundaryFlows(ioffset + nbf)
566  nn = 0
567  do m = 1, defn%npolyverts
568  if (.not. defn%ispv180(m)) then
569  nn = nn + 1
570  irectvert(nn) = m
571  end if
572  end do
573  irectvert(5) = irectvert(1)
574  m1 = irectvert(n)
575  m2 = irectvert(n + 1)
576  if (m2 .lt. m1) m2 = m2 + defn%npolyverts
577  mdiff = m2 - m1
578  if (mdiff .eq. 1) then
579  ! Assign BoundaryFlow to corresponding polygon face
580  defn%faceflow(m1) = defn%faceflow(m1) + qbf
581  else
582  ! Split BoundaryFlow between two faces on quad-refined edge
583  qbf = 5d-1 * qbf
584  defn%faceflow(m1) = defn%faceflow(m1) + qbf
585  defn%faceflow(m1 + 1) = defn%faceflow(m1 + 1) + qbf
586  end if
587  end do
588  ! Wrap around to 1 in position npolyverts+1
589  m = defn%npolyverts + 1
590  defn%faceflow(m) = defn%faceflow(1)
591  ! Bottom in position npolyverts+2
592  m = m + 1
593  defn%faceflow(m) = defn%faceflow(m) + &
594  this%fmi%BoundaryFlows(ioffset + 9)
595  ! Top in position npolyverts+3
596  m = m + 1
597  defn%faceflow(m) = defn%faceflow(m) + &
598  this%fmi%BoundaryFlows(ioffset + 10)
599 

◆ load_cell_defn()

subroutine methoddisvmodule::load_cell_defn ( class(methoddisvtype), intent(inout)  this,
integer(i4b), intent(in)  ic,
type(celldefntype), intent(inout), pointer  defn 
)
private

Definition at line 321 of file MethodDisv.f90.

322  ! dummy
323  class(MethodDisvType), intent(inout) :: this
324  integer(I4B), intent(in) :: ic
325  type(CellDefnType), pointer, intent(inout) :: defn
326 
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)

◆ load_disv()

subroutine methoddisvmodule::load_disv ( class(methoddisvtype), intent(inout)  this,
type(particletype), intent(inout), pointer  particle,
integer(i4b), intent(in)  next_level,
class(methodtype), intent(inout), pointer  submethod 
)
private

Definition at line 70 of file MethodDisv.f90.

71  use cellmodule
72  use cellrectmodule
74  use cellutilmodule
75  ! dummy
76  class(MethodDisvType), intent(inout) :: this
77  type(ParticleType), pointer, intent(inout) :: particle
78  integer(I4B), intent(in) :: next_level
79  class(MethodType), pointer, intent(inout) :: submethod
80  ! local
81  integer(I4B) :: ic
82  class(CellType), pointer :: base
83  type(CellRectType), pointer :: rect
84  type(CellRectQuadType), pointer :: quad
85 
86  select type (cell => this%cell)
87  type is (cellpolytype)
88  ic = particle%idomain(next_level)
89  call this%load_cell_defn(ic, cell%defn)
90  if (this%fmi%ibdgwfsat0(ic) == 0) then
91  call method_cell_ptb%init( &
92  fmi=this%fmi, &
93  cell=this%cell, &
94  trackctl=this%trackctl, &
95  tracktimes=this%tracktimes)
96  submethod => method_cell_ptb
97  else if (particle%ifrctrn > 0) then
98  call method_cell_tern%init( &
99  fmi=this%fmi, &
100  cell=this%cell, &
101  trackctl=this%trackctl, &
102  tracktimes=this%tracktimes)
103  submethod => method_cell_tern
104  else if (cell%defn%can_be_rect) then
105  call cell_poly_to_rect(cell, rect)
106  base => rect
107  call method_cell_plck%init( &
108  fmi=this%fmi, &
109  cell=base, &
110  trackctl=this%trackctl, &
111  tracktimes=this%tracktimes)
112  submethod => method_cell_plck
113  else if (cell%defn%can_be_quad) then
114  call cell_poly_to_quad(cell, quad)
115  base => quad
116  call method_cell_quad%init( &
117  fmi=this%fmi, &
118  cell=base, &
119  trackctl=this%trackctl, &
120  tracktimes=this%tracktimes)
121  submethod => method_cell_quad
122  else
123  call method_cell_tern%init( &
124  fmi=this%fmi, &
125  cell=this%cell, &
126  trackctl=this%trackctl, &
127  tracktimes=this%tracktimes)
128  submethod => method_cell_tern
129  end if
130  end select
subroutine, public cell_poly_to_rect(poly, rect)
Convert CellPoly representation to CellRect. Assumes the conversion is possible.
Definition: CellUtil.f90:14
subroutine, public cell_poly_to_quad(poly, quad)
Convert CellPoly representation to CellRectQuad. Assumes the conversion is possible.
Definition: CellUtil.f90:115
Here is the call graph for this function:

◆ load_face_flows_to_defn_poly()

subroutine methoddisvmodule::load_face_flows_to_defn_poly ( class(methoddisvtype), intent(inout)  this,
type(celldefntype), intent(inout)  defn 
)
private

Definition at line 490 of file MethodDisv.f90.

491  ! dummy
492  class(MethodDisvType), intent(inout) :: this
493  type(CellDefnType), intent(inout) :: defn
494  ! local
495  integer(I4B) :: m, n, nfaces
496  real(DP) :: q
497 
498  nfaces = defn%npolyverts + 3
499  do m = 1, nfaces
500  n = defn%facenbr(m)
501  if (n > 0) then
502  q = this%fmi%gwfflowja(this%fmi%dis%con%ia(defn%icell) + n)
503  defn%faceflow(m) = defn%faceflow(m) + q
504  end if
505  if (defn%faceflow(m) < dzero) defn%inoexitface = 0
506  end do

◆ load_flows()

subroutine methoddisvmodule::load_flows ( class(methoddisvtype), intent(inout)  this,
type(celldefntype), intent(inout)  defn 
)
private

Definition at line 455 of file MethodDisv.f90.

456  ! dummy
457  class(MethodDisvType), intent(inout) :: this
458  type(CellDefnType), intent(inout) :: defn
459  ! local
460  integer(I4B) :: nfaces
461  integer(I4B) :: nslots
462 
463  ! expand faceflow array if needed
464  nfaces = defn%npolyverts + 3
465  nslots = size(defn%faceflow)
466  if (nslots < nfaces) call expandarray(defn%faceflow, nfaces - nslots)
467 
468  ! Load face flows, including boundary flows. As with cell verts,
469  ! the face flow array wraps around. Top and bottom flows make up
470  ! the last two elements, respectively, for size npolyverts + 3.
471  ! If there is no flow through any face, set a no-exit-face flag.
472  defn%faceflow = dzero
473  defn%inoexitface = 1
474  call this%load_boundary_flows_to_defn_poly(defn)
475  call this%load_face_flows_to_defn_poly(defn)
476 
477  ! Add up net distributed flow
478  defn%distflow = this%fmi%SourceFlows(defn%icell) + &
479  this%fmi%SinkFlows(defn%icell) + &
480  this%fmi%StorageFlows(defn%icell)
481 
482  ! Set weak sink flag
483  if (this%fmi%SinkFlows(defn%icell) .ne. dzero) then
484  defn%iweaksink = 1
485  else
486  defn%iweaksink = 0
487  end if

◆ load_indicators()

subroutine methoddisvmodule::load_indicators ( class(methoddisvtype), intent(inout)  this,
type(celldefntype), intent(inout), pointer  defn 
)
private

Definition at line 637 of file MethodDisv.f90.

638  ! dummy
639  class(MethodDisvType), intent(inout) :: this
640  type(CellDefnType), pointer, intent(inout) :: defn
641  ! local
642  integer(I4B) :: npolyverts
643  integer(I4B) :: m
644  integer(I4B) :: m0
645  integer(I4B) :: m1
646  integer(I4B) :: m2
647  integer(I4B) :: ic
648  integer(I4B) :: num90
649  integer(I4B) :: num180
650  real(DP) :: x0
651  real(DP) :: y0
652  real(DP) :: x1
653  real(DP) :: y1
654  real(DP) :: x2
655  real(DP) :: y2
656  real(DP) :: epsang
657  real(DP) :: s0x
658  real(DP) :: s0y
659  real(DP) :: &
660  s0mag, s2x, s2y, s2mag, sinang, cosang, dotprod
661  logical(LGP) last180
662 
663  ic = defn%icell
664  npolyverts = defn%npolyverts
665 
666  if (size(defn%ispv180) < npolyverts + 3) &
667  call expandarray(defn%ispv180, npolyverts + 1)
668 
669  defn%ispv180(1:npolyverts + 1) = .false.
670  defn%can_be_rect = .false.
671  defn%can_be_quad = .false.
672  epsang = 1d-5
673  num90 = 0
674  num180 = 0
675  last180 = .false.
676 
677  ! assumes non-self-intersecting polygon
678  do m = 1, npolyverts
679  m1 = m
680  if (m1 .eq. 1) then
681  m0 = npolyverts
682  m2 = 2
683  else if (m1 .eq. npolyverts) then
684  m0 = npolyverts - 1
685  m2 = 1
686  else
687  m0 = m1 - 1
688  m2 = m1 + 1
689  end if
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)
696  s0x = x0 - x1
697  s0y = y0 - y1
698  s0mag = dsqrt(s0x * s0x + s0y * s0y)
699  s2x = x2 - x1
700  s2y = y2 - y1
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
707  num180 = num180 + 1
708  last180 = .true.
709  defn%ispv180(m) = .true.
710  end if
711  else
712  if (dabs(cosang) .lt. epsang) num90 = num90 + 1
713  last180 = .false.
714  end if
715  end do
716 
717  ! List of 180-degree indicators wraps around for convenience
718  defn%ispv180(npolyverts + 1) = defn%ispv180(1)
719 
720  ! Set rect/quad flags
721  if (num90 .eq. 4) then
722  if (num180 .eq. 0) then
723  defn%can_be_rect = .true.
724  else
725  defn%can_be_quad = .true.
726  end if
727  end if

◆ load_neighbors()

subroutine methoddisvmodule::load_neighbors ( class(methoddisvtype), intent(inout)  this,
type(celldefntype), intent(inout), pointer  defn 
)
private

Definition at line 381 of file MethodDisv.f90.

382  ! dummy
383  class(MethodDisvType), intent(inout) :: this
384  type(CellDefnType), pointer, intent(inout) :: defn
385  ! local
386  integer(I4B) :: ic1
387  integer(I4B) :: ic2
388  integer(I4B) :: icu1
389  integer(I4B) :: icu2
390  integer(I4B) :: j1
391  integer(I4B) :: j2
392  integer(I4B) :: k1
393  integer(I4B) :: k2
394  integer(I4B) :: iloc
395  integer(I4B) :: ipos
396  integer(I4B) :: istart1
397  integer(I4B) :: istart2
398  integer(I4B) :: istop1
399  integer(I4B) :: istop2
400  integer(I4B) :: isharedface
401  integer(I4B) :: ncpl
402  integer(I4B) :: nfaces
403  integer(I4B) :: nslots
404 
405  ! expand facenbr array if needed
406  nfaces = defn%npolyverts + 3
407  nslots = size(defn%facenbr)
408  if (nslots < nfaces) call expandarray(defn%facenbr, nfaces - nslots)
409 
410  select type (dis => this%fmi%dis)
411  type is (disvtype)
412  ! Load face neighbors.
413  defn%facenbr = 0
414  ic1 = defn%icell
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
422  ! mask could become relevant if PRT uses interface model
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
429  call shared_face(dis%javert(istart1:istop1), &
430  dis%javert(istart2:istop2), &
431  isharedface)
432  if (isharedface /= 0) then
433  ! Edge (polygon) face neighbor
434  defn%facenbr(isharedface) = int(iloc, 1)
435  else
436  if (k2 > k1) then
437  ! Bottom face neighbor
438  defn%facenbr(defn%npolyverts + 2) = int(iloc, 1)
439  else if (k2 < k1) then
440  ! Top face neighbor
441  defn%facenbr(defn%npolyverts + 3) = int(iloc, 1)
442  else
443  call pstop(1, "k2 should be <> k1, since no shared edge face")
444  end if
445  end if
446  end do
447  end select
448  ! List of edge (polygon) faces wraps around
449  defn%facenbr(defn%npolyverts + 1) = defn%facenbr(1)
Here is the call graph for this function:

◆ load_particle()

subroutine methoddisvmodule::load_particle ( class(methoddisvtype), intent(inout)  this,
type(cellpolytype), intent(inout), pointer  cell,
type(particletype), intent(inout), pointer  particle 
)

Definition at line 133 of file MethodDisv.f90.

134  ! modules
135  use disvmodule, only: disvtype
136  ! dummy
137  class(MethodDisvType), intent(inout) :: this
138  type(CellPolyType), pointer, intent(inout) :: cell
139  type(ParticleType), pointer, intent(inout) :: particle
140  ! local
141  integer(I4B) :: inface
142  integer(I4B) :: ipos
143  integer(I4B) :: ic
144  integer(I4B) :: icu
145  integer(I4B) :: inbr
146  integer(I4B) :: idiag
147  integer(I4B) :: icpl
148  integer(I4B) :: ilay
149  real(DP) :: z
150 
151  select type (dis => this%fmi%dis)
152  type is (disvtype)
153  inface = particle%iboundary(2)
154  idiag = dis%con%ia(cell%defn%icell)
155  inbr = cell%defn%facenbr(inface)
156  ipos = idiag + inbr
157  ic = dis%con%ja(ipos)
158  icu = dis%get_nodeuser(ic)
159  call get_jk(icu, dis%ncpl, dis%nlay, icpl, ilay)
160 
161  ! if returning to a cell through the bottom
162  ! face after previously leaving it through
163  ! that same face, we've entered a cycle
164  ! as can occur e.g. in wells. terminate
165  ! in the previous cell.
166  if (ic == particle%icp .and. inface == 7 .and. ilay < particle%ilay) then
167  particle%advancing = .false.
168  particle%idomain(2) = particle%icp
169  particle%istatus = 2
170  particle%izone = particle%izp
171  call this%save(particle, reason=3)
172  return
173  else
174  particle%icp = particle%idomain(2)
175  particle%izp = particle%izone
176  end if
177 
178  particle%idomain(2) = ic
179  particle%icu = icu
180  particle%ilay = ilay
181 
182  z = particle%z
183  call this%map_neighbor(cell%defn, inface, z)
184 
185  particle%iboundary(2) = inface
186  particle%idomain(3:) = 0
187  particle%iboundary(3:) = 0
188  particle%z = z
189  end select
190 
Vertex grid discretization.
Definition: Disv.f90:24
Here is the call graph for this function:

◆ load_polygon()

subroutine methoddisvmodule::load_polygon ( class(methoddisvtype), intent(inout)  this,
type(celldefntype), intent(inout), pointer  defn 
)
private

Definition at line 367 of file MethodDisv.f90.

368  ! dummy
369  class(MethodDisvType), intent(inout) :: this
370  type(CellDefnType), pointer, intent(inout) :: defn
371 
372  call this%fmi%dis%get_polyverts( &
373  defn%icell, &
374  defn%polyvert, &
375  closed=.true.)
376  defn%npolyverts = size(defn%polyvert, dim=2) - 1

◆ load_properties()

subroutine methoddisvmodule::load_properties ( class(methoddisvtype), intent(inout)  this,
integer(i4b), intent(in)  ic,
type(celldefntype), intent(inout), pointer  defn 
)
private

Definition at line 335 of file MethodDisv.f90.

336  ! dummy
337  class(MethodDisvType), intent(inout) :: this
338  integer(I4B), intent(in) :: ic
339  type(CellDefnType), pointer, intent(inout) :: defn
340  ! local
341  real(DP) :: top
342  real(DP) :: bot
343  real(DP) :: sat
344  integer(I4B) :: icu, icpl, ilay
345 
346  defn%icell = ic
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)
353  defn%top = top
354  defn%bot = bot
355  defn%sat = sat
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)
360  type is (disvtype)
361  icu = dis%get_nodeuser(ic)
362  call get_jk(icu, dis%ncpl, dis%nlay, icpl, ilay)
363  defn%ilay = ilay
364  end select
Here is the call graph for this function:

◆ map_neighbor()

subroutine methoddisvmodule::map_neighbor ( class(methoddisvtype), intent(inout)  this,
type(celldefntype), intent(inout), pointer  defn,
integer(i4b), intent(inout)  inface,
double precision, intent(inout)  z 
)
private

Definition at line 246 of file MethodDisv.f90.

247  ! dummy
248  class(MethodDisvType), intent(inout) :: this
249  type(CellDefnType), pointer, intent(inout) :: defn
250  integer(I4B), intent(inout) :: inface
251  double precision, intent(inout) :: z
252  ! local
253  integer(I4B) :: icin
254  integer(I4B) :: npolyvertsin
255  integer(I4B) :: ic
256  integer(I4B) :: npolyverts
257  integer(I4B) :: inbr
258  integer(I4B) :: inbrnbr
259  integer(I4B) :: j
260  integer(I4B) :: m
261  real(DP) :: zrel
262  real(DP) :: topfrom
263  real(DP) :: botfrom
264  real(DP) :: top
265  real(DP) :: bot
266  real(DP) :: sat
267 
268  ! Map to shared cell face of neighbor
269  inbr = defn%facenbr(inface)
270  if (inbr .eq. 0) then
271  ! Exterior face; no neighbor to map to
272  ! todo AMP: reconsider when multiple models allowed
273  inface = -1
274  else
275  ! Load definition for neighbor cell (neighbor with shared face)
276  icin = defn%icell
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)
280 
281  npolyvertsin = defn%npolyverts
282  npolyverts = this%neighbor%npolyverts
283  if (inface .eq. npolyvertsin + 2) then
284  ! Exits through bot, enters through top
285  inface = npolyverts + 3
286  else if (inface .eq. npolyvertsin + 3) then
287  ! Exits through top, enters through bot
288  inface = npolyverts + 2
289  else
290  ! Exits and enters through shared polygon face
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
295  inface = m
296  exit
297  end if
298  end do
299  ! Map z between cells
300  topfrom = defn%top
301  botfrom = defn%bot
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)
307  end if
308  end if

◆ pass_disv()

subroutine methoddisvmodule::pass_disv ( class(methoddisvtype), intent(inout)  this,
type(particletype), intent(inout), pointer  particle 
)
private

Definition at line 217 of file MethodDisv.f90.

218  ! dummy
219  class(MethodDisvType), intent(inout) :: this
220  type(ParticleType), pointer, intent(inout) :: particle
221  ! local
222  type(CellPolyType), pointer :: cell
223 
224  select type (c => this%cell)
225  type is (cellpolytype)
226  cell => c
227  ! If the entry face has no neighbors it's a
228  ! boundary face, so terminate the particle.
229  ! todo AMP: reconsider when multiple models supported
230  if (cell%defn%facenbr(particle%iboundary(2)) .eq. 0) then
231  particle%istatus = 2
232  particle%advancing = .false.
233  call this%save(particle, reason=3)
234  else
235  ! Update old to new cell properties
236  call this%load_particle(cell, particle)
237  if (.not. particle%advancing) return
238 
239  ! Update intercell mass flows
240  call this%update_flowja(cell, particle)
241  end if
242  end select

◆ update_flowja()

subroutine methoddisvmodule::update_flowja ( class(methoddisvtype), intent(inout)  this,
type(cellpolytype), intent(inout), pointer  cell,
type(particletype), intent(inout), pointer  particle 
)

Definition at line 193 of file MethodDisv.f90.

194  ! dummy
195  class(MethodDisvType), intent(inout) :: this
196  type(CellPolyType), pointer, intent(inout) :: cell
197  type(ParticleType), pointer, intent(inout) :: particle
198  ! local
199  integer(I4B) :: inbr
200  integer(I4B) :: idiag
201  integer(I4B) :: ipos
202 
203  idiag = this%fmi%dis%con%ia(cell%defn%icell)
204  inbr = cell%defn%facenbr(particle%iboundary(2))
205  ipos = idiag + inbr
206 
207  ! leaving old cell
208  this%flowja(ipos) = this%flowja(ipos) - done
209 
210  ! entering new cell
211  this%flowja(this%fmi%dis%con%isym(ipos)) &
212  = this%flowja(this%fmi%dis%con%isym(ipos)) + done
213