MODFLOW 6  version 6.7.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)
 Transform particle coordinates. More...
 
subroutine reset_transform (this)
 
subroutine get_model_coords (this, x, y, z)
 Return the particle's model 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 125 of file Particle.f90.

126  type(ParticleStoreType), pointer :: this !< store
127  integer(I4B), intent(in) :: np !< number of particles
128  character(*), intent(in) :: mempath !< path to memory
129 
130  allocate (this)
131  call mem_allocate(this%imdl, np, 'PLIMDL', mempath)
132  call mem_allocate(this%irpt, np, 'PLIRPT', mempath)
133  call mem_allocate(this%iprp, np, 'PLIPRP', mempath)
134  call mem_allocate(this%name, lenboundname, np, 'PLNAME', mempath)
135  call mem_allocate(this%icu, np, 'PLICU', mempath)
136  call mem_allocate(this%ilay, np, 'PLILAY', mempath)
137  call mem_allocate(this%izone, np, 'PLIZONE', mempath)
138  call mem_allocate(this%izp, np, 'PLIZP', mempath)
139  call mem_allocate(this%istatus, np, 'PLISTATUS', mempath)
140  call mem_allocate(this%x, np, 'PLX', mempath)
141  call mem_allocate(this%y, np, 'PLY', mempath)
142  call mem_allocate(this%z, np, 'PLZ', mempath)
143  call mem_allocate(this%trelease, np, 'PLTRELEASE', mempath)
144  call mem_allocate(this%tstop, np, 'PLTSTOP', mempath)
145  call mem_allocate(this%ttrack, np, 'PLTTRACK', mempath)
146  call mem_allocate(this%istopweaksink, np, 'PLISTOPWEAKSINK', mempath)
147  call mem_allocate(this%istopzone, np, 'PLISTOPZONE', mempath)
148  call mem_allocate(this%idrymeth, np, 'PLIDRYMETH', mempath)
149  call mem_allocate(this%ifrctrn, np, 'PLIFRCTRN', mempath)
150  call mem_allocate(this%iexmeth, np, 'PLIEXMETH', mempath)
151  call mem_allocate(this%extol, np, 'PLEXTOL', mempath)
152  call mem_allocate(this%extend, np, 'PLIEXTEND', mempath)
153  call mem_allocate(this%idomain, np, levelmax, 'PLIDOMAIN', mempath)
154  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 117 of file Particle.f90.

118  type(ParticleType), pointer :: particle !< particle
119  allocate (particle)
120  allocate (particle%idomain(levelmax))
121  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 158 of file Particle.f90.

159  class(ParticleStoreType), intent(inout) :: this !< store
160  character(*), intent(in) :: mempath !< path to memory
161 
162  call mem_deallocate(this%imdl, 'PLIMDL', mempath)
163  call mem_deallocate(this%iprp, 'PLIPRP', mempath)
164  call mem_deallocate(this%irpt, 'PLIRPT', mempath)
165  call mem_deallocate(this%name, 'PLNAME', mempath)
166  call mem_deallocate(this%icu, 'PLICU', mempath)
167  call mem_deallocate(this%ilay, 'PLILAY', mempath)
168  call mem_deallocate(this%izone, 'PLIZONE', mempath)
169  call mem_deallocate(this%izp, 'PLIZP', mempath)
170  call mem_deallocate(this%istatus, 'PLISTATUS', mempath)
171  call mem_deallocate(this%x, 'PLX', mempath)
172  call mem_deallocate(this%y, 'PLY', mempath)
173  call mem_deallocate(this%z, 'PLZ', mempath)
174  call mem_deallocate(this%trelease, 'PLTRELEASE', mempath)
175  call mem_deallocate(this%tstop, 'PLTSTOP', mempath)
176  call mem_deallocate(this%ttrack, 'PLTTRACK', mempath)
177  call mem_deallocate(this%istopweaksink, 'PLISTOPWEAKSINK', mempath)
178  call mem_deallocate(this%istopzone, 'PLISTOPZONE', mempath)
179  call mem_deallocate(this%idrymeth, 'PLIDRYMETH', mempath)
180  call mem_deallocate(this%ifrctrn, 'PLIFRCTRN', mempath)
181  call mem_deallocate(this%iexmeth, 'PLIEXMETH', mempath)
182  call mem_deallocate(this%extol, 'PLEXTOL', mempath)
183  call mem_deallocate(this%extend, 'PLIEXTEND', mempath)
184  call mem_deallocate(this%idomain, 'PLIDOMAIN', mempath)
185  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 
)
private
Parameters
[in,out]thisparticle
[out]xx coordinate
[out]yy coordinate
[out]zz coordinate

Definition at line 348 of file Particle.f90.

349  use geomutilmodule, only: transform, compose
350  class(ParticleType), intent(inout) :: this !< particle
351  real(DP), intent(out) :: x !< x coordinate
352  real(DP), intent(out) :: y !< y coordinate
353  real(DP), intent(out) :: z !< z coordinate
354 
355  if (this%transformed) then
356  ! Untransform coordinates
357  call transform(this%x, this%y, this%z, x, y, z, &
358  this%xorigin, this%yorigin, this%zorigin, &
359  this%sinrot, this%cosrot, invert=.true.)
360  else
361  x = this%x
362  y = this%y
363  z = this%z
364  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 227 of file Particle.f90.

228  class(ParticleType), intent(inout) :: this !< particle
229  type(ParticleStoreType), intent(in) :: store !< particle storage
230  integer(I4B), intent(in) :: imdl !< index of model particle originated in
231  integer(I4B), intent(in) :: iprp !< index of particle release package particle originated in
232  integer(I4B), intent(in) :: ip !< index into the particle list
233 
234  call this%reset_transform()
235  this%imdl = imdl
236  this%iprp = iprp
237  this%irpt = store%irpt(ip)
238  this%ip = ip
239  this%name = store%name(ip)
240  this%istopweaksink = store%istopweaksink(ip)
241  this%istopzone = store%istopzone(ip)
242  this%idrymeth = store%idrymeth(ip)
243  this%icp = 0
244  this%icu = store%icu(ip)
245  this%ilay = store%ilay(ip)
246  this%izone = store%izone(ip)
247  this%izp = store%izp(ip)
248  this%istatus = store%istatus(ip)
249  this%x = store%x(ip)
250  this%y = store%y(ip)
251  this%z = store%z(ip)
252  this%trelease = store%trelease(ip)
253  this%tstop = store%tstop(ip)
254  this%ttrack = store%ttrack(ip)
255  this%advancing = .true.
256  this%idomain(1:levelmax) = &
257  store%idomain(ip, 1:levelmax)
258  this%idomain(1) = imdl
259  this%iboundary(1:levelmax) = &
260  store%iboundary(ip, 1:levelmax)
261  this%ifrctrn = store%ifrctrn(ip)
262  this%iexmeth = store%iexmeth(ip)
263  this%extol = store%extol(ip)
264  this%iextend = store%extend(ip)

◆ num_stored()

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

Definition at line 367 of file Particle.f90.

368  class(ParticleStoreType) :: this
369  n = size(this%imdl)

◆ reset_transform()

subroutine particlemodule::reset_transform ( class(particletype), intent(inout)  this)
Parameters
[in,out]thisparticle

Definition at line 335 of file Particle.f90.

336  class(ParticleType), intent(inout) :: this !< particle
337 
338  this%xorigin = dzero
339  this%yorigin = dzero
340  this%zorigin = dzero
341  this%sinrot = dzero
342  this%cosrot = done
343  this%cosrot = done
344  this%transformed = .false.

◆ 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 189 of file Particle.f90.

190  ! dummy
191  class(ParticleStoreType), intent(inout) :: this !< particle store
192  integer(I4B), intent(in) :: np !< number of particles
193  character(*), intent(in) :: mempath !< path to memory
194 
195  ! resize arrays
196  call mem_reallocate(this%imdl, np, 'PLIMDL', mempath)
197  call mem_reallocate(this%iprp, np, 'PLIPRP', mempath)
198  call mem_reallocate(this%irpt, np, 'PLIRPT', mempath)
199  call mem_reallocate(this%name, lenboundname, np, 'PLNAME', mempath)
200  call mem_reallocate(this%icu, np, 'PLICU', mempath)
201  call mem_reallocate(this%ilay, np, 'PLILAY', mempath)
202  call mem_reallocate(this%izone, np, 'PLIZONE', mempath)
203  call mem_reallocate(this%izp, np, 'PLIZP', mempath)
204  call mem_reallocate(this%istatus, np, 'PLISTATUS', mempath)
205  call mem_reallocate(this%x, np, 'PLX', mempath)
206  call mem_reallocate(this%y, np, 'PLY', mempath)
207  call mem_reallocate(this%z, np, 'PLZ', mempath)
208  call mem_reallocate(this%trelease, np, 'PLTRELEASE', mempath)
209  call mem_reallocate(this%tstop, np, 'PLTSTOP', mempath)
210  call mem_reallocate(this%ttrack, np, 'PLTTRACK', mempath)
211  call mem_reallocate(this%istopweaksink, np, 'PLISTOPWEAKSINK', mempath)
212  call mem_reallocate(this%istopzone, np, 'PLISTOPZONE', mempath)
213  call mem_reallocate(this%idrymeth, np, 'PLIDRYMETH', mempath)
214  call mem_reallocate(this%ifrctrn, np, 'PLIFRCTRN', mempath)
215  call mem_reallocate(this%iexmeth, np, 'PLIEXMETH', mempath)
216  call mem_reallocate(this%extol, np, 'PLEXTOL', mempath)
217  call mem_reallocate(this%extend, np, 'PLIEXTEND', mempath)
218  call mem_reallocate(this%idomain, np, levelmax, 'PLIDOMAIN', mempath)
219  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 268 of file Particle.f90.

269  class(ParticleStoreType), intent(inout) :: this !< particle storage
270  type(ParticleType), intent(in) :: particle !< particle
271  integer(I4B), intent(in) :: ip !< particle index
272 
273  this%imdl(ip) = particle%imdl
274  this%iprp(ip) = particle%iprp
275  this%irpt(ip) = particle%irpt
276  this%name(ip) = particle%name
277  this%istopweaksink(ip) = particle%istopweaksink
278  this%istopzone(ip) = particle%istopzone
279  this%idrymeth(ip) = particle%idrymeth
280  this%icu(ip) = particle%icu
281  this%ilay(ip) = particle%ilay
282  this%izone(ip) = particle%izone
283  this%izp(ip) = particle%izp
284  this%istatus(ip) = particle%istatus
285  this%x(ip) = particle%x
286  this%y(ip) = particle%y
287  this%z(ip) = particle%z
288  this%trelease(ip) = particle%trelease
289  this%tstop(ip) = particle%tstop
290  this%ttrack(ip) = particle%ttrack
291  this%idomain( &
292  ip, &
293  1:levelmax) = &
294  particle%idomain(1:levelmax)
295  this%iboundary( &
296  ip, &
297  1:levelmax) = &
298  particle%iboundary(1:levelmax)
299  this%ifrctrn(ip) = particle%ifrctrn
300  this%iexmeth(ip) = particle%iexmeth
301  this%extol(ip) = particle%extol
302  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 
)
private

Apply a translation and/or rotation to particle coordinates. No rescaling. It's also possible to invert a transformation. Be sure to reset the transformation after using it.

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

Definition at line 311 of file Particle.f90.

313  use geomutilmodule, only: transform, compose
314  class(ParticleType), intent(inout) :: this !< particle
315  real(DP), intent(in), optional :: xorigin !< x coordinate of origin
316  real(DP), intent(in), optional :: yorigin !< y coordinate of origin
317  real(DP), intent(in), optional :: zorigin !< z coordinate of origin
318  real(DP), intent(in), optional :: sinrot !< sine of rotation angle
319  real(DP), intent(in), optional :: cosrot !< cosine of rotation angle
320  logical(LGP), intent(in), optional :: invert !< whether to invert
321 
322  call transform(this%x, this%y, this%z, &
323  this%x, this%y, this%z, &
324  xorigin, yorigin, zorigin, &
325  sinrot, cosrot, invert)
326 
327  call compose(this%xorigin, this%yorigin, this%zorigin, &
328  this%sinrot, this%cosrot, &
329  xorigin, yorigin, zorigin, &
330  sinrot, cosrot, invert)
331 
332  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