FMS  2024.03
Flexible Modeling System
horiz_interp_spherical.F90
1 !***********************************************************************
2 !* GNU Lesser General Public License
3 !*
4 !* This file is part of the GFDL Flexible Modeling System (FMS).
5 !*
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.
10 !*
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
14 !* for more details.
15 !*
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.
27 
28 !> @addtogroup horiz_interp_spherical_mod
29 !> @{
30 module horiz_interp_spherical_mod
31 
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, spherical
40 
41  implicit none
42  private
43 
45  module procedure horiz_interp_spherical_r4
46  module procedure horiz_interp_spherical_r8
47  end interface
48 
50  module procedure horiz_interp_spherical_new_r4
51  module procedure horiz_interp_spherical_new_r8
52  end interface
53 
55  module procedure horiz_interp_spherical_wght_r4
56  module procedure horiz_interp_spherical_wght_r8
57  end interface
58 
61 
62  !> private helper routines
63  interface full_search
64  module procedure full_search_r4
65  module procedure full_search_r8
66  end interface
67 
68  interface radial_search
69  module procedure radial_search_r4
70  module procedure radial_search_r8
71  end interface
72 
74  module procedure spherical_distance_r4
75  module procedure spherical_distance_r8
76  end interface
77 
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
83 
84  integer :: pe, root_pe
85 
86 
87  character(len=32) :: search_method = "radial_search" !< Namelist variable to indicate the searching
88  !! method to find the
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.
98 
99 !or "full_search"
100  namelist /horiz_interp_spherical_nml/ search_method
101 
102  !-----------------------------------------------------------------------
103  ! Include variable "version" to be written to log file.
104 #include<file_version.h>
105  logical :: module_is_initialized = .false.
106 
107 contains
108 
109  !#######################################################################
110 
111  !> Initializes module and writes version number to logfile.out
113  integer :: ierr, io
114 
115 
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')
120 
121  module_is_initialized = .true.
122 
123  end subroutine horiz_interp_spherical_init
124 
125  !#######################################################################
126 
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 )
130 
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.
135 
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)
140  endif
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)
144 
145  interp%horizInterpReals4_type%is_allocated = .false.
146  interp%horizInterpReals8_type%is_allocated = .false.
147 
148  end subroutine horiz_interp_spherical_del
149 
150  !#######################################################################
151 
152 #include "horiz_interp_spherical_r4.fh"
153 #include "horiz_interp_spherical_r8.fh"
154 
155 end module horiz_interp_spherical_mod
156 !> @}
157 ! close documentation grouping
integer function, public check_nml_error(IOSTAT, NML_NAME)
Checks the iostat argument that is returned after reading a namelist and determines if the error code...
Definition: fms.F90:580
subroutine, public write_version_number(version, tag, unit)
Prints to the log file (or a specified unit) the version id string and tag name.
Definition: fms.F90:758
subroutine, public horiz_interp_spherical_init
Initializes module and writes version number to logfile.out.
subroutine, public horiz_interp_spherical_del(Interp)
Deallocates memory used by "HI_KIND_TYPE" variables. Must be called before reinitializing with horiz_...
character(len=32) search_method
Namelist variable to indicate the searching method to find the nearest neighbor points....
Holds data pointers and metadata for horizontal interpolations, passed between the horiz_interp modul...
integer function stdout()
This function returns the current standard fortran unit numbers for output.
Definition: mpp_util.inc:43
integer function mpp_pe()
Returns processor ID.
Definition: mpp_util.inc:407
Error handler.
Definition: mpp.F90:382