FMS
2024.03
Flexible Modeling System
|
Performs spatial interpolation between grids using inverse-distance-weighted scheme. This module can interpolate data from rectangular/tripolar grid to rectangular/tripolar grid. The interpolation scheme is inverse-distance-weighted scheme. 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. More...
Data Types | |
interface | full_search |
private helper routines More... | |
interface | horiz_interp_spherical_new |
interface | horiz_interp_spherical_wght |
interface | radial_search |
interface | spherical_distance |
Functions/Subroutines | |
subroutine | full_search_ (theta_src, phi_src, theta_dst, phi_dst, map_src_add, map_src_dist, num_found, num_neighbors, max_src_dist) |
full_search_r4 | |
full_search_r8 | |
subroutine | horiz_interp_spherical_ (Interp, data_in, data_out, verbose, mask_in, mask_out, missing_value) |
Subroutine for performing the horizontal interpolation between two grids. HORIZ_INTERP_SPHERICAL_NEW_ must be called before calling this routine. More... | |
subroutine, public | horiz_interp_spherical_del (Interp) |
Deallocates memory used by "HI_KIND_TYPE" variables. Must be called before reinitializing with horiz_interp_spherical_new. More... | |
real(fms_hi_kind_) function | horiz_interp_spherical_distance_ (theta1, phi1, theta2, phi2) |
subroutine, public | horiz_interp_spherical_init |
Initializes module and writes version number to logfile.out. | |
subroutine | horiz_interp_spherical_new_ (Interp, lon_in, lat_in, lon_out, lat_out, num_nbrs, max_dist, src_modulo) |
horiz_interp_spherical_new_r4 | |
horiz_interp_spherical_new_r8 | |
horiz_interp_spherical_r4 | |
horiz_interp_spherical_r8 | |
subroutine | horiz_interp_spherical_wght_ (Interp, wt, verbose, mask_in, mask_out, missing_value) |
This routine isn't used internally it's similar to the routine above, just gets the weights as an out variable. | |
horiz_interp_spherical_wght_r4 | |
horiz_interp_spherical_wght_r8 | |
subroutine | radial_search_ (theta_src, phi_src, theta_dst, phi_dst, map_src_xsize, map_src_ysize, map_src_add, map_src_dist, num_found, num_neighbors, max_src_dist, src_is_modulo) |
radial_search_r4 | |
radial_search_r8 | |
spherical_distance_r4 | |
spherical_distance_r8 | |
logical function | update_dest_neighbors_ (map_src_add, map_src_dist, src_add, d, num_found, min_nbrs) |
Variables | |
real(r8_kind), parameter | epsln =1.e-10_r8_kind |
real(r8_kind), parameter | large =1.e20_r8_kind |
real(r8_kind), parameter | max_dist_default = 0.1_r8_kind |
integer, parameter | max_neighbors = 400 |
logical | module_is_initialized = .FALSE. |
integer, parameter | num_nbrs_default = 4 |
integer | pe |
integer | root_pe |
character(len=32) | search_method = "radial_search" |
Namelist variable to indicate the searching method to find the nearest neighbor points. Its value can be "radial_search" and "full_search", with default value "radial_search". when search_method is "radial_search", the search may be not quite accurate for some cases. Normally the search will be ok if you chose suitable max_dist. When search_method is "full_search", it will be always accurate, but will be slower comparing to "radial_search". Normally these two search algorithm will produce same results other than order of operation. "radial_search" are recommended to use. The purpose to add "full_search" is in case you think you interpolation results is not right, you have other option to verify. | |
Performs spatial interpolation between grids using inverse-distance-weighted scheme. This module can interpolate data from rectangular/tripolar grid to rectangular/tripolar grid. The interpolation scheme is inverse-distance-weighted scheme. 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.
Initialization routine.
Allocates space and initializes a derived-type variable that contains pre-computed interpolation indices and weights.
interface horiz_interp_spherical_mod::full_search |
private helper routines
Definition at line 63 of file horiz_interp_spherical.F90.
Private Member Functions | |
full_search_r4 | |
full_search_r8 | |
interface horiz_interp_spherical_mod::horiz_interp_spherical_new |
Definition at line 49 of file horiz_interp_spherical.F90.
Private Member Functions | |
horiz_interp_spherical_new_r4 | |
horiz_interp_spherical_new_r8 | |
interface horiz_interp_spherical_mod::horiz_interp_spherical_wght |
Definition at line 54 of file horiz_interp_spherical.F90.
Private Member Functions | |
horiz_interp_spherical_wght_r4 | |
horiz_interp_spherical_wght_r8 | |
interface horiz_interp_spherical_mod::radial_search |
Definition at line 68 of file horiz_interp_spherical.F90.
Private Member Functions | |
radial_search_r4 | |
radial_search_r8 | |
interface horiz_interp_spherical_mod::spherical_distance |
Definition at line 73 of file horiz_interp_spherical.F90.
Private Member Functions | |
spherical_distance_r4 | |
spherical_distance_r8 | |
subroutine horiz_interp_spherical_ | ( | type(horiz_interp_type), intent(in) | Interp, |
real(fms_hi_kind_), dimension(:,:), intent(in) | data_in, | ||
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, | ||
real(fms_hi_kind_), intent(in), optional | missing_value | ||
) |
Subroutine for performing the horizontal interpolation between two grids. HORIZ_INTERP_SPHERICAL_NEW_ must be called before calling this routine.
[in] | interp | A derived type variable containing indices and weights for subsequent interpolations. Returned by a previous call to HORIZ_INTERP_SPHERICAL_NEW_ |
[in] | data_in | Input data on source grid |
[out] | data_out | Output data on destination grid |
[in] | verbose | 0 = no output; 1 = min,max,means; 2 = most output |
[in] | mask_in | Input mask, must be the same size as the input data. The real(FMS_HI_KIND_) 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 |
[out] | mask_out | Output mask that specifies whether data was computed. |
[in] | missing_value | Used to indicate missing data |
Definition at line 202 of file horiz_interp_spherical.inc.
subroutine, public horiz_interp_spherical_mod::horiz_interp_spherical_del | ( | type (horiz_interp_type), intent(inout) | Interp | ) |
Deallocates memory used by "HI_KIND_TYPE" variables. Must be called before reinitializing with horiz_interp_spherical_new.
[in,out] | interp | A derived-type variable returned by previous call to horiz_interp_spherical_new. The input variable must have allocated arrays. The returned variable will contain deallocated arrays. |
Definition at line 129 of file horiz_interp_spherical.F90.
subroutine horiz_interp_spherical_new_ | ( | 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 | num_nbrs, | ||
real(fms_hi_kind_), intent(in), optional | max_dist, | ||
logical, intent(in), optional | src_modulo | ||
) |
[in,out] | interp | A derived type variable containing indices and weights for subsequent interpolations. To reinitialize for different grid-to-grid interpolation horiz_interp_del must be used first. |
[in] | lon_in | Latitude (radians) for source data grid |
[in] | lat_in | Longitude (radians) for source data grid |
[in] | lon_out | Longitude (radians) for output data grid |
[in] | lat_out | Latitude (radians) for output data grid |
[in] | src_modulo | indicates if the boundary condition along zonal boundary is cyclic or not. Cyclic when true |
[in] | num_nbrs | Number of nearest neighbors for regridding When number of neighbors within the radius max_dist is less than num_nbrs, All the neighbors will be used to interpolate onto destination grid. when number of neighbors within the radius max_dist is greater than num_nbrs, at least "num_nbrs" |
[in] | max_dist | Maximum region of influence around destination grid points |
Definition at line 25 of file horiz_interp_spherical.inc.