Wrapper to rgrd2. Allocates the work arrays internally.
Type | Intent | Optional | 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(:,:) | :: | p |
original p(x,y) |
|
real(kind=wp), | intent(in), | dimension(:) | :: | xx |
regridded xx |
|
real(kind=wp), | intent(in), | dimension(:) | :: | yy |
regridded yy |
|
real(kind=wp), | intent(out), | dimension(:,:) | :: | q |
regridded q(xx,yy) |
|
integer, | intent(in), | dimension(2) | :: | intpol | ||
integer, | intent(out) | :: | ier |
|
subroutine rgrd2_wrapper(x,y,p,xx,yy,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) :: p !! original p(x,y) real(wp),dimension(:),intent(in) :: xx !! regridded xx real(wp),dimension(:),intent(in) :: yy !! regridded yy real(wp),dimension(:,:),intent(out) :: q !! regridded q(xx,yy) integer,dimension(2),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 :: lw, liw integer :: nx, ny, mx, my integer,dimension(2) :: np, nq integer :: lwx, lwy real(wp),dimension(:),allocatable :: w integer,dimension(:),allocatable :: iw integer :: ierr1, ierr2 !get array sizes: nx = size(x) ny = size(y) np(1) = size(p,1) np(2) = size(p,2) mx = size(xx) my = size(yy) nq(1) = size(q,1) nq(2) = size(q,2) if (nx/=np(1) .or. ny/=np(2) .or. mx/=nq(1) .or. my/=nq(2)) 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 lw = lwx + lwy liw = mx + my allocate(w(lw), stat=ierr1) allocate(iw(liw), stat=ierr2) if (ierr1==0 .and. ierr2==0) then !call the main routine: call rgrd2(nx,ny,x,y,p,mx,my,xx,yy,q,intpol,w,lw,iw,liw,ier) else !error: out of memory ier = 100 end if !clean up: deallocate(w) deallocate(iw) end subroutine rgrd2_wrapper