get_distribution Subroutine

public subroutine get_distribution(fbeam, denf, pos, ind, coeffs)

Gets Guiding Center distribution at position pos or beam_grid indices ind

Arguments

Type IntentOptional AttributesName
real(kind=Float64), intent(out), dimension(:,:):: fbeam

Guiding Center Fast-ion distribution at pos/ind: F(E,p)

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

Guiding Center Fast-ion density at pos/ind [fast-ions/cm^3]

real(kind=Float64), intent(in), optional dimension(3):: pos

Position in beam grid coordinates

integer(kind=Int32), intent(in), optional dimension(3):: ind

beam_grid indices

type(InterpolCoeffs3D), intent(in), optional :: coeffs

Precomputed interpolation coefficients


Calls

proc~~get_distribution~~CallsGraph proc~get_distribution get_distribution interface~interpol interpol proc~get_distribution->interface~interpol proc~xyz_to_uvw xyz_to_uvw proc~get_distribution->proc~xyz_to_uvw proc~interpol3d_arr interpol3D_arr interface~interpol->proc~interpol3d_arr proc~interpol2d_arr interpol2D_arr interface~interpol->proc~interpol2d_arr proc~interpol1d_arr interpol1D_arr interface~interpol->proc~interpol1d_arr proc~interpol2d_2d_arr interpol2D_2D_arr interface~interpol->proc~interpol2d_2d_arr proc~interpol3d_2d_arr interpol3D_2D_arr interface~interpol->proc~interpol3d_2d_arr interface~interpol_coeff interpol_coeff proc~interpol3d_arr->interface~interpol_coeff proc~interpol2d_arr->interface~interpol_coeff proc~interpol1d_arr->interface~interpol_coeff proc~interpol2d_2d_arr->interface~interpol_coeff proc~interpol3d_2d_arr->interface~interpol_coeff 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~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~~get_distribution~~CalledByGraph proc~get_distribution get_distribution proc~mc_fastion mc_fastion proc~mc_fastion->proc~get_distribution proc~mc_fastion_pass_grid mc_fastion_pass_grid proc~mc_fastion_pass_grid->proc~get_distribution proc~fida_f fida_f proc~fida_f->proc~mc_fastion proc~pnpa_f pnpa_f proc~pnpa_f->proc~mc_fastion_pass_grid proc~pfida_f pfida_f proc~pfida_f->proc~mc_fastion_pass_grid proc~npa_f npa_f proc~npa_f->proc~mc_fastion program~fidasim fidasim program~fidasim->proc~fida_f program~fidasim->proc~pnpa_f program~fidasim->proc~pfida_f program~fidasim->proc~npa_f

Contents

Source Code


Source Code

subroutine get_distribution(fbeam, denf, pos, ind, coeffs)
    !+ Gets Guiding Center distribution at position `pos` or [[libfida:beam_grid]] indices `ind`
    real(Float64), dimension(:,:), intent(out)         :: fbeam
        !+ Guiding Center Fast-ion distribution at `pos`/`ind`: F(E,p)
    real(Float64), intent(out)                         :: denf
        !+ Guiding Center Fast-ion density at `pos`/`ind` [fast-ions/cm^3]
    real(Float64), dimension(3), intent(in), optional  :: pos
        !+ Position in beam grid coordinates
    integer(Int32), dimension(3), intent(in), optional :: ind
        !+ [[libfida:beam_grid]] indices
    type(InterpolCoeffs3D), intent(in), optional       :: coeffs
        !+ Precomputed interpolation coefficients

    real(Float64), dimension(3) :: xyz, uvw
    real(Float64) :: R, Z, Phi
    integer :: err

    if(present(coeffs)) then
        call interpol(fbm%r, fbm%z, fbm%phi, fbm%f, R, Z, Phi, fbeam, err, coeffs)
        call interpol(fbm%r, fbm%z, fbm%phi, fbm%denf, R, Z, Phi, denf, err, coeffs)
    else
        if(present(ind)) call get_position(ind,xyz)
        if(present(pos)) xyz = pos

        !! Convert to machine coordinates
        call xyz_to_uvw(xyz,uvw)
        R = sqrt(uvw(1)*uvw(1) + uvw(2)*uvw(2))
        Z = uvw(3)
        Phi = atan2(uvw(2),uvw(1))

        call interpol(fbm%r, fbm%z, fbm%phi, fbm%f, R, Z, Phi, fbeam, err)
        call interpol(fbm%r, fbm%z, fbm%phi, fbm%denf, R, Z, Phi, denf, err)
    endif

end subroutine get_distribution