MODFLOW 6  version 6.7.0.dev1
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  end type methodtype
56 
57  abstract interface
58  subroutine apply(this, particle, tmax)
59  import dp
60  import methodtype
61  import particletype
62  class(methodtype), intent(inout) :: this
63  type(particletype), pointer, intent(inout) :: particle
64  real(DP), intent(in) :: tmax
65  end subroutine apply
66  subroutine deallocate (this)
67  import methodtype
68  class(methodtype), intent(inout) :: this
69  end subroutine deallocate
70  end interface
71 
72 contains
73 
74  subroutine init(this, fmi, cell, subcell, trackctl, tracktimes, &
75  izone, flowja, porosity, retfactor)
76  class(methodtype), intent(inout) :: this
77  type(prtfmitype), intent(in), pointer, optional :: fmi
78  class(celltype), intent(in), pointer, optional :: cell
79  class(subcelltype), intent(in), pointer, optional :: subcell
80  type(trackcontroltype), intent(in), pointer, optional :: trackctl
81  type(timeselecttype), intent(in), pointer, optional :: tracktimes
82  integer(I4B), intent(in), pointer, optional :: izone(:)
83  real(DP), intent(in), pointer, optional :: flowja(:)
84  real(DP), intent(in), pointer, optional :: porosity(:)
85  real(DP), intent(in), pointer, optional :: retfactor(:)
86 
87  if (present(fmi)) this%fmi => fmi
88  if (present(cell)) this%cell => cell
89  if (present(subcell)) this%subcell => subcell
90  if (present(trackctl)) this%trackctl => trackctl
91  if (present(tracktimes)) this%tracktimes => tracktimes
92  if (present(izone)) this%izone => izone
93  if (present(flowja)) this%flowja => flowja
94  if (present(porosity)) this%porosity => porosity
95  if (present(retfactor)) this%retfactor => retfactor
96  end subroutine init
97 
98  !> @brief Track the particle over domains of the given
99  ! level until the particle terminates or tmax elapses.
100  recursive subroutine track(this, particle, level, tmax)
101  ! dummy
102  class(methodtype), intent(inout) :: this
103  type(particletype), pointer, intent(inout) :: particle
104  integer(I4B) :: level
105  real(dp), intent(in) :: tmax
106  ! local
107  logical(LGP) :: advancing
108  integer(I4B) :: nextlevel
109  class(methodtype), pointer :: submethod
110 
111  ! Advance the particle over subdomains
112  advancing = .true.
113  nextlevel = level + 1
114  do while (advancing)
115  call this%load(particle, nextlevel, submethod)
116  call submethod%apply(particle, tmax)
117  call this%try_pass(particle, nextlevel, advancing)
118  end do
119  end subroutine track
120 
121  !> @brief Try passing the particle to the next subdomain.
122  subroutine try_pass(this, particle, nextlevel, advancing)
123  class(methodtype), intent(inout) :: this
124  type(particletype), pointer, intent(inout) :: particle
125  integer(I4B) :: nextlevel
126  logical(LGP) :: advancing
127 
128  ! if the particle is done advancing, reset the domain boundary flag.
129  if (.not. particle%advancing) then
130  particle%iboundary = 0
131  advancing = .false.
132  else
133  ! otherwise pass the particle to the next subdomain.
134  ! if that leaves it on a boundary, stop advancing.
135  call this%pass(particle)
136  if (particle%iboundary(nextlevel - 1) .ne. 0) &
137  advancing = .false.
138  end if
139  end subroutine try_pass
140 
141  !> @brief Load the subdomain tracking method (submethod).
142  subroutine load(this, particle, next_level, submethod)
143  class(methodtype), intent(inout) :: this
144  type(particletype), pointer, intent(inout) :: particle
145  integer, intent(in) :: next_level
146  class(methodtype), pointer, intent(inout) :: submethod
147  call pstop(1, "load must be overridden")
148  end subroutine load
149 
150  !> @brief Pass the particle to the next subdomain.
151  subroutine pass(this, particle)
152  class(methodtype), intent(inout) :: this
153  type(particletype), pointer, intent(inout) :: particle
154  call pstop(1, "pass must be overridden")
155  end subroutine pass
156 
157  !> @brief Save the particle's state to output files.
158  subroutine save(this, particle, reason)
159  use tdismodule, only: kper, kstp, totimc
160  ! dummy
161  class(methodtype), intent(inout) :: this
162  type(particletype), pointer, intent(inout) :: particle
163  integer(I4B), intent(in) :: reason
164  ! local
165  integer(I4B) :: per, stp
166 
167  per = kper
168  stp = kstp
169 
170  ! If tracking time falls exactly on a boundary between time steps,
171  ! report the previous time step for this datum. This is to follow
172  ! MP7's behavior, and because the particle will have been tracked
173  ! up to this instant under the previous time step's conditions, so
174  ! the time step we're about to start shouldn't get "credit" for it.
175  if (particle%ttrack == totimc .and. (per > 1 .or. stp > 1)) then
176  if (stp > 1) then
177  stp = stp - 1
178  else if (per > 1) then
179  per = per - 1
180  stp = 1
181  end if
182  end if
183 
184  call this%trackctl%save(particle, kper=per, kstp=stp, reason=reason)
185  end subroutine save
186 
187 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:143
subroutine save(this, particle, reason)
Save the particle's state to output files.
Definition: Method.f90:159
recursive subroutine track(this, particle, level, tmax)
Track the particle over domains of the given.
Definition: Method.f90:101
subroutine try_pass(this, particle, nextlevel, advancing)
Try passing the particle to the next subdomain.
Definition: Method.f90:123
subroutine pass(this, particle)
Pass the particle to the next subdomain.
Definition: Method.f90:152
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:31
Particle tracked by the PRT model.
Definition: Particle.f90:78
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.