PUBLIC INTERFACE ~ PUBLIC DATA ~ PUBLIC ROUTINES ~ NAMELIST ~ DIAGNOSTIC FIELDS ~ ERROR MESSAGES ~ REFERENCES ~ NOTES

module horiz_interp_mod

Contact:  Bruce Wyman
Reviewers: 
Change History:  WebCVS Log
Last Modified: 2002/06/06 17:58:37


OVERVIEW

Performs spatial interpolation between rectangular latitude/longitude grids.

The interpolation algorithm uses a scheme that conserves the area-weighed integral of the input field. The domain of the output field must lie within the domain of the input field. Latitudes must be range between -¦Ð/2,+¦Ð/2, and longitudes of the input and output grids must be within +/-2¦Ð of each other. If the input and output fields both completely cover the sphere, then the global integrals of both fields will be the same to within machine precision.

There is an optional mask field for missing input data. An optional output mask field may be used in conjunction with the input mask to show where output data exists. A derived-type variable that contains interpolation indices and weights. It is allocated and initialized by calling horiz_interp_init, and it is deallocated by calling horiz_interp_end. The use of this type is recommended for performance when computing multiple interpolations between the same grids. The contents of this data type are private.


OTHER MODULES USED

      fms_mod
constants_mod

PUBLIC INTERFACE

use horiz_interp_mod [, only:  horiz_interp,
horiz_interp_init,
horiz_interp_end ]
horiz_interp:
Subroutine for performing the horizontal interpolation between two grids.
horiz_interp_init:
Allocates space and initializes a derived-type variable that contains pre-computed interpolation indices and weights.
horiz_interp_end:
Deallocates memory used by "horiz_interp_type" variables. Must be called before reinitializing with horiz_interp_init.


PUBLIC DATA

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


PUBLIC ROUTINES

  1. horiz_interp

    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 )
    DESCRIPTION
    Subroutine for performing the horizontal interpolation between two grids. There are two forms of this interface. Form A requires first calling horiz_interp_init, while Form B requires no initialization.


    INPUT
    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(:,:,:)]

    OUTPUT
    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(:,:,:)]

  2. horiz_interp_init

    call horiz_interp_init ( Interp, blon_in, blat_in, blon_out, blat_out, verbose )
    DESCRIPTION
    Allocates space and initializes a derived-type variable that contains pre-computed interpolation indices and weights for improved performance of multiple interpolations between the same grids. This routine does not need to be called if you are doing a single grid-to-grid interpolation.


    INPUT
    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]

    INPUT/OUTPUT
    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.

  3. horiz_interp_end

    call horiz_interp_end ( Interp )
    DESCRIPTION
    Deallocates memory used by "horiz_interp_type" variables. Must be called before reinitializing with horiz_interp_init.


    INPUT/OUTPUT
    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]


DATA SETS

None.


ERROR MESSAGES

FATAL in horiz_interp
size of input array incorrect
The input data array does not match the size of the input grid edges specified. If you are using the initialization interface make sure you have the correct grid size.
FATAL in horiz_interp
size of output array incorrect
The output data array does not match the size of the input grid edges specified. If you are using the initialization interface make sure you have the correct grid size.


REFERENCES

None.


COMPILER SPECIFICS

None.


PRECOMPILER OPTIONS

None.


LOADER OPTIONS

None.


TEST PROGRAM

INTERNAL
       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


KNOWN BUGS

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).



NOTES

Has not been checked with grids that do not cover the sphere.

Has not been checked with the optional mask arguments.

If a latitude or longitude index cannot be found the tolerance used for making this determination may need to be increased. This can be done by increasing the value of module variable num_iters (default 4).


FUTURE PLANS

None.


top