Performs spatial interpolation between rectangular latitude/longitude grids.
fms_mod
constants_mod
use horiz_interp_mod [, only: horiz_interp,
horiz_interp_init,
horiz_interp_end ]
| Name | Type | Value | Units | Description |
|---|---|---|---|---|
| dlon_in | real, dimension(:,:) | --- | --- | delta long |
| dlon_out | real, dimension(:,:) | --- | --- | delta long |
| dsph_in | real, dimension(:,:) | --- | --- | delta sin(lat) |
| dsph_out | real, dimension(:,:) | --- | --- | delta sin(lat) |
| faci | real, dimension(:,:) | --- | --- | weights |
| facj | real, dimension(:,:) | --- | --- | weights |
| ilon | integer, dimension(:,:) | --- | --- | indices |
| jlat | integer, dimension(:,:) | --- | --- | indices |
| wti | real, dimension(:,:,:) | --- | --- | weights for global output grids for bilinear interplation |
| wtj | real, dimension(:,:,:) | --- | --- | weights for global output grids for bilinear interplation |
| i_lon | integer, dimension(:,:,:) | --- | --- | indices |
| j_lat | integer, dimension(:,:,:) | --- | --- | indices |
call horiz_interp ( Interp, data_in, data_out, verbose, mask_in, mask_out )
call horiz_interp ( data_in, blon_in, blat_in, blon_out, blat_out, data_out, verbose, mask_in, mask_out )
| Interp | Derived-type variable containing interpolation indices and weights.
Returned by a previous call to horiz_interp_init. [Derived-type] [Derived-type] [Derived-type] |
| data_in | Input data on input grid defined by grid box edges
initialized in variable Interp. [real, dimension(:,:)] [real, dimension(:,:)] [real, dimension(:,:)] [real, dimension(:,:)] [real, dimension(:,:)] [real, dimension(:,:,:)] |
| verbose | flag for the amount of print output
verbose = 0, no output; = 1, min,max,means; = 2, still more [integer] [integer] [integer] [integer] [integer] |
| mask_in | Input mask, must be the same size as the input data.
The real value of mask_in must be in the range (0.,1.).
Set mask_in=0.0 for data points that should not be used
or have missing data. [real, dimension(:,:)] [real, dimension(:,:)] [real, dimension(:,:)] [real, dimension(:,:)] [real, dimension(:,:)] [real, dimension(:,:,:)] |
| data_out | Output data on output grid defined by grid box edges
initialized in variable Interp. [real, dimension(:,:)] [real, dimension(:,:)] [real, dimension(:,:)] [real, dimension(:,:)] [real, dimension(:,:)] [real, dimension(:,:,:)] |
| mask_out | Output mask that specifies whether data was computed. [real, dimension(:,:)] [real, dimension(:,:)] [real, dimension(:,:)] [real, dimension(:,:)] [real, dimension(:,:)] [real, dimension(:,:,:)] |
call horiz_interp_init ( Interp, blon_in, blat_in, blon_out, blat_out, verbose )
| blon_in | The longitude edges (in radians) for input data grid boxes.
The values are for adjacent grid boxes and must increase in
value. If there are M longitude grid boxes there must be
M+1 edge values. [real, dimension(:)] [real, dimension(:)] [real, dimension(:)] |
| blat_in | The latitude edges (in radians) for input data grid boxes.
The values are for adjacent grid boxes and may increase or
decrease in value. If there are N latitude grid boxes there
must be N+1 edge values. [real, dimension(:)] [real, dimension(:)] [real, dimension(:)] |
| blon_out | The longitude edges (in radians) for output data grid boxes.
The edge values may be stored as adjacent values or as pairs
for each grid box. The pairs do not have to be adjacent.
If there are MLON grid boxes in the output grid, then blon_out
is dimensioned by MLON+1 or (MLON,2). [real, dimension(:,:)] [real, dimension(:,:)] [real, dimension(:)] |
| blat_out | The latitude edges (in radians) for output data grid boxes.
The edge values may be stored as adjacent values or as pairs
for each grid box. The pairs do not have to be adjacent.
If there are NLAT grid boxes in the output grid, then blat_out
is dimensioned by NLAT+1 or (NLAT,2). [real, dimension(:,:)] [real, dimension(:,:)] [real, dimension(:)] |
| verbose | Integer flag that controls the amount of printed output.
verbose = 0, no output; = 1, min,max,means; = 2, still more [integer] |
| Interp | A derived-type variable containing indices and weights used for subsequent interpolations. To reinitialize this variable for a different grid-to-grid interpolation you must first use the "horiz_interp_end" interface. |
call horiz_interp_end ( Interp )
| Interp | A derived-type variable returned by previous call
to horiz_interp_init. The input variable must have
allocated arrays. The returned variable will contain
deallocated arrays. [horiz_interp_type] |
program test
use horiz_interp_mod
implicit none
integer, parameter :: nxi=177, nyi=91, nxo=133, nyo=77 ! resolution
real :: zi(nxi,nyi), zo(nxo,nyo) ! data
real :: xi(nxi+1), yi(nyi+1), xo(nxo+1), yo(nyo+1) ! grid edges
real :: pi, tpi, hpi, dx, dy
constants
hpi = acos(0.0)
pi = hpi*2.0
tpi = hpi*4.0
grid setup: west to east, south to north
dx = tpi/real(nxi); call setaxis (0.,dx,xi); xi(nxi+1) = xi(1)+tpi
dx = tpi/real(nxo); call setaxis (0.,dx,xo); xo(nxo+1) = xo(1)+tpi
dy = pi/real(nyi); call setaxis (-hpi,dy,yi); yi(nyi+1) = hpi
dy = pi/real(nyo); call setaxis (-hpi,dy,yo); yo(nyo+1) = hpi
random data on the input grid
call random_number (zi)
interpolate (flipping y-axis)
call horiz_interp (zi(:,1:nyi:+1), xi, yi(1:nyi+1:+1), xo, yo(1:nyo+1:+1), zo, verbose=2)
call horiz_interp (zi(:,nyi:1:-1), xi, yi(nyi+1:1:-1), xo, yo(1:nyo+1:+1), zo, verbose=2)
call horiz_interp (zi(:,nyi:1:-1), xi, yi(nyi+1:1:-1), xo, yo(nyo+1:1:-1), zo, verbose=2)
call horiz_interp (zi(:,1:nyi:+1), xi, yi(1:nyi+1:+1), xo, yo(nyo+1:1:-1), zo, verbose=2)
contains
set up a sequence of numbers
subroutine setaxis (xo,dx,x)
real, intent(in) :: xo, dx
real, intent(out) :: x(:)
integer :: i
x(1) = xo
do i=2,size(x)
x(i) = x(i-1)+dx
enddo
end subroutine setaxis
end program test
The verbose option for printing out min, max, and mean does so for the data local to each processor (instead on for the entire global field).