Performs spatial interpolation between grids using bilinear interpolation
mpp_mod
fms_mod
constants_mod
horiz_interp_type_mod
call horiz_interp_bilinear_init ( Interp, lon_in, lat_in, lon_out, lat_out, verbose, src_modulo )
lon_in | Longitude (in radians) for source data grid. [real, dimension(:,:)] |
lat_in | Latitude (in radians) for source data grid. [real, dimension(:,:)] |
lon_out | Longitude (in radians) for source data grid. [real, dimension(:,:)] |
lat_out | Latitude (in radians) for source data grid. [real, dimension(:,:)] |
src_modulo | logical variable to indicate if the boundary condition along zonal boundary
is cyclic or not. When true, the zonal boundary condition is cyclic. [logical, optional] |
verbose | flag for the amount of print output. [integer, optional] |
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. [type(horiz_interp_type)] |
call horiz_interp_bilinear ( Interp, data_in, data_out, verbose, mask_in,mask_out, missing_value, missing_permit)
Interp | Derived-type variable containing interpolation indices and weights.
Returned by a previous call to horiz_interp_init. [type(horiz_interp_type)] |
data_in | Input data on source grid. [real, dimension(:,:)] |
verbose | flag for the amount of print output.
verbose = 0, no output; = 1, min,max,means; = 2, still more [integer, optional] |
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(:,:),optional] |
missing_value | Use the missing_value to indicate missing data. [real, optional] |
missing_permit | numbers of points allowed to miss for the bilinear interpolation. The value should be between 0 and 3. |
data_out | Output data on destination grid. [real, dimension(:,:)] |
mask_out | Output mask that specifies whether data was computed. [real, dimension(:,:),optional] |
call horiz_interp_bilinear_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] |