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