MODFLOW 6  version 6.8.0.dev0
USGS Modular Hydrologic Model
MethodDisv.f90
Go to the documentation of this file.
2 
3  use kindmodule, only: dp, i4b, lgp
4  use errorutilmodule, only: pstop
5  use constantsmodule, only: done, dzero
16  use celldefnmodule
17  use cellpolymodule
18  use particlemodule
19  use prtfmimodule, only: prtfmitype
20  use disvmodule, only: disvtype
23  implicit none
24 
25  private
26  public :: methoddisvtype
27  public :: create_method_disv
28 
29  type, extends(methodmodeltype) :: methoddisvtype
30  private
31  type(celldefntype), pointer :: neighbor => null() !< ptr to a neighbor defn
32  type(methodcellpollocktype), pointer :: method_cell_plck => null()
33  type(methodcellpollockquadtype), pointer :: method_cell_quad => null()
34  type(methodcellternarytype), pointer :: method_cell_tern => null()
35  type(methodcellpasstobottype), pointer :: method_cell_ptb => null()
36  contains
37  procedure, public :: apply => apply_disv !< apply the DISV tracking method
38  procedure, public :: deallocate !< deallocate arrays and scalars
39  procedure, public :: load => load_disv !< load the cell method
40  procedure, public :: load_cell_defn !< load cell definition from the grid
41  procedure, public :: pass => pass_disv !< pass the particle to the next cell
42  procedure :: map_neighbor !< map a location on the cell face to the shared face of a neighbor
43  procedure :: update_flowja !< update intercell mass flows
44  procedure :: load_particle !< load particle properties
45  procedure :: load_cell_properties !< load cell properties
46  procedure :: load_cell_polygon !< load cell polygon
47  procedure :: load_cell_neighbors !< load cell face neighbors
48  procedure :: load_cell_flags !< load cell 180-degree vertex indicator
49  procedure :: load_cell_flows !< load the cell's flows
50  procedure :: load_cell_boundary_flows !< load boundary flows to a polygonal cell definition
51  procedure :: load_cell_face_flows !< load face flows to a polygonal cell definition
52  end type methoddisvtype
53 
54 contains
55 
56  !> @brief Create a new vertex grid (DISV) tracking method
57  subroutine create_method_disv(method)
58  ! dummy
59  type(methoddisvtype), pointer :: method
60  ! local
61  type(cellpolytype), pointer :: cell
62 
63  allocate (method)
64  allocate (method%name)
65  call create_cell_poly(cell)
66  method%cell => cell
67  method%name = "disv"
68  method%delegates = .true.
69  call create_defn(method%neighbor)
70 
71  ! Create cell methods this method can delegate to
72  call create_method_cell_pollock(method%method_cell_plck)
73  call create_method_cell_quad(method%method_cell_quad)
74  call create_method_cell_ternary(method%method_cell_tern)
75  call create_method_cell_ptb(method%method_cell_ptb)
76  end subroutine create_method_disv
77 
78  !> @brief Destroy the tracking method
79  subroutine deallocate (this)
80  class(methoddisvtype), intent(inout) :: this
81  deallocate (this%name)
82 
83  ! Deallocate owned cell methods
84  call this%method_cell_plck%deallocate()
85  deallocate (this%method_cell_plck)
86  call this%method_cell_quad%deallocate()
87  deallocate (this%method_cell_quad)
88  call this%method_cell_tern%deallocate()
89  deallocate (this%method_cell_tern)
90  call this%method_cell_ptb%deallocate()
91  deallocate (this%method_cell_ptb)
92  end subroutine deallocate
93 
94  !> @brief Load the cell and the tracking method
95  subroutine load_disv(this, particle, next_level, submethod)
96  use cellmodule
97  use cellrectmodule
99  use cellutilmodule
100  ! dummy
101  class(methoddisvtype), intent(inout) :: this
102  type(particletype), pointer, intent(inout) :: particle
103  integer(I4B), intent(in) :: next_level
104  class(methodtype), pointer, intent(inout) :: submethod
105  ! local
106  integer(I4B) :: ic
107  class(celltype), pointer :: base
108  type(cellrecttype), pointer :: rect
109  type(cellrectquadtype), pointer :: quad
110 
111  select type (cell => this%cell)
112  type is (cellpolytype)
113  ic = particle%itrdomain(next_level)
114  call this%load_cell_defn(ic, cell%defn)
115  if (this%fmi%ibdgwfsat0(ic) == 0) then
116  ! Cell is active but dry, so select and initialize pass-to-bottom
117  ! cell method and set cell method pointer
118  call this%method_cell_ptb%init( &
119  fmi=this%fmi, &
120  cell=this%cell, &
121  events=this%events, &
122  tracktimes=this%tracktimes)
123  submethod => this%method_cell_ptb
124  else if (particle%frctrn) then
125  ! Force the ternary method
126  call this%method_cell_tern%init( &
127  fmi=this%fmi, &
128  cell=this%cell, &
129  events=this%events, &
130  tracktimes=this%tracktimes)
131  submethod => this%method_cell_tern
132  else if (cell%defn%can_be_rect) then
133  ! Cell is a rectangle, convert it to a rectangular cell type and
134  ! initialize Pollock's method
135  call cell_poly_to_rect(cell, rect)
136  base => rect
137  call this%method_cell_plck%init( &
138  fmi=this%fmi, &
139  cell=base, &
140  events=this%events, &
141  tracktimes=this%tracktimes)
142  submethod => this%method_cell_plck
143  else if (cell%defn%can_be_quad) then
144  ! Cell is quad-refined, convert to a quad rect cell type and
145  ! initialize the corresponding method
146  call cell_poly_to_quad(cell, quad)
147  base => quad
148  call this%method_cell_quad%init( &
149  fmi=this%fmi, &
150  cell=base, &
151  events=this%events, &
152  tracktimes=this%tracktimes)
153  submethod => this%method_cell_quad
154  else
155  ! Default to the ternary method
156  call this%method_cell_tern%init( &
157  fmi=this%fmi, &
158  cell=this%cell, &
159  events=this%events, &
160  tracktimes=this%tracktimes)
161  submethod => this%method_cell_tern
162  end if
163  end select
164  end subroutine load_disv
165 
166  subroutine load_particle(this, cell, particle)
167  ! modules
168  use disvmodule, only: disvtype
169  use particlemodule, only: term_boundary
170  use particleeventmodule, only: terminate
171  ! dummy
172  class(methoddisvtype), intent(inout) :: this
173  type(cellpolytype), pointer, intent(inout) :: cell
174  type(particletype), pointer, intent(inout) :: particle
175  ! local
176  integer(I4B) :: inface
177  integer(I4B) :: ipos
178  integer(I4B) :: ic
179  integer(I4B) :: icu
180  integer(I4B) :: inbr
181  integer(I4B) :: idiag
182  integer(I4B) :: icpl
183  integer(I4B) :: ilay
184  real(DP) :: z
185 
186  select type (dis => this%fmi%dis)
187  type is (disvtype)
188  inface = particle%iboundary(level_feature)
189  idiag = dis%con%ia(cell%defn%icell)
190  inbr = cell%defn%facenbr(inface)
191  ipos = idiag + inbr
192  ic = dis%con%ja(ipos)
193  icu = dis%get_nodeuser(ic)
194  call get_jk(icu, dis%ncpl, dis%nlay, icpl, ilay)
195 
196  particle%itrdomain(level_feature) = ic
197  particle%icu = icu
198  particle%ilay = ilay
199 
200  z = particle%z
201  call this%map_neighbor(cell%defn, inface, z)
202 
203  particle%iboundary(level_feature) = inface
204  particle%itrdomain(level_subfeature:) = 0
205  particle%iboundary(level_subfeature:) = 0
206  particle%z = z
207  end select
208 
209  end subroutine
210 
211  subroutine update_flowja(this, cell, particle)
212  ! dummy
213  class(methoddisvtype), intent(inout) :: this
214  type(cellpolytype), pointer, intent(inout) :: cell
215  type(particletype), pointer, intent(inout) :: particle
216  ! local
217  integer(I4B) :: inbr
218  integer(I4B) :: idiag
219  integer(I4B) :: ipos
220 
221  idiag = this%fmi%dis%con%ia(cell%defn%icell)
222  inbr = cell%defn%facenbr(particle%iboundary(level_feature))
223  ipos = idiag + inbr
224 
225  ! leaving old cell
226  this%flowja(ipos) = this%flowja(ipos) - done
227 
228  ! entering new cell
229  this%flowja(this%fmi%dis%con%isym(ipos)) &
230  = this%flowja(this%fmi%dis%con%isym(ipos)) + done
231 
232  end subroutine update_flowja
233 
234  !> @brief Pass a particle to the next cell, if there is one
235  subroutine pass_disv(this, particle)
236  use particlemodule, only: term_boundary
237  use particleeventmodule, only: terminate
238  ! dummy
239  class(methoddisvtype), intent(inout) :: this
240  type(particletype), pointer, intent(inout) :: particle
241  ! local
242  type(cellpolytype), pointer :: cell
243  integer(I4B) :: iface
244  logical(LGP) :: no_neighbors
245 
246  select type (c => this%cell)
247  type is (cellpolytype)
248  cell => c
249  iface = particle%iboundary(level_feature)
250  no_neighbors = cell%defn%facenbr(iface) == 0
251 
252  ! todo AMP: reconsider when multiple models supported
253  if (no_neighbors) then
254  call this%terminate(particle, status=term_boundary)
255  return
256  end if
257 
258  call this%load_particle(cell, particle)
259  call this%update_flowja(cell, particle)
260  end select
261  end subroutine pass_disv
262 
263  !> @brief Map location on cell face to shared cell face of neighbor
264  subroutine map_neighbor(this, defn, inface, z)
265  ! dummy
266  class(methoddisvtype), intent(inout) :: this
267  type(celldefntype), pointer, intent(inout) :: defn
268  integer(I4B), intent(inout) :: inface
269  double precision, intent(inout) :: z
270  ! local
271  integer(I4B) :: icin
272  integer(I4B) :: npolyvertsin
273  integer(I4B) :: ic
274  integer(I4B) :: npolyverts
275  integer(I4B) :: inbr
276  integer(I4B) :: inbrnbr
277  integer(I4B) :: j
278  integer(I4B) :: m
279  real(DP) :: zrel
280  real(DP) :: topfrom
281  real(DP) :: botfrom
282  real(DP) :: top
283  real(DP) :: bot
284  real(DP) :: sat
285 
286  ! Map to shared cell face of neighbor
287  inbr = defn%facenbr(inface)
288  if (inbr .eq. 0) then
289  ! Exterior face; no neighbor to map to
290  ! todo AMP: reconsider when multiple models allowed
291  inface = -1
292  else
293  ! Load definition for neighbor cell (neighbor with shared face)
294  icin = defn%icell
295  j = this%fmi%dis%con%ia(icin)
296  ic = this%fmi%dis%con%ja(j + inbr)
297  call this%load_cell_defn(ic, this%neighbor)
298 
299  npolyvertsin = defn%npolyverts
300  npolyverts = this%neighbor%npolyverts
301  if (inface .eq. npolyvertsin + 2) then
302  ! Exits through bot, enters through top
303  inface = npolyverts + 3
304  else if (inface .eq. npolyvertsin + 3) then
305  ! Exits through top, enters through bot
306  inface = npolyverts + 2
307  else
308  ! Exits and enters through shared polygon face
309  j = this%fmi%dis%con%ia(ic)
310  do m = 1, npolyverts + 3
311  inbrnbr = this%neighbor%facenbr(m)
312  if (this%fmi%dis%con%ja(j + inbrnbr) .eq. icin) then
313  inface = m
314  exit
315  end if
316  end do
317  ! Map z between cells
318  topfrom = defn%top
319  botfrom = defn%bot
320  zrel = (z - botfrom) / (topfrom - botfrom)
321  top = this%fmi%dis%top(ic)
322  bot = this%fmi%dis%bot(ic)
323  sat = this%fmi%gwfsat(ic)
324  z = bot + zrel * sat * (top - bot)
325  end if
326  end if
327  end subroutine map_neighbor
328 
329  !> @brief Apply the DISV tracking method to a particle.
330  subroutine apply_disv(this, particle, tmax)
331  class(methoddisvtype), intent(inout) :: this
332  type(particletype), pointer, intent(inout) :: particle
333  real(DP), intent(in) :: tmax
334 
335  call this%track(particle, 1, tmax)
336  end subroutine apply_disv
337 
338  !> @brief Load cell definition from the grid
339  subroutine load_cell_defn(this, ic, defn)
340  ! dummy
341  class(methoddisvtype), intent(inout) :: this
342  integer(I4B), intent(in) :: ic
343  type(celldefntype), pointer, intent(inout) :: defn
344 
345  call this%load_cell_properties(ic, defn)
346  call this%load_cell_polygon(defn)
347  call this%load_cell_neighbors(defn)
348  call this%load_cell_saturation_status(defn)
349  call this%load_cell_flags(defn)
350  call this%load_cell_flows(defn)
351  end subroutine load_cell_defn
352 
353  !> @brief Loads cell properties to cell definition from the grid.
354  subroutine load_cell_properties(this, ic, defn)
355  ! dummy
356  class(methoddisvtype), intent(inout) :: this
357  integer(I4B), intent(in) :: ic
358  type(celldefntype), pointer, intent(inout) :: defn
359  ! local
360  real(DP) :: top
361  real(DP) :: bot
362  real(DP) :: sat
363  integer(I4B) :: icu, icpl, ilay
364 
365  defn%icell = ic
366  defn%iatop = get_iatop(this%fmi%dis%get_ncpl(), &
367  this%fmi%dis%get_nodeuser(ic))
368  top = this%fmi%dis%top(ic)
369  bot = this%fmi%dis%bot(ic)
370  sat = this%fmi%gwfsat(ic)
371  top = bot + sat * (top - bot)
372  defn%top = top
373  defn%bot = bot
374  defn%sat = sat
375  defn%porosity = this%porosity(ic)
376  defn%retfactor = this%retfactor(ic)
377  defn%izone = this%izone(ic)
378  select type (dis => this%fmi%dis)
379  type is (disvtype)
380  icu = dis%get_nodeuser(ic)
381  call get_jk(icu, dis%ncpl, dis%nlay, icpl, ilay)
382  defn%ilay = ilay
383  end select
384 
385  end subroutine load_cell_properties
386 
387  subroutine load_cell_polygon(this, defn)
388  ! dummy
389  class(methoddisvtype), intent(inout) :: this
390  type(celldefntype), pointer, intent(inout) :: defn
391 
392  call this%fmi%dis%get_polyverts( &
393  defn%icell, &
394  defn%polyvert, &
395  closed=.true.)
396  defn%npolyverts = size(defn%polyvert, dim=2) - 1
397  end subroutine load_cell_polygon
398 
399  !> @brief Loads face neighbors to cell definition from the grid
400  !! Assumes cell index and number of vertices are already loaded.
401  subroutine load_cell_neighbors(this, defn)
402  ! dummy
403  class(methoddisvtype), intent(inout) :: this
404  type(celldefntype), pointer, intent(inout) :: defn
405  ! local
406  integer(I4B) :: ic1
407  integer(I4B) :: ic2
408  integer(I4B) :: icu1
409  integer(I4B) :: icu2
410  integer(I4B) :: j1
411  integer(I4B) :: j2
412  integer(I4B) :: k1
413  integer(I4B) :: k2
414  integer(I4B) :: iloc
415  integer(I4B) :: ipos
416  integer(I4B) :: istart1
417  integer(I4B) :: istart2
418  integer(I4B) :: istop1
419  integer(I4B) :: istop2
420  integer(I4B) :: isharedface
421  integer(I4B) :: ncpl
422  integer(I4B) :: nfaces
423  integer(I4B) :: nslots
424 
425  ! expand facenbr array if needed
426  nfaces = defn%npolyverts + 3
427  nslots = size(defn%facenbr)
428  if (nslots < nfaces) call expandarray(defn%facenbr, nfaces - nslots)
429 
430  select type (dis => this%fmi%dis)
431  type is (disvtype)
432  ! Load face neighbors.
433  defn%facenbr = 0
434  ic1 = defn%icell
435  icu1 = dis%get_nodeuser(ic1)
436  ncpl = dis%get_ncpl()
437  call get_jk(icu1, ncpl, dis%nlay, j1, k1)
438  istart1 = dis%iavert(j1)
439  istop1 = dis%iavert(j1 + 1) - 1
440  do iloc = 1, dis%con%ia(ic1 + 1) - dis%con%ia(ic1) - 1
441  ipos = dis%con%ia(ic1) + iloc
442  ! mask could become relevant if PRT uses interface model
443  if (dis%con%mask(ipos) == 0) cycle
444  ic2 = dis%con%ja(ipos)
445  icu2 = dis%get_nodeuser(ic2)
446  call get_jk(icu2, ncpl, dis%nlay, j2, k2)
447  istart2 = dis%iavert(j2)
448  istop2 = dis%iavert(j2 + 1) - 1
449  call shared_face(dis%javert(istart1:istop1), &
450  dis%javert(istart2:istop2), &
451  isharedface)
452  if (isharedface /= 0) then
453  ! Edge (polygon) face neighbor
454  defn%facenbr(isharedface) = int(iloc, 1)
455  else
456  if (k2 > k1) then
457  ! Bottom face neighbor
458  defn%facenbr(defn%npolyverts + 2) = int(iloc, 1)
459  else if (k2 < k1) then
460  ! Top face neighbor
461  defn%facenbr(defn%npolyverts + 3) = int(iloc, 1)
462  else
463  call pstop(1, "k2 should be <> k1, since no shared edge face")
464  end if
465  end if
466  end do
467  end select
468  ! List of edge (polygon) faces wraps around
469  defn%facenbr(defn%npolyverts + 1) = defn%facenbr(1)
470  end subroutine load_cell_neighbors
471 
472  !> @brief Load flows into the cell definition.
473  !! These include face, boundary and net distributed flows.
474  !! Assumes cell index and number of vertices are already loaded.
475  subroutine load_cell_flows(this, defn)
476  ! dummy
477  class(methoddisvtype), intent(inout) :: this
478  type(celldefntype), pointer, intent(inout) :: defn
479  ! local
480  integer(I4B) :: nfaces, nslots
481 
482  ! expand faceflow array if needed
483  nfaces = defn%npolyverts + 3
484  nslots = size(defn%faceflow)
485  if (nslots < nfaces) call expandarray(defn%faceflow, nfaces - nslots)
486 
487  ! Load face flows, including boundary flows. As with cell verts,
488  ! the face flow array wraps around. Top and bottom flows make up
489  ! the last two elements, respectively, for size npolyverts + 3.
490  ! If there is no flow through any face, set a no-exit-face flag.
491  defn%faceflow = dzero
492  call this%load_cell_boundary_flows(defn)
493  call this%load_cell_face_flows(defn)
494  call this%cap_cell_wt_flow(defn)
495  call this%load_cell_no_exit_face(defn)
496 
497  ! Add up net distributed flow
498  defn%distflow = this%fmi%SourceFlows(defn%icell) + &
499  this%fmi%SinkFlows(defn%icell) + &
500  this%fmi%StorageFlows(defn%icell)
501 
502  ! Set weak sink flag
503  if (this%fmi%SinkFlows(defn%icell) .ne. dzero) then
504  defn%iweaksink = 1
505  else
506  defn%iweaksink = 0
507  end if
508  end subroutine load_cell_flows
509 
510  subroutine load_cell_face_flows(this, defn)
511  ! dummy
512  class(methoddisvtype), intent(inout) :: this
513  type(celldefntype), pointer, intent(inout) :: defn
514  ! local
515  integer(I4B) :: m, n, nfaces
516  real(DP) :: q
517 
518  nfaces = defn%npolyverts + 3
519  do m = 1, nfaces
520  n = defn%facenbr(m)
521  if (n > 0) then
522  q = this%fmi%gwfflowja(this%fmi%dis%con%ia(defn%icell) + n)
523  defn%faceflow(m) = defn%faceflow(m) + q
524  end if
525  end do
526  end subroutine load_cell_face_flows
527 
528  !> @brief Load boundary flows from the grid into a polygonal cell.
529  !! Assumes cell index and number of vertices are already loaded.
530  subroutine load_cell_boundary_flows(this, defn)
531  ! dummy
532  class(methoddisvtype), intent(inout) :: this
533  type(celldefntype), pointer, intent(inout) :: defn
534 
535  ! local
536  integer(I4B) :: ic, iv, npolyverts
537 
538  ic = defn%icell
539  npolyverts = defn%npolyverts
540  do iv = 1, npolyverts
541  defn%faceflow(iv) = &
542  defn%faceflow(iv) + &
543  this%fmi%BoundaryFlows(ic, iv)
544  end do
545  defn%faceflow(npolyverts + 1) = defn%faceflow(1)
546  defn%faceflow(npolyverts + 2) = &
547  defn%faceflow(npolyverts + 2) + &
548  this%fmi%BoundaryFlows(ic, this%fmi%max_faces - 1)
549  defn%faceflow(npolyverts + 3) = &
550  defn%faceflow(npolyverts + 3) + &
551  this%fmi%BoundaryFlows(ic, this%fmi%max_faces)
552 
553  end subroutine load_cell_boundary_flows
554 
555  !> @brief Load 180-degree vertex indicator array and set flags
556  !! indicating how cell can be represented. Assumes cell index
557  !! and number of vertices are already loaded.
558  !<
559  subroutine load_cell_flags(this, defn)
560  ! dummy
561  class(methoddisvtype), intent(inout) :: this
562  type(celldefntype), pointer, intent(inout) :: defn
563  ! local
564  integer(I4B) :: npolyverts
565  integer(I4B) :: m
566  integer(I4B) :: m0
567  integer(I4B) :: m1
568  integer(I4B) :: m2
569  integer(I4B) :: ic
570  integer(I4B) :: num90
571  integer(I4B) :: num180
572  real(DP) :: x0
573  real(DP) :: y0
574  real(DP) :: x1
575  real(DP) :: y1
576  real(DP) :: x2
577  real(DP) :: y2
578  real(DP) :: epsang
579  real(DP) :: s0x
580  real(DP) :: s0y
581  real(DP) :: &
582  s0mag, s2x, s2y, s2mag, sinang, cosang, dotprod
583  logical(LGP) last180
584 
585  ic = defn%icell
586  npolyverts = defn%npolyverts
587 
588  if (size(defn%ispv180) < npolyverts + 3) &
589  call expandarray(defn%ispv180, npolyverts + 1)
590 
591  defn%ispv180(1:npolyverts + 1) = .false.
592  defn%can_be_rect = .false.
593  defn%can_be_quad = .false.
594  epsang = 1d-5
595  num90 = 0
596  num180 = 0
597  last180 = .false.
598 
599  ! assumes non-self-intersecting polygon
600  do m = 1, npolyverts
601  m1 = m
602  if (m1 .eq. 1) then
603  m0 = npolyverts
604  m2 = 2
605  else if (m1 .eq. npolyverts) then
606  m0 = npolyverts - 1
607  m2 = 1
608  else
609  m0 = m1 - 1
610  m2 = m1 + 1
611  end if
612  x0 = defn%polyvert(1, m0)
613  y0 = defn%polyvert(2, m0)
614  x1 = defn%polyvert(1, m1)
615  y1 = defn%polyvert(2, m1)
616  x2 = defn%polyvert(1, m2)
617  y2 = defn%polyvert(2, m2)
618  s0x = x0 - x1
619  s0y = y0 - y1
620  s0mag = dsqrt(s0x * s0x + s0y * s0y)
621  s2x = x2 - x1
622  s2y = y2 - y1
623  s2mag = dsqrt(s2x * s2x + s2y * s2y)
624  sinang = (s0x * s2y - s0y * s2x) / (s0mag * s2mag)
625  cosang = dsqrt(dabs(done - (sinang * sinang)))
626  if (dabs(sinang) .lt. epsang) then
627  dotprod = s0x * s2x + s0y * s2y
628  if (dotprod .lt. dzero) then
629  num180 = num180 + 1
630  last180 = .true.
631  defn%ispv180(m) = .true.
632  end if
633  else
634  if (dabs(cosang) .lt. epsang) num90 = num90 + 1
635  last180 = .false.
636  end if
637  end do
638 
639  ! List of 180-degree indicators wraps around for convenience
640  defn%ispv180(npolyverts + 1) = defn%ispv180(1)
641 
642  ! Set rect/quad flags
643  if (num90 .eq. 4) then
644  if (num180 .eq. 0) then
645  defn%can_be_rect = .true.
646  else
647  defn%can_be_quad = .true.
648  end if
649  end if
650  end subroutine load_cell_flags
651 
652 end module methoddisvmodule
subroutine, public create_defn(cellDefn)
Create a new cell definition object.
Definition: CellDefn.f90:58
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_poly(cell)
Create a new polygonal cell.
Definition: CellPoly.f90:20
subroutine, public cell_poly_to_rect(poly, rect)
Convert CellPoly representation to CellRect. Assumes the conversion is possible.
Definition: CellUtil.f90:14
subroutine, public cell_poly_to_quad(poly, quad)
Convert CellPoly representation to CellRectQuad. Assumes the conversion is possible.
Definition: CellUtil.f90:115
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
subroutine pstop(status, message)
Stop the program, optionally specifying an error status code.
Definition: ErrorUtil.f90:24
subroutine, public shared_face(iverts1, iverts2, iface)
Find the lateral face shared by two cells.
Definition: GeomUtil.f90:408
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
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.
procedure subroutine, public create_method_cell_quad(method)
Create a new Pollock quad-refined cell method.
subroutine, public create_method_cell_ternary(method)
Create a tracking method.
subroutine, public create_method_disv(method)
Create a new vertex grid (DISV) tracking method.
Definition: MethodDisv.f90:58
subroutine load_cell_boundary_flows(this, defn)
Load boundary flows from the grid into a polygonal cell. Assumes cell index and number of vertices ar...
Definition: MethodDisv.f90:531
subroutine load_cell_polygon(this, defn)
Definition: MethodDisv.f90:388
subroutine load_cell_flows(this, defn)
Load flows into the cell definition. These include face, boundary and net distributed flows....
Definition: MethodDisv.f90:476
subroutine load_cell_flags(this, defn)
Load 180-degree vertex indicator array and set flags indicating how cell can be represented....
Definition: MethodDisv.f90:560
subroutine map_neighbor(this, defn, inface, z)
Map location on cell face to shared cell face of neighbor.
Definition: MethodDisv.f90:265
subroutine load_cell_properties(this, ic, defn)
Loads cell properties to cell definition from the grid.
Definition: MethodDisv.f90:355
subroutine update_flowja(this, cell, particle)
Definition: MethodDisv.f90:212
subroutine load_cell_face_flows(this, defn)
Definition: MethodDisv.f90:511
subroutine load_cell_neighbors(this, defn)
Loads face neighbors to cell definition from the grid Assumes cell index and number of vertices are a...
Definition: MethodDisv.f90:402
subroutine load_disv(this, particle, next_level, submethod)
Load the cell and the tracking method.
Definition: MethodDisv.f90:96
subroutine load_particle(this, cell, particle)
Definition: MethodDisv.f90:167
subroutine load_cell_defn(this, ic, defn)
Load cell definition from the grid.
Definition: MethodDisv.f90:340
subroutine apply_disv(this, particle, tmax)
Apply the DISV tracking method to a particle.
Definition: MethodDisv.f90:331
subroutine pass_disv(this, particle)
Pass a particle to the next cell, if there is one.
Definition: MethodDisv.f90:236
Particle tracking strategies.
Definition: Method.f90:2
@, public level_feature
Definition: Method.f90:41
@, public level_subfeature
Definition: Method.f90:42
@, public terminate
particle terminated
@ term_boundary
terminated at a boundary face
Definition: Particle.f90:31
Base grid cell definition.
Definition: CellDefn.f90:25
Base type for grid cells of a concrete type. Contains a cell-definition which is information shared b...
Definition: Cell.f90:12
Vertex grid discretization.
Definition: Disv.f90:24
Base type for particle tracking methods.
Definition: Method.f90:58
Particle tracked by the PRT model.
Definition: Particle.f90:62