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 321 of file MethodDisv.f90.

322  class(MethodDisvType), intent(inout) :: this
323  type(ParticleType), pointer, intent(inout) :: particle
324  real(DP), intent(in) :: tmax
325 
326  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 613 of file MethodDisv.f90.

614  ! dummy
615  class(MethodDisvType), intent(inout) :: this
616  type(CellDefnType), intent(inout) :: defn
617  ! local
618  integer(I4B) :: ic
619  integer(I4B) :: npolyverts
620  integer(I4B) :: ioffset
621  integer(I4B) :: iv
622 
623  ic = defn%icell
624  npolyverts = defn%npolyverts
625 
626  ioffset = (ic - 1) * max_poly_cells
627  do iv = 1, npolyverts
628  defn%faceflow(iv) = &
629  defn%faceflow(iv) + &
630  this%fmi%BoundaryFlows(ioffset + iv)
631  end do
632  defn%faceflow(npolyverts + 1) = defn%faceflow(1)
633  defn%faceflow(npolyverts + 2) = &
634  defn%faceflow(npolyverts + 2) + &
635  this%fmi%BoundaryFlows(ioffset + max_poly_cells - 1)
636  defn%faceflow(npolyverts + 3) = &
637  defn%faceflow(npolyverts + 3) + &
638  this%fmi%BoundaryFlows(ioffset + max_poly_cells)
639 

◆ 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 520 of file MethodDisv.f90.

521  ! dummy
522  class(MethodDisvType), intent(inout) :: this
523  type(CellDefnType), intent(inout) :: defn
524  ! local
525  integer(I4B) :: ioffset
526 
527  ! assignment of BoundaryFlows to faceflow below assumes clockwise
528  ! ordering of faces, with face 1 being the "western" face
529  ioffset = (defn%icell - 1) * 10
530  defn%faceflow(1) = defn%faceflow(1) + &
531  this%fmi%BoundaryFlows(ioffset + 4)
532  defn%faceflow(2) = defn%faceflow(2) + &
533  this%fmi%BoundaryFlows(ioffset + 2)
534  defn%faceflow(3) = defn%faceflow(3) + &
535  this%fmi%BoundaryFlows(ioffset + 3)
536  defn%faceflow(4) = defn%faceflow(4) + &
537  this%fmi%BoundaryFlows(ioffset + 1)
538  defn%faceflow(5) = defn%faceflow(1)
539  defn%faceflow(6) = defn%faceflow(6) + &
540  this%fmi%BoundaryFlows(ioffset + 9)
541  defn%faceflow(7) = defn%faceflow(7) + &
542  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 547 of file MethodDisv.f90.

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

◆ 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 330 of file MethodDisv.f90.

331  ! dummy
332  class(MethodDisvType), intent(inout) :: this
333  integer(I4B), intent(in) :: ic
334  type(CellDefnType), pointer, intent(inout) :: defn
335 
336  call this%load_properties(ic, defn)
337  call this%load_polygon(defn)
338  call this%load_neighbors(defn)
339  call this%load_indicators(defn)
340  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  ! Cell is active but dry, so select and initialize pass-to-bottom
92  ! cell method and set cell method pointer
93  call method_cell_ptb%init( &
94  fmi=this%fmi, &
95  cell=this%cell, &
96  trackctl=this%trackctl, &
97  tracktimes=this%tracktimes)
98  submethod => method_cell_ptb
99  else if (particle%ifrctrn > 0) then
100  ! Force the ternary method
101  call method_cell_tern%init( &
102  fmi=this%fmi, &
103  cell=this%cell, &
104  trackctl=this%trackctl, &
105  tracktimes=this%tracktimes)
106  submethod => method_cell_tern
107  else if (cell%defn%can_be_rect) then
108  ! Cell is a rectangle, convert it to a rectangular cell type and
109  ! initialize Pollock's method
110  call cell_poly_to_rect(cell, rect)
111  base => rect
112  call method_cell_plck%init( &
113  fmi=this%fmi, &
114  cell=base, &
115  trackctl=this%trackctl, &
116  tracktimes=this%tracktimes)
117  submethod => method_cell_plck
118  else if (cell%defn%can_be_quad) then
119  ! Cell is quad-refined, convert to a quad rect cell type and
120  ! initialize the corresponding method
121  call cell_poly_to_quad(cell, quad)
122  base => quad
123  call method_cell_quad%init( &
124  fmi=this%fmi, &
125  cell=base, &
126  trackctl=this%trackctl, &
127  tracktimes=this%tracktimes)
128  submethod => method_cell_quad
129  else
130  ! Default to the ternary method
131  call method_cell_tern%init( &
132  fmi=this%fmi, &
133  cell=this%cell, &
134  trackctl=this%trackctl, &
135  tracktimes=this%tracktimes)
136  submethod => method_cell_tern
137  end if
138  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 499 of file MethodDisv.f90.

500  ! dummy
501  class(MethodDisvType), intent(inout) :: this
502  type(CellDefnType), intent(inout) :: defn
503  ! local
504  integer(I4B) :: m, n, nfaces
505  real(DP) :: q
506 
507  nfaces = defn%npolyverts + 3
508  do m = 1, nfaces
509  n = defn%facenbr(m)
510  if (n > 0) then
511  q = this%fmi%gwfflowja(this%fmi%dis%con%ia(defn%icell) + n)
512  defn%faceflow(m) = defn%faceflow(m) + q
513  end if
514  if (defn%faceflow(m) < dzero) defn%inoexitface = 0
515  end do

◆ load_flows()

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

Definition at line 464 of file MethodDisv.f90.

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

◆ load_indicators()

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

Definition at line 646 of file MethodDisv.f90.

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

◆ load_neighbors()

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

Definition at line 390 of file MethodDisv.f90.

391  ! dummy
392  class(MethodDisvType), intent(inout) :: this
393  type(CellDefnType), pointer, intent(inout) :: defn
394  ! local
395  integer(I4B) :: ic1
396  integer(I4B) :: ic2
397  integer(I4B) :: icu1
398  integer(I4B) :: icu2
399  integer(I4B) :: j1
400  integer(I4B) :: j2
401  integer(I4B) :: k1
402  integer(I4B) :: k2
403  integer(I4B) :: iloc
404  integer(I4B) :: ipos
405  integer(I4B) :: istart1
406  integer(I4B) :: istart2
407  integer(I4B) :: istop1
408  integer(I4B) :: istop2
409  integer(I4B) :: isharedface
410  integer(I4B) :: ncpl
411  integer(I4B) :: nfaces
412  integer(I4B) :: nslots
413 
414  ! expand facenbr array if needed
415  nfaces = defn%npolyverts + 3
416  nslots = size(defn%facenbr)
417  if (nslots < nfaces) call expandarray(defn%facenbr, nfaces - nslots)
418 
419  select type (dis => this%fmi%dis)
420  type is (disvtype)
421  ! Load face neighbors.
422  defn%facenbr = 0
423  ic1 = defn%icell
424  icu1 = dis%get_nodeuser(ic1)
425  ncpl = dis%get_ncpl()
426  call get_jk(icu1, ncpl, dis%nlay, j1, k1)
427  istart1 = dis%iavert(j1)
428  istop1 = dis%iavert(j1 + 1) - 1
429  do iloc = 1, dis%con%ia(ic1 + 1) - dis%con%ia(ic1) - 1
430  ipos = dis%con%ia(ic1) + iloc
431  ! mask could become relevant if PRT uses interface model
432  if (dis%con%mask(ipos) == 0) cycle
433  ic2 = dis%con%ja(ipos)
434  icu2 = dis%get_nodeuser(ic2)
435  call get_jk(icu2, ncpl, dis%nlay, j2, k2)
436  istart2 = dis%iavert(j2)
437  istop2 = dis%iavert(j2 + 1) - 1
438  call shared_face(dis%javert(istart1:istop1), &
439  dis%javert(istart2:istop2), &
440  isharedface)
441  if (isharedface /= 0) then
442  ! Edge (polygon) face neighbor
443  defn%facenbr(isharedface) = int(iloc, 1)
444  else
445  if (k2 > k1) then
446  ! Bottom face neighbor
447  defn%facenbr(defn%npolyverts + 2) = int(iloc, 1)
448  else if (k2 < k1) then
449  ! Top face neighbor
450  defn%facenbr(defn%npolyverts + 3) = int(iloc, 1)
451  else
452  call pstop(1, "k2 should be <> k1, since no shared edge face")
453  end if
454  end if
455  end do
456  end select
457  ! List of edge (polygon) faces wraps around
458  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 141 of file MethodDisv.f90.

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

377  ! dummy
378  class(MethodDisvType), intent(inout) :: this
379  type(CellDefnType), pointer, intent(inout) :: defn
380 
381  call this%fmi%dis%get_polyverts( &
382  defn%icell, &
383  defn%polyvert, &
384  closed=.true.)
385  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 344 of file MethodDisv.f90.

345  ! dummy
346  class(MethodDisvType), intent(inout) :: this
347  integer(I4B), intent(in) :: ic
348  type(CellDefnType), pointer, intent(inout) :: defn
349  ! local
350  real(DP) :: top
351  real(DP) :: bot
352  real(DP) :: sat
353  integer(I4B) :: icu, icpl, ilay
354 
355  defn%icell = ic
356  defn%iatop = get_iatop(this%fmi%dis%get_ncpl(), &
357  this%fmi%dis%get_nodeuser(ic))
358  top = this%fmi%dis%top(ic)
359  bot = this%fmi%dis%bot(ic)
360  sat = this%fmi%gwfsat(ic)
361  top = bot + sat * (top - bot)
362  defn%top = top
363  defn%bot = bot
364  defn%sat = sat
365  defn%porosity = this%porosity(ic)
366  defn%retfactor = this%retfactor(ic)
367  defn%izone = this%izone(ic)
368  select type (dis => this%fmi%dis)
369  type is (disvtype)
370  icu = dis%get_nodeuser(ic)
371  call get_jk(icu, dis%ncpl, dis%nlay, icpl, ilay)
372  defn%ilay = ilay
373  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 255 of file MethodDisv.f90.

256  ! dummy
257  class(MethodDisvType), intent(inout) :: this
258  type(CellDefnType), pointer, intent(inout) :: defn
259  integer(I4B), intent(inout) :: inface
260  double precision, intent(inout) :: z
261  ! local
262  integer(I4B) :: icin
263  integer(I4B) :: npolyvertsin
264  integer(I4B) :: ic
265  integer(I4B) :: npolyverts
266  integer(I4B) :: inbr
267  integer(I4B) :: inbrnbr
268  integer(I4B) :: j
269  integer(I4B) :: m
270  real(DP) :: zrel
271  real(DP) :: topfrom
272  real(DP) :: botfrom
273  real(DP) :: top
274  real(DP) :: bot
275  real(DP) :: sat
276 
277  ! Map to shared cell face of neighbor
278  inbr = defn%facenbr(inface)
279  if (inbr .eq. 0) then
280  ! Exterior face; no neighbor to map to
281  ! todo AMP: reconsider when multiple models allowed
282  inface = -1
283  else
284  ! Load definition for neighbor cell (neighbor with shared face)
285  icin = defn%icell
286  j = this%fmi%dis%con%ia(icin)
287  ic = this%fmi%dis%con%ja(j + inbr)
288  call this%load_cell_defn(ic, this%neighbor)
289 
290  npolyvertsin = defn%npolyverts
291  npolyverts = this%neighbor%npolyverts
292  if (inface .eq. npolyvertsin + 2) then
293  ! Exits through bot, enters through top
294  inface = npolyverts + 3
295  else if (inface .eq. npolyvertsin + 3) then
296  ! Exits through top, enters through bot
297  inface = npolyverts + 2
298  else
299  ! Exits and enters through shared polygon face
300  j = this%fmi%dis%con%ia(ic)
301  do m = 1, npolyverts + 3
302  inbrnbr = this%neighbor%facenbr(m)
303  if (this%fmi%dis%con%ja(j + inbrnbr) .eq. icin) then
304  inface = m
305  exit
306  end if
307  end do
308  ! Map z between cells
309  topfrom = defn%top
310  botfrom = defn%bot
311  zrel = (z - botfrom) / (topfrom - botfrom)
312  top = this%fmi%dis%top(ic)
313  bot = this%fmi%dis%bot(ic)
314  sat = this%fmi%gwfsat(ic)
315  z = bot + zrel * sat * (top - bot)
316  end if
317  end if

◆ pass_disv()

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

Definition at line 225 of file MethodDisv.f90.

226  ! dummy
227  class(MethodDisvType), intent(inout) :: this
228  type(ParticleType), pointer, intent(inout) :: particle
229  ! local
230  type(CellPolyType), pointer :: cell
231 
232  select type (c => this%cell)
233  type is (cellpolytype)
234  cell => c
235  ! If the entry face has no neighbors it's a
236  ! boundary face, so terminate the particle.
237  ! todo AMP: reconsider when multiple models supported
238  if (cell%defn%facenbr(particle%iboundary(2)) .eq. 0) then
239  particle%istatus = 2
240  particle%advancing = .false.
241  call this%save(particle, reason=3)
242  else
243  ! Otherwise, load cell properties into the
244  ! particle. It may be marked to terminate.
245  call this%load_particle(cell, particle)
246  if (.not. particle%advancing) return
247 
248  ! Update intercell mass flows
249  call this%update_flowja(cell, particle)
250  end if
251  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 201 of file MethodDisv.f90.

202  ! dummy
203  class(MethodDisvType), intent(inout) :: this
204  type(CellPolyType), pointer, intent(inout) :: cell
205  type(ParticleType), pointer, intent(inout) :: particle
206  ! local
207  integer(I4B) :: inbr
208  integer(I4B) :: idiag
209  integer(I4B) :: ipos
210 
211  idiag = this%fmi%dis%con%ia(cell%defn%icell)
212  inbr = cell%defn%facenbr(particle%iboundary(2))
213  ipos = idiag + inbr
214 
215  ! leaving old cell
216  this%flowja(ipos) = this%flowja(ipos) - done
217 
218  ! entering new cell
219  this%flowja(this%fmi%dis%con%isym(ipos)) &
220  = this%flowja(this%fmi%dis%con%isym(ipos)) + done
221