rgrd3_wrapper Subroutine

private subroutine rgrd3_wrapper(x, y, z, p, xx, yy, zz, q, intpol, ier)

Wrapper to rgrd3. Allocates the work arrays internally.

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in), dimension(:) :: x

original x

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

original y

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

original z

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

original p(x,y,z)

real(kind=wp), intent(in), dimension(:) :: xx

regridded xx

real(kind=wp), intent(in), dimension(:) :: yy

regridded yy

real(kind=wp), intent(in), dimension(:) :: zz

regridded zz

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

regridded q(xx,yy,zz)

integer, intent(in), dimension(3) :: intpol
integer, intent(out) :: ier
  • 0 : no errors
  • 1-6 : error [see original code]
  • 10 : input vectors are the wrong size
  • 100 : out of memory

Calls

proc~~rgrd3_wrapper~~CallsGraph proc~rgrd3_wrapper regridpack_module::rgrd3_wrapper proc~rgrd3 regridpack_module::rgrd3 proc~rgrd3_wrapper->proc~rgrd3 proc~cubnmx regridpack_module::cubnmx proc~rgrd3->proc~cubnmx proc~cubt3 regridpack_module::cubt3 proc~rgrd3->proc~cubt3 proc~linmx regridpack_module::linmx proc~rgrd3->proc~linmx proc~lint3 regridpack_module::lint3 proc~rgrd3->proc~lint3 proc~cubt2 regridpack_module::cubt2 proc~cubt3->proc~cubt2 proc~lint2 regridpack_module::lint2 proc~cubt3->proc~lint2 proc~lint3->proc~cubt2 proc~lint3->proc~lint2 proc~cubt1 regridpack_module::cubt1 proc~cubt2->proc~cubt1 proc~lint1 regridpack_module::lint1 proc~cubt2->proc~lint1 proc~lint2->proc~cubt1 proc~lint2->proc~lint1

Called by

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

Source Code

    subroutine rgrd3_wrapper(x,y,z,p,xx,yy,zz,q,intpol,ier)

    implicit none

    real(wp),dimension(:),intent(in)              :: x            !! original x
    real(wp),dimension(:),intent(in)              :: y            !! original y
    real(wp),dimension(:),intent(in)              :: z            !! original z
    real(wp),dimension(:,:,:),intent(in)          :: p            !! original p(x,y,z)
    real(wp),dimension(:),intent(in)              :: xx           !! regridded xx
    real(wp),dimension(:),intent(in)              :: yy           !! regridded yy
    real(wp),dimension(:),intent(in)              :: zz           !! regridded zz
    real(wp),dimension(:,:,:),intent(out)         :: q            !! regridded q(xx,yy,zz)
    integer,dimension(3),intent(in)               :: intpol
    integer,intent(out)                           :: ier          !! * 0   : no errors
                                                                  !! * 1-6 : error [see original code]
                                                                  !! * 10  : input vectors are the wrong size
                                                                  !! * 100 : out of memory

    integer :: nx, ny, nz, mx, my, mz
    integer,dimension(3) :: np, nq
    integer :: lw, liw
    integer :: lwx, lwy, lwz
    real(wp),dimension(:),allocatable :: w
    integer,dimension(:),allocatable :: iw
    integer :: ierr1, ierr2

    !get array sizes:

    nx = size(x)
    ny = size(y)
    nz = size(z)
    np(1) = size(p,1)
    np(2) = size(p,2)
    np(3) = size(p,3)

    mx = size(xx)
    my = size(yy)
    mz = size(zz)
    nq(1) = size(q,1)
    nq(2) = size(q,2)
    nq(3) = size(q,3)

    if (nx/=np(1) .or. ny/=np(2) .or. nz/=np(3) .or. mx/=nq(1) .or. my/=nq(2) .or. mz/=nq(3)) then

        !Error: vectors are the wrong size
        ier = 10
        return

    end if

    !allocate work matrices:

    select case(intpol(1))
    case(1)
        lwx = mx
    case(3)
        lwx = 4*mx
        case default
        ier = 6     !Error: invalid intpol value
        return
    end select

    select case(intpol(2))
    case(1)
        lwy = my+2*mx
    case(3)
        lwy = 4*(mx+my)
    end select
    select case(intpol(3))
    case(1)
        lwz = 2*mx*my+mz
    case(3)
        lwz = 4*(mx*my+mz)
    end select

    lw  = lwx + lwy + lwz
    liw = mx + my + mz

    allocate(w(lw),   stat=ierr1)
    allocate(iw(liw), stat=ierr2)

    if (ierr1==0 .and. ierr2==0) then

        !call the main routine:
        call rgrd3(nx,ny,nz,x,y,z,p,mx,my,mz,xx,yy,zz,q,intpol,w,lw,iw,liw,ier)

    else

        !error: out of memory
        ier = 100

    end if

    !clean up:

    deallocate(w)
    deallocate(iw)

    end subroutine rgrd3_wrapper