MODFLOW 6  version 6.6.0.dev0
USGS Modular Hydrologic Model
particlemodule Module Reference

Data Types

type  particletype
 Particle tracked by the PRT model. More...
 
type  particlestoretype
 Structure of arrays to store particles. More...
 

Functions/Subroutines

subroutine, public create_particle (particle)
 Create a new particle. More...
 
subroutine, public allocate_particle_store (this, np, mempath)
 Create a new particle store. More...
 
subroutine deallocate (this, mempath)
 Deallocate particle arrays. More...
 
subroutine resize (this, np, mempath)
 Reallocate particle arrays. More...
 
subroutine load_particle (this, store, imdl, iprp, ip)
 Load a particle from the particle store. More...
 
subroutine save_particle (this, particle, ip)
 Save a particle's state to the particle store. More...
 
subroutine transform_coords (this, xorigin, yorigin, zorigin, sinrot, cosrot, invert, reset)
 Apply the given global-to-local transformation to the particle. More...
 
subroutine get_model_coords (this, x, y, z)
 Return the particle's model (global) coordinates. More...
 
integer function num_stored (this)
 

Variables

integer, parameter, public levelmax = 4
 

Function/Subroutine Documentation

◆ allocate_particle_store()

subroutine, public particlemodule::allocate_particle_store ( type(particlestoretype), pointer  this,
integer(i4b), intent(in)  np,
character(*), intent(in)  mempath 
)
Parameters
thisstore
[in]npnumber of particles
[in]mempathpath to memory

Definition at line 119 of file Particle.f90.

120  type(ParticleStoreType), pointer :: this !< store
121  integer(I4B), intent(in) :: np !< number of particles
122  character(*), intent(in) :: mempath !< path to memory
123 
124  allocate (this)
125  call mem_allocate(this%imdl, np, 'PLIMDL', mempath)
126  call mem_allocate(this%irpt, np, 'PLIRPT', mempath)
127  call mem_allocate(this%iprp, np, 'PLIPRP', mempath)
128  call mem_allocate(this%name, lenboundname, np, 'PLNAME', mempath)
129  call mem_allocate(this%icu, np, 'PLICU', mempath)
130  call mem_allocate(this%ilay, np, 'PLILAY', mempath)
131  call mem_allocate(this%izone, np, 'PLIZONE', mempath)
132  call mem_allocate(this%istatus, np, 'PLISTATUS', mempath)
133  call mem_allocate(this%x, np, 'PLX', mempath)
134  call mem_allocate(this%y, np, 'PLY', mempath)
135  call mem_allocate(this%z, np, 'PLZ', mempath)
136  call mem_allocate(this%trelease, np, 'PLTRELEASE', mempath)
137  call mem_allocate(this%tstop, np, 'PLTSTOP', mempath)
138  call mem_allocate(this%ttrack, np, 'PLTTRACK', mempath)
139  call mem_allocate(this%istopweaksink, np, 'PLISTOPWEAKSINK', mempath)
140  call mem_allocate(this%istopzone, np, 'PLISTOPZONE', mempath)
141  call mem_allocate(this%ifrctrn, np, 'PLIFRCTRN', mempath)
142  call mem_allocate(this%iexmeth, np, 'PLIEXMETH', mempath)
143  call mem_allocate(this%extol, np, 'PLEXTOL', mempath)
144  call mem_allocate(this%extend, np, 'PLIEXTEND', mempath)
145  call mem_allocate(this%idomain, np, levelmax, 'PLIDOMAIN', mempath)
146  call mem_allocate(this%iboundary, np, levelmax, 'PLIBOUNDARY', mempath)
Here is the caller graph for this function:

◆ create_particle()

subroutine, public particlemodule::create_particle ( type(particletype), pointer  particle)

Definition at line 111 of file Particle.f90.

112  type(ParticleType), pointer :: particle !< particle
113  allocate (particle)
114  allocate (particle%idomain(levelmax))
115  allocate (particle%iboundary(levelmax))
Here is the caller graph for this function:

◆ deallocate()

subroutine particlemodule::deallocate ( class(particlestoretype), intent(inout)  this,
character(*), intent(in)  mempath 
)
private
Parameters
[in,out]thisstore
[in]mempathpath to memory

Definition at line 150 of file Particle.f90.

151  class(ParticleStoreType), intent(inout) :: this !< store
152  character(*), intent(in) :: mempath !< path to memory
153 
154  call mem_deallocate(this%imdl, 'PLIMDL', mempath)
155  call mem_deallocate(this%iprp, 'PLIPRP', mempath)
156  call mem_deallocate(this%irpt, 'PLIRPT', mempath)
157  call mem_deallocate(this%name, 'PLNAME', mempath)
158  call mem_deallocate(this%icu, 'PLICU', mempath)
159  call mem_deallocate(this%ilay, 'PLILAY', mempath)
160  call mem_deallocate(this%izone, 'PLIZONE', mempath)
161  call mem_deallocate(this%istatus, 'PLISTATUS', mempath)
162  call mem_deallocate(this%x, 'PLX', mempath)
163  call mem_deallocate(this%y, 'PLY', mempath)
164  call mem_deallocate(this%z, 'PLZ', mempath)
165  call mem_deallocate(this%trelease, 'PLTRELEASE', mempath)
166  call mem_deallocate(this%tstop, 'PLTSTOP', mempath)
167  call mem_deallocate(this%ttrack, 'PLTTRACK', mempath)
168  call mem_deallocate(this%istopweaksink, 'PLISTOPWEAKSINK', mempath)
169  call mem_deallocate(this%istopzone, 'PLISTOPZONE', mempath)
170  call mem_deallocate(this%ifrctrn, 'PLIFRCTRN', mempath)
171  call mem_deallocate(this%iexmeth, 'PLIEXMETH', mempath)
172  call mem_deallocate(this%extol, 'PLEXTOL', mempath)
173  call mem_deallocate(this%extend, 'PLIEXTEND', mempath)
174  call mem_deallocate(this%idomain, 'PLIDOMAIN', mempath)
175  call mem_deallocate(this%iboundary, 'PLIBOUNDARY', mempath)

◆ get_model_coords()

subroutine particlemodule::get_model_coords ( class(particletype), intent(inout)  this,
real(dp), intent(out)  x,
real(dp), intent(out)  y,
real(dp), intent(out)  z 
)
Parameters
[in,out]thisparticle
[out]xx coordinate
[out]yy coordinate
[out]zz coordinate

Definition at line 337 of file Particle.f90.

338  use geomutilmodule, only: transform, compose
339  class(ParticleType), intent(inout) :: this !< particle
340  real(DP), intent(out) :: x !< x coordinate
341  real(DP), intent(out) :: y !< y coordinate
342  real(DP), intent(out) :: z !< z coordinate
343 
344  if (this%transformed) then
345  ! Transform back from local to model coordinates
346  call transform(this%x, this%y, this%z, x, y, z, &
347  this%xorigin, this%yorigin, this%zorigin, &
348  this%sinrot, this%cosrot, .true.)
349  else
350  ! Already in model coordinates
351  x = this%x
352  y = this%y
353  z = this%z
354  end if
subroutine, public transform(xin, yin, zin, xout, yout, zout, xorigin, yorigin, zorigin, sinrot, cosrot, invert)
Apply a 3D translation and optional 2D rotation to coordinates.
Definition: GeomUtil.f90:183
subroutine, public compose(xorigin, yorigin, zorigin, sinrot, cosrot, xorigin_new, yorigin_new, zorigin_new, sinrot_new, cosrot_new, invert)
Apply a 3D translation and 2D rotation to an existing transformation.
Definition: GeomUtil.f90:243
Here is the call graph for this function:

◆ load_particle()

subroutine particlemodule::load_particle ( class(particletype), intent(inout)  this,
type(particlestoretype), intent(in)  store,
integer(i4b), intent(in)  imdl,
integer(i4b), intent(in)  iprp,
integer(i4b), intent(in)  ip 
)
private

This routine is used to initialize a particle for tracking. The advancing flag and coordinate transformation are reset.

Parameters
[in,out]thisparticle
[in]storeparticle storage
[in]imdlindex of model particle originated in
[in]iprpindex of particle release package particle originated in
[in]ipindex into the particle list

Definition at line 215 of file Particle.f90.

216  class(ParticleType), intent(inout) :: this !< particle
217  type(ParticleStoreType), intent(in) :: store !< particle storage
218  integer(I4B), intent(in) :: imdl !< index of model particle originated in
219  integer(I4B), intent(in) :: iprp !< index of particle release package particle originated in
220  integer(I4B), intent(in) :: ip !< index into the particle list
221 
222  call this%transform(reset=.true.)
223  this%imdl = imdl
224  this%iprp = iprp
225  this%irpt = store%irpt(ip)
226  this%ip = ip
227  this%name = store%name(ip)
228  this%istopweaksink = store%istopweaksink(ip)
229  this%istopzone = store%istopzone(ip)
230  this%icu = store%icu(ip)
231  this%ilay = store%ilay(ip)
232  this%izone = store%izone(ip)
233  this%istatus = store%istatus(ip)
234  this%x = store%x(ip)
235  this%y = store%y(ip)
236  this%z = store%z(ip)
237  this%trelease = store%trelease(ip)
238  this%tstop = store%tstop(ip)
239  this%ttrack = store%ttrack(ip)
240  this%advancing = .true.
241  this%idomain(1:levelmax) = &
242  store%idomain(ip, 1:levelmax)
243  this%idomain(1) = imdl
244  this%iboundary(1:levelmax) = &
245  store%iboundary(ip, 1:levelmax)
246  this%ifrctrn = store%ifrctrn(ip)
247  this%iexmeth = store%iexmeth(ip)
248  this%extol = store%extol(ip)
249  this%iextend = store%extend(ip)

◆ num_stored()

integer function particlemodule::num_stored ( class(particlestoretype this)

Definition at line 357 of file Particle.f90.

358  class(ParticleStoreType) :: this
359  n = size(this%imdl)

◆ resize()

subroutine particlemodule::resize ( class(particlestoretype), intent(inout)  this,
integer(i4b), intent(in)  np,
character(*), intent(in)  mempath 
)
private
Parameters
[in,out]thisparticle store
[in]npnumber of particles
[in]mempathpath to memory

Definition at line 179 of file Particle.f90.

180  ! -- dummy
181  class(ParticleStoreType), intent(inout) :: this !< particle store
182  integer(I4B), intent(in) :: np !< number of particles
183  character(*), intent(in) :: mempath !< path to memory
184 
185  ! resize arrays
186  call mem_reallocate(this%imdl, np, 'PLIMDL', mempath)
187  call mem_reallocate(this%iprp, np, 'PLIPRP', mempath)
188  call mem_reallocate(this%irpt, np, 'PLIRPT', mempath)
189  call mem_reallocate(this%name, lenboundname, np, 'PLNAME', mempath)
190  call mem_reallocate(this%icu, np, 'PLICU', mempath)
191  call mem_reallocate(this%ilay, np, 'PLILAY', mempath)
192  call mem_reallocate(this%izone, np, 'PLIZONE', mempath)
193  call mem_reallocate(this%istatus, np, 'PLISTATUS', mempath)
194  call mem_reallocate(this%x, np, 'PLX', mempath)
195  call mem_reallocate(this%y, np, 'PLY', mempath)
196  call mem_reallocate(this%z, np, 'PLZ', mempath)
197  call mem_reallocate(this%trelease, np, 'PLTRELEASE', mempath)
198  call mem_reallocate(this%tstop, np, 'PLTSTOP', mempath)
199  call mem_reallocate(this%ttrack, np, 'PLTTRACK', mempath)
200  call mem_reallocate(this%istopweaksink, np, 'PLISTOPWEAKSINK', mempath)
201  call mem_reallocate(this%istopzone, np, 'PLISTOPZONE', mempath)
202  call mem_reallocate(this%ifrctrn, np, 'PLIFRCTRN', mempath)
203  call mem_reallocate(this%iexmeth, np, 'PLIEXMETH', mempath)
204  call mem_reallocate(this%extol, np, 'PLEXTOL', mempath)
205  call mem_reallocate(this%extend, np, 'PLIEXTEND', mempath)
206  call mem_reallocate(this%idomain, np, levelmax, 'PLIDOMAIN', mempath)
207  call mem_reallocate(this%iboundary, np, levelmax, 'PLIBOUNDARY', mempath)

◆ save_particle()

subroutine particlemodule::save_particle ( class(particlestoretype), intent(inout)  this,
type(particletype), intent(in)  particle,
integer(i4b), intent(in)  ip 
)
private
Parameters
[in,out]thisparticle storage
[in]ipparticle index

Definition at line 253 of file Particle.f90.

254  class(ParticleStoreType), intent(inout) :: this !< particle storage
255  type(ParticleType), intent(in) :: particle !< particle
256  integer(I4B), intent(in) :: ip !< particle index
257 
258  this%imdl(ip) = particle%imdl
259  this%iprp(ip) = particle%iprp
260  this%irpt(ip) = particle%irpt
261  this%name(ip) = particle%name
262  this%istopweaksink(ip) = particle%istopweaksink
263  this%istopzone(ip) = particle%istopzone
264  this%icu(ip) = particle%icu
265  this%ilay(ip) = particle%ilay
266  this%izone(ip) = particle%izone
267  this%istatus(ip) = particle%istatus
268  this%x(ip) = particle%x
269  this%y(ip) = particle%y
270  this%z(ip) = particle%z
271  this%trelease(ip) = particle%trelease
272  this%tstop(ip) = particle%tstop
273  this%ttrack(ip) = particle%ttrack
274  this%idomain( &
275  ip, &
276  1:levelmax) = &
277  particle%idomain(1:levelmax)
278  this%iboundary( &
279  ip, &
280  1:levelmax) = &
281  particle%iboundary(1:levelmax)
282  this%ifrctrn(ip) = particle%ifrctrn
283  this%iexmeth(ip) = particle%iexmeth
284  this%extol(ip) = particle%extol
285  this%extend(ip) = particle%iextend

◆ transform_coords()

subroutine particlemodule::transform_coords ( class(particletype), intent(inout)  this,
real(dp), intent(in), optional  xorigin,
real(dp), intent(in), optional  yorigin,
real(dp), intent(in), optional  zorigin,
real(dp), intent(in), optional  sinrot,
real(dp), intent(in), optional  cosrot,
logical(lgp), intent(in), optional  invert,
logical(lgp), intent(in), optional  reset 
)
private
Parameters
[in,out]thisparticle
[in]xoriginx coordinate of origin
[in]yoriginy coordinate of origin
[in]zoriginz coordinate of origin
[in]sinrotsine of rotation angle
[in]cosrotcosine of rotation angle
[in]invertwhether to invert
[in]resetwhether to reset

Definition at line 289 of file Particle.f90.

291  use geomutilmodule, only: transform, compose
292  class(ParticleType), intent(inout) :: this !< particle
293  real(DP), intent(in), optional :: xorigin !< x coordinate of origin
294  real(DP), intent(in), optional :: yorigin !< y coordinate of origin
295  real(DP), intent(in), optional :: zorigin !< z coordinate of origin
296  real(DP), intent(in), optional :: sinrot !< sine of rotation angle
297  real(DP), intent(in), optional :: cosrot !< cosine of rotation angle
298  logical(LGP), intent(in), optional :: invert !< whether to invert
299  logical(LGP), intent(in), optional :: reset !< whether to reset
300 
301  ! Reset if requested
302  if (present(reset)) then
303  if (reset) then
304  this%xorigin = dzero
305  this%yorigin = dzero
306  this%zorigin = dzero
307  this%sinrot = dzero
308  this%cosrot = done
309  this%cosrot = done
310  this%transformed = .false.
311  return
312  end if
313  end if
314 
315  ! Otherwise, transform coordinates
316  call transform(this%x, this%y, this%z, &
317  this%x, this%y, this%z, &
318  xorigin, yorigin, zorigin, &
319  sinrot, cosrot, invert)
320 
321  ! Modify transformation from model coordinates to particle's new
322  ! local coordinates by incorporating this latest transformation
323  call compose(this%xorigin, this%yorigin, this%zorigin, &
324  this%sinrot, this%cosrot, &
325  xorigin, yorigin, zorigin, &
326  sinrot, cosrot, invert)
327 
328  ! Set isTransformed flag to true. Note that there is no check
329  ! to see whether the modification brings the coordinates back
330  ! to model coordinates (in which case the origin would be very
331  ! close to zero and sinrot and cosrot would be very close to 0.
332  ! and 1., respectively, allowing for roundoff error).
333  this%transformed = .true.
Here is the call graph for this function:

Variable Documentation

◆ levelmax

integer, parameter, public particlemodule::levelmax = 4

Definition at line 14 of file Particle.f90.

14  integer, parameter, public :: levelmax = 4