rgrd1u Subroutine

private subroutine rgrd1u(nx, p, mx, q, intpol, w, lw, iw, liw, ier)

subroutine rgrd1u interpolates the nx vector p onto the mx vector q. it is assumed that p and q are values on uniform nx and mx grids which subdivide the same interval (INCLUDING END POINTS). if p and q are values on nonuniform grids and/or if q is defined on a grid which lies within the p grid then subroutine rgrd1 should be used.

method

linear or cubic interpolation (see intpol) is used when the mx uniform grid is not a subgrid of the nx uniform grid (i.e., whenever mx-1 does not divide nx-1). q is set directly from p in the subgrid case.

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: nx

the integer dimension of p. nx > 1 if intpol = 1 or nx > 3 if intpol = 3 is required (see ier = 2).

real(kind=wp), intent(in) :: p(nx)

a real(wp) nx dimensioned vector of given values

integer, intent(in) :: mx

the integer dimension of q. mx > 1 is required (see ier = 1)

real(kind=wp), intent(out), dimension(mx) :: q

a real(wp) mx dimensioned vector of values which are interpolated from p.

integer, intent(in) :: intpol

an integer which sets linear or cubic interpolation as follows:

  • intpol = 1 sets linear interpolation
  • intpol = 3 sets cubic interpolation

values other than 1 or 3 in intpol are not allowed (ier = 4).

real(kind=wp), intent(inout) :: w(lw)

a real(wp) work space of length lw.

integer, intent(in) :: lw

the integer length of the work space w in the routine calling rgrd1u. if mx-1 divides nx-1 then the mx uniform grid is a subgrid of the nx uniform grid. in this case let lwmin = 1. otherwise let lwmin = mx if intpol = 1 or lwmin = mx if intpol = 3. then lw must be greater than or equal to lwmin (see ier=4).

integer, intent(inout) :: iw(liw)

an integer work space of length liw

integer, intent(in) :: liw

the integer length of the integer work space iw in the routine calling rgrd1u. liw must be greater than or equal to mx.

integer, intent(out) :: ier

an integer error flag set as follows:

  • ier = 0 if no errors in input arguments are detected
  • ier = 1 if mx < 2
  • ier = 2 if nx < 2 when intpol=1 or nx < 4 when intpol=3.
  • ier = 3 if intpol is not equal to 1 or 3
  • ier = 4 if lw or liw is too small (insufficient work space)

Calls

proc~~rgrd1u~~CallsGraph proc~rgrd1u regridpack_module::rgrd1u proc~cubnmxu regridpack_module::cubnmxu proc~rgrd1u->proc~cubnmxu proc~cubt1u regridpack_module::cubt1u proc~rgrd1u->proc~cubt1u proc~linmxu regridpack_module::linmxu proc~rgrd1u->proc~linmxu proc~lint1u regridpack_module::lint1u proc~rgrd1u->proc~lint1u

Called by

proc~~rgrd1u~~CalledByGraph proc~rgrd1u regridpack_module::rgrd1u interface~regrid regridpack_module::regrid interface~regrid->proc~rgrd1u

Source Code

    subroutine rgrd1u(nx,p,mx,q,intpol,w,lw,iw,liw,ier)

    implicit none

    integer,intent(in)  :: intpol       !! an integer which sets linear or cubic interpolation as follows:
                                        !!
                                        !! * intpol = 1 sets linear interpolation
                                        !! * intpol = 3 sets cubic interpolation
                                        !!
                                        !! values other than 1 or 3 in intpol are not allowed (ier = 4).
    integer,intent(in)  :: liw          !! the integer length of the integer work space iw in the routine calling rgrd1u.
                                        !! liw must be greater than or equal to mx.
    integer,intent(inout)  :: iw(liw)   !! an integer work space of length liw
    integer,intent(in)  :: lw           !! the integer length of the work space w in the routine calling rgrd1u.
                                        !! if mx-1 divides nx-1 then the mx uniform grid is a subgrid of
                                        !! the nx uniform grid.  in this case let lwmin = 1.  otherwise
                                        !! let lwmin = mx if intpol = 1 or lwmin = mx if intpol = 3.
                                        !! then lw must be greater than or equal to lwmin (see ier=4).
    integer,intent(in) :: nx            !! the integer dimension of p.  nx > 1 if intpol = 1 or
                                        !! nx > 3 if intpol = 3 is required (see ier = 2).
    integer,intent(in)  :: mx           !! the integer dimension of q.  mx > 1 is required (see ier = 1)
    real(wp),intent(in) :: p(nx)        !! a real(wp) nx dimensioned vector of given values
    real(wp),intent(inout) :: w(lw)     !! a real(wp) work space of length lw.
    real(wp),dimension(mx),intent(out) :: q  !! a real(wp) mx dimensioned vector of values which are interpolated from p.
    integer,intent(out) :: ier  !! an integer error flag set as follows:
                                !!
                                !! * ier = 0 if no errors in input arguments are detected
                                !! * ier = 1 if  mx < 2
                                !! * ier = 2 if nx < 2 when intpol=1 or nx < 4 when intpol=3.
                                !! * ier = 3 if intpol is not equal to 1 or 3
                                !! * ier = 4 if lw or liw is too small (insufficient work space)

    integer  :: inmx,isubx,i2,i3,i4,i5,lwmin

    ! check input arguments

    ! check mx
    ier = 1
    if (mx < 2) return

    ! check intpol
    ier = 3
    if (intpol/=1 .and. intpol/=3) return

    ! check nx
    ier = 2
    if (intpol==1 .and. nx<2) return
    if (intpol==3 .and. nx<4) return

    ! set subgrid integer indicator
    inmx = (nx-1)/(mx-1)
    isubx = nx - inmx*(mx-1)

    ! set minimum and check work space
    ier = 4
    if (isubx/=1) then
        if (intpol==1) lwmin = mx
        if (intpol==3) lwmin = 4*mx
    else
        lwmin = 1
    end if
    if (lw < lwmin) return
    if (liw < mx) return

    ! input arguments o.k.

    ier = 0

    ! preset pointers

    i2 = 1
    i3 = 1
    i4 = 1
    i5 = 1
    if (intpol == 1) then
        ! linear interpolation in x
        if (isubx /= 1) then
            call linmxu(nx,mx,iw,w)
        end if
        call lint1u(nx,p,mx,q,iw,w,inmx,isubx)
    else
        ! cubic interpolation in x
        if (isubx /= 1) then
            i2 = 1
            i3 = i2+mx
            i4 = i3+mx
            i5 = i4+mx
            call cubnmxu(nx,mx,iw,w(i2),w(i3),w(i4),w(i5))
        end if
        call cubt1u(nx,p,mx,q,iw,w(i2),w(i3),w(i4),w(i5),inmx,isubx)
    end if

    end subroutine rgrd1u