1 !***********************************************************************
2 !* GNU Lesser General Public License
4 !* This file is part of the GFDL Flexible Modeling System (FMS).
6 !* FMS is free software: you can redistribute it and/or modify it under
7 !* the terms of the GNU Lesser General Public License as published by
8 !* the Free Software Foundation, either version 3 of the License, or (at
9 !* your option) any later version.
11 !* FMS is distributed in the hope that it will be useful, but WITHOUT
12 !* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
13 !* FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
16 !* You should have received a copy of the GNU Lesser General Public
17 !* License along with FMS. If not, see <http://www.gnu.org/licenses/>.
18 !***********************************************************************
19 !> @defgroup horiz_interp_spherical_mod horiz_interp_spherical_mod
20 !> @ingroup horiz_interp
21 !> @brief Performs spatial interpolation between grids using inverse-distance-weighted scheme.
22 !> This module can interpolate data from rectangular/tripolar grid
23 !! to rectangular/tripolar grid. The interpolation scheme is inverse-distance-weighted
24 !! scheme. There is an optional mask field for missing input data.
25 !! An optional output mask field may be used in conjunction with
26 !! the input mask to show where output data exists.
28 !> @addtogroup horiz_interp_spherical_mod
30 module horiz_interp_spherical_mod
32 use platform_mod, only : r4_kind, r8_kind
33 use mpp_mod, only : mpp_error, FATAL, WARNING, stdout
34 use mpp_mod, only : mpp_root_pe, mpp_pe
35 use mpp_mod, only : input_nml_file
36 use fms_mod, only : write_version_number
37 use fms_mod, only : check_nml_error
38 use constants_mod, only : pi
39 use horiz_interp_type_mod, only : horiz_interp_type, stats
44 interface horiz_interp_spherical
45 module procedure horiz_interp_spherical_r4
46 module procedure horiz_interp_spherical_r8
49 interface horiz_interp_spherical_new
50 module procedure horiz_interp_spherical_new_r4
51 module procedure horiz_interp_spherical_new_r8
54 interface horiz_interp_spherical_wght
55 module procedure horiz_interp_spherical_wght_r4
56 module procedure horiz_interp_spherical_wght_r8
59 public :: horiz_interp_spherical_new, horiz_interp_spherical, horiz_interp_spherical_del
60 public :: horiz_interp_spherical_init, horiz_interp_spherical_wght
62 !> private helper routines
64 module procedure full_search_r4
65 module procedure full_search_r8
68 interface radial_search
69 module procedure radial_search_r4
70 module procedure radial_search_r8
73 interface spherical_distance
74 module procedure spherical_distance_r4
75 module procedure spherical_distance_r8
78 integer, parameter :: max_neighbors = 400
79 real(R8_KIND), parameter :: max_dist_default = 0.1_r8_kind ! radians
80 integer, parameter :: num_nbrs_default = 4
81 real(R8_KIND), parameter :: large=1.e20_r8_kind
82 real(R8_KIND), parameter :: epsln=1.e-10_r8_kind
84 integer :: pe, root_pe
87 character(len=32) :: search_method = "radial_search" !< Namelist variable to indicate the searching
89 !! nearest neighbor points. Its value can be "radial_search" and "full_search",
90 !! with default value "radial_search". when search_method is "radial_search",
91 !! the search may be not quite accurate for some cases. Normally the search will
92 !! be ok if you chose suitable max_dist. When search_method is "full_search",
93 !! it will be always accurate, but will be slower comparing to "radial_search".
94 !! Normally these two search algorithm will produce same results other than
95 !! order of operation. "radial_search" are recommended to use. The purpose to
96 !! add "full_search" is in case you think you interpolation results is
97 !! not right, you have other option to verify.
100 namelist /horiz_interp_spherical_nml/ search_method
102 !-----------------------------------------------------------------------
103 ! Include variable "version" to be written to log file.
104 #include<file_version.h>
105 logical :: module_is_initialized = .FALSE.
109 !#######################################################################
111 !> Initializes module and writes version number to logfile.out
112 subroutine horiz_interp_spherical_init
116 if(module_is_initialized) return
117 call write_version_number("horiz_interp_spherical_mod", version)
118 read (input_nml_file, horiz_interp_spherical_nml, iostat=io)
119 ierr = check_nml_error(io,'horiz_interp_spherical_nml')
121 module_is_initialized = .true.
123 end subroutine horiz_interp_spherical_init
125 !#######################################################################
127 !> Deallocates memory used by "HI_KIND_TYPE" variables.
128 !! Must be called before reinitializing with horiz_interp_spherical_new.
129 subroutine horiz_interp_spherical_del( Interp )
131 type (horiz_interp_type), intent(inout) :: Interp !< A derived-type variable returned by previous
132 !! call to horiz_interp_spherical_new. The input variable
133 !! must have allocated arrays. The returned variable will
134 !! contain deallocated arrays.
136 if(Interp%horizInterpReals4_type%is_allocated) then
137 if(allocated(Interp%horizInterpReals4_type%src_dist)) deallocate(Interp%horizInterpReals4_type%src_dist)
138 else if (Interp%horizInterpReals8_type%is_allocated) then
139 if(allocated(Interp%horizInterpReals8_type%src_dist)) deallocate(Interp%horizInterpReals8_type%src_dist)
141 if(allocated(Interp%num_found)) deallocate(Interp%num_found)
142 if(allocated(Interp%i_lon)) deallocate(Interp%i_lon)
143 if(allocated(Interp%j_lat)) deallocate(Interp%j_lat)
145 Interp%horizInterpReals4_type%is_allocated = .false.
146 Interp%horizInterpReals8_type%is_allocated = .false.
148 end subroutine horiz_interp_spherical_del
150 !#######################################################################
152 #include "horiz_interp_spherical_r4.fh"
153 #include "horiz_interp_spherical_r8.fh"
155 end module horiz_interp_spherical_mod
157 ! close documentation grouping