40 double precision :: w(ixg^t,1:nw)
43 integer :: n, idir, igrid, ipe_particle, nparticles_local
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
98 call sample_update_payload(igrid,x(:,n),v(:,n),q(n),m(n),defpayload,
ndefpayload,0.d0)
117 subroutine sample_fill_gridvars
120 double precision,
dimension(ixG^T,1:nw) :: w
121 integer :: igrid, iigrid, idir
123 do iigrid=1,igridstail; igrid=igrids(iigrid);
126 w(ixg^t,1:nw) = ps(igrid)%w(ixg^t,1:nw)
127 call phys_to_primitive(ixg^
ll,ixg^
ll,w,ps(igrid)%x)
133 end subroutine sample_fill_gridvars
135 subroutine sample_integrate_particles(end_time)
139 double precision,
intent(in) :: end_time
141 double precision,
dimension(1:ndir) :: v, x
144 double precision :: tloc, tlocnew, dt_p, h1
145 double precision,
parameter :: eps=1.0
d-6, hmin=1.0
d-8
146 integer :: ipart, iipart, igrid
147 integer :: nok, nbad, ierror
151 dt_p = sample_get_particle_dt(
particle(ipart), end_time)
169 call sample_update_payload(igrid,x,v,0.d0,0.d0,defpayload,
ndefpayload,tlocnew)
178 end subroutine sample_integrate_particles
181 subroutine sample_update_payload(igrid,xpart,upart,qpart,mpart,mypayload,mynpayload,particle_time)
183 integer,
intent(in) :: igrid,mynpayload
184 double precision,
intent(in) :: xpart(1:
ndir),upart(1:
ndir),qpart,mpart,particle_time
185 double precision,
intent(out) :: mypayload(mynpayload)
186 double precision :: myw(ixg^t,1:nw),mywold(ixg^t,1:nw)
187 double precision :: wp, wpold, td
192 myw(ixg^t,1:nw) = gridvars(igrid)%w(ixg^t,1:nw)
193 if (
time_advance) mywold(ixg^t,1:nw) = gridvars(igrid)%wold(ixg^t,1:nw)
196 call phys_to_conserved(ixg^
ll,ixg^
ll,myw,ps(igrid)%x)
205 wp = wpold * (1.0d0 - td) + wp * td
210 end subroutine sample_update_payload
212 function sample_get_particle_dt(partp, end_time)
result(dt_p)
217 double precision,
intent(in) :: end_time
218 double precision :: dt_p
219 double precision :: tout, dt_cfl
220 double precision :: v(1:
ndir), xp(3), told, tnew
221 integer :: ipart, iipart, nout, id
233 told = partp%self%time
247 end function sample_get_particle_dt
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.
double precision, dimension(:), allocatable, parameter d
double precision dt
global time step
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.
logical time_advance
do time evolving
double precision, dimension(:,:), allocatable rnode
Corner coordinates.
double precision, dimension(^nd) dxlevel
store unstretched cell size of current level
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.
subroutine write_particle_output()
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
double precision t_next_output
Time to write next particle output.
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)
subroutine, public sample_init()
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