MODFLOW 6  version 6.7.0.dev2
USGS Modular Hydrologic Model
methodcellmodule Module Reference

Data Types

type  methodcelltype
 

Functions/Subroutines

subroutine try_pass (this, particle, nextlevel, advancing)
 Try passing the particle to the next subdomain. More...
 
subroutine assess (this, particle, cell_defn, tmax)
 Check reporting/terminating conditions before tracking the particle across the cell. More...
 
subroutine cellexit (this, particle)
 Particle exits a cell. More...
 
logical(lgp) function forms_cycle (this, particle, event)
 Check if the event forms a cycle in the particle path. More...
 
subroutine store_event (this, particle, event)
 Save the event in the particle's history. Acts like a queue, the oldest event is removed when the event count exceeds the maximum size. More...
 
integer(i4b) function get_level (this)
 Get the cell method's level. More...
 

Function/Subroutine Documentation

◆ assess()

subroutine methodcellmodule::assess ( class(methodcelltype), intent(inout)  this,
type(particletype), intent(inout), pointer  particle,
type(celldefntype), intent(inout), pointer  cell_defn,
real(dp), intent(in)  tmax 
)
private

Check a number of conditions determining whether to continue tracking the particle or terminate it, as well as whether to record any output data as per selected reporting conditions.

Definition at line 59 of file MethodCell.f90.

60  ! modules
65  ! dummy
66  class(MethodCellType), intent(inout) :: this
67  type(ParticleType), pointer, intent(inout) :: particle
68  type(CellDefnType), pointer, intent(inout) :: cell_defn
69  real(DP), intent(in) :: tmax
70  ! local
71  logical(LGP) :: dry_cell, dry_particle, no_exit_face, stop_zone, weak_sink
72  integer(I4B) :: i
73  real(DP) :: t, ttrackmax
74 
75  dry_cell = this%fmi%ibdgwfsat0(cell_defn%icell) == 0
76  dry_particle = particle%z > cell_defn%top
77  no_exit_face = cell_defn%inoexitface > 0
78  stop_zone = cell_defn%izone > 0 .and. particle%istopzone == cell_defn%izone
79  weak_sink = cell_defn%iweaksink > 0
80 
81  particle%izone = cell_defn%izone
82  if (stop_zone) then
83  call this%terminate(particle, status=term_stopzone)
84  return
85  end if
86 
87  if (no_exit_face .and. .not. dry_cell) then
88  call this%terminate(particle, status=term_no_exits)
89  return
90  end if
91 
92  if (weak_sink) then
93  if (particle%istopweaksink > 0) then
94  call this%terminate(particle, status=term_weaksink)
95  return
96  else
97  call this%weaksink(particle)
98  end if
99  end if
100 
101  if (dry_cell) then
102  if (particle%idrymeth == 0) then
103  ! drop to cell bottom. handled by pass
104  ! to bottom method, nothing to do here
105  no_exit_face = .false.
106  else if (particle%idrymeth == 1) then
107  ! stop
108  call this%terminate(particle, status=term_inactive)
109  return
110  else if (particle%idrymeth == 2) then
111  ! stay
112  particle%advancing = .false.
113  no_exit_face = .false.
114 
115  ! we might report tracking times
116  ! out of order here, but we want
117  ! the particle termination event
118  ! (if this is the last time step)
119  ! to have the maximum tracking t,
120  ! so we need to keep tabs on it.
121  ttrackmax = totim
122 
123  ! update tracking time to time
124  ! step end time and save record
125  particle%ttrack = totim
126  call this%timestep(particle)
127 
128  ! record user tracking times
129  call this%tracktimes%advance()
130  if (this%tracktimes%any()) then
131  do i = this%tracktimes%selection(1), this%tracktimes%selection(2)
132  t = this%tracktimes%times(i)
133  if (t < totimc) cycle
134  if (t >= tmax) exit
135  particle%ttrack = t
136  call this%usertime(particle)
137  if (t > ttrackmax) ttrackmax = t
138  end do
139  end if
140 
141  ! terminate if last period/step
142  if (endofsimulation) then
143  particle%ttrack = ttrackmax
144  call this%terminate(particle, status=term_no_exits)
145  return
146  end if
147  end if
148  else if (dry_particle .and. this%name /= "passtobottom") then
149  if (particle%idrymeth == 0) then
150  ! drop to water table
151  particle%z = cell_defn%top
152  else if (particle%idrymeth == 1) then
153  ! stop
154  call this%terminate(particle, status=term_inactive)
155  return
156  else if (particle%idrymeth == 2) then
157  ! stay
158  particle%advancing = .false.
159  no_exit_face = .false.
160 
161  ! we might report tracking times
162  ! out of order here, but we want
163  ! the particle termination event
164  ! (if this is the last time step)
165  ! to have the maximum tracking t,
166  ! so we need to keep tabs on it.
167  ttrackmax = totim
168 
169  ! update tracking time to time
170  ! step end time and save record
171  particle%ttrack = totim
172  call this%timestep(particle)
173 
174  ! record user tracking times
175  call this%tracktimes%advance()
176  if (this%tracktimes%any()) then
177  do i = this%tracktimes%selection(1), this%tracktimes%selection(2)
178  t = this%tracktimes%times(i)
179  if (t < totimc) cycle
180  if (t >= tmax) exit
181  particle%ttrack = t
182  call this%usertime(particle)
183  if (t > ttrackmax) ttrackmax = t
184  end do
185  end if
186  end if
187  end if
188 
189  if (no_exit_face) then
190  particle%advancing = .false.
191  particle%istatus = term_no_exits
192  call this%terminate(particle)
193  return
194  end if
195 
@, public weaksink
particle entered a weak sink
@, public usertime
user-specified tracking time
@, public terminate
particle terminated
@, public timestep
time step ended
@ term_weaksink
terminated in a weak sink cell
Definition: Particle.f90:31
@ term_inactive
terminated in an inactive cell
Definition: Particle.f90:34
@ term_no_exits
terminated in a cell with no exit face
Definition: Particle.f90:32
@ term_stopzone
terminated in a cell with a stop zone number
Definition: Particle.f90:33
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

◆ cellexit()

subroutine methodcellmodule::cellexit ( class(methodcelltype), intent(inout)  this,
type(particletype), intent(inout), pointer  particle 
)

Definition at line 199 of file MethodCell.f90.

200  class(MethodCellType), intent(inout) :: this
201  type(ParticleType), pointer, intent(inout) :: particle
202  class(ParticleEventType), pointer :: event
203  ! local
204  integer(I4B) :: i, nhist
205  class(*), pointer :: prev
206 
207  allocate (cellexiteventtype :: event)
208  select type (event)
209  type is (cellexiteventtype)
210  event%exit_face = particle%iboundary(level_feature)
211  end select
212  call this%events%dispatch(particle, event)
213  if (particle%icycwin == 0) then
214  deallocate (event)
215  return
216  end if
217  if (this%forms_cycle(particle, event)) then
218  ! print event history
219  print *, "Cyclic pathline detected"
220  nhist = particle%history%Count()
221  do i = 1, nhist
222  prev => particle%history%GetItem(i)
223  select type (prev)
224  class is (particleeventtype)
225  print *, "Back ", nhist - i + 1, ": ", prev%get_text()
226  end select
227  end do
228  print *, "Current :", event%get_text()
229  call pstop(1, 'Cyclic pathline detected, aborting')
230  else
231  call this%store_event(particle, event)
232  end if
Here is the call graph for this function:

◆ forms_cycle()

logical(lgp) function methodcellmodule::forms_cycle ( class(methodcelltype), intent(inout)  this,
type(particletype), intent(inout), pointer  particle,
class(particleeventtype), intent(in), pointer  event 
)
private

Definition at line 236 of file MethodCell.f90.

237  ! dummy
238  class(MethodCellType), intent(inout) :: this
239  type(ParticleType), pointer, intent(inout) :: particle
240  class(ParticleEventType), pointer, intent(in) :: event
241  ! local
242  class(IteratorType), allocatable :: itr
243  logical(LGP) :: found_cycle
244 
245  found_cycle = .false.
246  select type (event)
247  type is (cellexiteventtype)
248  itr = particle%history%Iterator()
249  do while (itr%has_next())
250  call itr%next()
251  select type (prev => itr%value())
252  class is (cellexiteventtype)
253  if (event%icu == prev%icu .and. &
254  event%ilay == prev%ilay .and. &
255  event%izone == prev%izone .and. &
256  event%exit_face == prev%exit_face .and. &
257  event%exit_face /= 0) then
258  found_cycle = .true.
259  exit
260  end if
261  end select
262  end do
263  end select

◆ get_level()

integer(i4b) function methodcellmodule::get_level ( class(methodcelltype), intent(in)  this)
private

Definition at line 287 of file MethodCell.f90.

288  class(MethodCellType), intent(in) :: this
289  integer(I4B) :: level
290  level = level_feature

◆ store_event()

subroutine methodcellmodule::store_event ( class(methodcelltype), intent(inout)  this,
type(particletype), intent(inout), pointer  particle,
class(particleeventtype), intent(in), pointer  event 
)
private

Definition at line 269 of file MethodCell.f90.

270  ! dummy
271  class(MethodCellType), intent(inout) :: this
272  type(ParticleType), pointer, intent(inout) :: particle
273  class(ParticleEventType), pointer, intent(in) :: event
274  ! local
275  class(*), pointer :: p
276 
277  select type (event)
278  type is (cellexiteventtype)
279  p => event
280  call particle%history%Add(p)
281  if (particle%history%Count() > particle%icycwin) &
282  call particle%history%RemoveNode(1, .true.)
283  end select

◆ try_pass()

subroutine methodcellmodule::try_pass ( class(methodcelltype), intent(inout)  this,
type(particletype), intent(inout), pointer  particle,
integer(i4b)  nextlevel,
logical(lgp)  advancing 
)
private

Definition at line 29 of file MethodCell.f90.

30  class(MethodCellType), intent(inout) :: this
31  type(ParticleType), pointer, intent(inout) :: particle
32  integer(I4B) :: nextlevel
33  logical(LGP) :: advancing
34 
35  if (particle%advancing) then
36  ! if still advancing, pass to the next subdomain.
37  ! if that puts us on a boundary, then we're done.
38  ! raise a cell exit event.
39  call this%pass(particle)
40  if (particle%iboundary(nextlevel - 1) .ne. 0) then
41  advancing = .false.
42  call this%cellexit(particle)
43  end if
44  else
45  ! otherwise we're already done so
46  ! reset the domain boundary value.
47  advancing = .false.
48  particle%iboundary = 0
49  end if