MODFLOW 6  version 6.6.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-grid method. 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 295 of file MethodDisv.f90.

296  class(MethodDisvType), intent(inout) :: this
297  type(ParticleType), pointer, intent(inout) :: particle
298  real(DP), intent(in) :: tmax
299  call this%track(particle, 1, tmax)

◆ create_method_disv()

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

Definition at line 49 of file MethodDisv.f90.

50  ! -- dummy
51  type(MethodDisvType), pointer :: method
52  ! -- local
53  type(CellPolyType), pointer :: cell
54 
55  allocate (method)
56  allocate (method%type)
57  call create_cell_poly(cell)
58  method%cell => cell
59  method%type = "disv"
60  method%delegates = .true.
61  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 65 of file MethodDisv.f90.

66  class(MethodDisvType), intent(inout) :: this
67  deallocate (this%type)

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

594  ! -- dummy
595  class(MethodDisvType), intent(inout) :: this
596  type(CellDefnType), intent(inout) :: defn
597  ! -- local
598  integer(I4B) :: ic
599  integer(I4B) :: npolyverts
600  integer(I4B) :: ioffset
601  integer(I4B) :: iv
602 
603  ic = defn%icell
604  npolyverts = defn%npolyverts
605 
606  ioffset = (ic - 1) * max_poly_cells
607  do iv = 1, npolyverts
608  defn%faceflow(iv) = &
609  defn%faceflow(iv) + &
610  this%fmi%BoundaryFlows(ioffset + iv)
611  end do
612  defn%faceflow(npolyverts + 1) = defn%faceflow(1)
613  defn%faceflow(npolyverts + 2) = &
614  defn%faceflow(npolyverts + 2) + &
615  this%fmi%BoundaryFlows(ioffset + max_poly_cells - 1)
616  defn%faceflow(npolyverts + 3) = &
617  defn%faceflow(npolyverts + 3) + &
618  this%fmi%BoundaryFlows(ioffset + max_poly_cells)
619 

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

500  ! -- dummy
501  class(MethodDisvType), intent(inout) :: this
502  type(CellDefnType), intent(inout) :: defn
503  ! -- local
504  integer(I4B) :: ioffset
505 
506  ! assignment of BoundaryFlows to faceflow below assumes clockwise
507  ! ordering of faces, with face 1 being the "western" face
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)
522 

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

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

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

304  ! -- dummy
305  class(MethodDisvType), intent(inout) :: this
306  integer(I4B), intent(in) :: ic
307  type(CellDefnType), pointer, intent(inout) :: defn
308 
309  ! -- Load basic cell properties
310  call this%load_properties(ic, defn)
311 
312  ! -- Load polygon vertices
313  call this%load_polygon(defn)
314 
315  ! -- Load face neighbors
316  call this%load_neighbors(defn)
317 
318  ! -- Load 180-degree indicator
319  call this%load_indicators(defn)
320 
321  ! -- Load flows (assumes face neighbors already loaded)
322  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 71 of file MethodDisv.f90.

72  use cellmodule
73  use cellrectmodule
75  use cellutilmodule
76  ! -- dummy
77  class(MethodDisvType), intent(inout) :: this
78  type(ParticleType), pointer, intent(inout) :: particle
79  integer(I4B), intent(in) :: next_level
80  class(MethodType), pointer, intent(inout) :: submethod
81  ! -- local
82  integer(I4B) :: ic
83  class(CellType), pointer :: base
84  type(CellRectType), pointer :: rect
85  type(CellRectQuadType), pointer :: quad
86 
87  select type (cell => this%cell)
88  type is (cellpolytype)
89  ! load cell definition
90  ic = particle%idomain(next_level)
91  call this%load_cell_defn(ic, cell%defn)
92  if (this%fmi%ibdgwfsat0(ic) == 0) then
93  ! -- Cell is active but dry, so select and initialize pass-to-bottom
94  ! -- cell method and set cell method pointer
95  call method_cell_ptb%init( &
96  cell=this%cell, &
97  trackfilectl=this%trackfilectl, &
98  tracktimes=this%tracktimes)
99  submethod => method_cell_ptb
100  else
101  ! -- Select and initialize cell method and set cell method pointer
102  if (particle%ifrctrn > 0) then
103  call method_cell_tern%init( &
104  cell=this%cell, &
105  trackfilectl=this%trackfilectl, &
106  tracktimes=this%tracktimes)
107  submethod => method_cell_tern
108  else if (cell%defn%can_be_rect) then
109  call cell_poly_to_rect(cell, rect)
110  base => rect
111  call method_cell_plck%init( &
112  cell=base, &
113  trackfilectl=this%trackfilectl, &
114  tracktimes=this%tracktimes)
115  submethod => method_cell_plck
116  else if (cell%defn%can_be_quad) then
117  call cell_poly_to_quad(cell, quad)
118  base => quad
119  call method_cell_quad%init( &
120  cell=base, &
121  trackfilectl=this%trackfilectl, &
122  tracktimes=this%tracktimes)
123  submethod => method_cell_quad
124  else
125  call method_cell_tern%init( &
126  cell=this%cell, &
127  trackfilectl=this%trackfilectl, &
128  tracktimes=this%tracktimes)
129  submethod => method_cell_tern
130  end if
131  end if
132  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 478 of file MethodDisv.f90.

479  ! dummy
480  class(MethodDisvType), intent(inout) :: this
481  type(CellDefnType), intent(inout) :: defn
482  ! local
483  integer(I4B) :: m, n, nfaces
484  real(DP) :: q
485 
486  nfaces = defn%npolyverts + 3
487  do m = 1, nfaces
488  n = defn%facenbr(m)
489  if (n > 0) then
490  q = this%fmi%gwfflowja(this%fmi%dis%con%ia(defn%icell) + n)
491  defn%faceflow(m) = defn%faceflow(m) + q
492  end if
493  if (defn%faceflow(m) < dzero) defn%inoexitface = 0
494  end do

◆ load_flows()

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

Definition at line 442 of file MethodDisv.f90.

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

◆ load_indicators()

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

Definition at line 625 of file MethodDisv.f90.

626  ! -- dummy
627  class(MethodDisvType), intent(inout) :: this
628  type(CellDefnType), pointer, intent(inout) :: defn
629  ! -- local
630  integer(I4B) :: npolyverts
631  integer(I4B) :: m
632  integer(I4B) :: m0
633  integer(I4B) :: m1
634  integer(I4B) :: m2
635  integer(I4B) :: ic
636  integer(I4B) :: num90
637  integer(I4B) :: num180
638  real(DP) :: x0
639  real(DP) :: y0
640  real(DP) :: x1
641  real(DP) :: y1
642  real(DP) :: x2
643  real(DP) :: y2
644  real(DP) :: epsang
645  real(DP) :: s0x
646  real(DP) :: s0y
647  real(DP) :: &
648  s0mag, s2x, s2y, s2mag, sinang, cosang, dotprod
649  logical(LGP) last180
650 
651  ic = defn%icell
652  npolyverts = defn%npolyverts
653 
654  ! -- expand ispv180 array if needed
655  if (size(defn%ispv180) < npolyverts + 3) &
656  call expandarray(defn%ispv180, npolyverts + 1)
657 
658  ! -- Load 180-degree indicator.
659  ! -- Also, set flags that indicate how cell can be represented.
660  defn%ispv180(1:npolyverts + 1) = .false.
661  defn%can_be_rect = .false.
662  defn%can_be_quad = .false.
663  epsang = 1d-5
664  num90 = 0
665  num180 = 0
666  last180 = .false.
667  ! assumes non-self-intersecting polygon
668  do m = 1, npolyverts
669  m1 = m
670  if (m1 .eq. 1) then
671  m0 = npolyverts
672  m2 = 2
673  else if (m1 .eq. npolyverts) then
674  m0 = npolyverts - 1
675  m2 = 1
676  else
677  m0 = m1 - 1
678  m2 = m1 + 1
679  end if
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)
686  s0x = x0 - x1
687  s0y = y0 - y1
688  s0mag = dsqrt(s0x * s0x + s0y * s0y)
689  s2x = x2 - x1
690  s2y = y2 - y1
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
697  num180 = num180 + 1
698  last180 = .true.
699  defn%ispv180(m) = .true.
700  end if
701  else
702  if (dabs(cosang) .lt. epsang) num90 = num90 + 1
703  last180 = .false.
704  end if
705  end do
706 
707  ! -- List of 180-degree indicators wraps around for convenience
708  defn%ispv180(npolyverts + 1) = defn%ispv180(1)
709 
710  if (num90 .eq. 4) then
711  if (num180 .eq. 0) then
712  defn%can_be_rect = .true.
713  else
714  defn%can_be_quad = .true.
715  end if
716  end if
717 

◆ load_neighbors()

subroutine methoddisvmodule::load_neighbors ( 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  ! -- local
372  integer(I4B) :: ic1
373  integer(I4B) :: ic2
374  integer(I4B) :: icu1
375  integer(I4B) :: icu2
376  integer(I4B) :: j1
377  integer(I4B) :: j2
378  integer(I4B) :: k1
379  integer(I4B) :: k2
380  integer(I4B) :: iloc
381  integer(I4B) :: ipos
382  integer(I4B) :: istart1
383  integer(I4B) :: istart2
384  integer(I4B) :: istop1
385  integer(I4B) :: istop2
386  integer(I4B) :: isharedface
387  integer(I4B) :: ncpl
388  integer(I4B) :: nfaces
389  integer(I4B) :: nslots
390 
391  ! -- expand facenbr array if needed
392  nfaces = defn%npolyverts + 3
393  nslots = size(defn%facenbr)
394  if (nslots < nfaces) call expandarray(defn%facenbr, nfaces - nslots)
395 
396  select type (dis => this%fmi%dis)
397  type is (disvtype)
398  ! -- Load face neighbors.
399  defn%facenbr = 0
400  ic1 = defn%icell
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
408  ! mask could become relevant if PRT uses interface model
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
415  call shared_face(dis%javert(istart1:istop1), &
416  dis%javert(istart2:istop2), &
417  isharedface)
418  if (isharedface /= 0) then
419  ! -- Edge (polygon) face neighbor
420  defn%facenbr(isharedface) = int(iloc, 1)
421  else
422  if (k2 > k1) then
423  ! -- Bottom face neighbor
424  defn%facenbr(defn%npolyverts + 2) = int(iloc, 1)
425  else if (k2 < k1) then
426  ! -- Top face neighbor
427  defn%facenbr(defn%npolyverts + 3) = int(iloc, 1)
428  else
429  call pstop(1, "k2 should be <> k1, since no shared edge face")
430  end if
431  end if
432  end do
433  end select
434  ! -- List of edge (polygon) faces wraps around
435  defn%facenbr(defn%npolyverts + 1) = defn%facenbr(1)
436 
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 135 of file MethodDisv.f90.

136  ! modules
137  use disvmodule, only: disvtype
138  ! dummy
139  class(MethodDisvType), intent(inout) :: this
140  type(CellPolyType), pointer, intent(inout) :: cell
141  type(ParticleType), pointer, intent(inout) :: particle
142  ! local
143  integer(I4B) :: inface
144  integer(I4B) :: ipos
145  integer(I4B) :: ic
146  integer(I4B) :: icu
147  integer(I4B) :: inbr
148  integer(I4B) :: idiag
149  integer(I4B) :: icpl
150  integer(I4B) :: ilay
151  real(DP) :: z
152 
153  select type (dis => this%fmi%dis)
154  type is (disvtype)
155  ! compute and set reduced/user node numbers and layer
156  inface = particle%iboundary(2)
157  idiag = dis%con%ia(cell%defn%icell)
158  inbr = cell%defn%facenbr(inface)
159  ipos = idiag + inbr
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
164  particle%icu = icu
165  particle%ilay = ilay
166 
167  ! map/set particle entry face and z coordinate
168  z = particle%z
169  call this%map_neighbor(cell%defn, inface, z)
170  particle%iboundary(2) = inface
171  particle%idomain(3:) = 0
172  particle%iboundary(3:) = 0
173  particle%z = z
174  end select
175 
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 352 of file MethodDisv.f90.

353  ! dummy
354  class(MethodDisvType), intent(inout) :: this
355  type(CellDefnType), pointer, intent(inout) :: defn
356 
357  call this%fmi%dis%get_polyverts( &
358  defn%icell, &
359  defn%polyvert, &
360  closed=.true.)
361  defn%npolyverts = size(defn%polyvert, dim=2) - 1
362 

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

327  ! dummy
328  class(MethodDisvType), intent(inout) :: this
329  integer(I4B), intent(in) :: ic
330  type(CellDefnType), pointer, intent(inout) :: defn
331  ! local
332  real(DP) :: top
333  real(DP) :: bot
334  real(DP) :: sat
335 
336  defn%icell = 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)
343  defn%top = top
344  defn%bot = bot
345  defn%sat = sat
346  defn%porosity = this%porosity(ic)
347  defn%retfactor = this%retfactor(ic)
348  defn%izone = this%izone(ic)
349 
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 229 of file MethodDisv.f90.

230  ! dummy
231  class(MethodDisvType), intent(inout) :: this
232  type(CellDefnType), pointer, intent(inout) :: defn
233  integer(I4B), intent(inout) :: inface
234  double precision, intent(inout) :: z
235  ! local
236  integer(I4B) :: icin
237  integer(I4B) :: npolyvertsin
238  integer(I4B) :: ic
239  integer(I4B) :: npolyverts
240  integer(I4B) :: inbr
241  integer(I4B) :: inbrnbr
242  integer(I4B) :: j
243  integer(I4B) :: m
244  real(DP) :: zrel
245  real(DP) :: topfrom
246  real(DP) :: botfrom
247  real(DP) :: top
248  real(DP) :: bot
249  real(DP) :: sat
250 
251  ! -- Map to shared cell face of neighbor
252  inbr = defn%facenbr(inface)
253  if (inbr .eq. 0) then
254  ! -- Exterior face; no neighbor to map to
255  ! todo AMP: reconsider when multiple models allowed
256  inface = -1
257  else
258  ! -- Load definition for neighbor cell (neighbor with shared face)
259  icin = defn%icell
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)
263 
264  npolyvertsin = defn%npolyverts
265  npolyverts = this%neighbor%npolyverts
266  if (inface .eq. npolyvertsin + 2) then
267  ! -- Exits through bot, enters through top
268  inface = npolyverts + 3
269  else if (inface .eq. npolyvertsin + 3) then
270  ! -- Exits through top, enters through bot
271  inface = npolyverts + 2
272  else
273  ! -- Exits and enters through shared polygon face
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
278  inface = m
279  exit
280  end if
281  end do
282  ! -- Map z between cells
283  topfrom = defn%top
284  botfrom = defn%bot
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)
290  end if
291  end if

◆ pass_disv()

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

Definition at line 202 of file MethodDisv.f90.

203  ! dummy
204  class(MethodDisvType), intent(inout) :: this
205  type(ParticleType), pointer, intent(inout) :: particle
206  ! local
207  type(CellPolyType), pointer :: cell
208 
209  select type (c => this%cell)
210  type is (cellpolytype)
211  cell => c
212  ! If the entry face has no neighbors it's a
213  ! boundary face, so terminate the particle.
214  ! todo AMP: reconsider when multiple models supported
215  if (cell%defn%facenbr(particle%iboundary(2)) .eq. 0) then
216  particle%istatus = 2
217  particle%advancing = .false.
218  call this%save(particle, reason=3)
219  else
220  ! Otherwise, load cell properties into the
221  ! particle and update intercell mass flows
222  call this%load_particle(cell, particle)
223  call this%update_flowja(cell, particle)
224  end if
225  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 178 of file MethodDisv.f90.

179  ! dummy
180  class(MethodDisvType), intent(inout) :: this
181  type(CellPolyType), pointer, intent(inout) :: cell
182  type(ParticleType), pointer, intent(inout) :: particle
183  ! local
184  integer(I4B) :: inbr
185  integer(I4B) :: idiag
186  integer(I4B) :: ipos
187 
188  idiag = this%fmi%dis%con%ia(cell%defn%icell)
189  inbr = cell%defn%facenbr(particle%iboundary(2))
190  ipos = idiag + inbr
191 
192  ! leaving old cell
193  this%flowja(ipos) = this%flowja(ipos) - done
194 
195  ! entering new cell
196  this%flowja(this%fmi%dis%con%isym(ipos)) &
197  = this%flowja(this%fmi%dis%con%isym(ipos)) + done
198