interpol2D_2D_arr Subroutine

public subroutine interpol2D_2D_arr(x, y, z, xout, yout, zout, err, coeffs)

Performs bilinear interpolation on a 2D grid of 2D arrays z(:,:,x,y)

Arguments

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

The abscissa values of z

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

The ordinate values of z

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

Values at the abscissa/ordinates: z(:,:,x,y)

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

The abscissa value to interpolate

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

The ordinate value to interpolate

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

Interpolant: z(:,:,xout,yout)

integer, intent(out), optional :: err

Error code

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

Precomputed Linear Interpolation Coefficients


Calls

proc~~interpol2d_2d_arr~~CallsGraph proc~interpol2d_2d_arr interpol2D_2D_arr interface~interpol_coeff interpol_coeff proc~interpol2d_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~~interpol2d_2d_arr~~CalledByGraph proc~interpol2d_2d_arr interpol2D_2D_arr interface~interpol interpol interface~interpol->proc~interpol2d_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 interpol2D_2D_arr(x, y, z, xout, yout, zout, err, coeffs)
    !+ Performs bilinear interpolation on a 2D grid of 2D arrays z(:,:,x,y)
    real(Float64), dimension(:), intent(in)       :: x
        !+ The abscissa values of `z`
    real(Float64), dimension(:), intent(in)       :: y
        !+ The ordinate values of `z`
    real(Float64), dimension(:,:,:,:), intent(in) :: z
        !+ Values at the abscissa/ordinates: z(:,:,x,y)
    real(Float64), intent(in)                     :: xout
        !+ The abscissa value to interpolate
    real(Float64), intent(in)                     :: yout
        !+ The ordinate value to interpolate
    real(Float64), dimension(:,:), intent(out)    :: zout
        !+ Interpolant: z(:,:,xout,yout)
    integer, intent(out), optional                :: err
        !+ Error code
    type(InterpolCoeffs2D), intent(in), optional :: coeffs
        !+ Precomputed Linear Interpolation Coefficients

    type(InterpolCoeffs2D) :: c
    integer :: i, j, err_status

    err_status = 1
    if(present(coeffs)) then
        c = coeffs
        err_status = 0
    else
        call interpol_coeff(x,y,xout,yout,c,err_status)
    endif

    if(err_status.eq.0) then
        i = c%i
        j = c%j
        zout = c%b11*z(:,:,i,j) + c%b12*z(:,:,i,j+1) + c%b21*z(:,:,i+1,j) + c%b22*z(:,:,i+1,j+1)
    else
        zout = 0.0
    endif

    if(present(err)) err = err_status

end subroutine interpol2D_2D_arr