get_ep_denf Subroutine

public subroutine get_ep_denf(energy, pitch, denf, pos, ind, coeffs)

Get guiding center fast-ion density at given energy and pitch at position pos or beam_grid indices ind

Arguments

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

Energy [keV]

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

Pitch

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

Fast-ion density [fast-ions/(cm^3dEdp)]

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_ep_denf~~CallsGraph proc~get_ep_denf get_ep_denf interface~interpol interpol proc~get_ep_denf->interface~interpol proc~xyz_to_uvw xyz_to_uvw proc~get_ep_denf->proc~xyz_to_uvw proc~get_position get_position proc~get_ep_denf->proc~get_position proc~cyl_to_xyz cyl_to_xyz proc~get_position->proc~cyl_to_xyz proc~cyl_to_uvw cyl_to_uvw proc~cyl_to_xyz->proc~cyl_to_uvw proc~uvw_to_xyz uvw_to_xyz proc~cyl_to_xyz->proc~uvw_to_xyz

Called by

proc~~get_ep_denf~~CalledByGraph proc~get_ep_denf get_ep_denf proc~cfpd_f cfpd_f proc~cfpd_f->proc~get_ep_denf proc~npa_weights npa_weights proc~npa_weights->proc~get_ep_denf proc~fida_weights_mc fida_weights_mc proc~fida_weights_mc->proc~get_ep_denf proc~fida_weights_los fida_weights_los proc~fida_weights_los->proc~get_ep_denf program~fidasim fidasim program~fidasim->proc~cfpd_f program~fidasim->proc~npa_weights program~fidasim->proc~fida_weights_mc program~fidasim->proc~fida_weights_los

Contents

Source Code


Source Code

subroutine get_ep_denf(energy, pitch, denf, pos, ind, coeffs)
    !+ Get guiding center fast-ion density at given energy and pitch
    !+ at position `pos` or [[libfida:beam_grid]] indices `ind`
    real(Float64), intent(in)                          :: energy
        !+ Energy [keV]
    real(Float64), intent(in)                          :: pitch
        !+ Pitch
    real(Float64), intent(out)                         :: denf
        !+ Fast-ion density [fast-ions/(cm^3*dE*dp)]
    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), dimension(fbm%nenergy,fbm%npitch)  :: fbeam
    integer(Int32), dimension(2) :: epi
    real(Float64) :: R, Phi, Z
    real(Float64) :: dE, dp
    integer :: err

    epi(1) = minloc(abs(fbm%energy - energy),1)
    epi(2) = minloc(abs(fbm%pitch - pitch),1)
    dE = abs(fbm%energy(epi(1)) - energy)
    dp = abs(fbm%pitch(epi(2)) - pitch)

    if((dE.le.fbm%dE).and.(dp.le.fbm%dp)) then
        if(present(coeffs)) then
            call interpol(inter_grid%r, inter_grid%z, inter_grid%phi, fbm%f, R, Z, Phi, fbeam, 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(inter_grid%r, inter_grid%z, inter_grid%phi, fbm%f, R, Z, Phi, fbeam, err)
        endif
        denf = fbeam(epi(1),epi(2))
    else
        denf = 0.0
    endif

end subroutine get_ep_denf