10 integer,
parameter :: Boris=1, vay=2, hc=3, lm=4
21 if(physics_type/=
'mhd')
call mpistop(
"Lorentz particles need magnetic field!")
22 if(
ndir/=3)
call mpistop(
"Lorentz particles need ndir=3!")
66 case(
'HC',
'hc',
'higueracary')
68 case(
'LM',
'lm',
'lapentamarkidis')
83 integer :: n, idir, igrid, ipe_particle, nparticles_local
84 double precision :: lfac
102 rrd(n,idir) =
rng%unif_01()
106 {^
d&x(^
d,n) = xprobmin^
d + rrd(n+1,^
d) * (xprobmax^
d - xprobmin^
d)\}
128 if (ipe_particle ==
mype)
then
139 nparticles_local = nparticles_local + 1
173 integer :: igrid, iigrid, idir, idim
174 double precision,
dimension(ixG^T,1:ndir) :: beta
175 double precision,
dimension(ixG^T,1:nw) :: w,wold
176 double precision :: current(ixG^T,7-2*ndir:3)
178 double precision,
dimension(ixG^T,1:ndir) :: vE, bhat
179 double precision,
dimension(ixG^T) :: kappa, kappa_B, absB, tmp
191 double precision,
intent(in) :: end_time
192 double precision :: defpayload(ndefpayload)
193 double precision :: usrpayload(nusrpayload)
194 integer :: ipart, iipart, igrid
195 double precision :: lfac, q, m, dt_p
196 double precision :: xp(ndir), xpm(ndir), xpc(ndir), xpcm(ndir)
197 double precision :: up(ndir), upc(ndir), tp
198 double precision,
dimension(ndir) :: b, e, bc, ec, vfluid, current
199 double precision :: rho, rhoold, td
218 xpcm = xpc + 0.5d0 * dt_p * upc/lfac
268 xpc = xpcm + 0.5d0 * dt_p * upc/lfac
284 particle(ipart)%payload(1:ndefpayload) = defpayload
298 double precision,
intent(in) :: e(ndir), b(ndir), q, m, dtp
299 double precision,
intent(inout) :: upart(ndir)
300 double precision :: lfac, cosphi, sinphi, phi1, phi2, r, re, sigma, td
301 double precision,
dimension(ndir) :: emom, uprime, tau, s, tmp, uplus, xcart1, xcart2, ucart2, radmom
302 double precision,
dimension(ndir) :: upartk, vbar, Fk, C1, C2, dupartk
303 double precision :: abserr, tol, lfack, J11, J12, J13, J21, J22, J23, J31, J32, J33
304 double precision :: iJ11, iJ12, iJ13, iJ21, iJ22, iJ23, iJ31, iJ32, iJ33, Det
313 emom = q * e * dtp /(2.0d0 * m)
314 uprime = upart + emom
328 tau = q * b * dtp / (2.0d0 * lfac * m)
329 s = 2.0d0 * tau / (1.0d0+sum(tau(:)**2))
331 call cross(uprime,tau,tmp)
332 call cross(uprime+tmp,s,tmp)
349 emom = q * e * dtp /(2.0d0 * m)
350 tau = q * b * dtp / (2.0d0 * m)
352 call cross(upart,tau,tmp)
353 uprime = upart + 2.d0*emom + tmp/lfac
356 sigma = lfac**2 - sum(tau(:)*tau(:))
357 lfac = sqrt((sigma + sqrt(sigma**2 &
358 + 4.d0 * (sum(tau(:)*tau(:)) + sum(uprime(:)*tau(:)/
c_norm)**2))) / 2.d0)
360 call cross(uprime,tau,tmp)
361 upart = (uprime + sum(uprime(:)*tau(:))*tau/lfac**2 + tmp/lfac) / (1.d0+sum(tau(:)*tau(:))/lfac**2)
366 emom = q * e * dtp /(2.0d0 * m)
367 tau = q * b * dtp / (2.0d0 * m)
368 uprime = upart + emom
371 sigma = lfac**2 - sum(tau(:)*tau(:))
372 lfac = sqrt((sigma + sqrt(sigma**2 &
373 + 4.d0 * (sum(tau(:)*tau(:)) + sum(uprime(:)*tau(:)/
c_norm)**2))) / 2.d0)
375 call cross(uprime,tau,tmp)
376 upart = (uprime + sum(uprime(:)*tau(:))*tau/lfac**2 + tmp/lfac) / (1.d0+sum(tau(:)*tau(:))/lfac**2) &
391 do while(abserr > tol .and. nk < nkmax)
396 vbar = (upart + upartk) / (lfac + lfack)
397 call cross(vbar,tau,tmp)
400 fk = upartk - upart - q*dtp/m * (e + tmp)
403 c1 = (lfack + lfac - upartk(1:
ndim) / lfack /
c_norm**2 * (upartk + upart)) / (lfack + lfac)**2
404 c2 = - upartk / lfack /
c_norm**2 / (lfack + lfac)**2
407 j11 = 1. - q*dtp/m * (c2(1) * (upartk(2) + upart(2)) * tau(3) - c2(1) * (upartk(3) + upart(3)) * tau(2))
408 j12 = - q*dtp/m * (c1(2) * tau(3) - c2(2) * (upartk(3) + upart(3)) * tau(2))
409 j13 = - q*dtp/m * (c2(3) * (upartk(2) + upart(2)) * tau(3) - c1(3) * tau(2))
410 j21 = - q*dtp/m * (- c1(1) * tau(3) + c2(1) * (upartk(3) + upart(3)) * tau(1))
411 j22 = 1. - q*dtp/m * (- c2(2) * (upartk(1) + upart(1)) * tau(3) + c2(2) * (upartk(3) + upart(3)) * tau(1))
412 j23 = - q*dtp/m * (- c2(3) * (upartk(1) + upart(1)) * tau(3) + c1(3) * tau(1))
413 j31 = - q*dtp/m * (c1(1) * tau(2) - c2(1) * (upartk(2) + upart(2)) * tau(1))
414 j32 = - q*dtp/m * (c2(2) * (upartk(1) + upart(1)) * tau(2) - c1(2) * tau(1))
415 j33 = 1. - q*dtp/m * (c2(3) * (upartk(1) + upart(1)) * tau(2) - c2(3) * (upartk(2) + upart(2)) * tau(1))
418 det = j11*j22*j33 + j21*j32*j13 + j31*j12*j23 - j11*j32*j23 - j31*j22*j13 - j21*j12*j33
419 ij11 = (j22*j33 - j23*j32) / det
420 ij12 = (j13*j32 - j12*j33) / det
421 ij13 = (j12*j23 - j13*j22) / det
422 ij21 = (j23*j31 - j21*j33) / det
423 ij22 = (j11*j33 - j13*j31) / det
424 ij23 = (j13*j21 - j11*j23) / det
425 ij31 = (j21*j32 - j22*j31) / det
426 ij32 = (j12*j31 - j11*j32) / det
427 ij33 = (j11*j22 - j12*j21) / det
430 dupartk(1) = - (ij11 * fk(1) + ij12 * fk(2) + ij13 * fk(3))
431 dupartk(2) = - (ij21 * fk(1) + ij22 * fk(2) + ij23 * fk(3))
432 dupartk(3) = - (ij31 * fk(1) + ij32 * fk(2) + ij33 * fk(3))
435 upartk=upartk+dupartk
436 abserr=sqrt(sum(dupartk(:)*dupartk(:)))
452 integer,
intent(in) :: igrid,mynpayload
453 double precision,
intent(in) :: xpart(1:ndir),upart(1:ndir),qpart,mpart,particle_time
454 double precision,
intent(out) :: mypayload(mynpayload)
455 double precision :: b(3), e(3), tmp(3), lfac, vfluid(3), current(3), rho, rhoold, td
457 call get_bfield(igrid, xpart, particle_time, b)
458 call get_efield(igrid, xpart, particle_time, e)
494 if (mynpayload > 0)
then
500 if (mynpayload > 1)
then
501 call cross(upart,b,tmp)
502 tmp = tmp / sqrt(sum(b(:)**2))
503 mypayload(2) = mpart/abs(qpart)*sqrt(sum(tmp(:)**2)) / sqrt(sum(b(:)**2))
507 if (mynpayload > 2)
then
508 mypayload(3) = mpart*sum(tmp(:)**2)/(2.0d0*sqrt(sum(b(:)**2)))
512 if (mynpayload > 3)
then
513 mypayload(4) = sum(e(:)*b(:))
522 double precision,
intent(in) :: end_time
523 double precision :: dt_p
524 integer :: ipart, iipart, nout
525 double precision,
dimension(ndir) :: b,v
526 double precision :: lfac,absb,dt_cfl
527 double precision :: tout
528 double precision,
parameter :: cfl = 0.5d0
535 call get_vec(
bp, partp%igrid,partp%self%x,partp%self%time,b)
536 absb = sqrt(sum(b(:)**2))
541 v(:) = abs(partp%self%u(:) / lfac)
548 dt_cfl = min(bigdouble, &
553 if(
phi_ .gt.
ndim) dt_cfl = min(dt_cfl, &
554 sqrt(
rnode(rpdx1_,partp%igrid)/partp%self%x(
r_)) &
557 dt_cfl = min(dt_cfl,0.1d0/max(v(
phi_), smalldouble))
559 dt_cfl = min(dt_cfl,(partp%self%x(
r_)+smalldouble)/max(v(
r_), smalldouble))
562 dt_cfl = dt_cfl * cfl
565 dt_p = abs(
dtheta * partp%self%m * lfac &
566 / (partp%self%q * absb) )
568 dt_p = min(dt_p, dt_cfl)
Module with geometry-related routines (e.g., divergence, curl)
integer, parameter cylindrical
This module contains definitions of global parameters and variables and some generic functions/subrou...
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 ierrmpi
A global MPI error return code.
double precision c_norm
Normalised speed of light.
double precision, dimension(:,:), allocatable rnode
Corner coordinates.
integer r_
Indices for cylindrical coordinates FOR TESTS, negative value when not used:
Module with shared functionality for all the particle movers.
subroutine partcoord_from_cartesian(xp, xpcart)
Convert to noncartesian coordinates.
pure subroutine limit_dt_endtime(t_left, dt_p)
integer num_particles
Number of particles.
subroutine fix_phi_crossing(xp, igrid)
Fix particle position when crossing the 0,2pi boundary in noncartesian coordinates.
double precision const_dt_particles
If positive, a constant time step for the particles.
integer ngridvars
Number of variables for grid field.
subroutine get_efield(igrid, x, tloc, e)
double precision particles_etah
integer nusrpayload
Number of user-defined payload variables for a particle.
pure subroutine get_lfac(u, lfac)
Get Lorentz factor from relativistic momentum.
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 partvec_from_cartesian(xp, up, upcart)
Convert vector from Cartesian coordinates.
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.
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
subroutine partcoord_to_cartesian(xp, xpcart)
Convert position to Cartesian coordinates.
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
subroutine partvec_to_cartesian(xp, up, upcart)
Convert vector to Cartesian coordinates.
integer ndefpayload
Number of default payload variables for a particle.
integer, dimension(:), allocatable jp
Variable index for current.
integer rhop
Variable index for density.
Particle mover with Newtonian/relativistic Boris scheme for Lorentz dynamics By Jannis Teunissen,...
double precision function lorentz_get_particle_dt(partp, end_time)
subroutine lorentz_update_payload(igrid, xpart, upart, qpart, mpart, mypayload, mynpayload, particle_time)
Update payload subroutine.
subroutine lorentz_fill_gridvars
subroutine, public lorentz_init()
subroutine lorentz_integrate_particles(end_time)
Relativistic particle integrator.
subroutine lorentz_kick(upart, e, b, q, m, dtp)
Momentum update subroutine for full Lorentz dynamics.
subroutine, public lorentz_create_particles()
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