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_conserve_mod horiz_interp_conserve_mod
20 !> @ingroup horiz_interp
21 !> @brief Performs spatial interpolation between grids using conservative interpolation
23 !> @author Bruce Wyman, Zhi Liang
25 !> This module can conservatively interpolate data from any logically rectangular grid
26 !! to any rectangular grid. The interpolation scheme is area-averaging
27 !! conservative scheme. There is an optional mask field for missing input data in both
28 !! horiz_interp__conserveinit and horiz_interp_conserve. For efficiency purpose, mask should only be
29 !! kept in horiz_interp_init (will remove the mask in horiz_interp in the future).
30 !! There are 1-D and 2-D version of horiz_interp_conserve_init for 1-D and 2-D grid.
31 !! There is a optional argument mask in horiz_interp_conserve_init_2d and no mask should
32 !! to passed into horiz_interp_conserv. optional argument mask will not be passed into
33 !! horiz_interp_conserve_init_1d and optional argument mask may be passed into
34 !! horiz_interp_conserve (For the purpose of reproduce Memphis??? results).
35 !! An optional output mask field may be used in conjunction with the input mask to show
36 !! where output data exists.
38 module horiz_interp_conserve_mod
40 use platform_mod, only: r4_kind, r8_kind
41 use mpp_mod, only: mpp_send, mpp_recv, mpp_pe, mpp_root_pe, mpp_npes
42 use mpp_mod, only: mpp_error, FATAL, mpp_sync_self
43 use mpp_mod, only: COMM_TAG_1, COMM_TAG_2
44 use fms_mod, only: write_version_number
45 use grid2_mod, only: get_great_circle_algorithm
46 use constants_mod, only: PI
47 use horiz_interp_type_mod, only: horiz_interp_type, CONSERVE
56 !> @brief Allocates space and initializes a derived-type variable
57 !! that contains pre-computed interpolation indices and weights.
59 !> Allocates space and initializes a derived-type variable
60 !! that contains pre-computed interpolation indices and weights
61 !! for improved performance of multiple interpolations between
64 !! Longitude (in radians) for source data grid.
67 !! Latitude (in radians) for source data grid.
70 !! Longitude (in radians) for destination data grid.
73 !! Latitude (in radians) for destination data grid.
76 !! flag for the amount of print output.
79 !! Input mask. must be the size (size(lon_in)-1, size(lon. The real value of
80 !! mask_in must be in the range (0.,1.). Set mask_in=0.0 for data points
81 !! that should not be used or have missing data.
84 !! Output mask that specifies whether data was computed.
87 !! A derived-type variable containing indices and weights used for subsequent
88 !! interpolations. To reinitialize this variable for a different grid-to-grid
89 !! interpolation you must first use the "horiz_interp_del" interface.
91 !> @ingroup horiz_interp_conserve_mod
92 interface horiz_interp_conserve_new
93 module procedure horiz_interp_conserve_new_1dx1d_r4
94 module procedure horiz_interp_conserve_new_1dx2d_r4
95 module procedure horiz_interp_conserve_new_2dx1d_r4
96 module procedure horiz_interp_conserve_new_2dx2d_r4
97 module procedure horiz_interp_conserve_new_1dx1d_r8
98 module procedure horiz_interp_conserve_new_1dx2d_r8
99 module procedure horiz_interp_conserve_new_2dx1d_r8
100 module procedure horiz_interp_conserve_new_2dx2d_r8
103 interface horiz_interp_conserve
104 module procedure horiz_interp_conserve_r4
105 module procedure horiz_interp_conserve_r8
108 !> private helper routines
110 module procedure data_sum_r4
111 module procedure data_sum_r8
115 module procedure stats_r4
116 module procedure stats_r8
119 interface horiz_interp_conserve_version1
120 module procedure horiz_interp_conserve_version1_r8
121 module procedure horiz_interp_conserve_version1_r4
124 interface horiz_interp_conserve_version2
125 module procedure horiz_interp_conserve_version2_r8
126 module procedure horiz_interp_conserve_version2_r4
131 !> @addtogroup horiz_interp_conserve_mod
133 public :: horiz_interp_conserve_init
134 public :: horiz_interp_conserve_new, horiz_interp_conserve, horiz_interp_conserve_del
136 integer :: pe, root_pe
137 !-----------------------------------------------------------------------
138 ! Include variable "version" to be written to log file.
139 #include<file_version.h>
140 logical :: module_is_initialized = .FALSE.
142 logical :: great_circle_algorithm = .false.
146 !> Writes version number to logfile.
147 subroutine horiz_interp_conserve_init
149 if(module_is_initialized) return
150 call write_version_number("HORIZ_INTERP_CONSERVE_MOD", version)
152 great_circle_algorithm = get_great_circle_algorithm()
154 module_is_initialized = .true.
156 end subroutine horiz_interp_conserve_init
158 !> Deallocates memory used by "HI_KIND_TYPE" variables.
159 !! Must be called before reinitializing with horiz_interp_new.
160 subroutine horiz_interp_conserve_del ( Interp )
162 type (horiz_interp_type), intent(inout) :: Interp !< A derived-type variable returned by
163 !! previous call to horiz_interp_new. The input variable must have
164 !! allocated arrays. The returned variable will contain deallocated arrays.
166 select case(Interp%version)
168 if( Interp%horizInterpReals8_type%is_allocated) then
169 if(allocated(Interp%horizInterpReals8_type%area_src)) deallocate(Interp%horizInterpReals8_type%area_src)
170 if(allocated(Interp%horizInterpReals8_type%area_dst)) deallocate(Interp%horizInterpReals8_type%area_dst)
171 if(allocated(Interp%horizInterpReals8_type%facj)) deallocate(Interp%horizInterpReals8_type%facj)
172 if(allocated(Interp%jlat)) deallocate(Interp%jlat)
173 if(allocated(Interp%horizInterpReals8_type%faci)) deallocate(Interp%horizInterpReals8_type%faci)
174 if(allocated(Interp%ilon)) deallocate(Interp%ilon)
175 else if( Interp%horizInterpReals4_type%is_allocated) then
176 if(allocated(Interp%horizInterpReals4_type%area_src)) deallocate(Interp%horizInterpReals4_type%area_src)
177 if(allocated(Interp%horizInterpReals4_type%area_dst)) deallocate(Interp%horizInterpReals4_type%area_dst)
178 if(allocated(Interp%horizInterpReals4_type%facj)) deallocate(Interp%horizInterpReals4_type%facj)
179 if(allocated(Interp%jlat)) deallocate(Interp%jlat)
180 if(allocated(Interp%horizInterpReals4_type%faci)) deallocate(Interp%horizInterpReals4_type%faci)
181 if(allocated(Interp%ilon)) deallocate(Interp%ilon)
184 if( Interp%horizInterpReals8_type%is_allocated) then
185 if(allocated(Interp%i_src)) deallocate(Interp%i_src)
186 if(allocated(Interp%j_src)) deallocate(Interp%j_src)
187 if(allocated(Interp%i_dst)) deallocate(Interp%i_dst)
188 if(allocated(Interp%j_dst)) deallocate(Interp%j_dst)
189 if(allocated(Interp%horizInterpReals8_type%area_frac_dst)) &
190 deallocate(Interp%horizInterpReals8_type%area_frac_dst)
191 else if( Interp%horizInterpReals4_type%is_allocated ) then
192 if(allocated(Interp%i_src)) deallocate(Interp%i_src)
193 if(allocated(Interp%j_src)) deallocate(Interp%j_src)
194 if(allocated(Interp%i_dst)) deallocate(Interp%i_dst)
195 if(allocated(Interp%j_dst)) deallocate(Interp%j_dst)
196 if(allocated(Interp%horizInterpReals4_type%area_frac_dst)) &
197 deallocate(Interp%horizInterpReals4_type%area_frac_dst)
200 Interp%horizInterpReals4_type%is_allocated = .false.
201 Interp%horizInterpReals8_type%is_allocated = .false.
203 end subroutine horiz_interp_conserve_del
205 #include "horiz_interp_conserve_r4.fh"
206 #include "horiz_interp_conserve_r8.fh"
208 end module horiz_interp_conserve_mod
210 ! close documentation grouping