gyro_correction Subroutine

public subroutine gyro_correction(fields, energy, pitch, rp, vp, theta_in)

Calculates gyro correction for Guiding Center MC distribution calculation

Arguments

Type IntentOptional AttributesName
type(LocalEMFields), intent(in) :: fields

Electromagnetic fields at guiding center

real(kind=Float64), intent(in) :: energy

Energy of particle

real(kind=Float64), intent(in) :: pitch

Particle pitch w.r.t the magnetic field

real(kind=Float64), intent(out), dimension(3):: rp

Particle position

real(kind=Float64), intent(out), dimension(3):: vp

Particle velocity

real(kind=Float64), intent(in), optional :: theta_in

Gyro-angle


Calls

proc~~gyro_correction~~CallsGraph proc~gyro_correction gyro_correction interface~randu randu proc~gyro_correction->interface~randu proc~gyro_step gyro_step proc~gyro_correction->proc~gyro_step proc~pitch_to_vec pitch_to_vec proc~gyro_correction->proc~pitch_to_vec proc~cross_product cross_product proc~gyro_step->proc~cross_product

Called by

proc~~gyro_correction~~CalledByGraph proc~gyro_correction gyro_correction proc~fida_weights_mc fida_weights_mc proc~fida_weights_mc->proc~gyro_correction proc~fida_f fida_f proc~fida_f->proc~gyro_correction proc~neutron_mc neutron_mc proc~neutron_mc->proc~gyro_correction proc~pfida_mc pfida_mc proc~pfida_mc->proc~gyro_correction proc~pfida_f pfida_f proc~pfida_f->proc~gyro_correction proc~neutron_f neutron_f proc~neutron_f->proc~gyro_correction proc~fida_mc fida_mc proc~fida_mc->proc~gyro_correction program~fidasim fidasim program~fidasim->proc~fida_weights_mc program~fidasim->proc~fida_f program~fidasim->proc~neutron_mc program~fidasim->proc~pfida_mc program~fidasim->proc~pfida_f program~fidasim->proc~neutron_f program~fidasim->proc~fida_mc

Contents

Source Code


Source Code

subroutine gyro_correction(fields, energy, pitch, rp, vp, theta_in)
    !+ Calculates gyro correction for Guiding Center MC distribution calculation
    type(LocalEMFields), intent(in)          :: fields
        !+ Electromagnetic fields at guiding center
    real(Float64), intent(in)                :: energy
        !+ Energy of particle
    real(Float64), intent(in)                :: pitch
        !+ Particle pitch w.r.t the magnetic field
    real(Float64), dimension(3), intent(out) :: rp
        !+ Particle position
    real(Float64), dimension(3), intent(out) :: vp
        !+ Particle velocity
    real(Float64), intent(in), optional      :: theta_in
        !+ Gyro-angle
    real(Float64), dimension(3) :: vi_norm, r_step
    real(Float64), dimension(1) :: randomu
    real(Float64) :: vabs, theta

    vabs  = sqrt(energy/(v2_to_E_per_amu*inputs%ab))

    if(present(theta_in)) then
        theta = theta_in
    else
        !! Sample gyroangle
        call randu(randomu)
        theta = 2*pi*randomu(1)
    endif

    !! Calculate velocity vector
    call pitch_to_vec(pitch, theta, fields, vi_norm)
    vp = vabs*vi_norm

    !! Move to particle location
    call gyro_step(vp, fields, r_step)
    rp = fields%pos - r_step

end subroutine gyro_correction