MPI-AMRVAC  3.0
The MPI - Adaptive Mesh Refinement - Versatile Advection Code
mod_particle_sample.t
Go to the documentation of this file.
1 !> Scattered sampling based on fixed-particle interpolation
2 !> By Fabio Bacchini (2020)
5 
6  private
7 
8  public :: sample_init
10 
11 contains
12 
13  subroutine sample_init()
15  integer :: idir
16 
17  allocate(vp(ndir))
18  do idir = 1, ndir
19  vp(idir) = idir
20  end do
21  ngridvars=nw
22 
24 
25  if (associated(particles_define_additional_gridvars)) then
27  end if
28 
30 
31  end subroutine sample_init
32 
34  ! initialise the particles (=fixed interpolation points)
38 
39  integer :: n, idir, igrid, ipe_particle, nparticles_local
40  double precision :: x(3, num_particles)
41  double precision :: v(3, num_particles)
42  double precision :: q(num_particles)
43  double precision :: m(num_particles)
44  double precision :: rrd(num_particles,ndir)
45  double precision :: w(ixg^t,1:nw)
46  double precision :: defpayload(ndefpayload)
47  double precision :: usrpayload(nusrpayload)
48  logical :: follow(num_particles), check
49 
50  follow = .false.
51  x = 0.0d0
52 
53  if (mype==0) then
54  if (.not. associated(usr_create_particles)) then
55  ! Randomly distributed
56  do idir=1,ndir
57  do n = 1, num_particles
58  rrd(n,idir) = rng%unif_01()
59  end do
60  end do
61  do n=1, num_particles
62  {^d&x(^d,n) = xprobmin^d + rrd(n+1,^d) * (xprobmax^d - xprobmin^d)\}
63  end do
64  else
65  call usr_create_particles(num_particles, x, v, q, m, follow)
66  end if
67  end if
68 
69  call mpi_bcast(x,3*num_particles,mpi_double_precision,0,icomm,ierrmpi)
70  call mpi_bcast(follow,num_particles,mpi_logical,0,icomm,ierrmpi)
71 
72  nparticles_local = 0
73 
74  do n=1,num_particles
75  call find_particle_ipe(x(:,n),igrid,ipe_particle)
76  particle(n)%igrid = igrid
77  particle(n)%ipe = ipe_particle
78 
79  if(ipe_particle == mype) then
80  check = .true.
81 
82  ! Check for user-defined modifications or rejection conditions
83  if (associated(usr_check_particle)) call usr_check_particle(igrid, x(:,n), v(:,n), q(n), m(n), follow(n), check)
84  if (check) then
86  else
87  cycle
88  end if
89 
90  nparticles_local = nparticles_local + 1
91 
92  allocate(particle(n)%self)
93  particle(n)%self%follow = follow(n)
94  particle(n)%self%index = n
95  particle(n)%self%time = global_time
96  particle(n)%self%dt = 0.0d0
97  particle(n)%self%x = 0.d0
98  particle(n)%self%x(:) = x(:,n)
99  particle(n)%self%u(:) = 0.d0
100  allocate(particle(n)%payload(npayload))
101  call sample_update_payload(igrid,ps(igrid)%w,pso(igrid)%w,ps(igrid)%x,x(:,n),v(:,n),q(n),m(n),defpayload,ndefpayload,0.d0)
102  particle(n)%payload(1:ndefpayload) = defpayload
103  if (associated(usr_update_payload)) then
104  call usr_update_payload(igrid,ps(igrid)%w,pso(igrid)%w,ps(igrid)%x,x(:,n),v(:,n),q(n),m(n),usrpayload,nusrpayload,0.d0)
105  particle(n)%payload(ndefpayload+1:npayload)=usrpayload
106  end if
107  end if
108 
109  end do
110 
111  call mpi_allreduce(nparticles_local,nparticles,1,mpi_integer,mpi_sum,icomm,ierrmpi)
112 
113  end subroutine sample_create_particles
114 
117 
118  integer :: igrid, iigrid, idir
119  double precision, dimension(ixG^T,1:nw) :: w
120 
121  do iigrid=1,igridstail; igrid=igrids(iigrid);
122 
123  ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
124 
125  gridvars(igrid)%w(ixg^t,1:ngridvars) = 0.0d0
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)
128  ! fill all variables:
129  gridvars(igrid)%w(ixg^t,1:ngridvars) = w(ixg^t,1:ngridvars)
130 
131  if(time_advance) then
132  gridvars(igrid)%wold(ixg^t,1:ngridvars) = 0.0d0
133  w(ixg^t,1:nw) = pso(igrid)%w(ixg^t,1:nw)
134  call phys_to_primitive(ixg^ll,ixg^ll,w,ps(igrid)%x)
135  gridvars(igrid)%wold(ixg^t,1:ngridvars) = w(ixg^t,1:ngridvars)
136  end if
137 
138  end do
139 
140  end subroutine sample_fill_gridvars
141 
142  subroutine sample_integrate_particles(end_time)
143  ! this interpolates the HD/MHD quantities at the particle positions
146  double precision, intent(in) :: end_time
147 
148  double precision, dimension(1:ndir) :: v, x
149  double precision :: defpayload(ndefpayload)
150  double precision :: usrpayload(nusrpayload)
151  double precision :: tloc, tlocnew, dt_p, h1
152  double precision,parameter :: eps=1.0d-6, hmin=1.0d-8
153  integer :: ipart, iipart, igrid
154  integer :: nok, nbad, ierror
155 
156  do iipart=1,nparticles_active_on_mype
157  ipart = particles_active_on_mype(iipart);
158  dt_p = sample_get_particle_dt(particle(ipart), end_time)
159  particle(ipart)%self%dt = dt_p
160 
161  igrid = particle(ipart)%igrid
162  igrid_working = igrid
163  tloc = particle(ipart)%self%time
164  x(1:ndir) = particle(ipart)%self%x(1:ndir)
165  tlocnew = tloc+dt_p
166 
167  ! Position update (if defined)
168  ! TODO: this may create problems with interpolation out of boundaries
169  if (associated(usr_particle_position)) call usr_particle_position(x,ipart,tloc,tlocnew)
170  particle(ipart)%self%x(1:ndir) = x
171 
172  ! Time update
173  particle(ipart)%self%time = tlocnew
174 
175  ! Update payload
176  call sample_update_payload(igrid,ps(igrid)%w,pso(igrid)%w,ps(igrid)%x,x,v,0.d0,0.d0,defpayload,ndefpayload,tlocnew)
177  particle(ipart)%payload(1:ndefpayload) = defpayload
178  if (associated(usr_update_payload)) then
179  call usr_update_payload(igrid,ps(igrid)%w,pso(igrid)%w,ps(igrid)%x,x,v,0.d0,0.d0,usrpayload,nusrpayload,tlocnew)
180  particle(ipart)%payload(ndefpayload+1:npayload) = usrpayload
181  end if
182 
183  end do
184 
185  end subroutine sample_integrate_particles
186 
187  !> Payload update
188  subroutine sample_update_payload(igrid,w,wold,xgrid,xpart,upart,qpart,mpart,mypayload,mynpayload,particle_time)
190  integer, intent(in) :: igrid,mynpayload
191  double precision, intent(in) :: w(ixG^T,1:nw),wold(ixG^T,1:nw)
192  double precision, intent(in) :: xgrid(ixG^T,1:ndim),xpart(1:ndir),upart(1:ndir),qpart,mpart,particle_time
193  double precision, intent(out) :: mypayload(mynpayload)
194  double precision :: wp, wp1, wp2, td
195  integer :: ii
196 
197  td = (particle_time - global_time) / dt
198 
199  ! There are npayload=nw payloads, one for each primitive fluid quantity
200  do ii=1,mynpayload
201  if (.not.time_advance) then
202  call interpolate_var(igrid,ixg^ll,ixm^ll,w(ixg^t,ii),xgrid,xpart,wp)
203  else
204  call interpolate_var(igrid,ixg^ll,ixm^ll,wold(ixg^t,ii),xgrid,xpart,wp1)
205  call interpolate_var(igrid,ixg^ll,ixm^ll,w(ixg^t,ii),xgrid,xpart,wp2)
206  wp = wp1 * (1.0d0 - td) + wp2 * td
207  end if
208  mypayload(ii) = wp
209  end do
210 
211  end subroutine sample_update_payload
212 
213  pure function sample_get_particle_dt(partp, end_time) result(dt_p)
215  use mod_geometry
216  type(particle_ptr), intent(in) :: partp
217  double precision, intent(in) :: end_time
218  double precision :: dt_p
219  integer :: ipart, iipart, nout
220  double precision :: tout, dt_cfl
221  double precision :: v(1:ndir)
222 
223  if (const_dt_particles > 0) then
224  dt_p = const_dt_particles
225  return
226  end if
227 
228  dt_p = dtsave_particles
229 
230  ! Make sure we don't advance beyond end_time
231  call limit_dt_endtime(end_time - partp%self%time, dt_p)
232 
233  end function sample_get_particle_dt
234 
235 end module mod_particle_sample
Module with geometry-related routines (e.g., divergence, curl)
Definition: mod_geometry.t:2
This module contains definitions of global parameters and variables and some generic functions/subrou...
double precision global_time
The global simulation time.
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 (within a block with 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
integer, dimension(:), allocatable vp
Variable index for fluid velocity.
procedure(sub_integrate), pointer particles_integrate
procedure(sub_fill_gridvars), pointer particles_fill_gridvars
integer ndefpayload
Number of default payload variables for a particle.
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-particle interpolation By Fabio Bacchini (2020)
subroutine, public sample_init()
pure double precision function sample_get_particle_dt(partp, end_time)
subroutine sample_integrate_particles(end_time)
subroutine sample_update_payload(igrid, w, wold, xgrid, 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