FMS
2024.03
Flexible Modeling System
|
Performs spatial interpolation between grids. More...
Data Types | |
interface | horiz_interp |
Subroutine for performing the horizontal interpolation between two grids. More... | |
Functions/Subroutines | |
subroutine | horiz_interp_base_2d_ (Interp, data_in, data_out, verbose, mask_in, mask_out, missing_value, missing_permit, err_msg, new_missing_handle) |
subroutine | horiz_interp_base_3d_ (Interp, data_in, data_out, verbose, mask_in, mask_out, missing_value, missing_permit, err_msg) |
Overload of interface HORIZ_INTERP_BASE_2D_ uses 3d arrays for data and mask this allows for multiple interpolations with one call. | |
subroutine, public | horiz_interp_del (Interp) |
Deallocates memory used by "horiz_interp_type" variables. Must be called before reinitializing with horiz_interp_new. More... | |
subroutine, public | horiz_interp_end |
Dummy routine. | |
subroutine, public | horiz_interp_init |
Initialize module and writes version number to logfile.out. | |
subroutine | horiz_interp_new_1d_ (Interp, lon_in, lat_in, lon_out, lat_out, verbose, interp_method, num_nbrs, max_dist, src_modulo, grid_at_center, mask_in, mask_out) |
subroutine | horiz_interp_new_1d_dst_ (Interp, lon_in, lat_in, lon_out, lat_out, verbose, interp_method, num_nbrs, max_dist, src_modulo, mask_in, mask_out, is_latlon_in) |
subroutine | horiz_interp_new_1d_src_ (Interp, lon_in, lat_in, lon_out, lat_out, verbose, interp_method, num_nbrs, max_dist, src_modulo, grid_at_center, mask_in, mask_out, is_latlon_out) |
subroutine | horiz_interp_new_2d_ (Interp, lon_in, lat_in, lon_out, lat_out, verbose, interp_method, num_nbrs, max_dist, src_modulo, mask_in, mask_out, is_latlon_in, is_latlon_out) |
subroutine | horiz_interp_read_weights_ (Interp, weight_filename, lon_out, lat_out, lon_in, lat_in, weight_file_source, interp_method, isw, iew, jsw, jew, nglon, nglat) |
Subroutine for reading a weight file and use it to fill in the horiz interp type for the bilinear interpolation method. More... | |
subroutine | horiz_interp_solo_1d_ (data_in, lon_in, lat_in, lon_out, lat_out, data_out, verbose, mask_in, mask_out, interp_method, missing_value, missing_permit, num_nbrs, max_dist, src_modulo, grid_at_center) |
Interpolates from a rectangular grid to rectangular grid. interp_method can be the value conservative, bilinear or spherical. horiz_interp_new don't need to be called before calling this routine. | |
subroutine | horiz_interp_solo_1d_dst_ (data_in, lon_in, lat_in, lon_out, lat_out, data_out, verbose, mask_in, mask_out, interp_method, missing_value, missing_permit, num_nbrs, max_dist, src_modulo) |
interpolates from any grid to rectangular longitude/latitude grid. interp_method should be "spherical". horiz_interp_new don't need to be called before calling this routine. | |
subroutine | horiz_interp_solo_1d_src_ (data_in, lon_in, lat_in, lon_out, lat_out, data_out, verbose, mask_in, mask_out, interp_method, missing_value, missing_permit, num_nbrs, max_dist, src_modulo, grid_at_center) |
Interpolates from a uniformly spaced grid to any output grid. interp_method can be the value "onservative","bilinear" or "spherical". horiz_interp_new don't need to be called before calling this routine. | |
subroutine | horiz_interp_solo_2d_ (data_in, lon_in, lat_in, lon_out, lat_out, data_out, verbose, mask_in, mask_out, interp_method, missing_value, missing_permit, num_nbrs, max_dist, src_modulo) |
Interpolates from any grid to any grid. interp_method should be "spherical" horiz_interp_new don't need to be called before calling this routine. | |
subroutine | horiz_interp_solo_old_ (data_in, wb, sb, dx, dy, lon_out, lat_out, data_out, verbose, mask_in, mask_out) |
Overloaded version of interface horiz_interp_solo_2. More... | |
logical function | is_lat_lon_ (lon, lat) |
Variables | |
logical | module_is_initialized = .FALSE. |
logical | reproduce_siena = .false. |
Set reproduce_siena = .true. to reproduce siena results. Set reproduce_siena = .false. to decrease truncation error in routine poly_area in file mosaic_util.c. The truncation error of second order conservative remapping might be big for high resolution grid. | |
Performs spatial interpolation between grids.
Creates a 1D horiz_interp_type with the given parameters.
This module can interpolate data from any logically rectangular grid to any logically rectangular grid. Four interpolation schems are used here: conservative, bilinear, bicubic and inverse of square distance weighted. The four interpolation schemes are implemented seperately in horiz_interp_conserver_mod, horiz_interp_blinear_mod, horiz_interp_bicubic_mod and horiz_interp_spherical_mod. bicubic interpolation requires the source grid is regular lon/lat grid. User can choose the interpolation method in the public interface horiz_interp_new through optional argument interp_method, with acceptable value "conservative", "bilinear", "bicubic" and "spherical". The default value is "conservative". 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.
interface horiz_interp_mod::horiz_interp |
Subroutine for performing the horizontal interpolation between two grids.
Subroutine for performing the horizontal interpolation between two grids. There are two forms of this interface. Form A requires first calling horiz_interp_new, while Form B requires no initialization.
Interp | Derived-type variable containing interpolation indices and weights. Returned by a previous call to horiz_interp_new. |
data_in | Input data on source grid. |
verbose | flag for the amount of print output. verbose = 0, no output; = 1, min,max,means; = 2, still more |
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. It is Not needed for spherical regrid. |
missing_value | Use the missing_value to indicate missing data. |
missing_permit | numbers of points allowed to miss for the bilinear interpolation. The value should be between 0 and 3. |
lon_in,lat_in | longitude and latitude (in radians) of source grid. More explanation can be found in the documentation of horiz_interp_new. |
lon_out,lat_out | longitude and latitude (in radians) of destination grid. More explanation can be found in the documentation of horiz_interp_new. |
data_out | Output data on destination grid. |
mask_out | Output mask that specifies whether data was computed. |
FATAL,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,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. |
Definition at line 204 of file horiz_interp.F90.
subroutine, public horiz_interp_mod::horiz_interp_del | ( | type (horiz_interp_type), intent(inout) | Interp | ) |
Deallocates memory used by "horiz_interp_type" variables. Must be called before reinitializing with horiz_interp_new.
[in,out] | interp | A derived-type variable returned by previous call to horiz_interp_new. The input variable must have allocated arrays. The returned variable will contain deallocated arrays |
Definition at line 286 of file horiz_interp.F90.
subroutine horiz_interp_new_1d_ | ( | type(horiz_interp_type), intent(inout) | Interp, |
real(fms_hi_kind_), dimension(:), intent(in) | lon_in, | ||
real(fms_hi_kind_), dimension(:), intent(in) | lat_in, | ||
real(fms_hi_kind_), dimension(:), intent(in) | lon_out, | ||
real(fms_hi_kind_), dimension(:), intent(in) | lat_out, | ||
integer, intent(in), optional | verbose, | ||
character(len=*), intent(in), optional | interp_method, | ||
integer, intent(in), optional | num_nbrs, | ||
real(fms_hi_kind_), intent(in), optional | max_dist, | ||
logical, intent(in), optional | src_modulo, | ||
logical, intent(in), optional | grid_at_center, | ||
real(fms_hi_kind_), dimension(:,:), intent(in), optional | mask_in, | ||
real(fms_hi_kind_), dimension(:,:), intent(inout), optional | mask_out | ||
) |
[in] | mask_in | dummy variable |
[in,out] | mask_out | dummy variable |
Definition at line 22 of file horiz_interp.inc.
subroutine horiz_interp_new_1d_src_ | ( | type(horiz_interp_type), intent(inout) | Interp, |
real(fms_hi_kind_), dimension(:), intent(in) | lon_in, | ||
real(fms_hi_kind_), dimension(:), intent(in) | lat_in, | ||
real(fms_hi_kind_), dimension(:,:), intent(in) | lon_out, | ||
real(fms_hi_kind_), dimension(:,:), intent(in) | lat_out, | ||
integer, intent(in), optional | verbose, | ||
character(len=*), intent(in), optional | interp_method, | ||
integer, intent(in), optional | num_nbrs, | ||
real(fms_hi_kind_), intent(in), optional | max_dist, | ||
logical, intent(in), optional | src_modulo, | ||
logical, intent(in), optional | grid_at_center, | ||
real(fms_hi_kind_), dimension(:,:), intent(in), optional | mask_in, | ||
real(fms_hi_kind_), dimension(:,:), intent(out), optional | mask_out, | ||
logical, intent(in), optional | is_latlon_out | ||
) |
[in] | num_nbrs | minimum number of neighbors |
Definition at line 155 of file horiz_interp.inc.
subroutine horiz_interp_read_weights_ | ( | type(horiz_interp_type), intent(inout) | Interp, |
character(len=*), intent(in) | weight_filename, | ||
real(fms_hi_kind_), dimension(:,:), intent(in) | lon_out, | ||
real(fms_hi_kind_), dimension(:,:), intent(in) | lat_out, | ||
real(fms_hi_kind_), dimension(:), intent(in) | lon_in, | ||
real(fms_hi_kind_), dimension(:), intent(in) | lat_in, | ||
character(len=*), intent(in) | weight_file_source, | ||
character(len=*), intent(in) | interp_method, | ||
integer, intent(in) | isw, | ||
integer, intent(in) | iew, | ||
integer, intent(in) | jsw, | ||
integer, intent(in) | jew, | ||
integer, intent(in) | nglon, | ||
integer, intent(in) | nglat | ||
) |
Subroutine for reading a weight file and use it to fill in the horiz interp type for the bilinear interpolation method.
[in,out] | interp | Horiz interp time to fill |
[in] | weight_filename | Name of the weight file |
[in] | lat_out | Output (model) latitude |
[in] | lon_out | Output (model) longitude |
[in] | lat_in | Input (data) latitude |
[in] | lon_in | Input (data) longitude |
[in] | weight_file_source | Source of the weight file |
[in] | interp_method | The interp method to use |
[in] | jew | Starting and ending indices of the compute domain |
[in] | nglon | Number of longitudes in the global domain |
[in] | nglat | Number of latitudes in the globl domain |
Definition at line 846 of file horiz_interp.inc.
subroutine horiz_interp_solo_old_ | ( | real(fms_hi_kind_), dimension(:,:), intent(in) | data_in, |
real(fms_hi_kind_), intent(in) | wb, | ||
real(fms_hi_kind_), intent(in) | sb, | ||
real(fms_hi_kind_), intent(in) | dx, | ||
real(fms_hi_kind_), intent(in) | dy, | ||
real(fms_hi_kind_), dimension(:), intent(in) | lon_out, | ||
real(fms_hi_kind_), dimension(:), intent(in) | lat_out, | ||
real(fms_hi_kind_), dimension(:,:), intent(out) | data_out, | ||
integer, intent(in), optional | verbose, | ||
real(fms_hi_kind_), dimension(:,:), intent(in), optional | mask_in, | ||
real(fms_hi_kind_), dimension(:,:), intent(out), optional | mask_out | ||
) |
Overloaded version of interface horiz_interp_solo_2.
[in] | data_in | Global input data stored from west to east (1st dimension), south to north (2nd dimension) |
[in] | wb | Longitude (radians) that correspond to western-most boundary of grid box j=1 in array data_in |
[in] | sb | Latitude (radians) that correspond to western-most boundary of grid box j=1 in array data_in |
[in] | dx | Grid spacing (in radians) for the longitude axis (first dimension) for the input data |
[in] | dy | Grid spacing (in radians) for the latitude axis (first dimension) for the input data |
[in] | lon_out | The longitude edges (in radians) for output data grid boxes. The values are for adjacent grid boxes and must increase in value. If there are MLON grid boxes there must be MLON+1 edge values |
[in] | lat_out | The latitude edges (in radians) for output data grid boxes. The values are for adjacent grid boxes and may increase or decrease in value. If there are NLAT grid boxes there must be NLAT+1 edge values |
[out] | data_out | Output data on the output grid defined by grid box |
Definition at line 739 of file horiz_interp.inc.