mc_fastion Subroutine

public subroutine mc_fastion(ind, fields, eb, ptch, denf)

Samples a Guiding Center Fast-ion distribution function at a given beam_grid index

Arguments

Type IntentOptional AttributesName
integer, intent(in), dimension(3):: ind

beam_grid index

type(LocalEMFields), intent(out) :: fields

Electromagnetic fields at the guiding center

real(kind=Float64), intent(out) :: eb

Energy of the fast ion

real(kind=Float64), intent(out) :: ptch

Pitch of the fast ion

real(kind=Float64), intent(out) :: denf

Fast-ion density at guiding center


Calls

proc~~mc_fastion~~CallsGraph proc~mc_fastion mc_fastion proc~get_fields get_fields proc~mc_fastion->proc~get_fields interface~randu randu proc~mc_fastion->interface~randu proc~get_distribution get_distribution proc~mc_fastion->proc~get_distribution interface~randind randind proc~mc_fastion->interface~randind proc~calc_perp_vectors calc_perp_vectors proc~get_fields->proc~calc_perp_vectors proc~uvw_to_xyz uvw_to_xyz proc~get_fields->proc~uvw_to_xyz proc~in_plasma in_plasma proc~get_fields->proc~in_plasma proc~xyz_to_uvw xyz_to_uvw proc~get_fields->proc~xyz_to_uvw interface~interpol interpol proc~get_distribution->interface~interpol proc~get_distribution->proc~xyz_to_uvw proc~interpol3d_arr interpol3D_arr interface~interpol->proc~interpol3d_arr proc~interpol2d_2d_arr interpol2D_2D_arr interface~interpol->proc~interpol2d_2d_arr proc~interpol2d_arr interpol2D_arr interface~interpol->proc~interpol2d_arr proc~interpol1d_arr interpol1D_arr interface~interpol->proc~interpol1d_arr proc~interpol3d_2d_arr interpol3D_2D_arr interface~interpol->proc~interpol3d_2d_arr proc~in_plasma->proc~xyz_to_uvw interface~interpol_coeff interpol_coeff proc~in_plasma->interface~interpol_coeff proc~cyl_to_uvw cyl_to_uvw proc~in_plasma->proc~cyl_to_uvw proc~cyl_interpol3d_coeff cyl_interpol3D_coeff interface~interpol_coeff->proc~cyl_interpol3d_coeff proc~interpol1d_coeff interpol1D_coeff interface~interpol_coeff->proc~interpol1d_coeff proc~interpol2d_coeff_arr interpol2D_coeff_arr interface~interpol_coeff->proc~interpol2d_coeff_arr proc~interpol2d_coeff interpol2D_coeff interface~interpol_coeff->proc~interpol2d_coeff proc~cyl_interpol3d_coeff_arr cyl_interpol3D_coeff_arr interface~interpol_coeff->proc~cyl_interpol3d_coeff_arr proc~interpol1d_coeff_arr interpol1D_coeff_arr interface~interpol_coeff->proc~interpol1d_coeff_arr proc~interpol3d_arr->interface~interpol_coeff proc~interpol2d_2d_arr->interface~interpol_coeff proc~interpol2d_arr->interface~interpol_coeff proc~interpol1d_arr->interface~interpol_coeff proc~interpol3d_2d_arr->interface~interpol_coeff proc~cyl_interpol3d_coeff->proc~interpol2d_coeff proc~interpol2d_coeff_arr->proc~interpol2d_coeff proc~cyl_interpol3d_coeff_arr->proc~cyl_interpol3d_coeff proc~cyl_interpol3d_coeff_arr->proc~interpol2d_coeff proc~interpol1d_coeff_arr->proc~interpol1d_coeff

Called by

proc~~mc_fastion~~CalledByGraph proc~mc_fastion mc_fastion proc~fida_f fida_f proc~fida_f->proc~mc_fastion proc~npa_f npa_f proc~npa_f->proc~mc_fastion program~fidasim fidasim program~fidasim->proc~fida_f program~fidasim->proc~npa_f

Contents

Source Code


Source Code

subroutine mc_fastion(ind,fields,eb,ptch,denf)
    !+ Samples a Guiding Center Fast-ion distribution function at a given [[libfida:beam_grid]] index
    integer, dimension(3), intent(in) :: ind
        !+ [[libfida:beam_grid]] index
    type(LocalEMFields), intent(out)       :: fields
        !+ Electromagnetic fields at the guiding center
    real(Float64), intent(out)             :: eb
        !+ Energy of the fast ion
    real(Float64), intent(out)             :: ptch
        !+ Pitch of the fast ion
    real(Float64), intent(out)             :: denf
        !+ Fast-ion density at guiding center

    real(Float64), dimension(fbm%nenergy,fbm%npitch) :: fbeam
    real(Float64), dimension(3) :: rg
    real(Float64), dimension(3) :: randomu3
    integer, dimension(2,1) :: ep_ind

    call randu(randomu3)
    rg(1) = beam_grid%xc(ind(1)) + beam_grid%dr(1)*(randomu3(1) - 0.5)
    rg(2) = beam_grid%yc(ind(2)) + beam_grid%dr(2)*(randomu3(2) - 0.5)
    rg(3) = beam_grid%zc(ind(3)) + beam_grid%dr(3)*(randomu3(3) - 0.5)
    denf=0.d0

    call get_fields(fields,pos=rg)
    if(.not.fields%in_plasma) return

    call get_distribution(fbeam,denf,pos=rg, coeffs=fields%b)
    call randind(fbeam,ep_ind)
    call randu(randomu3)
    eb = fbm%energy(ep_ind(1,1)) + fbm%dE*(randomu3(1)-0.5)
    ptch = fbm%pitch(ep_ind(2,1)) + fbm%dp*(randomu3(2)-0.5)

end subroutine mc_fastion