MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
methoddismodule Module Reference

Data Types

type  methoddistype
 

Functions/Subroutines

subroutine, public create_method_dis (method)
 Create a new structured grid (DIS) tracking method. More...
 
subroutine destroy (this)
 Destructor the tracking method. More...
 
subroutine load_dis (this, particle, next_level, submethod)
 Load the cell geometry and method (tracking strategy) More...
 
subroutine pass_dis (this, particle)
 Pass a particle to the next cell, if there is one. More...
 
subroutine apply_dis (this, particle, tmax)
 Apply the method to a particle. More...
 
double precision function get_top (this, iatop)
 Returns a top elevation based on index iatop. More...
 
subroutine load_cell_defn (this, ic, defn)
 Loads cell definition from the grid. More...
 
subroutine load_nbrs_to_defn (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_to_defn (this, defn)
 Load flows into the cell definition. These include face flows and net distributed flows. Assumes cell index and number of vertices are already loaded. More...
 
subroutine load_boundary_flows_to_defn (this, defn)
 Add boundary flows to the cell definition faceflow array. Assumes cell index and number of vertices are already loaded. More...
 

Function/Subroutine Documentation

◆ apply_dis()

subroutine methoddismodule::apply_dis ( class(methoddistype), intent(inout)  this,
type(particletype), intent(inout), pointer  particle,
real(dp), intent(in)  tmax 
)
private

Definition at line 229 of file MethodDis.f90.

230  class(MethodDisType), intent(inout) :: this
231  type(ParticleType), pointer, intent(inout) :: particle
232  real(DP), intent(in) :: tmax
233 
234  call this%track(particle, 1, tmax) ! kluge, hardwired to level 1

◆ create_method_dis()

subroutine, public methoddismodule::create_method_dis ( type(methoddistype), pointer  method)

Definition at line 36 of file MethodDis.f90.

37  ! -- dummy
38  type(MethodDisType), pointer :: method
39  ! -- local
40  type(CellRectType), pointer :: cell
41 
42  allocate (method)
43  allocate (method%type)
44  call create_cell_rect(cell)
45  method%cell => cell
46  method%type = "dis"
47  method%delegates = .true.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ destroy()

subroutine methoddismodule::destroy ( class(methoddistype), intent(inout)  this)
private

Definition at line 51 of file MethodDis.f90.

52  class(MethodDisType), intent(inout) :: this
53  deallocate (this%type)

◆ get_top()

double precision function methoddismodule::get_top ( class(methoddistype), intent(inout)  this,
integer, intent(in)  iatop 
)
private

Definition at line 238 of file MethodDis.f90.

239  class(MethodDisType), intent(inout) :: this
240  integer, intent(in) :: iatop
241  double precision :: top
242 
243  if (iatop .lt. 0) then
244  top = this%fmi%dis%top(-iatop)
245  else
246  top = this%fmi%dis%bot(iatop)
247  end if

◆ load_boundary_flows_to_defn()

subroutine methoddismodule::load_boundary_flows_to_defn ( class(methoddistype), intent(inout)  this,
type(celldefntype), intent(inout)  defn 
)
private

Definition at line 409 of file MethodDis.f90.

410  ! -- dummy
411  class(MethodDisType), intent(inout) :: this
412  type(CellDefnType), intent(inout) :: defn
413  ! -- local
414  integer(I4B) :: ic
415  integer(I4B) :: ioffset
416 
417  ic = defn%icell
418  ioffset = (ic - 1) * 10
419  defn%faceflow(1) = defn%faceflow(1) + &
420  this%fmi%BoundaryFlows(ioffset + 1) ! kluge note: should these be additive (seems so)???
421  defn%faceflow(2) = defn%faceflow(2) + &
422  this%fmi%BoundaryFlows(ioffset + 2)
423  defn%faceflow(3) = defn%faceflow(3) + &
424  this%fmi%BoundaryFlows(ioffset + 3)
425  defn%faceflow(4) = defn%faceflow(4) + &
426  this%fmi%BoundaryFlows(ioffset + 4)
427  defn%faceflow(5) = defn%faceflow(1)
428  defn%faceflow(6) = defn%faceflow(6) + &
429  this%fmi%BoundaryFlows(ioffset + 9)
430  defn%faceflow(7) = defn%faceflow(7) + &
431  this%fmi%BoundaryFlows(ioffset + 10)

◆ load_cell_defn()

subroutine methoddismodule::load_cell_defn ( class(methoddistype), intent(inout)  this,
integer(i4b), intent(in)  ic,
type(celldefntype), intent(inout), pointer  defn 
)
private

Definition at line 251 of file MethodDis.f90.

252  ! -- dummy
253  class(MethodDisType), intent(inout) :: this
254  integer(I4B), intent(in) :: ic
255  type(CellDefnType), pointer, intent(inout) :: defn
256 
257  select type (dis => this%fmi%dis)
258  type is (distype)
259  ! -- Set basic cell properties
260  defn%icell = ic
261  defn%npolyverts = 4 ! rectangular cell always has 4 vertices
262  defn%iatop = get_iatop(dis%get_ncpl(), &
263  dis%get_nodeuser(ic))
264  defn%top = dis%bot(ic) + &
265  this%fmi%gwfsat(ic) * (dis%top(ic) - dis%bot(ic))
266  defn%bot = dis%bot(ic)
267  defn%sat = this%fmi%gwfsat(ic)
268  defn%porosity = this%porosity(ic)
269  defn%retfactor = this%retfactor(ic)
270  defn%izone = this%izone(ic)
271  defn%can_be_rect = .true.
272  defn%can_be_quad = .false.
273 
274  ! -- Load cell polygon vertices
275  call dis%get_polyverts( &
276  defn%icell, &
277  defn%polyvert, &
278  closed=.true.)
279 
280  ! -- Load face neighbors
281  call this%load_nbrs_to_defn(defn)
282 
283  ! -- Load 180 degree face indicators
284  defn%ispv180(1:defn%npolyverts + 1) = .false.
285 
286  ! -- Load flows (assumes face neighbors already loaded)
287  call this%load_flows_to_defn(defn)
288  end select
Here is the call graph for this function:

◆ load_dis()

subroutine methoddismodule::load_dis ( class(methoddistype), 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 57 of file MethodDis.f90.

58  ! -- dummy
59  class(MethodDisType), intent(inout) :: this
60  type(ParticleType), pointer, intent(inout) :: particle
61  integer(I4B), intent(in) :: next_level
62  class(MethodType), pointer, intent(inout) :: submethod
63  ! -- local
64  integer(I4B) :: ic
65  integer(I4B) :: icu
66  integer(I4B) :: irow
67  integer(I4B) :: jcol
68  integer(I4B) :: klay
69  real(DP) :: areax
70  real(DP) :: areay
71  real(DP) :: areaz
72  real(DP) :: dx
73  real(DP) :: dy
74  real(DP) :: dz
75  real(DP) :: factor
76  real(DP) :: term
77 
78  select type (cell => this%cell)
79  type is (cellrecttype)
80  select type (dis => this%fmi%dis)
81  type is (distype)
82  ic = particle%idomain(next_level)
83  call this%load_cell_defn(ic, cell%defn)
84 
85  ! -- If cell is active but dry, select and initialize
86  ! -- pass-to-bottom method and set cell method pointer
87  if (this%fmi%ibdgwfsat0(ic) == 0) then ! kluge note: use cellDefn%sat == DZERO here instead?
88  call method_cell_ptb%init( &
89  cell=this%cell, &
90  trackfilectl=this%trackfilectl, &
91  tracktimes=this%tracktimes)
92  submethod => method_cell_ptb
93  else
94  ! -- load rectangular cell (todo: refactor into separate routine)
95  icu = dis%get_nodeuser(ic)
96  call get_ijk(icu, dis%nrow, dis%ncol, dis%nlay, &
97  irow, jcol, klay)
98  dx = dis%delr(jcol)
99  dy = dis%delc(irow)
100  dz = cell%defn%top - cell%defn%bot
101  cell%dx = dx
102  cell%dy = dy
103  cell%dz = dz
104  cell%sinrot = dzero
105  cell%cosrot = done
106  cell%xOrigin = cell%defn%polyvert(1, 1) ! kluge note: could avoid using polyvert here
107  cell%yOrigin = cell%defn%polyvert(2, 1)
108  cell%zOrigin = cell%defn%bot
109  cell%ipvOrigin = 1
110  areax = dy * dz
111  areay = dx * dz
112  areaz = dx * dy
113  factor = done / cell%defn%retfactor
114  factor = factor / cell%defn%porosity
115  term = factor / areax
116  cell%vx1 = cell%defn%faceflow(1) * term
117  cell%vx2 = -cell%defn%faceflow(3) * term
118  term = factor / areay
119  cell%vy1 = cell%defn%faceflow(4) * term
120  cell%vy2 = -cell%defn%faceflow(2) * term
121  term = factor / areaz
122  cell%vz1 = cell%defn%faceflow(6) * term
123  cell%vz2 = -cell%defn%faceflow(7) * term
124 
125  ! -- Select and initialize Pollock's method and set method pointer
126  call method_cell_plck%init( &
127  cell=this%cell, &
128  trackfilectl=this%trackfilectl, &
129  tracktimes=this%tracktimes)
130  submethod => method_cell_plck
131  end if
132  end select
133  end select
Here is the call graph for this function:

◆ load_flows_to_defn()

subroutine methoddismodule::load_flows_to_defn ( class(methoddistype), intent(inout)  this,
type(celldefntype), intent(inout)  defn 
)
private

Definition at line 362 of file MethodDis.f90.

363  ! -- dummy
364  class(MethodDisType), intent(inout) :: this
365  type(CellDefnType), intent(inout) :: defn
366  ! -- local
367  integer(I4B) :: ic
368  integer(I4B) :: m
369  integer(I4B) :: n
370  integer(I4B) :: npolyverts
371 
372  ic = defn%icell
373  npolyverts = defn%npolyverts
374 
375  ! -- Load face flows.
376  defn%faceflow = 0d0 ! kluge note: eventually use DZERO for 0d0 throughout
377  ! -- As with polygon nbrs, polygon face flows wrap around for
378  ! -- convenience at position npolyverts+1, and bot and top flows
379  ! -- are tacked on the end of the list
380  do m = 1, npolyverts + 3
381  n = defn%facenbr(m)
382  if (n > 0) &
383  defn%faceflow(m) = this%fmi%gwfflowja(this%fmi%dis%con%ia(ic) + n)
384  end do
385 
386  ! -- Add BoundaryFlows to face flows
387  call this%load_boundary_flows_to_defn(defn)
388  ! -- Set inoexitface flag
389  defn%inoexitface = 1
390  do m = 1, npolyverts + 3 ! kluge note: can be streamlined with above code
391  if (defn%faceflow(m) < 0d0) defn%inoexitface = 0
392  end do
393 
394  ! -- Add up net distributed flow
395  defn%distflow = this%fmi%SourceFlows(ic) + this%fmi%SinkFlows(ic) + &
396  this%fmi%StorageFlows(ic)
397 
398  ! -- Set weak sink flag
399  if (this%fmi%SinkFlows(ic) .ne. 0d0) then
400  defn%iweaksink = 1
401  else
402  defn%iweaksink = 0
403  end if
404 

◆ load_nbrs_to_defn()

subroutine methoddismodule::load_nbrs_to_defn ( class(methoddistype), intent(inout)  this,
type(celldefntype), intent(inout), pointer  defn 
)
private

Definition at line 293 of file MethodDis.f90.

294  ! -- dummy
295  class(MethodDisType), intent(inout) :: this
296  type(CellDefnType), pointer, intent(inout) :: defn
297  ! -- local
298  integer(I4B) :: ic1
299  integer(I4B) :: ic2
300  integer(I4B) :: icu1
301  integer(I4B) :: icu2
302  integer(I4B) :: j1
303  integer(I4B) :: iloc
304  integer(I4B) :: ipos
305  integer(I4B) :: irow1
306  integer(I4B) :: irow2
307  integer(I4B) :: jcol1
308  integer(I4B) :: jcol2
309  integer(I4B) :: klay1
310  integer(I4B) :: klay2
311  integer(I4B) :: iedgeface
312 
313  select type (dis => this%fmi%dis)
314  type is (distype)
315  ! -- Load face neighbors
316  defn%facenbr = 0
317  ic1 = defn%icell
318  icu1 = dis%get_nodeuser(ic1)
319  call get_ijk(icu1, dis%nrow, dis%ncol, dis%nlay, &
320  irow1, jcol1, klay1)
321  call get_jk(icu1, dis%get_ncpl(), dis%nlay, j1, klay1)
322  do iloc = 1, dis%con%ia(ic1 + 1) - dis%con%ia(ic1) - 1
323  ipos = dis%con%ia(ic1) + iloc
324  if (dis%con%mask(ipos) == 0) cycle ! kluge note: need mask here???
325  ic2 = dis%con%ja(ipos)
326  icu2 = dis%get_nodeuser(ic2)
327  call get_ijk(icu2, dis%nrow, dis%ncol, dis%nlay, &
328  irow2, jcol2, klay2)
329  if (klay2 == klay1) then
330  ! -- Edge (polygon) face neighbor
331  if (irow2 > irow1) then
332  ! Neighbor to the S
333  iedgeface = 4 ! kluge note: make sure this numbering is consistent with numbering in cell method
334  else if (jcol2 > jcol1) then
335  ! Neighbor to the E
336  iedgeface = 3
337  else if (irow2 < irow1) then
338  ! Neighbor to the N
339  iedgeface = 2
340  else
341  ! Neighbor to the W
342  iedgeface = 1
343  end if
344  defn%facenbr(iedgeface) = int(iloc, 1)
345  else if (klay2 > klay1) then
346  ! -- Bottom face neighbor
347  defn%facenbr(defn%npolyverts + 2) = int(iloc, 1)
348  else
349  ! -- Top face neighbor
350  defn%facenbr(defn%npolyverts + 3) = int(iloc, 1)
351  end if
352  end do
353  end select
354  ! -- List of edge (polygon) faces wraps around
355  ! todo: why need to wrap around? no analog to "closing" a polygon?
356  defn%facenbr(defn%npolyverts + 1) = defn%facenbr(1)
Here is the call graph for this function:

◆ pass_dis()

subroutine methoddismodule::pass_dis ( class(methoddistype), intent(inout)  this,
type(particletype), intent(inout), pointer  particle 
)
private

Definition at line 137 of file MethodDis.f90.

138  ! -- dummy
139  class(MethodDisType), intent(inout) :: this
140  type(ParticleType), pointer, intent(inout) :: particle
141  ! -- local
142  integer(I4B) :: inface
143  integer(I4B) :: ipos
144  integer(I4B) :: ic
145  integer(I4B) :: icu
146  integer(I4B) :: inbr
147  integer(I4B) :: idiag
148  integer(I4B) :: ilay
149  integer(I4B) :: irow
150  integer(I4B) :: icol
151  real(DP) :: z
152  real(DP) :: zrel
153  real(DP) :: topfrom
154  real(DP) :: botfrom
155  real(DP) :: top
156  real(DP) :: bot
157  real(DP) :: sat
158 
159  inface = particle%iboundary(2)
160  z = particle%z
161 
162  select type (cell => this%cell)
163  type is (cellrecttype)
164  select type (dis => this%fmi%dis)
165  type is (distype)
166  inbr = cell%defn%facenbr(inface)
167  if (inbr .eq. 0) then
168  ! -- Exterior face; no neighbor to map to
169  ! particle%idomain(1) = 0
170  ! particle%idomain(2) = 0 ! kluge note: set a "has_exited" attribute instead???
171  ! particle%idomain(1) = -abs(particle%idomain(1)) ! kluge???
172  ! particle%idomain(2) = -abs(particle%idomain(2)) ! kluge???
173  particle%istatus = 2 ! kluge note: use -2 to allow check for transfer to another model???
174  particle%advancing = .false.
175  call this%save(particle, reason=3) ! reason=3: termination
176  ! particle%iboundary(2) = -1
177  else
178  idiag = dis%con%ia(cell%defn%icell)
179  ipos = idiag + inbr
180  ic = dis%con%ja(ipos) ! kluge note: use PRT model's DIS instead of fmi's???
181  particle%idomain(2) = ic
182 
183  ! compute and set user node number and layer on particle
184  icu = dis%get_nodeuser(ic)
185  call get_ijk(icu, dis%nrow, dis%ncol, dis%nlay, &
186  irow, icol, ilay)
187  particle%icu = icu
188  particle%ilay = ilay
189 
190  ! call this%mapToNbrCell(cellRect%cellDefn,inface,z)
191  if (inface .eq. 1) then
192  inface = 3
193  else if (inface .eq. 2) then
194  inface = 4
195  else if (inface .eq. 3) then
196  inface = 1
197  else if (inface .eq. 4) then
198  inface = 2
199  else if (inface .eq. 6) then
200  inface = 7
201  else if (inface .eq. 7) then
202  inface = 6
203  end if
204  particle%iboundary(2) = inface
205  if (inface < 5) then
206  ! -- Map z between cells
207  topfrom = cell%defn%top
208  botfrom = cell%defn%bot
209  zrel = (z - botfrom) / (topfrom - botfrom)
210  top = dis%top(ic) ! kluge note: use PRT model's DIS instead of fmi's???
211  bot = dis%bot(ic)
212  sat = this%fmi%gwfsat(ic)
213  z = bot + zrel * sat * (top - bot)
214  end if
215  particle%z = z
216  ! -- Update cell-cell flows of particle mass.
217  ! Every particle is currently assigned unit mass.
218  ! -- leaving old cell
219  this%flowja(ipos) = this%flowja(ipos) - done
220  ! -- entering new cell
221  this%flowja(dis%con%isym(ipos)) &
222  = this%flowja(dis%con%isym(ipos)) + done
223  end if
224  end select
225  end select
Here is the call graph for this function: