interpol3D_2D_arr Subroutine

public subroutine interpol3D_2D_arr(r, z, phi, f, rout, zout, phiout, fout, err, coeffs)

Performs cylindrical interpolation on a 3D grid of 2D arrays f(:,:,r,z,phi)

Arguments

Type IntentOptional AttributesName
real(kind=Float64), intent(in), dimension(:):: r

R values

real(kind=Float64), intent(in), dimension(:):: z

Z values

real(kind=Float64), intent(in), dimension(:):: phi

Phi values

real(kind=Float64), intent(in), dimension(:,:,:,:,:):: f

Values at r,z,phi: f(:,:,r,z,phi)

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

R value to interpolate

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

Z value to interpolate

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

Phi value to interpolate

real(kind=Float64), intent(out), dimension(:,:):: fout

Interpolant: f(:,:,rout,zout,phiout)

integer, intent(out), optional :: err

Error code

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

Precomputed Interpolation Coefficients


Calls

proc~~interpol3d_2d_arr~~CallsGraph proc~interpol3d_2d_arr interpol3D_2D_arr interface~interpol_coeff 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~~interpol3d_2d_arr~~CalledByGraph proc~interpol3d_2d_arr interpol3D_2D_arr interface~interpol interpol interface~interpol->proc~interpol3d_2d_arr proc~get_distribution get_distribution proc~get_distribution->interface~interpol proc~get_ep_denf get_ep_denf proc~get_ep_denf->interface~interpol 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_weights_mc fida_weights_mc proc~fida_weights_mc->proc~get_ep_denf 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_weights_mc 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 interpol3D_2D_arr(r, z, phi, f, rout, zout, phiout, fout, err, coeffs)
    !+ Performs cylindrical interpolation on a 3D grid of 2D arrays
    !+ f(:,:,r,z,phi)
    real(Float64), dimension(:), intent(in) :: r
        !+ R values
    real(Float64), dimension(:), intent(in) :: z
        !+ Z values
    real(Float64), dimension(:), intent(in) :: phi
        !+ Phi values
    real(Float64), dimension(:,:,:,:,:), intent(in) :: f
        !+ Values at r,z,phi: f(:,:,r,z,phi)
    real(Float64), intent(in)               :: rout
        !+ R value to interpolate
    real(Float64), intent(in)               :: zout
        !+ Z value to interpolate
    real(Float64), intent(in)               :: phiout
        !+ Phi value to interpolate
    real(Float64), dimension(:,:), intent(out)    :: fout
        !+ Interpolant: f(:,:,rout,zout,phiout)
    integer, intent(out), optional          :: err
        !+ Error code
    type(InterpolCoeffs3D), intent(in), optional :: coeffs
        !+ Precomputed Interpolation Coefficients

    type(InterpolCoeffs3D) :: b
    integer :: i, j, k, k2, err_status
    integer :: nphi

    err_status = 1

    nphi = size(phi)

    if(present(coeffs)) then
        b = coeffs
        if(nphi .eq. 1) then
            b%b212 = 0
            b%b222 = 0
            b%b122 = 0
            b%b112 = 0
            b%k = 1
        endif
        err_status = 0
    else
        call interpol_coeff(r,z,phi,rout,zout,phiout,b,err_status)
    endif

    if(err_status.eq.0) then
        i = b%i
        j = b%j
        k = b%k
        if(nphi .eq. 1) then
            k2 = min(k+1,nphi)
        else
            k2 = k+1
        endif
        fout = b%b111*f(:,:,i,j,k)    + b%b121*f(:,:,i,j+1,k) +   &
               b%b112*f(:,:,i,j,k2)   + b%b122*f(:,:,i,j+1,k2) +  &
               b%b211*f(:,:,i+1,j,k)  + b%b221*f(:,:,i+1,j+1,k) + &
               b%b212*f(:,:,i+1,j,k2) + b%b222*f(:,:,i+1,j+1,k2)
    else
        fout = 0.0
    endif

    if(present(err)) err = err_status

end subroutine interpol3D_2D_arr