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

300  class(MethodDisvType), intent(inout) :: this
301  type(ParticleType), pointer, intent(inout) :: particle
302  real(DP), intent(in) :: tmax
303  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%type)
56  call create_cell_poly(cell)
57  method%cell => cell
58  method%type = "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%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 597 of file MethodDisv.f90.

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

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

504  ! -- dummy
505  class(MethodDisvType), intent(inout) :: this
506  type(CellDefnType), intent(inout) :: defn
507  ! -- local
508  integer(I4B) :: ioffset
509 
510  ! assignment of BoundaryFlows to faceflow below assumes clockwise
511  ! ordering of faces, with face 1 being the "western" face
512  ioffset = (defn%icell - 1) * 10
513  defn%faceflow(1) = defn%faceflow(1) + &
514  this%fmi%BoundaryFlows(ioffset + 4)
515  defn%faceflow(2) = defn%faceflow(2) + &
516  this%fmi%BoundaryFlows(ioffset + 2)
517  defn%faceflow(3) = defn%faceflow(3) + &
518  this%fmi%BoundaryFlows(ioffset + 3)
519  defn%faceflow(4) = defn%faceflow(4) + &
520  this%fmi%BoundaryFlows(ioffset + 1)
521  defn%faceflow(5) = defn%faceflow(1)
522  defn%faceflow(6) = defn%faceflow(6) + &
523  this%fmi%BoundaryFlows(ioffset + 9)
524  defn%faceflow(7) = defn%faceflow(7) + &
525  this%fmi%BoundaryFlows(ioffset + 10)
526 

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

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

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

308  ! -- dummy
309  class(MethodDisvType), intent(inout) :: this
310  integer(I4B), intent(in) :: ic
311  type(CellDefnType), pointer, intent(inout) :: defn
312 
313  ! -- Load basic cell properties
314  call this%load_properties(ic, defn)
315 
316  ! -- Load polygon vertices
317  call this%load_polygon(defn)
318 
319  ! -- Load face neighbors
320  call this%load_neighbors(defn)
321 
322  ! -- Load 180-degree indicator
323  call this%load_indicators(defn)
324 
325  ! -- Load flows (assumes face neighbors already loaded)
326  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  ! load cell definition
89  ic = particle%idomain(next_level)
90  call this%load_cell_defn(ic, cell%defn)
91  if (this%fmi%ibdgwfsat0(ic) == 0) then
92  ! -- Cell is active but dry, so select and initialize pass-to-bottom
93  ! -- cell method and set cell method pointer
94  call method_cell_ptb%init( &
95  fmi=this%fmi, &
96  cell=this%cell, &
97  trackctl=this%trackctl, &
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  fmi=this%fmi, &
105  cell=this%cell, &
106  trackctl=this%trackctl, &
107  tracktimes=this%tracktimes)
108  submethod => method_cell_tern
109  else if (cell%defn%can_be_rect) then
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  call cell_poly_to_quad(cell, quad)
120  base => quad
121  call method_cell_quad%init( &
122  fmi=this%fmi, &
123  cell=base, &
124  trackctl=this%trackctl, &
125  tracktimes=this%tracktimes)
126  submethod => method_cell_quad
127  else
128  call method_cell_tern%init( &
129  fmi=this%fmi, &
130  cell=this%cell, &
131  trackctl=this%trackctl, &
132  tracktimes=this%tracktimes)
133  submethod => method_cell_tern
134  end if
135  end if
136  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 482 of file MethodDisv.f90.

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

◆ load_flows()

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

Definition at line 446 of file MethodDisv.f90.

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

◆ load_indicators()

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

Definition at line 629 of file MethodDisv.f90.

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

◆ load_neighbors()

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

Definition at line 371 of file MethodDisv.f90.

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

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

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

◆ 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 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  ! local
336  real(DP) :: top
337  real(DP) :: bot
338  real(DP) :: sat
339 
340  defn%icell = ic
341  defn%iatop = get_iatop(this%fmi%dis%get_ncpl(), &
342  this%fmi%dis%get_nodeuser(ic))
343  top = this%fmi%dis%top(ic)
344  bot = this%fmi%dis%bot(ic)
345  sat = this%fmi%gwfsat(ic)
346  top = bot + sat * (top - bot)
347  defn%top = top
348  defn%bot = bot
349  defn%sat = sat
350  defn%porosity = this%porosity(ic)
351  defn%retfactor = this%retfactor(ic)
352  defn%izone = this%izone(ic)
353 
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 233 of file MethodDisv.f90.

234  ! dummy
235  class(MethodDisvType), intent(inout) :: this
236  type(CellDefnType), pointer, intent(inout) :: defn
237  integer(I4B), intent(inout) :: inface
238  double precision, intent(inout) :: z
239  ! local
240  integer(I4B) :: icin
241  integer(I4B) :: npolyvertsin
242  integer(I4B) :: ic
243  integer(I4B) :: npolyverts
244  integer(I4B) :: inbr
245  integer(I4B) :: inbrnbr
246  integer(I4B) :: j
247  integer(I4B) :: m
248  real(DP) :: zrel
249  real(DP) :: topfrom
250  real(DP) :: botfrom
251  real(DP) :: top
252  real(DP) :: bot
253  real(DP) :: sat
254 
255  ! -- Map to shared cell face of neighbor
256  inbr = defn%facenbr(inface)
257  if (inbr .eq. 0) then
258  ! -- Exterior face; no neighbor to map to
259  ! todo AMP: reconsider when multiple models allowed
260  inface = -1
261  else
262  ! -- Load definition for neighbor cell (neighbor with shared face)
263  icin = defn%icell
264  j = this%fmi%dis%con%ia(icin)
265  ic = this%fmi%dis%con%ja(j + inbr)
266  call this%load_cell_defn(ic, this%neighbor)
267 
268  npolyvertsin = defn%npolyverts
269  npolyverts = this%neighbor%npolyverts
270  if (inface .eq. npolyvertsin + 2) then
271  ! -- Exits through bot, enters through top
272  inface = npolyverts + 3
273  else if (inface .eq. npolyvertsin + 3) then
274  ! -- Exits through top, enters through bot
275  inface = npolyverts + 2
276  else
277  ! -- Exits and enters through shared polygon face
278  j = this%fmi%dis%con%ia(ic)
279  do m = 1, npolyverts + 3
280  inbrnbr = this%neighbor%facenbr(m)
281  if (this%fmi%dis%con%ja(j + inbrnbr) .eq. icin) then
282  inface = m
283  exit
284  end if
285  end do
286  ! -- Map z between cells
287  topfrom = defn%top
288  botfrom = defn%bot
289  zrel = (z - botfrom) / (topfrom - botfrom)
290  top = this%fmi%dis%top(ic)
291  bot = this%fmi%dis%bot(ic)
292  sat = this%fmi%gwfsat(ic)
293  z = bot + zrel * sat * (top - bot)
294  end if
295  end if

◆ pass_disv()

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

Definition at line 206 of file MethodDisv.f90.

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

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