MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
particlemodule Module Reference

Data Types

type  particletype
 A 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 destroy_particle (this)
 Destroy a particle. More...
 
subroutine, public create_particle_store (this, np, mempath)
 Create a new particle store. More...
 
subroutine destroy_store (this, mempath)
 Deallocate particle arrays. More...
 
subroutine resize_store (this, np, mempath)
 Reallocate particle arrays. More...
 
subroutine load_from_store (this, store, imdl, iprp, ip)
 Initialize particle from particle list. More...
 
subroutine load_from_particle (this, particle, ip)
 Update particle store from particle. 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...
 
pure character(len=lenmempath) function, public get_particle_id (particle)
 Return the particle's composite ID. More...
 

Variables

integer, parameter, public levelmin = 0
 
integer, parameter, public levelmax = 4
 

Function/Subroutine Documentation

◆ create_particle()

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

Definition at line 102 of file Particle.f90.

103  type(ParticleType), pointer :: particle !< particle
104  allocate (particle)
105  allocate (particle%idomain(levelmin:levelmax))
106  allocate (particle%iboundary(levelmin:levelmax))
Here is the caller graph for this function:

◆ create_particle_store()

subroutine, public particlemodule::create_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 117 of file Particle.f90.

118  type(ParticleStoreType), pointer :: this !< store
119  integer(I4B), intent(in) :: np !< number of particles
120  character(*), intent(in) :: mempath !< path to memory
121 
122  allocate (this)
123  call mem_allocate(this%imdl, np, 'PLIMDL', mempath)
124  call mem_allocate(this%irpt, np, 'PLIRPT', mempath)
125  call mem_allocate(this%iprp, np, 'PLIPRP', mempath)
126  call mem_allocate(this%name, lenboundname, np, 'PLNAME', mempath)
127  ! -- kluge todo: update mem_allocate to allow custom range of indices?
128  ! e.g. here we want to allocate 0-4 for trackdomain levels, not 1-5
129  allocate (this%idomain(np, levelmin:levelmax))
130  allocate (this%iboundary(np, levelmin:levelmax))
131  call mem_allocate(this%icu, np, 'PLICU', mempath)
132  call mem_allocate(this%ilay, np, 'PLILAY', mempath)
133  call mem_allocate(this%izone, np, 'PLIZONE', mempath)
134  call mem_allocate(this%istatus, np, 'PLISTATUS', mempath)
135  call mem_allocate(this%x, np, 'PLX', mempath)
136  call mem_allocate(this%y, np, 'PLY', mempath)
137  call mem_allocate(this%z, np, 'PLZ', mempath)
138  call mem_allocate(this%trelease, np, 'PLTRELEASE', mempath)
139  call mem_allocate(this%tstop, np, 'PLTSTOP', mempath)
140  call mem_allocate(this%ttrack, np, 'PLTTRACK', mempath)
141  call mem_allocate(this%istopweaksink, np, 'PLISTOPWEAKSINK', mempath)
142  call mem_allocate(this%istopzone, np, 'PLISTOPZONE', mempath)
Here is the caller graph for this function:

◆ destroy_particle()

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

Definition at line 110 of file Particle.f90.

111  class(ParticleType), intent(inout) :: this !< particle
112  deallocate (this%idomain)
113  deallocate (this%iboundary)

◆ destroy_store()

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

Definition at line 146 of file Particle.f90.

147  class(ParticleStoreType), intent(inout) :: this !< store
148  character(*), intent(in) :: mempath !< path to memory
149 
150  call mem_deallocate(this%imdl, 'PLIMDL', mempath)
151  call mem_deallocate(this%iprp, 'PLIPRP', mempath)
152  call mem_deallocate(this%irpt, 'PLIRPT', mempath)
153  call mem_deallocate(this%name, 'PLNAME', mempath)
154  deallocate (this%idomain)
155  deallocate (this%iboundary)
156  call mem_deallocate(this%icu, 'PLICU', mempath)
157  call mem_deallocate(this%ilay, 'PLILAY', mempath)
158  call mem_deallocate(this%izone, 'PLIZONE', mempath)
159  call mem_deallocate(this%istatus, 'PLISTATUS', mempath)
160  call mem_deallocate(this%x, 'PLX', mempath)
161  call mem_deallocate(this%y, 'PLY', mempath)
162  call mem_deallocate(this%z, 'PLZ', mempath)
163  call mem_deallocate(this%trelease, 'PLTRELEASE', mempath)
164  call mem_deallocate(this%tstop, 'PLTSTOP', mempath)
165  call mem_deallocate(this%ttrack, 'PLTTRACK', mempath)
166  call mem_deallocate(this%istopweaksink, 'PLISTOPWEAKSINK', mempath)
167  call mem_deallocate(this%istopzone, 'PLISTOPZONE', 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 328 of file Particle.f90.

329  use geomutilmodule, only: transform, compose
330  class(ParticleType), intent(inout) :: this !< particle
331  real(DP), intent(out) :: x !< x coordinate
332  real(DP), intent(out) :: y !< y coordinate
333  real(DP), intent(out) :: z !< z coordinate
334 
335  if (this%transformed) then
336  ! -- Transform back from local to model coordinates
337  call transform(this%x, this%y, this%z, x, y, z, &
338  this%xorigin, this%yorigin, this%zorigin, &
339  this%sinrot, this%cosrot, .true.)
340  else
341  ! -- Already in model coordinates
342  x = this%x
343  y = this%y
344  z = this%z
345  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:180
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:240
Here is the call graph for this function:

◆ get_particle_id()

pure character(len=lenmempath) function, public particlemodule::get_particle_id ( class(particletype), intent(in)  particle)

Particles are uniquely identified by model index, PRP index, location index, and release time.

Returns
particle id

Definition at line 353 of file Particle.f90.

354  class(ParticleType), intent(in) :: particle !< particle
355  character(len=LENMEMPATH) :: id !< particle id
356 
357  write (id, '(I0,"-",I0,"-",I0,"-",F0.0)') &
358  particle%imdl, particle%iprp, particle%irpt, particle%trelease
Here is the caller graph for this function:

◆ load_from_particle()

subroutine particlemodule::load_from_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 248 of file Particle.f90.

249  class(ParticleStoreType), intent(inout) :: this !< particle storage
250  type(ParticleType), intent(in) :: particle !< particle
251  integer(I4B), intent(in) :: ip !< particle index
252 
253  this%imdl(ip) = particle%imdl
254  this%iprp(ip) = particle%iprp
255  this%irpt(ip) = particle%irpt
256  this%name(ip) = particle%name
257  this%istopweaksink(ip) = particle%istopweaksink
258  this%istopzone(ip) = particle%istopzone
259  this%icu(ip) = particle%icu
260  this%ilay(ip) = particle%ilay
261  this%izone(ip) = particle%izone
262  this%istatus(ip) = particle%istatus
263  this%x(ip) = particle%x
264  this%y(ip) = particle%y
265  this%z(ip) = particle%z
266  this%trelease(ip) = particle%trelease
267  this%tstop(ip) = particle%tstop
268  this%ttrack(ip) = particle%ttrack
269  this%idomain( &
270  ip, &
271  levelmin:levelmax) = &
272  particle%idomain(levelmin:levelmax)
273  this%iboundary( &
274  ip, &
275  levelmin:levelmax) = &
276  particle%iboundary(levelmin:levelmax)

◆ load_from_store()

subroutine particlemodule::load_from_store ( 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 
)

This routine is used to initialize a particle from the list so it can be tracked by prt_solve. The particle's advancing flag is set and local coordinate transformations 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 214 of file Particle.f90.

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

◆ resize_store()

subroutine particlemodule::resize_store ( 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 171 of file Particle.f90.

172  ! -- modules
174  ! -- dummy
175  class(ParticleStoreType), intent(inout) :: this !< particle store
176  integer(I4B), intent(in) :: np !< number of particles
177  character(*), intent(in) :: mempath !< path to memory
178 
179  ! resize 1D arrays
180  call mem_reallocate(this%imdl, np, 'PLIMDL', mempath)
181  call mem_reallocate(this%iprp, np, 'PLIPRP', mempath)
182  call mem_reallocate(this%irpt, np, 'PLIRPT', mempath)
183  call mem_reallocate(this%name, lenboundname, np, 'PLNAME', mempath)
184  call mem_reallocate(this%icu, np, 'PLICU', mempath)
185  call mem_reallocate(this%ilay, np, 'PLILAY', mempath)
186  call mem_reallocate(this%izone, np, 'PLIZONE', mempath)
187  call mem_reallocate(this%istatus, np, 'PLISTATUS', mempath)
188  call mem_reallocate(this%x, np, 'PLX', mempath)
189  call mem_reallocate(this%y, np, 'PLY', mempath)
190  call mem_reallocate(this%z, np, 'PLZ', mempath)
191  call mem_reallocate(this%trelease, np, 'PLTRELEASE', mempath)
192  call mem_reallocate(this%tstop, np, 'PLTSTOP', mempath)
193  call mem_reallocate(this%ttrack, np, 'PLTTRACK', mempath)
194  call mem_reallocate(this%istopweaksink, np, 'PLISTOPWEAKSINK', mempath)
195  call mem_reallocate(this%istopzone, np, 'PLISTOPZONE', mempath)
196  ! resize first dimension of 2D arrays
197  ! todo: memory manager support?
198  call expandarray2d( &
199  this%idomain, &
200  np - size(this%idomain, 1), &
201  0)
202  call expandarray2d( &
203  this%iboundary, &
204  np - size(this%iboundary, 1), &
205  0)

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

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

Variable Documentation

◆ levelmax

integer, parameter, public particlemodule::levelmax = 4

Definition at line 15 of file Particle.f90.

◆ levelmin

integer, parameter, public particlemodule::levelmin = 0

Definition at line 15 of file Particle.f90.

15  integer, parameter, public :: levelmin = 0, levelmax = 4