MODFLOW 6  version 6.6.0.dev0
USGS Modular Hydrologic Model
MethodDis.f90
Go to the documentation of this file.
2 
3  use kindmodule, only: dp, i4b, lgp
4  use constantsmodule, only: done, dzero
5  use methodmodule, only: methodtype
7  use cellmodule
10  use particlemodule
11  use prtfmimodule, only: prtfmitype
12  use dismodule, only: distype
13  use geomutilmodule, only: get_ijk, get_jk
14  use mathutilmodule, only: is_close
15  implicit none
16 
17  private
18  public :: methoddistype
19  public :: create_method_dis
20 
21  type, extends(methodtype) :: methoddistype
22  private
23  contains
24  procedure, public :: apply => apply_dis !< apply the DIS tracking method
25  procedure, public :: deallocate !< deallocate arrays and scalars
26  procedure, public :: load => load_dis !< load the method
27  procedure, public :: pass => pass_dis !< pass the particle to the next domain
28  procedure :: get_top !< get cell top elevation
29  procedure :: update_flowja !< load intercell mass flows
30  procedure :: load_particle !< load particle properties
31  procedure :: load_properties !< load cell properties
32  procedure :: load_neighbors !< load cell face neighbors
33  procedure :: load_flows !< load cell face flows
34  procedure :: load_boundary_flows_to_defn !< load boundary flows to the cell definition
35  procedure :: load_face_flows_to_defn !< load face flows to the cell definition
36  procedure :: load_celldefn !< load cell definition from the grid
37  procedure :: load_cell !< load cell geometry and flows
38  end type methoddistype
39 
40 contains
41 
42  !> @brief Create a new structured grid (DIS) tracking method
43  subroutine create_method_dis(method)
44  ! dummy
45  type(methoddistype), pointer :: method
46  ! local
47  type(cellrecttype), pointer :: cell
48 
49  allocate (method)
50  allocate (method%name)
51  call create_cell_rect(cell)
52  method%cell => cell
53  method%name = "dis"
54  method%delegates = .true.
55  end subroutine create_method_dis
56 
57  !> @brief Destructor the tracking method
58  subroutine deallocate (this)
59  class(methoddistype), intent(inout) :: this
60  deallocate (this%name)
61  end subroutine deallocate
62 
63  subroutine load_cell(this, ic, cell)
64  ! dummy
65  class(methoddistype), intent(inout) :: this
66  integer(I4B), intent(in) :: ic
67  type(cellrecttype), intent(inout) :: cell
68  ! local
69  integer(I4B) :: icu
70  integer(I4B) :: irow
71  integer(I4B) :: jcol
72  integer(I4B) :: klay
73  real(DP) :: areax
74  real(DP) :: areay
75  real(DP) :: areaz
76  real(DP) :: dx
77  real(DP) :: dy
78  real(DP) :: dz
79  real(DP) :: factor
80  real(DP) :: term
81 
82  select type (dis => this%fmi%dis)
83  type is (distype)
84  icu = dis%get_nodeuser(ic)
85 
86  call get_ijk(icu, dis%nrow, dis%ncol, dis%nlay, &
87  irow, jcol, klay)
88 
89  dx = dis%delr(jcol)
90  dy = dis%delc(irow)
91  dz = cell%defn%top - cell%defn%bot
92 
93  cell%dx = dx
94  cell%dy = dy
95  cell%dz = dz
96  cell%sinrot = dzero
97  cell%cosrot = done
98  cell%xOrigin = cell%defn%polyvert(1, 1)
99  cell%yOrigin = cell%defn%polyvert(2, 1)
100  cell%zOrigin = cell%defn%bot
101  cell%ipvOrigin = 1
102 
103  factor = done / cell%defn%retfactor
104  factor = factor / cell%defn%porosity
105 
106  areaz = dx * dy
107  term = factor / areaz
108 
109  cell%vz1 = cell%defn%faceflow(6) * term
110  cell%vz2 = -cell%defn%faceflow(7) * term
111 
112  if (this%fmi%ibdgwfsat0(ic) == 0) then
113  cell%vx1 = dzero
114  cell%vx2 = dzero
115  cell%vy1 = dzero
116  cell%vy2 = dzero
117  cell%vz1 = dzero
118  cell%vz2 = dzero
119  return
120  end if
121 
122  areax = dy * dz
123  term = factor / areax
124  cell%vx1 = cell%defn%faceflow(1) * term
125  cell%vx2 = -cell%defn%faceflow(3) * term
126 
127  areay = dx * dz
128  term = factor / areay
129  cell%vy1 = cell%defn%faceflow(4) * term
130  cell%vy2 = -cell%defn%faceflow(2) * term
131 
132  end select
133  end subroutine load_cell
134 
135  !> @brief Load the cell geometry and tracking method
136  subroutine load_dis(this, particle, next_level, submethod)
137  ! dummy
138  class(methoddistype), intent(inout) :: this
139  type(particletype), pointer, intent(inout) :: particle
140  integer(I4B), intent(in) :: next_level
141  class(methodtype), pointer, intent(inout) :: submethod
142  ! local
143  integer(I4B) :: ic
144 
145  select type (cell => this%cell)
146  type is (cellrecttype)
147  ic = particle%idomain(next_level)
148  call this%load_celldefn(ic, cell%defn)
149  call this%load_cell(ic, cell)
150  if (this%fmi%ibdgwfsat0(ic) == 0) then
151  call method_cell_ptb%init( &
152  fmi=this%fmi, &
153  cell=this%cell, &
154  trackctl=this%trackctl, &
155  tracktimes=this%tracktimes)
156  submethod => method_cell_ptb
157  else
158  call method_cell_plck%init( &
159  fmi=this%fmi, &
160  cell=this%cell, &
161  trackctl=this%trackctl, &
162  tracktimes=this%tracktimes)
163  submethod => method_cell_plck
164  end if
165  end select
166  end subroutine load_dis
167 
168  !> @brief Load cell properties into the particle, including
169  ! the z coordinate, entry face, and node and layer numbers.
170  subroutine load_particle(this, cell, particle)
171  ! dummy
172  class(methoddistype), intent(inout) :: this
173  type(cellrecttype), pointer, intent(inout) :: cell
174  type(particletype), pointer, intent(inout) :: particle
175  ! local
176  integer(I4B) :: ic
177  integer(I4B) :: icu
178  integer(I4B) :: ilay
179  integer(I4B) :: irow
180  integer(I4B) :: icol
181  integer(I4B) :: inface
182  integer(I4B) :: idiag
183  integer(I4B) :: inbr
184  integer(I4B) :: ipos
185  real(DP) :: z
186  real(DP) :: zrel
187  real(DP) :: topfrom
188  real(DP) :: botfrom
189  real(DP) :: top
190  real(DP) :: bot
191  real(DP) :: sat
192 
193  select type (dis => this%fmi%dis)
194  type is (distype)
195  ! compute reduced/user node numbers and layer
196  inface = particle%iboundary(2)
197  inbr = cell%defn%facenbr(inface)
198  idiag = this%fmi%dis%con%ia(cell%defn%icell)
199  ipos = idiag + inbr
200  ic = dis%con%ja(ipos)
201  icu = dis%get_nodeuser(ic)
202  call get_ijk(icu, dis%nrow, dis%ncol, dis%nlay, &
203  irow, icol, ilay)
204 
205  ! if returning to a cell through the bottom
206  ! face after previously leaving it through
207  ! that same face, we've entered a cycle
208  ! as can occur e.g. in wells. terminate
209  ! in the previous cell.
210  if (ic == particle%icp .and. inface == 7 .and. ilay < particle%ilay) then
211  particle%advancing = .false.
212  particle%idomain(2) = particle%icp
213  particle%istatus = 2
214  particle%izone = particle%izp
215  call this%save(particle, reason=3)
216  return
217  else
218  particle%icp = particle%idomain(2)
219  particle%izp = particle%izone
220  end if
221 
222  ! update node numbers and layer
223  particle%idomain(2) = ic
224  particle%icu = icu
225  particle%ilay = ilay
226 
227  ! map/set particle entry face
228  if (inface .eq. 1) then
229  inface = 3
230  else if (inface .eq. 2) then
231  inface = 4
232  else if (inface .eq. 3) then
233  inface = 1
234  else if (inface .eq. 4) then
235  inface = 2
236  else if (inface .eq. 6) then
237  inface = 7
238  else if (inface .eq. 7) then
239  inface = 6
240  end if
241  particle%iboundary(2) = inface
242 
243  ! map z between cells
244  z = particle%z
245  if (inface < 5) then
246  topfrom = cell%defn%top
247  botfrom = cell%defn%bot
248  zrel = (z - botfrom) / (topfrom - botfrom)
249  top = dis%top(ic)
250  bot = dis%bot(ic)
251  sat = this%fmi%gwfsat(ic)
252  z = bot + zrel * sat * (top - bot)
253  end if
254  particle%z = z
255  end select
256  end subroutine load_particle
257 
258  !> @brief Update cell-cell flows of particle mass.
259  !! Every particle is currently assigned unit mass.
260  subroutine update_flowja(this, cell, particle)
261  ! dummy
262  class(methoddistype), intent(inout) :: this
263  type(cellrecttype), pointer, intent(inout) :: cell
264  type(particletype), pointer, intent(inout) :: particle
265  ! local
266  integer(I4B) :: idiag
267  integer(I4B) :: inbr
268  integer(I4B) :: inface
269  integer(I4B) :: ipos
270 
271  inface = particle%iboundary(2)
272  inbr = cell%defn%facenbr(inface)
273  idiag = this%fmi%dis%con%ia(cell%defn%icell)
274  ipos = idiag + inbr
275 
276  ! leaving old cell
277  this%flowja(ipos) = this%flowja(ipos) - done
278 
279  ! entering new cell
280  this%flowja(this%fmi%dis%con%isym(ipos)) &
281  = this%flowja(this%fmi%dis%con%isym(ipos)) + done
282  end subroutine update_flowja
283 
284  !> @brief Pass a particle to the next cell, if there is one
285  subroutine pass_dis(this, particle)
286  ! dummy
287  class(methoddistype), intent(inout) :: this
288  type(particletype), pointer, intent(inout) :: particle
289  ! local
290  type(cellrecttype), pointer :: cell
291 
292  select type (c => this%cell)
293  type is (cellrecttype)
294  cell => c
295  ! If the entry face has no neighbors it's a
296  ! boundary face, so terminate the particle.
297  ! todo AMP: reconsider when multiple models supported
298  if (cell%defn%facenbr(particle%iboundary(2)) .eq. 0) then
299  particle%istatus = 2
300  particle%advancing = .false.
301  call this%save(particle, reason=3)
302  else
303  ! Update old to new cell properties
304  call this%load_particle(cell, particle)
305  if (.not. particle%advancing) return
306 
307  ! Update intercell mass flows
308  call this%update_flowja(cell, particle)
309  end if
310  end select
311  end subroutine pass_dis
312 
313  !> @brief Apply the structured tracking method to a particle.
314  subroutine apply_dis(this, particle, tmax)
315  class(methoddistype), intent(inout) :: this
316  type(particletype), pointer, intent(inout) :: particle
317  real(DP), intent(in) :: tmax
318 
319  call this%track(particle, 1, tmax)
320  end subroutine apply_dis
321 
322  !> @brief Returns a top elevation based on index iatop
323  function get_top(this, iatop) result(top)
324  class(methoddistype), intent(inout) :: this
325  integer, intent(in) :: iatop
326  double precision :: top
327 
328  if (iatop .lt. 0) then
329  top = this%fmi%dis%top(-iatop)
330  else
331  top = this%fmi%dis%bot(iatop)
332  end if
333  end function get_top
334 
335  !> @brief Loads cell definition from the grid
336  subroutine load_celldefn(this, ic, defn)
337  ! dummy
338  class(methoddistype), intent(inout) :: this
339  integer(I4B), intent(in) :: ic
340  type(celldefntype), pointer, intent(inout) :: defn
341 
342  ! Load basic cell properties
343  call this%load_properties(ic, defn)
344 
345  ! Load cell polygon vertices
346  call this%fmi%dis%get_polyverts( &
347  defn%icell, &
348  defn%polyvert, &
349  closed=.true.)
350  call this%load_neighbors(defn)
351 
352  ! Load 180 degree face indicators
353  defn%ispv180(1:defn%npolyverts + 1) = .false.
354 
355  call this%load_flows(defn)
356 
357  end subroutine load_celldefn
358 
359  subroutine load_properties(this, ic, defn)
360  ! dummy
361  class(methoddistype), intent(inout) :: this
362  integer(I4B), intent(in) :: ic
363  type(celldefntype), pointer, intent(inout) :: defn
364  ! local
365  integer(I4B) :: irow, icol, ilay, icu
366 
367  defn%icell = ic
368  defn%npolyverts = 4 ! rectangular cell always has 4 vertices
369  defn%iatop = get_iatop(this%fmi%dis%get_ncpl(), &
370  this%fmi%dis%get_nodeuser(ic))
371  defn%top = this%fmi%dis%bot(ic) + &
372  this%fmi%gwfsat(ic) * &
373  (this%fmi%dis%top(ic) - this%fmi%dis%bot(ic))
374  defn%bot = this%fmi%dis%bot(ic)
375  defn%sat = this%fmi%gwfsat(ic)
376  defn%porosity = this%porosity(ic)
377  defn%retfactor = this%retfactor(ic)
378  select type (dis => this%fmi%dis)
379  type is (distype)
380  icu = dis%get_nodeuser(ic)
381  call get_ijk(icu, dis%nrow, dis%ncol, dis%nlay, irow, icol, ilay)
382  defn%ilay = ilay
383  end select
384  defn%izone = this%izone(ic)
385  defn%can_be_rect = .true.
386  defn%can_be_quad = .false.
387  end subroutine load_properties
388 
389  !> @brief Loads face neighbors to cell definition from the grid.
390  !! Assumes cell index and number of vertices are already loaded.
391  subroutine load_neighbors(this, defn)
392  ! dummy
393  class(methoddistype), intent(inout) :: this
394  type(celldefntype), pointer, intent(inout) :: defn
395  ! local
396  integer(I4B) :: ic1
397  integer(I4B) :: ic2
398  integer(I4B) :: icu1
399  integer(I4B) :: icu2
400  integer(I4B) :: j1
401  integer(I4B) :: iloc
402  integer(I4B) :: ipos
403  integer(I4B) :: irow1
404  integer(I4B) :: irow2
405  integer(I4B) :: jcol1
406  integer(I4B) :: jcol2
407  integer(I4B) :: klay1
408  integer(I4B) :: klay2
409  integer(I4B) :: iedgeface
410 
411  select type (dis => this%fmi%dis)
412  type is (distype)
413  ! Load face neighbors
414  defn%facenbr = 0
415  ic1 = defn%icell
416  icu1 = dis%get_nodeuser(ic1)
417  call get_ijk(icu1, dis%nrow, dis%ncol, dis%nlay, &
418  irow1, jcol1, klay1)
419  call get_jk(icu1, dis%get_ncpl(), dis%nlay, j1, klay1)
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_ijk(icu2, dis%nrow, dis%ncol, dis%nlay, &
427  irow2, jcol2, klay2)
428  if (klay2 == klay1) then
429  ! Edge (polygon) face neighbor
430  if (irow2 > irow1) then
431  ! Neighbor to the S
432  iedgeface = 4
433  else if (jcol2 > jcol1) then
434  ! Neighbor to the E
435  iedgeface = 3
436  else if (irow2 < irow1) then
437  ! Neighbor to the N
438  iedgeface = 2
439  else
440  ! Neighbor to the W
441  iedgeface = 1
442  end if
443  defn%facenbr(iedgeface) = int(iloc, 1)
444  else if (klay2 > klay1) then
445  ! Bottom face neighbor
446  defn%facenbr(defn%npolyverts + 2) = int(iloc, 1)
447  else
448  ! Top face neighbor
449  defn%facenbr(defn%npolyverts + 3) = int(iloc, 1)
450  end if
451  end do
452  end select
453  ! List of edge (polygon) faces wraps around
454  defn%facenbr(defn%npolyverts + 1) = defn%facenbr(1)
455  end subroutine load_neighbors
456 
457  !> @brief Load flows into the cell definition.
458  !! These include face, boundary and net distributed flows.
459  !! Assumes cell index and number of vertices are already loaded.
460  subroutine load_flows(this, defn)
461  class(methoddistype), intent(inout) :: this
462  type(celldefntype), intent(inout) :: defn
463 
464  ! Load face flows, including boundary flows. As with cell verts,
465  ! the face flow array wraps around. Top and bottom flows make up
466  ! the last two elements, respectively, for size npolyverts + 3.
467  ! If there is no flow through any face, set a no-exit-face flag.
468  defn%faceflow = dzero
469  defn%inoexitface = 1
470  call this%load_boundary_flows_to_defn(defn)
471  call this%load_face_flows_to_defn(defn)
472 
473  ! Add up net distributed flow
474  defn%distflow = this%fmi%SourceFlows(defn%icell) + &
475  this%fmi%SinkFlows(defn%icell) + &
476  this%fmi%StorageFlows(defn%icell)
477 
478  ! Set weak sink flag
479  if (this%fmi%SinkFlows(defn%icell) .ne. dzero) then
480  defn%iweaksink = 1
481  else
482  defn%iweaksink = 0
483  end if
484  end subroutine load_flows
485 
486  subroutine load_face_flows_to_defn(this, defn)
487  ! dummy
488  class(methoddistype), intent(inout) :: this
489  type(celldefntype), intent(inout) :: defn
490  ! local
491  integer(I4B) :: m, n, nfaces
492  real(DP) :: q
493 
494  nfaces = defn%npolyverts + 3
495  do m = 1, nfaces
496  n = defn%facenbr(m)
497  if (n > 0) then
498  q = this%fmi%gwfflowja(this%fmi%dis%con%ia(defn%icell) + n)
499  defn%faceflow(m) = defn%faceflow(m) + q
500  end if
501  if (defn%faceflow(m) < dzero) defn%inoexitface = 0
502  end do
503  end subroutine load_face_flows_to_defn
504 
505  !> @brief Add boundary flows to the cell definition faceflow array.
506  !! Assumes cell index and number of vertices are already loaded.
507  subroutine load_boundary_flows_to_defn(this, defn)
508  ! dummy
509  class(methoddistype), intent(inout) :: this
510  type(celldefntype), intent(inout) :: defn
511  ! local
512  integer(I4B) :: ioffset
513 
514  ioffset = (defn%icell - 1) * max_poly_cells
515  defn%faceflow(1) = defn%faceflow(1) + &
516  this%fmi%BoundaryFlows(ioffset + 1)
517  defn%faceflow(2) = defn%faceflow(2) + &
518  this%fmi%BoundaryFlows(ioffset + 2)
519  defn%faceflow(3) = defn%faceflow(3) + &
520  this%fmi%BoundaryFlows(ioffset + 3)
521  defn%faceflow(4) = defn%faceflow(4) + &
522  this%fmi%BoundaryFlows(ioffset + 4)
523  defn%faceflow(5) = defn%faceflow(1)
524  defn%faceflow(6) = defn%faceflow(6) + &
525  this%fmi%BoundaryFlows(ioffset + 9)
526  defn%faceflow(7) = defn%faceflow(7) + &
527  this%fmi%BoundaryFlows(ioffset + 10)
528  end subroutine load_boundary_flows_to_defn
529 
530 end module methoddismodule
integer(i4b) function, public get_iatop(ncpl, icu)
Get the index corresponding to top elevation of a cell in the grid. This is -1 if the cell is in the ...
Definition: CellDefn.f90:56
integer(i4b), parameter, public max_poly_cells
Definition: Cell.f90:9
subroutine, public create_cell_rect(cell)
Create a new rectangular cell.
Definition: CellRect.f90:40
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
real(dp), parameter done
real constant 1
Definition: Constants.f90:76
Definition: Dis.f90:1
subroutine, public get_ijk(nodenumber, nrow, ncol, nlay, irow, icol, ilay)
Get row, column and layer indices from node number and grid dimensions. If nodenumber is invalid,...
Definition: GeomUtil.f90:100
subroutine, public get_jk(nodenumber, ncpl, nlay, icpl, ilay)
Get layer index and within-layer node index from node number and grid dimensions. If nodenumber is in...
Definition: GeomUtil.f90:128
This module defines variable data types.
Definition: kind.f90:8
pure logical function, public is_close(a, b, rtol, atol, symmetric)
Check if a real value is approximately equal to another.
Definition: MathUtil.f90:46
Cell-level tracking methods.
type(methodcellpollocktype), pointer, public method_cell_plck
type(methodcellpasstobottype), pointer, public method_cell_ptb
subroutine load_flows(this, defn)
Load flows into the cell definition. These include face, boundary and net distributed flows....
Definition: MethodDis.f90:461
subroutine load_boundary_flows_to_defn(this, defn)
Add boundary flows to the cell definition faceflow array. Assumes cell index and number of vertices a...
Definition: MethodDis.f90:508
double precision function get_top(this, iatop)
Returns a top elevation based on index iatop.
Definition: MethodDis.f90:324
subroutine load_particle(this, cell, particle)
Load cell properties into the particle, including.
Definition: MethodDis.f90:171
subroutine load_properties(this, ic, defn)
Definition: MethodDis.f90:360
subroutine update_flowja(this, cell, particle)
Update cell-cell flows of particle mass. Every particle is currently assigned unit mass.
Definition: MethodDis.f90:261
subroutine, public create_method_dis(method)
Create a new structured grid (DIS) tracking method.
Definition: MethodDis.f90:44
subroutine load_celldefn(this, ic, defn)
Loads cell definition from the grid.
Definition: MethodDis.f90:337
subroutine apply_dis(this, particle, tmax)
Apply the structured tracking method to a particle.
Definition: MethodDis.f90:315
subroutine load_dis(this, particle, next_level, submethod)
Load the cell geometry and tracking method.
Definition: MethodDis.f90:137
subroutine load_cell(this, ic, cell)
Definition: MethodDis.f90:64
subroutine load_face_flows_to_defn(this, defn)
Definition: MethodDis.f90:487
subroutine load_neighbors(this, defn)
Loads face neighbors to cell definition from the grid. Assumes cell index and number of vertices are ...
Definition: MethodDis.f90:392
subroutine pass_dis(this, particle)
Pass a particle to the next cell, if there is one.
Definition: MethodDis.f90:286
Particle tracking strategies.
Definition: Method.f90:2
Base grid cell definition.
Definition: CellDefn.f90:12
Structured grid discretization.
Definition: Dis.f90:23
Base type for particle tracking methods.
Definition: Method.f90:31
Particle tracked by the PRT model.
Definition: Particle.f90:32