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 errorutilmodule, only: pstop
6  use subcellmodule, only: subcelltype
8  use basedismodule, only: disbasetype
9  use prtfmimodule, only: prtfmitype
10  use cellmodule, only: celltype
11  use celldefnmodule, only: celldefntype
14  implicit none
15 
16  private
17  public :: methodtype
18 
19  !> @brief Base type for particle tracking methods.
20  !!
21  !! The PRT tracking algorithm invokes a "tracking method" for each
22  !! domain. A domain can be a model, cell in a model, or subcell in
23  !! a cell. Tracking proceeds recursively, delegating to a possibly
24  !! arbitrary number of subdomains (currently, only the three above
25  !! are recognized). A tracking method is responsible for advancing
26  !! a particle through a domain, delegating to subdomains as needed
27  !! depending on cell geometry (implementing the strategy pattern).
28  !<
29  type, abstract :: methodtype
30  character(len=40), pointer, public :: type !< method name
31  logical(LGP), public :: delegates !< whether the method delegates
32  type(prtfmitype), pointer, public :: fmi => null() !< ptr to fmi
33  class(celltype), pointer, public :: cell => null() !< ptr to the current cell
34  class(subcelltype), pointer, public :: subcell => null() !< ptr to the current subcell
35  type(trackfilecontroltype), pointer, public :: trackfilectl => null() !< ptr to track file control
36  type(timeselecttype), pointer, public :: tracktimes => null() !< ptr to user-defined tracking times
37  integer(I4B), dimension(:), pointer, contiguous, public :: izone => null() !< pointer to zone numbers
38  real(dp), dimension(:), pointer, contiguous, public :: flowja => null() !< pointer to intercell flows
39  real(dp), dimension(:), pointer, contiguous, public :: porosity => null() !< pointer to aquifer porosity
40  real(dp), dimension(:), pointer, contiguous, public :: retfactor => null() !< pointer to retardation factor
41  contains
42  ! Implemented in all subtypes
43  procedure(apply), deferred :: apply
44  procedure(deallocate), deferred :: deallocate
45  ! Overridden in subtypes that delegate
46  procedure :: pass
47  procedure :: load
48  ! Implemented here
49  procedure :: init
50  procedure :: save
51  procedure :: track
52  procedure :: try_pass
53  procedure :: update
54  end type methodtype
55 
56  abstract interface
57  subroutine apply(this, particle, tmax)
58  import dp
59  import methodtype
60  import particletype
61  class(methodtype), intent(inout) :: this
62  type(particletype), pointer, intent(inout) :: particle
63  real(DP), intent(in) :: tmax
64  end subroutine apply
65  subroutine deallocate (this)
66  import methodtype
67  class(methodtype), intent(inout) :: this
68  end subroutine deallocate
69  end interface
70 
71 contains
72 
73  subroutine init(this, fmi, cell, subcell, trackfilectl, tracktimes, &
74  izone, flowja, porosity, retfactor)
75  class(methodtype), intent(inout) :: this
76  type(prtfmitype), intent(in), pointer, optional :: fmi
77  class(celltype), intent(in), pointer, optional :: cell
78  class(subcelltype), intent(in), pointer, optional :: subcell
79  type(trackfilecontroltype), intent(in), pointer, optional :: trackfilectl
80  type(timeselecttype), intent(in), pointer, optional :: tracktimes
81  integer(I4B), intent(in), pointer, optional :: izone(:)
82  real(DP), intent(in), pointer, optional :: flowja(:)
83  real(DP), intent(in), pointer, optional :: porosity(:)
84  real(DP), intent(in), pointer, optional :: retfactor(:)
85 
86  if (present(fmi)) this%fmi => fmi
87  if (present(cell)) this%cell => cell
88  if (present(subcell)) this%subcell => subcell
89  if (present(trackfilectl)) this%trackfilectl => trackfilectl
90  if (present(tracktimes)) this%tracktimes => tracktimes
91  if (present(izone)) this%izone => izone
92  if (present(flowja)) this%flowja => flowja
93  if (present(porosity)) this%porosity => porosity
94  if (present(retfactor)) this%retfactor => retfactor
95  end subroutine init
96 
97  !> @brief Track particle through subdomains
98  recursive subroutine track(this, particle, level, tmax)
99  ! dummy
100  class(methodtype), intent(inout) :: this
101  type(particletype), pointer, intent(inout) :: particle
102  integer(I4B) :: level
103  real(dp), intent(in) :: tmax
104  ! local
105  logical(LGP) :: advancing
106  integer(I4B) :: nextlevel
107  class(methodtype), pointer :: submethod
108 
109  advancing = .true.
110  nextlevel = level + 1
111  do while (advancing)
112  call this%load(particle, nextlevel, submethod)
113  call submethod%apply(particle, tmax)
114  call this%try_pass(particle, nextlevel, advancing)
115  end do
116  end subroutine track
117 
118  !> @brief Try passing the particle to the next subdomain
119  subroutine try_pass(this, particle, nextlevel, advancing)
120  class(methodtype), intent(inout) :: this
121  type(particletype), pointer, intent(inout) :: particle
122  integer(I4B) :: nextlevel
123  logical(LGP) :: advancing
124 
125  ! tracking submethod marked tracking complete?
126  ! reset domain boundary flag and don't advance
127  if (.not. particle%advancing) then
128  particle%iboundary = 0
129  advancing = .false.
130  else
131  ! otherwise pass particle to next subdomain
132  ! and if it's on a boundary, stop advancing
133  call this%pass(particle)
134  if (particle%iboundary(nextlevel - 1) .ne. 0) &
135  advancing = .false.
136  end if
137  end subroutine try_pass
138 
139  !> @brief Load subdomain tracking method (submethod)
140  subroutine load(this, particle, next_level, submethod)
141  class(methodtype), intent(inout) :: this
142  type(particletype), pointer, intent(inout) :: particle
143  integer, intent(in) :: next_level
144  class(methodtype), pointer, intent(inout) :: submethod
145  call pstop(1, "load must be overridden")
146  end subroutine load
147 
148  !> @brief Pass a particle to the next subdomain, internal use only
149  subroutine pass(this, particle)
150  class(methodtype), intent(inout) :: this
151  type(particletype), pointer, intent(inout) :: particle
152  call pstop(1, "pass must be overridden")
153  end subroutine pass
154 
155  !> @brief Save a particle's current state.
156  subroutine save(this, particle, reason)
157  use tdismodule, only: kper, kstp, totimc
158  ! dummy
159  class(methodtype), intent(inout) :: this
160  type(particletype), pointer, intent(inout) :: particle
161  integer(I4B), intent(in) :: reason
162  ! local
163  integer(I4B) :: per, stp
164 
165  per = kper
166  stp = kstp
167 
168  ! If tracking time falls exactly on a boundary between time steps,
169  ! report the previous time step for this datum. This is to follow
170  ! MP7's behavior, and because the particle will have been tracked
171  ! up to this instant under the previous time step's conditions, so
172  ! the time step we're about to start shouldn't get "credit" for it.
173  if (particle%ttrack == totimc .and. (per > 1 .or. stp > 1)) then
174  if (stp > 1) then
175  stp = stp - 1
176  else if (per > 1) then
177  per = per - 1
178  stp = 1
179  end if
180  end if
181 
182  ! Save the particle's state to any registered tracking output files
183  call this%trackfilectl%save(particle, kper=per, &
184  kstp=stp, reason=reason)
185  end subroutine save
186 
187  !> @brief Update particle state and check termination conditions
188  !!
189  !! Update the particle's properties (e.g. advancing flag, zone number,
190  !! status). If any termination conditions apply, the particle's status
191  !! will be set to the appropriate termination value. If any reporting
192  !! conditions apply, save particle state with the proper reason code.
193  !<
194  subroutine update(this, particle, cell_defn)
195  ! -- dummy
196  class(methodtype), intent(inout) :: this
197  type(particletype), pointer, intent(inout) :: particle
198  type(celldefntype), pointer, intent(inout) :: cell_defn
199 
200  particle%izone = cell_defn%izone
201  if (cell_defn%izone .ne. 0) then
202  if (particle%istopzone .eq. cell_defn%izone) then
203  particle%advancing = .false.
204  particle%istatus = 6
205  call this%save(particle, reason=3) ! particle terminated
206  return
207  end if
208  end if
209  if (cell_defn%inoexitface .ne. 0) then
210  particle%advancing = .false.
211  particle%istatus = 5
212  call this%save(particle, reason=3) ! particle terminated
213  return
214  end if
215  if (cell_defn%iweaksink .ne. 0) then
216  if (particle%istopweaksink .ne. 0) then
217  particle%advancing = .false.
218  particle%istatus = 3
219  call this%save(particle, reason=3) ! particle terminated
220  return
221  else
222  call this%save(particle, reason=4) ! particle exited weak sink
223  return
224  end if
225  end if
226  end subroutine update
227 
228 end module methodmodule
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
Particle tracking strategies.
Definition: Method.f90:2
subroutine load(this, particle, next_level, submethod)
Load subdomain tracking method (submethod)
Definition: Method.f90:141
subroutine save(this, particle, reason)
Save a particle's current state.
Definition: Method.f90:157
recursive subroutine track(this, particle, level, tmax)
Track particle through subdomains.
Definition: Method.f90:99
subroutine update(this, particle, cell_defn)
Update particle state and check termination conditions.
Definition: Method.f90:195
subroutine try_pass(this, particle, nextlevel, advancing)
Try passing the particle to the next subdomain.
Definition: Method.f90:120
subroutine pass(this, particle)
Pass a particle to the next subdomain, internal use only.
Definition: Method.f90:150
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:10
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:29
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.
Definition: TrackData.f90:38