MODFLOW 6  version 6.6.0.dev0
USGS Modular Hydrologic Model
Method.f90
Go to the documentation of this file.
1 !> @brief Particle tracking strategies
3 
4  use kindmodule, only: dp, i4b, lgp
5  use constantsmodule, only: dzero
6  use errorutilmodule, only: pstop
7  use subcellmodule, only: subcelltype
9  use basedismodule, only: disbasetype
10  use prtfmimodule, only: prtfmitype
11  use cellmodule, only: celltype
12  use celldefnmodule, only: celldefntype
15  use mathutilmodule, only: is_close
16  implicit none
17 
18  private
19  public :: methodtype
20 
21  !> @brief Base type for particle tracking methods.
22  !!
23  !! The PRT tracking algorithm invokes a "tracking method" for each
24  !! domain. A domain can be a model, cell in a model, or subcell in
25  !! a cell. Tracking proceeds recursively, delegating to a possibly
26  !! arbitrary number of subdomains (currently, only the three above
27  !! are recognized). A tracking method is responsible for advancing
28  !! a particle through a domain, delegating to subdomains as needed
29  !! depending on cell geometry (implementing the strategy pattern).
30  !<
31  type, abstract :: methodtype
32  character(len=40), pointer, public :: name !< method name
33  logical(LGP), public :: delegates !< whether the method delegates
34  type(prtfmitype), pointer, public :: fmi => null() !< ptr to fmi
35  class(celltype), pointer, public :: cell => null() !< ptr to the current cell
36  class(subcelltype), pointer, public :: subcell => null() !< ptr to the current subcell
37  type(trackcontroltype), pointer, public :: trackctl => null() !< ptr to track file control
38  type(timeselecttype), pointer, public :: tracktimes => null() !< ptr to user-defined tracking times
39  integer(I4B), dimension(:), pointer, contiguous, public :: izone => null() !< pointer to zone numbers
40  real(dp), dimension(:), pointer, contiguous, public :: flowja => null() !< pointer to intercell flows
41  real(dp), dimension(:), pointer, contiguous, public :: porosity => null() !< pointer to aquifer porosity
42  real(dp), dimension(:), pointer, contiguous, public :: retfactor => null() !< pointer to retardation factor
43  contains
44  ! Implemented in all subtypes
45  procedure(apply), deferred :: apply
46  procedure(deallocate), deferred :: deallocate
47  ! Overridden in subtypes that delegate
48  procedure :: pass
49  procedure :: load
50  ! Implemented here
51  procedure :: init
52  procedure :: save
53  procedure :: track
54  procedure :: try_pass
55  procedure :: check
56  end type methodtype
57 
58  abstract interface
59  subroutine apply(this, particle, tmax)
60  import dp
61  import methodtype
62  import particletype
63  class(methodtype), intent(inout) :: this
64  type(particletype), pointer, intent(inout) :: particle
65  real(DP), intent(in) :: tmax
66  end subroutine apply
67  subroutine deallocate (this)
68  import methodtype
69  class(methodtype), intent(inout) :: this
70  end subroutine deallocate
71  end interface
72 
73 contains
74 
75  subroutine init(this, fmi, cell, subcell, trackctl, tracktimes, &
76  izone, flowja, porosity, retfactor)
77  class(methodtype), intent(inout) :: this
78  type(prtfmitype), intent(in), pointer, optional :: fmi
79  class(celltype), intent(in), pointer, optional :: cell
80  class(subcelltype), intent(in), pointer, optional :: subcell
81  type(trackcontroltype), intent(in), pointer, optional :: trackctl
82  type(timeselecttype), intent(in), pointer, optional :: tracktimes
83  integer(I4B), intent(in), pointer, optional :: izone(:)
84  real(DP), intent(in), pointer, optional :: flowja(:)
85  real(DP), intent(in), pointer, optional :: porosity(:)
86  real(DP), intent(in), pointer, optional :: retfactor(:)
87 
88  if (present(fmi)) this%fmi => fmi
89  if (present(cell)) this%cell => cell
90  if (present(subcell)) this%subcell => subcell
91  if (present(trackctl)) this%trackctl => trackctl
92  if (present(tracktimes)) this%tracktimes => tracktimes
93  if (present(izone)) this%izone => izone
94  if (present(flowja)) this%flowja => flowja
95  if (present(porosity)) this%porosity => porosity
96  if (present(retfactor)) this%retfactor => retfactor
97  end subroutine init
98 
99  !> @brief Track the particle over domains of the given
100  ! level until the particle terminates or tmax elapses.
101  recursive subroutine track(this, particle, level, tmax)
102  ! dummy
103  class(methodtype), intent(inout) :: this
104  type(particletype), pointer, intent(inout) :: particle
105  integer(I4B) :: level
106  real(dp), intent(in) :: tmax
107  ! local
108  logical(LGP) :: advancing
109  integer(I4B) :: nextlevel
110  class(methodtype), pointer :: submethod
111 
112  ! Advance the particle over subdomains
113  advancing = .true.
114  nextlevel = level + 1
115  do while (advancing)
116  call this%load(particle, nextlevel, submethod)
117  call submethod%apply(particle, tmax)
118  call this%try_pass(particle, nextlevel, advancing)
119  end do
120  end subroutine track
121 
122  !> @brief Try passing the particle to the next subdomain.
123  subroutine try_pass(this, particle, nextlevel, advancing)
124  class(methodtype), intent(inout) :: this
125  type(particletype), pointer, intent(inout) :: particle
126  integer(I4B) :: nextlevel
127  logical(LGP) :: advancing
128 
129  ! if the particle is done advancing, reset the domain boundary flag.
130  if (.not. particle%advancing) then
131  particle%iboundary = 0
132  advancing = .false.
133  else
134  ! otherwise pass the particle to the next subdomain.
135  ! if that leaves it on a boundary, stop advancing.
136  call this%pass(particle)
137  if (particle%iboundary(nextlevel - 1) .ne. 0) then
138  advancing = .false.
139  end if
140  end if
141  end subroutine try_pass
142 
143  !> @brief Load the subdomain tracking method (submethod).
144  subroutine load(this, particle, next_level, submethod)
145  class(methodtype), intent(inout) :: this
146  type(particletype), pointer, intent(inout) :: particle
147  integer, intent(in) :: next_level
148  class(methodtype), pointer, intent(inout) :: submethod
149  call pstop(1, "load must be overridden")
150  end subroutine load
151 
152  !> @brief Pass the particle to the next subdomain.
153  subroutine pass(this, particle)
154  class(methodtype), intent(inout) :: this
155  type(particletype), pointer, intent(inout) :: particle
156  call pstop(1, "pass must be overridden")
157  end subroutine pass
158 
159  !> @brief Save the particle's state to output files.
160  subroutine save(this, particle, reason)
161  use tdismodule, only: kper, kstp, totimc
162  ! dummy
163  class(methodtype), intent(inout) :: this
164  type(particletype), pointer, intent(inout) :: particle
165  integer(I4B), intent(in) :: reason
166  ! local
167  integer(I4B) :: per, stp
168 
169  per = kper
170  stp = kstp
171 
172  ! If tracking time falls exactly on a boundary between time steps,
173  ! report the previous time step for this datum. This is to follow
174  ! MP7's behavior, and because the particle will have been tracked
175  ! up to this instant under the previous time step's conditions, so
176  ! the time step we're about to start shouldn't get "credit" for it.
177  if (particle%ttrack == totimc .and. (per > 1 .or. stp > 1)) then
178  if (stp > 1) then
179  stp = stp - 1
180  else if (per > 1) then
181  per = per - 1
182  stp = 1
183  end if
184  end if
185 
186  call this%trackctl%save(particle, kper=per, kstp=stp, reason=reason)
187  end subroutine save
188 
189  !> @brief Check reporting/terminating conditions before tracking.
190  !!
191  !! Check a number of conditions determining whether to continue
192  !! tracking the particle or terminate it, as well as whether to
193  !! record any output data as per selected reporting conditions.
194  !<
195  subroutine check(this, particle, cell_defn)
196  ! modules
197  use tdismodule, only: endofsimulation, totim
198  ! dummy
199  class(methodtype), intent(inout) :: this
200  type(particletype), pointer, intent(inout) :: particle
201  type(celldefntype), pointer, intent(inout) :: cell_defn
202  ! local
203  logical(LGP) :: dry_cell, dry_particle, no_exit_face, stop_zone, weak_sink
204 
205  dry_cell = this%fmi%ibdgwfsat0(cell_defn%icell) == 0
206  dry_particle = particle%z > cell_defn%top
207  no_exit_face = cell_defn%inoexitface > 0
208  stop_zone = cell_defn%izone > 0 .and. particle%istopzone == cell_defn%izone
209  weak_sink = cell_defn%iweaksink > 0
210 
211  particle%izone = cell_defn%izone
212  if (stop_zone) then
213  particle%advancing = .false.
214  particle%istatus = 6
215  call this%save(particle, reason=3)
216  return
217  end if
218 
219  if (no_exit_face .and. .not. dry_cell) then
220  particle%advancing = .false.
221  particle%istatus = 5
222  call this%save(particle, reason=3)
223  return
224  end if
225 
226  if (weak_sink) then
227  if (particle%istopweaksink > 0) then
228  particle%advancing = .false.
229  particle%istatus = 3
230  call this%save(particle, reason=3)
231  return
232  else
233  call this%save(particle, reason=4)
234  end if
235  end if
236 
237  if (dry_cell) then
238  if (particle%idrymeth == 0) then
239  no_exit_face = .false.
240  else if (particle%idrymeth == 1) then
241  ! stop
242  particle%advancing = .false.
243  particle%istatus = 7
244  call this%save(particle, reason=3)
245  return
246  else if (particle%idrymeth == 2) then
247  ! stay
248  no_exit_face = .false.
249  particle%advancing = .false.
250  particle%ttrack = totim
251  ! terminate if last period/step
252  if (endofsimulation) then
253  particle%istatus = 5
254  call this%save(particle, reason=3)
255  return
256  end if
257  call this%save(particle, reason=2)
258  end if
259  else if (dry_particle .and. this%name /= "passtobottom") then
260  ! dry particle
261  if (particle%idrymeth == 0) then
262  ! drop to water table
263  particle%z = cell_defn%top
264  call this%save(particle, reason=1)
265  else if (particle%idrymeth == 1) then
266  ! terminate
267  particle%advancing = .false.
268  particle%istatus = 7
269  call this%save(particle, reason=3)
270  return
271  else if (particle%idrymeth == 2) then
272  ! stay
273  no_exit_face = .false.
274  particle%advancing = .false.
275  return
276  end if
277  end if
278 
279  if (no_exit_face) then
280  particle%advancing = .false.
281  particle%istatus = 5
282  call this%save(particle, reason=3)
283  return
284  end if
285 
286  end subroutine check
287 
288 end module methodmodule
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
subroutine pstop(status, message)
Stop the program, optionally specifying an error status code.
Definition: ErrorUtil.f90:24
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
Particle tracking strategies.
Definition: Method.f90:2
subroutine load(this, particle, next_level, submethod)
Load the subdomain tracking method (submethod).
Definition: Method.f90:145
subroutine save(this, particle, reason)
Save the particle's state to output files.
Definition: Method.f90:161
subroutine check(this, particle, cell_defn)
Check reporting/terminating conditions before tracking.
Definition: Method.f90:196
recursive subroutine track(this, particle, level, tmax)
Track the particle over domains of the given.
Definition: Method.f90:102
subroutine try_pass(this, particle, nextlevel, advancing)
Try passing the particle to the next subdomain.
Definition: Method.f90:124
subroutine pass(this, particle)
Pass the particle to the next subdomain.
Definition: Method.f90:154
logical(lgp), pointer, public endofsimulation
flag indicating end of simulation
Definition: tdis.f90:28
real(dp), pointer, public totim
time relative to start of simulation
Definition: tdis.f90:32
real(dp), pointer, public totimc
simulation time at start of time step
Definition: tdis.f90:33
integer(i4b), pointer, public kstp
current time step number
Definition: tdis.f90:24
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:23
Specify times for some event to occur.
Definition: TimeSelect.f90:2
Base grid cell definition.
Definition: CellDefn.f90:12
Base type for grid cells of a concrete type. Contains a cell-definition which is information shared b...
Definition: Cell.f90:13
Base type for particle tracking methods.
Definition: Method.f90:31
Particle tracked by the PRT model.
Definition: Particle.f90:32
A subcell of a cell.
Definition: Subcell.f90:9
Represents a series of instants at which some event should occur.
Definition: TimeSelect.f90:30
Manages particle track (i.e. pathline) files.