MODFLOW 6  version 6.8.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
11  use cellmodule
12  use celldefnmodule
13  use cellrectmodule
14  use particlemodule
15  use prtfmimodule, only: prtfmitype
16  use dismodule, only: distype
17  use geomutilmodule, only: get_ijk, get_jk
18  use mathutilmodule, only: is_close
19  implicit none
20 
21  private
22  public :: methoddistype
23  public :: create_method_dis
24 
25  type, extends(methodmodeltype) :: methoddistype
26  private
27  type(methodcellpollocktype), pointer :: method_cell_plck => null()
28  type(methodcellpasstobottype), pointer :: method_cell_ptb => null()
29  contains
30  procedure, public :: apply => apply_dis !< apply the DIS tracking method
31  procedure, public :: deallocate !< deallocate arrays and scalars
32  procedure, public :: load => load_dis !< load the method
33  procedure, public :: pass => pass_dis !< pass the particle to the next domain
34  procedure :: get_top !< get cell top elevation
35  procedure :: update_flowja !< load intercell mass flows
36  procedure :: load_particle !< load particle properties
37  procedure :: load_cell_properties !< load basic cell properties
38  procedure :: load_cell_neighbors !< load cell face neighbors
39  procedure :: load_cell_flows !< load cell flows
40  procedure :: load_cell_boundary_flows !< load boundary flows to the cell definition
41  procedure :: load_cell_face_flows !< load face flows to the cell definition
42  procedure :: load_cell_defn !< load cell definition from the grid
43  procedure :: load_cell !< load cell geometry and flows
44  end type methoddistype
45 
46 contains
47 
48  !> @brief Create a new structured grid (DIS) tracking method
49  subroutine create_method_dis(method)
50  ! dummy
51  type(methoddistype), pointer :: method
52  ! local
53  type(cellrecttype), pointer :: cell
54 
55  allocate (method)
56  allocate (method%name)
57  call create_cell_rect(cell)
58  method%cell => cell
59  method%name = "dis"
60  method%delegates = .true.
61 
62  ! Create cell methods this method can delegate to
63  call create_method_cell_pollock(method%method_cell_plck)
64  call create_method_cell_ptb(method%method_cell_ptb)
65  end subroutine create_method_dis
66 
67  !> @brief Destructor the tracking method
68  subroutine deallocate (this)
69  class(methoddistype), intent(inout) :: this
70  deallocate (this%name)
71 
72  ! Deallocate owned cell methods
73  call this%method_cell_plck%deallocate()
74  deallocate (this%method_cell_plck)
75  call this%method_cell_ptb%deallocate()
76  deallocate (this%method_cell_ptb)
77  end subroutine deallocate
78 
79  subroutine load_cell(this, ic, cell)
80  ! dummy
81  class(methoddistype), intent(inout) :: this
82  integer(I4B), intent(in) :: ic
83  type(cellrecttype), intent(inout) :: cell
84  ! local
85  integer(I4B) :: icu
86  integer(I4B) :: irow
87  integer(I4B) :: jcol
88  integer(I4B) :: klay
89  real(DP) :: areax
90  real(DP) :: areay
91  real(DP) :: areaz
92  real(DP) :: dx
93  real(DP) :: dy
94  real(DP) :: dz
95  real(DP) :: factor
96  real(DP) :: term
97 
98  select type (dis => this%fmi%dis)
99  type is (distype)
100  icu = dis%get_nodeuser(ic)
101 
102  call get_ijk(icu, dis%nrow, dis%ncol, dis%nlay, &
103  irow, jcol, klay)
104 
105  dx = dis%delr(jcol)
106  dy = dis%delc(irow)
107  dz = cell%defn%top - cell%defn%bot
108 
109  cell%dx = dx
110  cell%dy = dy
111  cell%dz = dz
112  cell%sinrot = dzero
113  cell%cosrot = done
114  cell%xOrigin = cell%defn%polyvert(1, 1)
115  cell%yOrigin = cell%defn%polyvert(2, 1)
116  cell%zOrigin = cell%defn%bot
117  cell%ipvOrigin = 1
118 
119  factor = done / cell%defn%retfactor
120  factor = factor / cell%defn%porosity
121 
122  areaz = dx * dy
123  term = factor / areaz
124 
125  cell%vz1 = cell%defn%faceflow(6) * term
126  cell%vz2 = -cell%defn%faceflow(7) * term
127 
128  if (this%fmi%ibdgwfsat0(ic) == 0) then
129  cell%vx1 = dzero
130  cell%vx2 = dzero
131  cell%vy1 = dzero
132  cell%vy2 = dzero
133  cell%vz1 = dzero
134  cell%vz2 = dzero
135  return
136  end if
137 
138  areax = dy * dz
139  term = factor / areax
140  cell%vx1 = cell%defn%faceflow(1) * term
141  cell%vx2 = -cell%defn%faceflow(3) * term
142 
143  areay = dx * dz
144  term = factor / areay
145  cell%vy1 = cell%defn%faceflow(4) * term
146  cell%vy2 = -cell%defn%faceflow(2) * term
147 
148  end select
149  end subroutine load_cell
150 
151  !> @brief Load the cell geometry and tracking method
152  subroutine load_dis(this, particle, next_level, submethod)
153  ! dummy
154  class(methoddistype), intent(inout) :: this
155  type(particletype), pointer, intent(inout) :: particle
156  integer(I4B), intent(in) :: next_level
157  class(methodtype), pointer, intent(inout) :: submethod
158  ! local
159  integer(I4B) :: ic
160 
161  select type (cell => this%cell)
162  type is (cellrecttype)
163  ic = particle%itrdomain(next_level)
164  call this%load_cell_defn(ic, cell%defn)
165  call this%load_cell(ic, cell)
166  if (this%fmi%ibdgwfsat0(ic) == 0) then
167  call this%method_cell_ptb%init( &
168  fmi=this%fmi, &
169  cell=this%cell, &
170  events=this%events, &
171  tracktimes=this%tracktimes)
172  submethod => this%method_cell_ptb
173  else
174  call this%method_cell_plck%init( &
175  fmi=this%fmi, &
176  cell=this%cell, &
177  events=this%events, &
178  tracktimes=this%tracktimes)
179  submethod => this%method_cell_plck
180  end if
181  end select
182  end subroutine load_dis
183 
184  !> @brief Load cell properties into the particle, including
185  ! the z coordinate, entry face, and node and layer numbers.
186  subroutine load_particle(this, cell, particle)
187  use particlemodule, only: term_boundary
188  use particleeventmodule, only: terminate
189  ! dummy
190  class(methoddistype), intent(inout) :: this
191  type(cellrecttype), pointer, intent(inout) :: cell
192  type(particletype), pointer, intent(inout) :: particle
193  ! local
194  integer(I4B) :: ic
195  integer(I4B) :: icu
196  integer(I4B) :: ilay
197  integer(I4B) :: irow
198  integer(I4B) :: icol
199  integer(I4B) :: inface
200  integer(I4B) :: idiag
201  integer(I4B) :: inbr
202  integer(I4B) :: ipos
203  real(DP) :: z
204  real(DP) :: zrel
205  real(DP) :: topfrom
206  real(DP) :: botfrom
207  real(DP) :: top
208  real(DP) :: bot
209  real(DP) :: sat
210 
211  select type (dis => this%fmi%dis)
212  type is (distype)
213  ! compute reduced/user node numbers and layer
214  inface = particle%iboundary(level_feature)
215  inbr = cell%defn%facenbr(inface)
216  idiag = this%fmi%dis%con%ia(cell%defn%icell)
217  ipos = idiag + inbr
218  ic = dis%con%ja(ipos)
219  icu = dis%get_nodeuser(ic)
220  call get_ijk(icu, dis%nrow, dis%ncol, dis%nlay, irow, icol, ilay)
221 
222  ! update node numbers and layer
223  particle%itrdomain(level_feature) = 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(level_feature) = 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(level_feature)
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  use particlemodule, only: term_boundary
287  use particleeventmodule, only: terminate
288  ! dummy
289  class(methoddistype), intent(inout) :: this
290  type(particletype), pointer, intent(inout) :: particle
291  ! local
292  type(cellrecttype), pointer :: cell
293  integer(I4B) :: iface
294  logical(LGP) :: no_neighbors
295 
296  select type (c => this%cell)
297  type is (cellrecttype)
298  cell => c
299  iface = particle%iboundary(level_feature)
300  no_neighbors = cell%defn%facenbr(iface) == 0
301 
302  ! todo AMP: reconsider when multiple models supported
303  if (no_neighbors) then
304  call this%terminate(particle, status=term_boundary)
305  return
306  end if
307 
308  call this%load_particle(cell, particle)
309  call this%update_flowja(cell, particle)
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_cell_defn(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  call this%load_cell_properties(ic, defn)
343  call this%fmi%dis%get_polyverts( &
344  defn%icell, &
345  defn%polyvert, &
346  closed=.true.)
347  call this%load_cell_neighbors(defn)
348  call this%load_cell_saturation_status(defn)
349  defn%ispv180(1:defn%npolyverts + 1) = .false.
350  call this%load_cell_flows(defn)
351 
352  end subroutine load_cell_defn
353 
354  subroutine load_cell_properties(this, ic, defn)
355  ! dummy
356  class(methoddistype), intent(inout) :: this
357  integer(I4B), intent(in) :: ic
358  type(celldefntype), pointer, intent(inout) :: defn
359  ! local
360  integer(I4B) :: irow, icol, ilay, icu
361 
362  defn%icell = ic
363  defn%npolyverts = 4 ! rectangular cell always has 4 vertices
364  defn%iatop = get_iatop(this%fmi%dis%get_ncpl(), &
365  this%fmi%dis%get_nodeuser(ic))
366  defn%top = this%fmi%dis%bot(ic) + &
367  this%fmi%gwfsat(ic) * &
368  (this%fmi%dis%top(ic) - this%fmi%dis%bot(ic))
369  defn%bot = this%fmi%dis%bot(ic)
370  defn%sat = this%fmi%gwfsat(ic)
371  defn%porosity = this%porosity(ic)
372  defn%retfactor = this%retfactor(ic)
373  select type (dis => this%fmi%dis)
374  type is (distype)
375  icu = dis%get_nodeuser(ic)
376  call get_ijk(icu, dis%nrow, dis%ncol, dis%nlay, irow, icol, ilay)
377  defn%ilay = ilay
378  end select
379  defn%izone = this%izone(ic)
380  defn%can_be_rect = .true.
381  defn%can_be_quad = .false.
382 
383  end subroutine load_cell_properties
384 
385  !> @brief Loads face neighbors to cell definition from the grid.
386  !! Assumes cell index and number of vertices are already loaded.
387  subroutine load_cell_neighbors(this, defn)
388  ! dummy
389  class(methoddistype), intent(inout) :: this
390  type(celldefntype), pointer, intent(inout) :: defn
391  ! local
392  integer(I4B) :: ic1
393  integer(I4B) :: ic2
394  integer(I4B) :: icu1
395  integer(I4B) :: icu2
396  integer(I4B) :: j1
397  integer(I4B) :: iloc
398  integer(I4B) :: ipos
399  integer(I4B) :: irow1
400  integer(I4B) :: irow2
401  integer(I4B) :: jcol1
402  integer(I4B) :: jcol2
403  integer(I4B) :: klay1
404  integer(I4B) :: klay2
405  integer(I4B) :: iedgeface
406 
407  select type (dis => this%fmi%dis)
408  type is (distype)
409  ! Load face neighbors
410  defn%facenbr = 0
411  ic1 = defn%icell
412  icu1 = dis%get_nodeuser(ic1)
413  call get_ijk(icu1, dis%nrow, dis%ncol, dis%nlay, &
414  irow1, jcol1, klay1)
415  call get_jk(icu1, dis%get_ncpl(), dis%nlay, j1, klay1)
416  do iloc = 1, dis%con%ia(ic1 + 1) - dis%con%ia(ic1) - 1
417  ipos = dis%con%ia(ic1) + iloc
418  ! mask could become relevant if PRT uses interface model
419  if (dis%con%mask(ipos) == 0) cycle
420  ic2 = dis%con%ja(ipos)
421  icu2 = dis%get_nodeuser(ic2)
422  call get_ijk(icu2, dis%nrow, dis%ncol, dis%nlay, &
423  irow2, jcol2, klay2)
424  if (klay2 == klay1) then
425  ! Edge (polygon) face neighbor
426  if (irow2 > irow1) then
427  ! Neighbor to the S
428  iedgeface = 4
429  else if (jcol2 > jcol1) then
430  ! Neighbor to the E
431  iedgeface = 3
432  else if (irow2 < irow1) then
433  ! Neighbor to the N
434  iedgeface = 2
435  else
436  ! Neighbor to the W
437  iedgeface = 1
438  end if
439  defn%facenbr(iedgeface) = int(iloc, 1)
440  else if (klay2 > klay1) then
441  ! Bottom face neighbor
442  defn%facenbr(defn%npolyverts + 2) = int(iloc, 1)
443  else
444  ! Top face neighbor
445  defn%facenbr(defn%npolyverts + 3) = int(iloc, 1)
446  end if
447  end do
448  end select
449  ! List of edge (polygon) faces wraps around
450  defn%facenbr(defn%npolyverts + 1) = defn%facenbr(1)
451  end subroutine load_cell_neighbors
452 
453  !> @brief Load flows into the cell definition.
454  !! These include face, boundary and net distributed flows.
455  !! Assumes cell index and number of vertices are already loaded.
456  subroutine load_cell_flows(this, defn)
457  class(methoddistype), intent(inout) :: this
458  type(celldefntype), pointer, intent(inout) :: defn
459 
460  ! Load face flows, including boundary flows. As with cell verts,
461  ! the face flow array wraps around. Top and bottom flows make up
462  ! the last two elements, respectively, for size npolyverts + 3.
463  ! If there is no flow through any face, set a no-exit-face flag.
464  defn%faceflow = dzero
465  call this%load_cell_boundary_flows(defn)
466  call this%load_cell_face_flows(defn)
467  call this%cap_cell_wt_flow(defn)
468  call this%load_cell_no_exit_face(defn)
469 
470  ! Add up net distributed flow
471  defn%distflow = this%fmi%SourceFlows(defn%icell) + &
472  this%fmi%SinkFlows(defn%icell) + &
473  this%fmi%StorageFlows(defn%icell)
474 
475  ! Set weak sink flag
476  if (this%fmi%SinkFlows(defn%icell) .ne. dzero) then
477  defn%iweaksink = 1
478  else
479  defn%iweaksink = 0
480  end if
481  end subroutine load_cell_flows
482 
483  subroutine load_cell_face_flows(this, defn)
484  ! dummy
485  class(methoddistype), intent(inout) :: this
486  type(celldefntype), pointer, intent(inout) :: defn
487  ! local
488  integer(I4B) :: m, n, nfaces
489  real(DP) :: q
490 
491  nfaces = defn%npolyverts + 3
492  do m = 1, nfaces
493  n = defn%facenbr(m)
494  if (n > 0) then
495  q = this%fmi%gwfflowja(this%fmi%dis%con%ia(defn%icell) + n)
496  defn%faceflow(m) = defn%faceflow(m) + q
497  end if
498  end do
499  end subroutine load_cell_face_flows
500 
501  !> @brief Add boundary flows to the cell definition faceflow array.
502  !! Assumes cell index and number of vertices are already loaded.
503  subroutine load_cell_boundary_flows(this, defn)
504  ! dummy
505  class(methoddistype), intent(inout) :: this
506  type(celldefntype), pointer, intent(inout) :: defn
507  ! local
508 
509  defn%faceflow(1) = defn%faceflow(1) + &
510  this%fmi%BoundaryFlows(defn%icell, 1)
511  defn%faceflow(2) = defn%faceflow(2) + &
512  this%fmi%BoundaryFlows(defn%icell, 2)
513  defn%faceflow(3) = defn%faceflow(3) + &
514  this%fmi%BoundaryFlows(defn%icell, 3)
515  defn%faceflow(4) = defn%faceflow(4) + &
516  this%fmi%BoundaryFlows(defn%icell, 4)
517  defn%faceflow(5) = defn%faceflow(1)
518  defn%faceflow(6) = defn%faceflow(6) + &
519  this%fmi%BoundaryFlows(defn%icell, this%fmi%max_faces - 1)
520  defn%faceflow(7) = defn%faceflow(7) + &
521  this%fmi%BoundaryFlows(defn%icell, this%fmi%max_faces)
522  end subroutine load_cell_boundary_flows
523 
524 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:70
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
subroutine, public create_method_cell_ptb(method)
Create a new pass-to-bottom tracking method.
procedure subroutine, public create_method_cell_pollock(method)
Create a tracking method.
subroutine load_cell_face_flows(this, defn)
Definition: MethodDis.f90:484
subroutine load_cell_defn(this, ic, defn)
Loads cell definition from the grid.
Definition: MethodDis.f90:337
double precision function get_top(this, iatop)
Returns a top elevation based on index iatop.
Definition: MethodDis.f90:324
subroutine load_cell_properties(this, ic, defn)
Definition: MethodDis.f90:355
subroutine load_particle(this, cell, particle)
Load cell properties into the particle, including.
Definition: MethodDis.f90:187
subroutine load_cell_boundary_flows(this, defn)
Add boundary flows to the cell definition faceflow array. Assumes cell index and number of vertices a...
Definition: MethodDis.f90:504
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:50
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:153
subroutine load_cell_neighbors(this, defn)
Loads face neighbors to cell definition from the grid. Assumes cell index and number of vertices are ...
Definition: MethodDis.f90:388
subroutine load_cell(this, ic, cell)
Definition: MethodDis.f90:80
subroutine load_cell_flows(this, defn)
Load flows into the cell definition. These include face, boundary and net distributed flows....
Definition: MethodDis.f90:457
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
@, public level_feature
Definition: Method.f90:41
@, public terminate
particle terminated
@ term_boundary
terminated at a boundary face
Definition: Particle.f90:31
Base grid cell definition.
Definition: CellDefn.f90:25
Structured grid discretization.
Definition: Dis.f90:23
Base type for particle tracking methods.
Definition: Method.f90:58
Particle tracked by the PRT model.
Definition: Particle.f90:62