MPI-AMRVAC  3.0
The MPI - Adaptive Mesh Refinement - Versatile Advection Code
mod_particle_advect.t
Go to the documentation of this file.
1 !> Tracer for advected particles moving with fluid flows
2 !> By Jannis Teunissen, Bart Ripperda, Oliver Porth, and Fabio Bacchini (2017-2020)
5 
6  private
7 
8  public :: advect_init
10 
11 contains
12 
13  subroutine advect_init()
15  integer :: idir
16 
17  allocate(vp(ndir))
18  do idir = 1, ndir
19  vp(idir) = idir
20  end do
22 
24 
25  if (associated(particles_define_additional_gridvars)) then
27  end if
28 
30 
31  end subroutine advect_init
32 
34  ! initialise the particles
37 
38  integer :: n, idir, igrid, ipe_particle, nparticles_local
39  double precision :: x(3, num_particles)
40  double precision :: v(3, num_particles)
41  double precision :: q(num_particles)
42  double precision :: m(num_particles)
43  double precision :: rrd(num_particles,ndir)
44  double precision :: w(ixg^t,1:nw)
45  double precision :: defpayload(ndefpayload)
46  double precision :: usrpayload(nusrpayload)
47  logical :: follow(num_particles), check
48 
49  follow = .false.
50  x = 0.0d0
51 
52  if (mype==0) then
53  if (.not. associated(usr_create_particles)) then
54  ! Randomly distributed
55  do idir=1,ndir
56  do n = 1, num_particles
57  rrd(n,idir) = rng%unif_01()
58  end do
59  end do
60  do n=1, num_particles
61  {^d&x(^d,n) = xprobmin^d + rrd(n+1,^d) * (xprobmax^d - xprobmin^d)\}
62  end do
63  else
64  call usr_create_particles(num_particles, x, v, q, m, follow)
65  end if
66  end if
67 
68  call mpi_bcast(x,3*num_particles,mpi_double_precision,0,icomm,ierrmpi)
69  call mpi_bcast(follow,num_particles,mpi_logical,0,icomm,ierrmpi)
70 
71  nparticles_local = 0
72 
73  do n=1,num_particles
74  call find_particle_ipe(x(:,n),igrid,ipe_particle)
75  particle(n)%igrid = igrid
76  particle(n)%ipe = ipe_particle
77 
78  if(ipe_particle == mype) then
79  check = .true.
80 
81  ! Check for user-defined modifications or rejection conditions
82  if (associated(usr_check_particle)) call usr_check_particle(igrid, x(:,n), v(:,n), q(n), m(n), follow(n), check)
83  if (check) then
85  else
86  cycle
87  end if
88 
89  nparticles_local = nparticles_local + 1
90 
91  allocate(particle(n)%self)
92  particle(n)%self%follow = follow(n)
93  particle(n)%self%index = n
94  particle(n)%self%time = global_time
95  particle(n)%self%dt = 0.0d0
96  particle(n)%self%x = 0.d0
97  particle(n)%self%x(:) = x(:,n)
98  w=ps(igrid)%w
99  call phys_to_primitive(ixg^ll,ixg^ll,w,ps(igrid)%x)
100  do idir=1,ndir
101  call interpolate_var(igrid,ixg^ll,ixm^ll,&
102  w(ixg^t,iw_mom(idir)),ps(igrid)%x,x(:,n),v(idir,n))
103  end do
104  particle(n)%self%u(:) = 0.d0
105  particle(n)%self%u(1:ndir) = v(1:ndir,n)
106  allocate(particle(n)%payload(npayload))
107  ! Compute default and user-defined payloads
108  call advect_update_payload(igrid,ps(igrid)%w,pso(igrid)%w,ps(igrid)%x,x(:,n),v(:,n),q(n),m(n),defpayload,ndefpayload,0.d0)
109  particle(n)%payload(1:ndefpayload)=defpayload
110  if (associated(usr_update_payload)) then
111  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)
112  particle(n)%payload(ndefpayload+1:npayload)=usrpayload
113  end if
114  end if
115 
116  end do
117 
118  call mpi_allreduce(nparticles_local,nparticles,1,mpi_integer,mpi_sum,icomm,ierrmpi)
119 
120  end subroutine advect_create_particles
121 
124 
125  integer :: igrid, iigrid, idir
126  double precision, dimension(ixG^T,1:nw) :: w
127 
128  do iigrid=1,igridstail; igrid=igrids(iigrid);
129 
130  ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
131 
132  gridvars(igrid)%w(ixg^t,1:ngridvars) = 0.0d0
133  w(ixg^t,1:nw) = ps(igrid)%w(ixg^t,1:nw)
134  call phys_to_primitive(ixg^ll,ixg^ll,w,ps(igrid)%x)
135  ! fill with velocity:
136  gridvars(igrid)%w(ixg^t,vp(:)) = w(ixg^t,iw_mom(:))
137 
138  if(time_advance) then
139  gridvars(igrid)%wold(ixg^t,1:ngridvars) = 0.0d0
140  w(ixg^t,1:nw) = pso(igrid)%w(ixg^t,1:nw)
141  call phys_to_primitive(ixg^ll,ixg^ll,w,ps(igrid)%x)
142  gridvars(igrid)%wold(ixg^t,vp(:)) = w(ixg^t,iw_mom(:))
143  end if
144 
145  end do
146 
147  end subroutine advect_fill_gridvars
148 
149  subroutine advect_integrate_particles(end_time)
150  ! this solves dx/dt=v for particles
151  use mod_odeint
154  double precision, intent(in) :: end_time
155 
156  double precision, dimension(1:ndir) :: v, x
157  double precision :: defpayload(ndefpayload)
158  double precision :: usrpayload(nusrpayload)
159  double precision :: tloc, tlocnew, dt_p, h1
160  double precision,parameter :: eps=1.0d-6, hmin=1.0d-8
161  integer :: ipart, iipart, igrid
162  integer :: nok, nbad, ierror
163 
164  do iipart=1,nparticles_active_on_mype
165  ipart = particles_active_on_mype(iipart);
166  dt_p = advect_get_particle_dt(particle(ipart), end_time)
167  particle(ipart)%self%dt = dt_p
168 
169  igrid = particle(ipart)%igrid
170  igrid_working = igrid
171  tloc = particle(ipart)%self%time
172  x(1:ndir) = particle(ipart)%self%x(1:ndir)
173  tlocnew = tloc+dt_p
174 
175  ! Position update
176  ! Simple forward Euler start
177  !call get_vec_advect(igrid,x,tloc,v,vp(1),vp(ndir))
178  !particle(ipart)%self%u(1:ndir) = v(1:ndir)
179  !particle(ipart)%self%x(1:ndir) = particle(ipart)%self%x(1:ndir) &
180  ! + dt_p * v(1:ndir)
181  ! Simple forward Euler end
182 
183  ! Adaptive stepwidth RK4:
184  h1 = dt_p/2.0d0
185  call odeint(x,ndir,tloc,tlocnew,eps,h1,hmin,nok,nbad,derivs_advect,rkqs,ierror)
186 
187  if (ierror /= 0) then
188  print *, "odeint returned error code", ierror
189  print *, "1 means hmin too small, 2 means MAXSTP exceeded"
190  print *, "Having a problem with particle", iipart
191  end if
192 
193  particle(ipart)%self%x(1:ndir) = x(1:ndir)
194 
195  ! Velocity update
196  call get_vec_advect(igrid,x,tlocnew,v,vp(1),vp(ndir))
197  particle(ipart)%self%u(1:ndir) = v(1:ndir)
198 
199  ! Time update
200  particle(ipart)%self%time = tlocnew
201 
202  ! Update payload
203  call advect_update_payload(igrid,ps(igrid)%w,pso(igrid)%w,ps(igrid)%x,x,v,0.d0,0.d0,defpayload,ndefpayload,tlocnew)
204  particle(ipart)%payload(1:ndefpayload) = defpayload
205  if (associated(usr_update_payload)) then
206  call usr_update_payload(igrid,ps(igrid)%w,pso(igrid)%w,ps(igrid)%x,x,v,0.d0,0.d0,usrpayload,nusrpayload,tlocnew)
207  end if
208  particle(ipart)%payload(ndefpayload+1:npayload) = usrpayload
209 
210  end do
211 
212  end subroutine advect_integrate_particles
213 
214  !> Payload update
215  subroutine advect_update_payload(igrid,w,wold,xgrid,xpart,upart,qpart,mpart,mypayload,mynpayload,particle_time)
217  integer, intent(in) :: igrid,mynpayload
218  double precision, intent(in) :: w(ixG^T,1:nw),wold(ixG^T,1:nw)
219  double precision, intent(in) :: xgrid(ixG^T,1:ndim),xpart(1:ndir),upart(1:ndir),qpart,mpart,particle_time
220  double precision, intent(out) :: mypayload(mynpayload)
221  double precision :: rho, rho1, rho2, td
222 
223  td = (particle_time - global_time) / dt
224 
225  ! Payload 1 is density
226  if (mynpayload > 0 ) then
227  if (.not.time_advance) then
228  call interpolate_var(igrid,ixg^ll,ixm^ll,w(ixg^t,iw_rho),xgrid,xpart,rho)
229  else
230  call interpolate_var(igrid,ixg^ll,ixm^ll,wold(ixg^t,iw_rho),xgrid,xpart,rho1)
231  call interpolate_var(igrid,ixg^ll,ixm^ll,w(ixg^t,iw_rho),xgrid,xpart,rho2)
232  rho = rho1 * (1.0d0 - td) + rho2 * td
233  end if
234  mypayload(1) = rho * w_convert_factor(1)
235  end if
236 
237  end subroutine advect_update_payload
238 
239  subroutine derivs_advect(t_s,x,dxdt)
241  double precision :: t_s, x(ndir)
242  double precision :: dxdt(ndir)
243 
244  double precision :: v(ndir)
245 
246  call get_vec_advect(igrid_working,x,t_s,v,vp(1),vp(ndir))
247  dxdt(:) = v(:)
248 
249  end subroutine derivs_advect
250 
251  pure function advect_get_particle_dt(partp, end_time) result(dt_p)
253  use mod_geometry
254  type(particle_ptr), intent(in) :: partp
255  double precision, intent(in) :: end_time
256  double precision :: dt_p
257  integer :: ipart, iipart, nout
258  double precision :: tout, dt_cfl
259  double precision :: v(1:ndir)
260 
261  if (const_dt_particles > 0) then
262  dt_p = const_dt_particles
263  return
264  end if
265 
266  ! make sure we step only one cell at a time:
267  v(1:ndir)=abs(partp%self%u(1:ndir))
268 
269  ! convert to angular velocity:
270  if(coordinate ==cylindrical.and.phi_>0) v(phi_) = abs(v(phi_)/partp%self%x(r_))
271 
272  dt_cfl = min({rnode(rpdx^d_,partp%igrid)/v(^d)},bigdouble)
273 
274  if(coordinate ==cylindrical.and.phi_>0) then
275  ! phi-momentum leads to radial velocity:
276  if(phi_ .gt. ndim) dt_cfl = min(dt_cfl, &
277  sqrt(rnode(rpdx1_,partp%igrid)/partp%self%x(r_)) &
278  / v(phi_))
279  ! limit the delta phi of the orbit (just for aesthetic reasons):
280  dt_cfl = min(dt_cfl,0.1d0/v(phi_))
281  ! take some care at the axis:
282  dt_cfl = min(dt_cfl,(partp%self%x(r_)+smalldouble)/v(r_))
283  end if
284 
285  dt_p = dt_cfl
286 
287  ! Make sure we don't advance beyond end_time
288  call limit_dt_endtime(end_time - partp%self%time, dt_p)
289 
290  end function advect_get_particle_dt
291 
292  subroutine get_vec_advect(igrid,x,tloc,var,ibeg,iend)
294 
295  integer,intent(in) :: igrid, ibeg, iend
296  double precision,dimension(ndir), intent(in) :: x
297  double precision, intent(in) :: tloc
298  double precision,dimension(iend-ibeg+1), intent(out) :: var
299  double precision,dimension(iend-ibeg+1) :: e1, e2
300  integer :: ivar, iloc
301  double precision :: td
302 
303  if(.not.time_advance) then
304  do ivar=ibeg,iend
305  iloc = ivar-ibeg+1
306  call interpolate_var(igrid,ixg^ll,ixm^ll,gridvars(igrid)%w(ixg^t,ivar),ps(igrid)%x(ixg^t,1:ndim),x,var(iloc))
307  end do
308  else
309  td = (tloc - global_time) / dt
310  do ivar=ibeg,iend
311  iloc = ivar-ibeg+1
312  call interpolate_var(igrid,ixg^ll,ixm^ll,gridvars(igrid)%wold(ixg^t,ivar),ps(igrid)%x(ixg^t,1:ndim),x,e1(iloc))
313  call interpolate_var(igrid,ixg^ll,ixm^ll,gridvars(igrid)%w(ixg^t,ivar),ps(igrid)%x(ixg^t,1:ndim),x,e2(iloc))
314  var(iloc) = e1(iloc) * (1.0d0 - td) + e2(iloc) * td
315  end do
316  end if
317 
318  end subroutine get_vec_advect
319 
320 end module mod_particle_advect
Module with geometry-related routines (e.g., divergence, curl)
Definition: mod_geometry.t:2
integer coordinate
Definition: mod_geometry.t:6
integer, parameter cylindrical
Definition: mod_geometry.t:9
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.
integer, parameter ndim
Number of spatial dimensions for grid variables.
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.
integer r_
Indices for cylindrical coordinates FOR TESTS, negative value when not used:
double precision, dimension(ndim) dxlevel
This module packages odeint from numerical recipes.
Definition: mod_odeint.t:2
subroutine odeint(ystart, nvar, x1, x2, eps, h1, hmin, nok, nbad, derivs, rkqs, ierror)
Definition: mod_odeint.t:12
subroutine rkqs(y, dydx, n, x, htry, eps, yscal, hdid, hnext, derivs)
Definition: mod_odeint.t:115
Tracer for advected particles moving with fluid flows By Jannis Teunissen, Bart Ripperda,...
subroutine, public advect_init()
pure double precision function advect_get_particle_dt(partp, end_time)
subroutine advect_integrate_particles(end_time)
subroutine derivs_advect(t_s, x, dxdt)
subroutine advect_update_payload(igrid, w, wold, xgrid, xpart, upart, qpart, mpart, mypayload, mynpayload, particle_time)
Payload update.
subroutine get_vec_advect(igrid, x, tloc, var, ibeg, iend)
subroutine, public advect_create_particles()
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.
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)
Module with all the methods that users can customize in AMRVAC.
procedure(check_particle), pointer usr_check_particle
procedure(create_particles), pointer usr_create_particles
procedure(update_payload), pointer usr_update_payload