35 integer :: n, idir, igrid, ipe_particle, nparticles_local
41 double precision :: w(ixg^t,1:nw)
54 rrd(n,idir) =
rng%unif_01()
58 {^
d&x(^
d,n) = xprobmin^
d + rrd(n+1,^
d) * (xprobmax^
d - xprobmin^
d)\}
75 if(ipe_particle ==
mype)
then
86 nparticles_local = nparticles_local + 1
115 integer :: igrid, iigrid, idir
116 double precision,
dimension(ixG^T,1:nw) :: w
118 do iigrid=1,igridstail; igrid=igrids(iigrid);
121 w(ixg^t,1:nw) = ps(igrid)%w(ixg^t,1:nw)
122 call phys_to_primitive(ixg^
ll,ixg^
ll,w,ps(igrid)%x)
134 double precision,
intent(in) :: end_time
136 double precision,
dimension(1:ndir) :: v, x
137 double precision :: defpayload(ndefpayload)
138 double precision :: usrpayload(nusrpayload)
139 double precision :: tloc, tlocnew, dt_p, h1
140 double precision,
parameter :: eps=1.0
d-6, hmin=1.0
d-8
141 integer :: ipart, iipart, igrid
142 integer :: nok, nbad, ierror
165 particle(ipart)%payload(1:ndefpayload) = defpayload
178 integer,
intent(in) :: igrid,mynpayload
179 double precision,
intent(in) :: xpart(1:ndir),upart(1:ndir),qpart,mpart,particle_time
180 double precision,
intent(out) :: mypayload(mynpayload)
181 double precision :: myw(ixG^T,1:nw),mywold(ixG^T,1:nw)
182 double precision :: wp, wpold, td
187 myw(ixg^t,1:nw) = gridvars(igrid)%w(ixg^t,1:nw)
188 if (
time_advance) mywold(ixg^t,1:nw) = gridvars(igrid)%wold(ixg^t,1:nw)
191 call phys_to_conserved(ixg^
ll,ixg^
ll,myw,ps(igrid)%x)
200 wp = wpold * (1.0d0 - td) + wp * td
212 double precision,
intent(in) :: end_time
213 double precision :: dt_p
214 integer :: ipart, iipart, nout, id
215 double precision :: tout, dt_cfl
216 double precision :: v(1:
ndir), xp(3), told, tnew
228 told = partp%self%time
Module with geometry-related routines (e.g., divergence, curl)
This module contains definitions of global parameters and variables and some generic functions/subrou...
double precision, dimension(:), allocatable w_convert_factor
Conversion factors the primitive variables.
double precision global_time
The global simulation time.
logical saveprim
If true, convert from conservative to primitive variables in output.
integer icomm
The MPI communicator.
integer mype
The rank of the current MPI task.
integer, dimension(:), allocatable, parameter d
integer ndir
Number of spatial dimensions (components) for vector variables.
integer ixm
the mesh range of a physical block without ghost cells
integer ierrmpi
A global MPI error return code.
double precision, dimension(:,:), allocatable rnode
Corner coordinates.
double precision, dimension(ndim) dxlevel
Module with shared functionality for all the particle movers.
pure subroutine limit_dt_endtime(t_left, dt_p)
integer num_particles
Number of particles.
double precision const_dt_particles
If positive, a constant time step for the particles.
double precision dtsave_particles
Time interval to save particles.
integer ngridvars
Number of variables for grid field.
integer nusrpayload
Number of user-defined payload variables for a particle.
subroutine find_particle_ipe(x, igrid_particle, ipe_particle)
integer npayload
Number of total payload variables for a particle.
integer, dimension(:), allocatable particles_active_on_mype
Array of identity numbers of active particles in current processor.
subroutine push_particle_into_particles_on_mype(ipart)
integer nparticles_active_on_mype
Number of active particles in current processor.
procedure(sub_define_additional_gridvars), pointer particles_define_additional_gridvars
integer nparticles
Identity number and total number of particles.
type(particle_ptr), dimension(:), allocatable particle
procedure(sub_integrate), pointer particles_integrate
procedure(sub_fill_gridvars), pointer particles_fill_gridvars
integer ndefpayload
Number of default payload variables for a particle.
logical function point_in_igrid_ghostc(x, igrid, ngh)
Quick check if particle coordinate is inside igrid (ghost cells included, except the last ngh)
integer igrid_working
set the current igrid for the particle integrator:
subroutine interpolate_var(igrid, ixIL, ixOL, gf, x, xloc, gfloc)
Scattered sampling based on fixed- or moving-particle interpolation By Fabio Bacchini (2020)
double precision function sample_get_particle_dt(partp, end_time)
subroutine, public sample_init()
subroutine sample_fill_gridvars
subroutine sample_integrate_particles(end_time)
subroutine sample_update_payload(igrid, xpart, upart, qpart, mpart, mypayload, mynpayload, particle_time)
Payload update.
subroutine, public sample_create_particles()
Module with all the methods that users can customize in AMRVAC.
procedure(particle_position), pointer usr_particle_position
procedure(check_particle), pointer usr_check_particle
procedure(create_particles), pointer usr_create_particles
procedure(update_payload), pointer usr_update_payload