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