24 integer,
parameter :: rk4=1, ark4=2
37 if (physics_type/=
'mhd')
call mpistop(
"GCA particles need magnetic field!")
38 if (
ndir/=3)
call mpistop(
"GCA particles need ndir=3!")
98 case(
'RK4',
'Rk4',
'rk4')
100 case(
'ARK4',
'ARk4',
'Ark4',
'ark4')
115 double precision :: b(
ndir), u(
ndir), magmom
116 double precision :: bnorm, lfac, vnorm, vperp, vpar
117 integer :: igrid, ipe_particle
118 integer :: n, idir, nparticles_local
133 rrd(n,idir) =
rng%unif_01()
137 {^
d&x(^
d,n) = xprobmin^
d + rrd(n+1,^
d) * (xprobmax^
d - xprobmin^
d)\}
159 if(ipe_particle ==
mype)
then
170 nparticles_local = nparticles_local + 1
186 vnorm = norm2(v(:, n))
187 vpar = sum(v(:, n) * b/bnorm)
188 vperp = sqrt(vnorm**2 - vpar**2)
196 magmom = m(n) * (vperp * lfac)**2 / (2.0d0 * bnorm)
222 integer :: igrid, iigrid, idir, idim
223 double precision,
dimension(ixG^T,1:ndir) :: beta
224 double precision,
dimension(ixG^T,1:nw) :: w,wold
225 double precision :: current(ixG^T,7-2*ndir:3)
227 double precision,
dimension(ixG^T,1:ndir) :: vE, bhat
228 double precision,
dimension(ixG^T) :: kappa, kappa_B, absB, tmp
232 do iigrid=1,igridstail; igrid=igrids(iigrid);
233 w(ixg^t,1:nw) = ps(igrid)%w(ixg^t,1:nw)
234 call phys_to_primitive(ixg^
ll,ixg^
ll,w,ps(igrid)%x)
237 absb(ixg^t) = sqrt(sum(gridvars(igrid)%w(ixg^t,
bp(:))**2,dim=
ndim+1))
238 ve(ixg^t,1) = gridvars(igrid)%w(ixg^t,
ep(2)) * gridvars(igrid)%w(ixg^t,
bp(3)) &
239 - gridvars(igrid)%w(ixg^t,
ep(3)) * gridvars(igrid)%w(ixg^t,
bp(2))
240 ve(ixg^t,2) = gridvars(igrid)%w(ixg^t,
ep(3)) * gridvars(igrid)%w(ixg^t,
bp(1)) &
241 - gridvars(igrid)%w(ixg^t,
ep(1)) * gridvars(igrid)%w(ixg^t,
bp(3))
242 ve(ixg^t,3) = gridvars(igrid)%w(ixg^t,
ep(1)) * gridvars(igrid)%w(ixg^t,
bp(2)) &
243 - gridvars(igrid)%w(ixg^t,
ep(2)) * gridvars(igrid)%w(ixg^t,
bp(1))
245 where (absb(ixg^t) .gt. 0.d0)
246 ve(ixg^t,idir) = ve(ixg^t,idir) / absb(ixg^t)**2
248 ve(ixg^t,idir) = 0.d0
252 call mpistop(
"GCA FILL GRIDVARS: vE>c! ABORTING...")
254 if (any(ve .ne. ve))
then
255 call mpistop(
"GCA FILL GRIDVARS: NaNs IN vE! ABORTING...")
259 kappa(ixg^t) = 1.d0/sqrt(1.0d0 - sum(ve(ixg^t,:)**2,dim=
ndim+1)/
c_norm**2)
263 kappa_b(ixg^t) = absb(ixg^t) / kappa(ixg^t)
265 if (any(kappa_b .ne. kappa_b))
then
266 call mpistop(
"GCA FILL GRIDVARS: NaNs IN kappa_B! ABORTING...")
272 gridvars(igrid)%w(ixg^t,
grad_kappa_b(idim)) = tmp(ixg^t)
277 where (absb(ixg^t) .gt. 0.d0)
278 bhat(ixg^t,idir) = gridvars(igrid)%w(ixg^t,
bp(idir)) / absb(ixg^t)
280 bhat(ixg^t,idir) = 0.d0
283 if (any(bhat .ne. bhat))
then
284 call mpistop(
"GCA FILL GRIDVARS: NaNs IN bhat! ABORTING...")
290 call gradient(bhat(ixg^t,idir),ixg^
ll,ixg^
ll^lsub1,idim,tmp)
292 + bhat(ixg^t,idim) * tmp(ixg^t)
294 + ve(ixg^t,idim) * tmp(ixg^t)
295 call gradient(ve(ixg^t,idir),ixg^
ll,ixg^
ll^lsub1,idim,tmp)
297 + bhat(ixg^t,idim) * tmp(ixg^t)
299 + ve(ixg^t,idim) * tmp(ixg^t)
311 double precision,
intent(in) :: end_time
313 double precision :: lfac, absS
314 double precision :: defpayload(ndefpayload)
315 double precision :: usrpayload(nusrpayload)
316 double precision :: dt_p, tloc, y(ndir+2),dydt(ndir+2),ytmp(ndir+2), euler_cfl, int_factor
317 double precision :: tk, k1(ndir+2),k2(ndir+2),k3(ndir+2),k4(ndir+2)
318 double precision,
dimension(1:ndir) :: x, vE, e, b, bhat, x_new, vfluid, current
319 double precision,
dimension(1:ndir) :: drift1, drift2
320 double precision,
dimension(1:ndir) :: drift3, drift4, drift5, drift6, drift7
321 double precision,
dimension(1:ndir) :: bdotgradb, vEdotgradb, gradkappaB
322 double precision,
dimension(1:ndir) :: bdotgradvE, vEdotgradvE
323 double precision,
dimension(1:ndir) :: gradBdrift, reldrift, bdotgradbdrift
324 double precision,
dimension(1:ndir) :: vEdotgradbdrift, bdotgradvEdrift
325 double precision,
dimension(1:ndir) :: vEdotgradvEdrift
326 double precision :: kappa, Mr, upar, m, absb, gamma, q, mompar, vpar, vEabs
327 double precision :: gradBdrift_abs, reldrift_abs, epar
328 double precision :: bdotgradbdrift_abs, vEdotgradbdrift_abs
329 double precision :: bdotgradvEdrift_abs, vEdotgradvEdrift_abs
330 double precision :: momentumpar1, momentumpar2, momentumpar3, momentumpar4
332 double precision,
parameter :: eps=1.0
d-6
334 double precision :: h1, hmin, h_old
335 integer :: nok, nbad, ic1^D, ic2^D, ierror, nvar
336 integer :: ipart, iipart, seed, ic^D,igrid_particle, ipe_particle, ipe_working
337 logical :: int_choice
338 logical :: BC_applied
351 x(1:ndir) =
particle(ipart)%self%x(1:ndir)
357 y(1:ndir) = x(1:ndir)
358 y(ndir+1) =
particle(ipart)%self%u(1)
359 y(ndir+2) =
particle(ipart)%self%u(2)
361 tk = tloc + dt_p/2.d0
362 y(1:ndir) = x + dt_p*k1(1:ndir)/2.d0
363 y(ndir+1) =
particle(ipart)%self%u(1) + dt_p*k1(ndir+1)/2.d0
365 tk = tloc + dt_p/2.d0
366 y(1:ndir) = x + dt_p*k2(1:ndir)/2.d0
367 y(ndir+1) =
particle(ipart)%self%u(1) + dt_p*k2(ndir+1)/2.d0
370 y(1:ndir) = x + dt_p*k3(1:ndir)
371 y(ndir+1) =
particle(ipart)%self%u(1) + dt_p*k3(ndir+1)
373 y(1:ndir) = x(1:ndir)
374 y(ndir+1) =
particle(ipart)%self%u(1)
375 y = y + dt_p/6.d0*(k1 + 2.d0*k2 + 2.d0*k3 + k4)
376 particle(ipart)%self%x(1:ndir) = y(1:ndir)
377 particle(ipart)%self%u(1) = y(ndir+1)
382 y(1:ndir) = x(1:ndir)
383 y(ndir+1) =
particle(ipart)%self%u(1)
384 y(ndir+2) =
particle(ipart)%self%u(2)
399 y(1:ndir+2) = y(1:ndir+2) + euler_cfl * dt_p * dydt(1:ndir+2)
400 particle(ipart)%self%x(1:ndir) = y(1:ndir)
401 particle(ipart)%self%u(1) = y(ndir+1)
402 particle(ipart)%self%u(2) = y(ndir+2)
431 if (ic1^d .le.
ixglo^d-2 .or. ic2^d .ge.
ixghi^d+2)
then
437 if (ic1^d .eq.
ixglo^d-1 .or. ic2^d .eq.
ixghi^d+1)
then
446 y(1:ndir+2) = ytmp(1:ndir+2)
448 particle(ipart)%self%x(1:ndir) = ytmp(1:ndir)
449 particle(ipart)%self%u(1) = ytmp(ndir+1)
450 particle(ipart)%self%u(2) = ytmp(ndir+2)
454 h1 = dt_p/2.0d0; hmin=1.0d-9; h_old=dt_p/2.0d0
461 if (any(y .ne. y))
then
462 call mpistop(
"NaNs DETECTED IN GCA_INTEGRATE BEFORE ODEINT CALL! ABORTING...")
466 call odeint(y,nvar,tloc,tloc+dt_p,eps,h1,hmin,nok,nbad,
derivs_gca_rk,
rkqs,ierror)
468 if (ierror /= 0)
then
469 print *,
"odeint returned error code", ierror
470 print *,
"1 means hmin too small, 2 means MAXSTP exceeded"
471 print *,
"Having a problem with particle", iipart
474 if (any(y .ne. y))
then
475 call mpistop(
"NaNs DETECTED IN GCA_INTEGRATE AFTER ODEINT CALL! ABORTING...")
482 particle(ipart)%self%x(1:ndir) = y(1:ndir)
483 particle(ipart)%self%u(1) = y(ndir+1)
484 particle(ipart)%self%u(2) = y(ndir+2)
498 absb = sqrt(sum(b(:)**2))
499 bhat(1:ndir) = b(1:ndir) / absb
501 epar = sum(e(:)*bhat(:))
503 call cross(e,bhat,ve)
505 ve(1:ndir) = ve(1:ndir) / absb
506 veabs = sqrt(sum(ve(:)**2))
508 kappa = 1.d0/sqrt(1.0d0 - sum(ve(:)**2)/
c_norm**2)
512 mr = y(ndir+2); upar = y(ndir+1); m=
particle(ipart)%self%m; q=
particle(ipart)%self%q
514 gamma = sqrt(1.0d0+upar**2/
c_norm**2+2.0d0*mr*absb/m/
c_norm**2)*kappa
527 particle(ipart)%payload(1:ndefpayload) = defpayload
541 double precision :: t_s, y(ndir+2)
542 double precision :: dydt(ndir+2)
544 double precision,
dimension(ndir):: vE, b, e, x, bhat, bdotgradb, vEdotgradb, gradkappaB, vfluid, current
545 double precision,
dimension(ndir):: bdotgradvE, vEdotgradvE, u, utmp1, utmp2, utmp3
546 double precision :: upar, Mr, gamma, absb, q, m, epar, kappa
556 x(1:ndir) = y(1:ndir)
561 if (any(x .ne. x))
then
562 write(*,*)
"ERROR IN DERIVS_GCA_RK: NaNs IN X OR Y!"
564 write(*,*)
"y",y(ndir+1:ndir+2)
565 call mpistop(
"ABORTING...")
582 if (any(b .ne. b) .or. any(e .ne. e) &
583 .or. any(bdotgradb .ne. bdotgradb) .or. any(vedotgradb .ne. vedotgradb) &
584 .or. any(gradkappab .ne. gradkappab) .or. any(bdotgradve .ne. bdotgradve) &
585 .or. any(vedotgradve .ne. vedotgradve))
then
586 write(*,*)
"ERROR IN DERIVS_GCA_RK: NaNs IN FIELD QUANTITIES!"
589 write(*,*)
"bdotgradb",bdotgradb
590 write(*,*)
"vEdotgradb",vedotgradb
591 write(*,*)
"gradkappaB",gradkappab
592 write(*,*)
"bdotgradvE",bdotgradve
593 write(*,*)
"vEdotgradvE",vedotgradve
594 call mpistop(
"ABORTING...")
597 absb = sqrt(sum(b(:)**2))
598 if (absb .gt. 0.d0)
then
599 bhat(1:ndir) = b(1:ndir) / absb
603 epar = sum(e(:)*bhat(:))
605 call cross(e,bhat,ve)
606 if (absb .gt. 0.d0)
then
607 ve(1:ndir) = ve(1:ndir) / absb
613 kappa = 1.d0/sqrt(1.0d0 - sum(ve(:)**2)/
c_norm**2)
614 gamma = sqrt(1.0d0+upar**2/
c_norm**2+2.0d0*mr*absb/m/
c_norm**2)*kappa
620 if (absb .gt. 0.d0)
then
621 utmp1(1:ndir) = bhat(1:ndir)/(absb/kappa**2)
625 utmp2(1:ndir) = mr/(gamma*q)*gradkappab(1:ndir) &
626 + m/q* (upar**2/gamma*bdotgradb(1:ndir) + upar*vedotgradb(1:ndir) &
627 + upar*bdotgradve(1:ndir) + gamma*vedotgradve(1:ndir))
629 utmp2(1:ndir) = utmp2(1:ndir) + upar*epar/(gamma)*ve(1:ndir)
632 call cross(utmp1,utmp2,utmp3)
633 u(1:ndir) = ve(1:ndir) + utmp3(1:ndir)
636 dydt(1:ndir) = ( u(1:ndir) + upar/gamma * bhat(1:ndir) )
637 dydt(ndir+1) = q/m*epar - mr/(m*gamma) * sum(bhat(:)*gradkappab(:)) &
638 + sum(ve(:)*(upar*bdotgradb(:)+gamma*vedotgradb(:)))
646 double precision :: t_s, y(ndir+2)
647 double precision :: dydt(ndir+2)
649 double precision,
dimension(ndir):: vE, b, e, x, bhat, bdotgradb, vEdotgradb, gradkappaB, vfluid, current
650 double precision,
dimension(ndir):: bdotgradvE, vEdotgradvE, u, utmp1, utmp2, utmp3
651 double precision :: upar, Mr, gamma, absb, q, m, epar, kappa
658 x(1:ndir) = y(1:ndir)
663 if (any(x .ne. x))
then
664 call mpistop(
"ERROR IN DERIVS_GCA: NaNs IN X! ABORTING...")
681 absb = sqrt(sum(b(:)**2))
682 if (absb .gt. 0.d0)
then
683 bhat(1:ndir) = b(1:ndir) / absb
688 epar = sum(e(:)*bhat(:))
689 call cross(e,bhat,ve)
690 if (absb .gt. 0.d0) ve(1:ndir) = ve(1:ndir) / absb
693 kappa = 1.d0/sqrt(1.0d0 - sum(ve(:)**2)/
c_norm**2)
694 gamma = sqrt(1.0d0+upar**2/
c_norm**2+2.0d0*mr*absb/m/
c_norm**2)*kappa
699 if (absb .gt. 0.d0)
then
700 utmp1(1:ndir) = bhat(1:ndir)/(absb/kappa**2)
704 utmp2(1:ndir) = mr/(gamma*q)*gradkappab(1:ndir) &
705 + m/q* (upar**2/gamma*bdotgradb(1:ndir) + upar*vedotgradb(1:ndir) &
706 + upar*bdotgradve(1:ndir) + gamma*vedotgradve(1:ndir))
708 utmp2(1:ndir) = utmp2(1:ndir) + upar*epar/(gamma)*ve(1:ndir)
711 call cross(utmp1,utmp2,utmp3)
712 u(1:ndir) = ve(1:ndir) + utmp3(1:ndir)
715 dydt(1:ndir) = ( u(1:ndir) + upar/gamma * bhat(1:ndir) )
716 dydt(ndir+1) = q/m*epar - mr/(m*gamma) * sum(bhat(:)*gradkappab(:)) &
717 + sum(ve(:)*(upar*bdotgradb(:)+gamma*vedotgradb(:)))
725 integer,
intent(in) :: igrid,mynpayload
726 double precision,
intent(in) :: xpart(1:ndir),upart(1:ndir),qpart,mpart,particle_time
727 double precision,
intent(out) :: mypayload(mynpayload)
728 double precision,
dimension(1:ndir) :: vE, e, b, bhat, vfluid, current
729 double precision,
dimension(1:ndir) :: drift1, drift2
730 double precision,
dimension(1:ndir) :: drift3, drift4, drift5, drift6, drift7
731 double precision,
dimension(1:ndir) :: bdotgradb, vEdotgradb, gradkappaB
732 double precision,
dimension(1:ndir) :: bdotgradvE, vEdotgradvE
733 double precision,
dimension(1:ndir) :: gradBdrift, reldrift, bdotgradbdrift
734 double precision,
dimension(1:ndir) :: vEdotgradbdrift, bdotgradvEdrift
735 double precision,
dimension(1:ndir) :: vEdotgradvEdrift
736 double precision :: kappa, upar, absb, gamma, vpar, vEabs
737 double precision :: gradBdrift_abs, reldrift_abs, epar
738 double precision :: bdotgradbdrift_abs, vEdotgradbdrift_abs
739 double precision :: bdotgradvEdrift_abs, vEdotgradvEdrift_abs
740 double precision :: momentumpar1, momentumpar2, momentumpar3, momentumpar4
742 call get_bfield(igrid, xpart(1:ndir), particle_time, b)
743 call get_efield(igrid, xpart(1:ndir), particle_time, e)
751 absb = sqrt(sum(b(:)**2))
752 bhat(1:ndir) = b(1:ndir) / absb
753 epar = sum(e(:)*bhat(:))
754 call cross(e,bhat,ve)
755 ve(1:ndir) = ve(1:ndir) / absb
756 veabs = sqrt(sum(ve(:)**2))
758 kappa = 1.d0/sqrt(1.0d0 - sum(ve(:)**2)/
c_norm**2)
762 vpar = upart(1)/upart(3)
771 drift1(1:ndir) = bhat(1:ndir)/(absb/kappa**2)
772 drift2(1:ndir) = upart(2)/(upart(3)*qpart)*gradkappab(1:ndir)
774 call cross(drift1,drift2,gradbdrift)
775 gradbdrift_abs = sqrt(sum(gradbdrift(:)**2))
777 drift3(1:ndir) = upar*epar/upart(3)*ve(1:ndir)
778 call cross(drift1,drift3,reldrift)
779 reldrift_abs = sqrt(sum(reldrift(:)**2))
781 drift4(1:ndir) = mpart/qpart* ( upar**2/upart(3)*bdotgradb(1:ndir))
782 call cross(drift1,drift4,bdotgradbdrift)
783 bdotgradbdrift_abs = sqrt(sum(bdotgradbdrift(:)**2))
785 drift5(1:ndir) = mpart/qpart* ( upar*vedotgradb(1:ndir))
786 call cross(drift1,drift5,vedotgradbdrift)
787 vedotgradbdrift_abs = sqrt(sum(vedotgradbdrift(:)**2))
789 drift6(1:ndir) = mpart/qpart* ( upar*bdotgradve(1:ndir))
790 call cross(drift1,drift6,bdotgradvedrift)
791 bdotgradvedrift_abs = sqrt(sum(bdotgradvedrift(:)**2))
793 drift7(1:ndir) = mpart/qpart* (upart(3)*vedotgradve(1:ndir))
794 call cross(drift1,drift7,vedotgradvedrift)
795 vedotgradvedrift_abs = sqrt(sum(vedotgradvedrift(:)**2))
797 momentumpar1 = qpart/mpart*epar
798 momentumpar2 = -(upart(2)/mpart/upart(3))*sum(bhat(:)*gradkappab(:))
799 momentumpar3 = upar*sum(ve(:)*bdotgradb(:))
800 momentumpar4 = upart(3)*sum(ve(:)*vedotgradb(:))
803 if (mynpayload > 0)
then
805 mypayload(1) = sqrt(2.0d0*mpart*upart(2)*absb)/abs(qpart*absb)
807 if (mynpayload > 1)
then
809 mypayload(2) = atan2(sqrt((2.0d0*upart(2)*absb)/(mpart*upart(3)**2)),vpar)
811 if (mynpayload > 2)
then
813 mypayload(3) = sqrt((2.0d0*upart(2)*absb)/(mpart*upart(3)**2))
815 if (mynpayload > 3)
then
817 mypayload(4) = momentumpar1
819 if (mynpayload > 4)
then
821 mypayload(5) = momentumpar2
823 if (mynpayload > 5)
then
825 mypayload(6) = momentumpar3
827 if (mynpayload > 6)
then
829 mypayload(7) = momentumpar4
831 if (mynpayload > 7)
then
835 if (mynpayload > 8)
then
837 mypayload(9) = gradbdrift_abs
839 if (mynpayload > 9)
then
841 mypayload(10) = reldrift_abs
843 if (mynpayload > 10)
then
845 mypayload(11) = bdotgradbdrift_abs
847 if (mynpayload > 11)
then
849 mypayload(12) = vedotgradbdrift_abs
851 if (mynpayload > 12)
then
853 mypayload(13) = bdotgradvedrift_abs
855 if (mynpayload > 13)
then
857 mypayload(14) = vedotgradvedrift_abs
866 double precision,
intent(in) :: end_time
867 double precision :: dt_p
869 double precision :: tout, dt_particles_mype, dt_cfl0, dt_cfl1, dt_a
870 double precision :: dxmin,
vp, a, gammap
872 double precision :: ap0, ap1, dt_cfl_ap0, dt_cfl_ap1, dt_cfl_ap
873 double precision :: dt_euler, dt_tmp
875 double precision :: cfl, uparcfl
876 double precision,
parameter :: uparmin=1.0
d-6*const_c
877 integer :: ipart, iipart, nout, ic^
d, igrid_particle, ipe_particle, ipe
878 logical :: bc_applied
890 dt_tmp = (end_time - partp%self%time)
891 if(dt_tmp .le. 0.0d0) dt_tmp = smalldouble
899 y(
ndir+1) = partp%self%u(1)
900 y(
ndir+2) = partp%self%u(2)
910 vp = sqrt(sum(v(:)**2))
912 dt_cfl0 = dxmin / max(
vp, smalldouble)
913 dt_cfl_ap0 = uparcfl * abs(max(abs(y(
ndir+1)),uparmin) / max(abs(ap0), smalldouble))
918 dt_euler = min(dt_tmp,dt_cfl0,dt_cfl_ap0)
922 partp%self%u(1) = y(
ndir+1)
923 partp%self%u(2) = y(
ndir+2)
938 vp = sqrt(sum(v(:)**2))
940 dt_cfl1 = dxmin / max(
vp, smalldouble)
941 dt_cfl_ap1 = uparcfl * abs(max(abs(y(
ndir+1)),uparmin) / max(abs(ap1), smalldouble))
944 dt_tmp = min(dt_euler, dt_cfl1, dt_cfl_ap1)
947 partp%self%u(1) = ytmp(
ndir+1)
948 partp%self%u(2) = ytmp(
ndir+2)
Module with geometry-related routines (e.g., divergence, curl)
subroutine gradient(q, ixIL, ixOL, idir, gradq)
Calculate gradient of a scalar q within ixL in direction idir.
This module contains definitions of global parameters and variables and some generic functions/subrou...
integer ixghi
Upper index of grid block arrays.
double precision global_time
The global simulation time.
integer, parameter ndim
Number of spatial dimensions for grid variables.
integer, parameter rpxmin
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 ierrmpi
A global MPI error return code.
double precision c_norm
Normalised speed of light.
double precision, dimension(:,:), allocatable rnode
Corner coordinates.
integer nghostcells
Number of ghost cells surrounding a grid.
integer, parameter ixglo
Lower index of grid block arrays (always 1)
This module packages odeint from numerical recipes.
subroutine odeint(ystart, nvar, x1, x2, eps, h1, hmin, nok, nbad, derivs, rkqs, ierror)
subroutine rkqs(y, dydx, n, x, htry, eps, yscal, hdid, hnext, derivs)
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.
logical function particle_in_igrid(ipart, igrid)
Quick check if particle is still in igrid.
integer ngridvars
Number of variables for grid field.
subroutine get_efield(igrid, x, tloc, e)
integer nusrpayload
Number of user-defined payload variables for a particle.
subroutine find_particle_ipe(x, igrid_particle, ipe_particle)
pure subroutine get_lfac_from_velocity(v, lfac)
Get Lorentz factor from velocity.
integer integrator
Integrator to be used for particles.
integer, dimension(:), allocatable ep
Variable index for electric field.
subroutine fill_gridvars_default
This routine fills the particle grid variables with the default for mhd, i.e. only E and B.
subroutine get_bfield(igrid, x, tloc, b)
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.
logical relativistic
Use a relativistic particle mover?
character(len=name_len) integrator_type_particles
String describing the particle integrator type.
subroutine get_vec(ix, igrid, x, tloc, vec)
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.
integer, dimension(:), allocatable bp
Variable index for magnetic field.
type(particle_ptr), dimension(:), allocatable particle
integer, dimension(:), allocatable vp
Variable index for fluid velocity.
procedure(sub_integrate), pointer particles_integrate
subroutine cross(a, b, c)
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:
integer, dimension(:), allocatable jp
Variable index for current.
integer ipart_working
set the current ipart for the particle integrator:
double precision particles_cfl
Particle CFL safety factor.
integer rhop
Variable index for density.
Particle mover with Newtonian/relativistic Guiding Center Approximation (GCA) By Jannis Teunissen,...
integer, dimension(:), allocatable, public, protected b_dot_grad_b
Variable index for (B . grad)B (curvature B drift)
subroutine derivs_gca(t_s, y, dydt)
subroutine derivs_gca_rk(t_s, y, dydt)
subroutine, public gca_init()
double precision function gca_get_particle_dt(partp, end_time)
subroutine, public gca_create_particles()
subroutine gca_fill_gridvars
integer, dimension(:), allocatable, public, protected b_dot_grad_ve
Variable index for polarization drift.
integer, dimension(:), allocatable, public, protected ve_dot_grad_ve
Variable index for polarization drift.
subroutine gca_integrate_particles(end_time)
integer, dimension(:), allocatable, public, protected grad_kappa_b
Variable index for gradient B, with relativistic correction 1/kappa where kappa = 1/sqrt(1 - E_perp^2...
subroutine gca_update_payload(igrid, xpart, upart, qpart, mpart, mypayload, mynpayload, particle_time)
Update payload subroutine.
integer, dimension(:), allocatable, public, protected ve_dot_grad_b
Variable index for curvature drift.
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
procedure(particle_fields), pointer usr_particle_fields