MODFLOW 6  version 6.7.0.dev2
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
17  use basedismodule, only: disbasetype
18  use prtfmimodule, only: prtfmitype
19  use cellmodule, only: celltype
20  use celldefnmodule, only: celldefntype
22  use mathutilmodule, only: is_close
23  implicit none
24 
26 
27  !> @brief Tracking method level enumeration.
28  !!
29  !> Tracking levels: 1: model, 2: grid feature, 3: grid subfeature.
30  !! A tracking level identifies the domain through which a tracking
31  !! method is responsible for moving a particle. Methods operate on
32  !! a particular level and delegate to submethods for levels higher
33  !! than (i.e. below the scope of) their own.
34  !<
35  enum, bind(C)
36  enumerator :: level_model = 1
37  enumerator :: level_feature = 2
38  enumerator :: level_subfeature = 3
39  end enum
40 
41  private
42  public :: methodtype
43 
44  !> @brief Base type for particle tracking methods.
45  !!
46  !! The PRT tracking algorithm invokes a "tracking method" for each
47  !! domain. A domain can be a model, cell in a model, or subcell in
48  !! a cell. Tracking proceeds recursively, delegating to a possibly
49  !! arbitrary number of subdomains (currently, only the three above
50  !! are recognized). A tracking method is responsible for advancing
51  !! a particle through a domain, delegating to subdomains as needed
52  !! depending on cell geometry (implementing the strategy pattern).
53  !<
54  type, abstract :: methodtype
55  character(len=40), pointer, public :: name !< method name
56  logical(LGP), public :: delegates !< whether the method delegates
57  type(prtfmitype), pointer, public :: fmi => null() !< ptr to fmi
58  class(celltype), pointer, public :: cell => null() !< ptr to the current cell
59  class(subcelltype), pointer, public :: subcell => null() !< ptr to the current subcell
60  type(particleeventdispatchertype), pointer, public :: events => null() !< ptr to event dispatcher
61  type(timeselecttype), pointer, public :: tracktimes => null() !< ptr to user-defined tracking times
62  integer(I4B), dimension(:), pointer, contiguous, public :: izone => null() !< pointer to zone numbers
63  real(dp), dimension(:), pointer, contiguous, public :: flowja => null() !< pointer to intercell flows
64  real(dp), dimension(:), pointer, contiguous, public :: porosity => null() !< pointer to aquifer porosity
65  real(dp), dimension(:), pointer, contiguous, public :: retfactor => null() !< pointer to retardation factor
66  contains
67  ! Implemented in all subtypes
68  procedure(apply), deferred :: apply !< apply the method to the particle
69  procedure(assess), deferred :: assess !< assess conditions before tracking
70  procedure(deallocate), deferred :: deallocate !< deallocate the method object
71  ! Overridden in subtypes that delegate
72  procedure :: pass !< pass the particle to the next subdomain
73  procedure :: load !< load the subdomain tracking method
74  procedure :: get_level !< get the tracking method level
75  ! Implemented here
76  procedure :: init
77  procedure :: track
78  procedure :: try_pass
79  ! Event firing methods
80  procedure :: release
81  procedure :: terminate
82  procedure :: timestep
83  procedure :: weaksink
84  procedure :: usertime
85  end type methodtype
86 
87  abstract interface
88  subroutine apply(this, particle, tmax)
89  import dp
90  import methodtype
91  import particletype
92  class(methodtype), intent(inout) :: this
93  type(particletype), pointer, intent(inout) :: particle
94  real(DP), intent(in) :: tmax
95  end subroutine apply
96  subroutine assess(this, particle, cell_defn, tmax)
97  import dp
98  import methodtype
99  import particletype
100  import celldefntype
101  class(methodtype), intent(inout) :: this
102  type(particletype), pointer, intent(inout) :: particle
103  type(celldefntype), pointer, intent(inout) :: cell_defn
104  real(DP), intent(in) :: tmax
105  end subroutine assess
106  subroutine deallocate (this)
107  import methodtype
108  class(methodtype), intent(inout) :: this
109  end subroutine deallocate
110  end interface
111 
112 contains
113 
114  subroutine init(this, fmi, cell, subcell, events, tracktimes, &
115  izone, flowja, porosity, retfactor)
116  class(methodtype), intent(inout) :: this
117  type(prtfmitype), intent(in), pointer, optional :: fmi
118  class(celltype), intent(in), pointer, optional :: cell
119  class(subcelltype), intent(in), pointer, optional :: subcell
120  type(particleeventdispatchertype), intent(in), pointer, optional :: events
121  type(timeselecttype), intent(in), pointer, optional :: tracktimes
122  integer(I4B), intent(in), pointer, optional :: izone(:)
123  real(DP), intent(in), pointer, optional :: flowja(:)
124  real(DP), intent(in), pointer, optional :: porosity(:)
125  real(DP), intent(in), pointer, optional :: retfactor(:)
126 
127  if (present(fmi)) this%fmi => fmi
128  if (present(cell)) this%cell => cell
129  if (present(subcell)) this%subcell => subcell
130  if (present(events)) this%events => events
131  if (present(tracktimes)) this%tracktimes => tracktimes
132  if (present(izone)) this%izone => izone
133  if (present(flowja)) this%flowja => flowja
134  if (present(porosity)) this%porosity => porosity
135  if (present(retfactor)) this%retfactor => retfactor
136  end subroutine init
137 
138  !> @brief Track the particle over subdomains of the given
139  ! level until the particle terminates or tmax is reached.
140  recursive subroutine track(this, particle, level, tmax)
141  ! dummy
142  class(methodtype), intent(inout) :: this
143  type(particletype), pointer, intent(inout) :: particle
144  integer(I4B) :: level
145  real(dp), intent(in) :: tmax
146  ! local
147  logical(LGP) :: advancing
148  integer(I4B) :: nextlevel
149  class(methodtype), pointer :: submethod
150 
151  advancing = .true.
152  nextlevel = level + 1
153  do while (advancing)
154  call this%load(particle, nextlevel, submethod)
155  call submethod%apply(particle, tmax)
156  call this%try_pass(particle, nextlevel, advancing)
157  end do
158  end subroutine track
159 
160  !> @brief Try passing the particle to the next subdomain.
161  subroutine try_pass(this, particle, nextlevel, advancing)
162  class(methodtype), intent(inout) :: this
163  type(particletype), pointer, intent(inout) :: particle
164  integer(I4B) :: nextlevel
165  logical(LGP) :: advancing
166 
167  if (particle%advancing) then
168  ! if still advancing, pass to the next subdomain.
169  ! if that puts us on a boundary, then we're done.
170  call this%pass(particle)
171  if (particle%iboundary(nextlevel - 1) .ne. 0) &
172  advancing = .false.
173  else
174  ! otherwise we're already done so
175  ! reset the domain boundary value.
176  advancing = .false.
177  particle%iboundary = 0
178  end if
179  end subroutine try_pass
180 
181  !> @brief Load the subdomain tracking method (submethod).
182  subroutine load(this, particle, next_level, submethod)
183  class(methodtype), intent(inout) :: this
184  type(particletype), pointer, intent(inout) :: particle
185  integer, intent(in) :: next_level
186  class(methodtype), pointer, intent(inout) :: submethod
187  call pstop(1, "load must be overridden")
188  end subroutine load
189 
190  !> @brief Get the tracking method's level.
191  function get_level(this) result(level)
192  class(methodtype), intent(in) :: this
193  integer(I4B) :: level
194  level = -1 ! suppress compiler warning
195  call pstop(1, "get_level must be overridden")
196  end function get_level
197 
198  !> @brief Pass the particle to the next subdomain.
199  subroutine pass(this, particle)
200  class(methodtype), intent(inout) :: this
201  type(particletype), pointer, intent(inout) :: particle
202  call pstop(1, "pass must be overridden")
203  end subroutine pass
204 
205  !> @brief Particle is released.
206  subroutine release(this, particle)
207  class(methodtype), intent(inout) :: this
208  type(particletype), pointer, intent(inout) :: particle
209  class(particleeventtype), pointer :: event
210  allocate (releaseeventtype :: event)
211  call this%events%dispatch(particle, event)
212  deallocate (event)
213  end subroutine release
214 
215  !> @brief Particle terminates.
216  subroutine terminate(this, particle, status)
217  class(methodtype), intent(inout) :: this
218  type(particletype), pointer, intent(inout) :: particle
219  integer(I4B), intent(in), optional :: status
220  class(particleeventtype), pointer :: event
221  particle%advancing = .false.
222  if (present(status)) particle%istatus = status
223  allocate (terminationeventtype :: event)
224  call this%events%dispatch(particle, event)
225  deallocate (event)
226  end subroutine terminate
227 
228  !> @brief Time step ends.
229  subroutine timestep(this, particle)
230  class(methodtype), intent(inout) :: this
231  type(particletype), pointer, intent(inout) :: particle
232  class(particleeventtype), pointer :: event
233  allocate (timestepeventtype :: event)
234  call this%events%dispatch(particle, event)
235  deallocate (event)
236  end subroutine timestep
237 
238  !> @brief Particle leaves a weak sink.
239  subroutine weaksink(this, particle)
240  class(methodtype), intent(inout) :: this
241  type(particletype), pointer, intent(inout) :: particle
242  class(particleeventtype), pointer :: event
243  allocate (weaksinkeventtype :: event)
244  call this%events%dispatch(particle, event)
245  deallocate (event)
246  end subroutine weaksink
247 
248  !> @brief User-defined tracking time occurs.
249  subroutine usertime(this, particle)
250  class(methodtype), intent(inout) :: this
251  type(particletype), pointer, intent(inout) :: particle
252  class(particleeventtype), pointer :: event
253  allocate (usertimeeventtype :: event)
254  call this%events%dispatch(particle, event)
255  deallocate (event)
256  end subroutine usertime
257 
258 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:183
subroutine terminate(this, particle, status)
Particle terminates.
Definition: Method.f90:217
subroutine release(this, particle)
Particle is released.
Definition: Method.f90:207
recursive subroutine track(this, particle, level, tmax)
Track the particle over subdomains of the given.
Definition: Method.f90:141
subroutine weaksink(this, particle)
Particle leaves a weak sink.
Definition: Method.f90:240
subroutine try_pass(this, particle, nextlevel, advancing)
Try passing the particle to the next subdomain.
Definition: Method.f90:162
integer(i4b) function get_level(this)
Get the tracking method's level.
Definition: Method.f90:192
subroutine timestep(this, particle)
Time step ends.
Definition: Method.f90:230
@, public level_feature
Definition: Method.f90:37
subroutine usertime(this, particle)
User-defined tracking time occurs.
Definition: Method.f90:250
@, public level_subfeature
Definition: Method.f90:38
@, public level_model
Definition: Method.f90:36
subroutine pass(this, particle)
Pass the particle to the next subdomain.
Definition: Method.f90:200
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:11
Base type for particle tracking methods.
Definition: Method.f90:54
Base type for particle events.
Particle tracked by the PRT model.
Definition: Particle.f90:56
A subcell of a cell.
Definition: Subcell.f90:9
Represents a series of instants at which some event should occur.
Definition: TimeSelect.f90:30