chore: remove outdated am4 ci workflow (#1639)
[FMS.git] / horiz_interp / horiz_interp_bilinear.F90
blobd8db732b222bba0772d1e0a4b46619939195c938
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.
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.
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_bilinear_mod horiz_interp_bilinear_mod
20 !> @ingroup horiz_interp
21 !> @brief Performs spatial interpolation between grids using bilinear interpolation
23 !> @author Zhi Liang <Zhi.Liang@noaa.gov>
24 !> This module can interpolate data from regular rectangular grid
25 !! to rectangular/tripolar grid. The interpolation scheme is bilinear interpolation.
26 !! There is an optional mask field for missing input data.
27 !! An optional output mask field may be used in conjunction with
28 !! the input mask to show where output data exists.
30 module horiz_interp_bilinear_mod
32   use mpp_mod,               only: mpp_error, FATAL, stdout, mpp_pe, mpp_root_pe
33   use fms_mod,               only: write_version_number
34   use constants_mod,         only: PI
35   use horiz_interp_type_mod, only: horiz_interp_type, stats, BILINEAR
36   use platform_mod,          only: r4_kind, r8_kind
37   use axis_utils2_mod,       only: nearest_index
38   use fms2_io_mod,           only: open_file, close_file, read_data, FmsNetcdfFile_t, get_dimension_size
39   use fms_string_utils_mod,  only: string
41   implicit none
42   private
45   public :: horiz_interp_bilinear_new, horiz_interp_bilinear, horiz_interp_bilinear_del
46   public :: horiz_interp_bilinear_init, horiz_interp_read_weights_bilinear
48   !> Creates a @ref horiz_interp_type for bilinear interpolation.
49   !> @ingroup horiz_interp_bilinear_mod
50   interface horiz_interp_bilinear_new
51     module procedure horiz_interp_bilinear_new_1d_r4
52     module procedure horiz_interp_bilinear_new_1d_r8
53     module procedure horiz_interp_bilinear_new_2d_r4
54     module procedure horiz_interp_bilinear_new_2d_r8
55   end interface
57   !> Subroutines for reading in weight files and using that to fill in the horiz_interp type instead
58   !! calculating it
59   !> @ingroup horiz_interp_bilinear_mod
60   interface horiz_interp_read_weights_bilinear
61     module procedure horiz_interp_read_weights_bilinear_r4
62     module procedure horiz_interp_read_weights_bilinear_r8
63   end interface
65   interface horiz_interp_bilinear
66     module procedure horiz_interp_bilinear_r4
67     module procedure horiz_interp_bilinear_r8
68   end interface
70 !> @addtogroup horiz_interp_bilinear_mod
71 !> @{
73   real(r8_kind), parameter :: epsln=1.e-10_r8_kind
74   integer, parameter :: DUMMY = -999
76 !! Private helper routines, interfaces for mixed real precision support
77   interface intersect
78     module procedure intersect_r4
79     module procedure intersect_r8
80   end interface
82   !-----------------------------------------------------------------------
83 ! Include variable "version" to be written to log file.
84 #include<file_version.h>
85   logical            :: module_is_initialized = .FALSE.
87 contains
89   !> Initialize this module and writes version number to logfile.
90   subroutine horiz_interp_bilinear_init
92     if(module_is_initialized) return
93     call write_version_number("HORIZ_INTERP_BILINEAR_MOD", version)
94     module_is_initialized = .true.
96   end subroutine horiz_interp_bilinear_init
98   !> @brief Deallocates memory used by "horiz_interp_type" variables.
99   !!
100   !> Must be called before reinitializing with horiz_interp_bilinear_new.
101   subroutine horiz_interp_bilinear_del( Interp )
103     type (horiz_interp_type), intent(inout) :: Interp!< A derived-type variable returned by previous
104                                    !! call to horiz_interp_bilinear_new. The input variable must
105                                    !! have allocated arrays. The returned variable will contain
106                                    !! deallocated arrays
108     if( Interp%horizInterpReals4_type%is_allocated) then
109       if(allocated(Interp%horizInterpReals4_type%wti))   deallocate(Interp%horizInterpReals4_type%wti)
110       if(allocated(Interp%horizInterpReals4_type%wtj))   deallocate(Interp%horizInterpReals4_type%wtj)
111     else if (Interp%horizInterpReals8_type%is_allocated) then
112       if(allocated(Interp%horizInterpReals8_type%wti))   deallocate(Interp%horizInterpReals8_type%wti)
113       if(allocated(Interp%horizInterpReals8_type%wtj))   deallocate(Interp%horizInterpReals8_type%wtj)
114     endif
115     if(allocated(Interp%i_lon)) deallocate(Interp%i_lon)
116     if(allocated(Interp%j_lat)) deallocate(Interp%j_lat)
118     Interp%horizInterpReals4_type%is_allocated = .false.
119     Interp%horizInterpReals8_type%is_allocated = .false.
121   end subroutine horiz_interp_bilinear_del
123 #include "horiz_interp_bilinear_r4.fh"
124 #include "horiz_interp_bilinear_r8.fh"
126 end module horiz_interp_bilinear_mod
127 !> @}
128 ! close documentation grouping