FMS  2025.04
Flexible Modeling System
horiz_interp_type.F90
1 !***********************************************************************
2 !* Apache License 2.0
3 !*
4 !* This file is part of the GFDL Flexible Modeling System (FMS).
5 !*
6 !* Licensed under the Apache License, Version 2.0 (the "License");
7 !* you may not use this file except in compliance with the License.
8 !* You may obtain a copy of the License at
9 !*
10 !* http://www.apache.org/licenses/LICENSE-2.0
11 !*
12 !* FMS is distributed in the hope that it will be useful, but WITHOUT
13 !* WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied;
14 !* without even the implied warranty of MERCHANTABILITY or FITNESS FOR A
15 !* PARTICULAR PURPOSE. See the License for the specific language
16 !* governing permissions and limitations under the License.
17 !***********************************************************************
18 !> @defgroup horiz_interp_type_mod horiz_interp_type_mod
19 !> @ingroup horiz_interp
20 !> @brief define derived data type that contains indices and weights used for subsequent
21 !! interpolations.
22 !> @author Zhi Liang
23 
24 !> @addtogroup
25 !> @{
26 module horiz_interp_type_mod
27 
28 use mpp_mod, only : mpp_send, mpp_recv, mpp_sync_self, mpp_error, fatal
29 use mpp_mod, only : mpp_pe, mpp_root_pe, mpp_npes
30 use mpp_mod, only : comm_tag_1, comm_tag_2
31 use platform_mod, only: r4_kind, r8_kind
32 
33 implicit none
34 private
35 
36 
37 ! parameter to determine interpolation method
38  integer, parameter :: CONSERVE = 1
39  integer, parameter :: BILINEAR = 2
40  integer, parameter :: SPHERICAL = 3
41  integer, parameter :: BICUBIC = 4
42 
43 public :: conserve, bilinear, spherical, bicubic
44 public :: horiz_interp_type, stats, assignment(=)
45 
46 !> @}
47 
48 !> @ingroup horiz_interp_type_mod
49 interface assignment(=)
50  module procedure horiz_interp_type_eq
51 end interface
52 
53 !> @ingroup horiz_interp_type_mod
54 interface stats
55  module procedure stats_r4
56  module procedure stats_r8
57 end interface
58 
59 
60 !> real(8) pointers for use in horiz_interp_type
62  real(kind=r8_kind), dimension(:,:), allocatable :: faci !< weights for conservative scheme
63  real(kind=r8_kind), dimension(:,:), allocatable :: facj !< weights for conservative scheme
64  real(kind=r8_kind), dimension(:,:), allocatable :: area_src !< area of the source grid
65  real(kind=r8_kind), dimension(:,:), allocatable :: area_dst !< area of the destination grid
66  real(kind=r8_kind), dimension(:,:,:), allocatable :: wti !< weights for bilinear interpolation
67  !! wti ist used for derivative "weights" in bicubic
68  real(kind=r8_kind), dimension(:,:,:), allocatable :: wtj !< weights for bilinear interpolation
69  !! wti ist used for derivative "weights" in bicubic
70  real(kind=r8_kind), dimension(:,:,:), allocatable :: src_dist !< distance between destination grid and
71  !! neighbor source grid.
72  real(kind=r8_kind), dimension(:,:), allocatable :: rat_x !< the ratio of coordinates of the dest grid
73  !! (x_dest -x_src_r)/(x_src_l -x_src_r)
74  !! and (y_dest -y_src_r)/(y_src_l -y_src_r)
75  real(kind=r8_kind), dimension(:,:), allocatable :: rat_y !< the ratio of coordinates of the dest grid
76  !! (x_dest -x_src_r)/(x_src_l -x_src_r)
77  !! and (y_dest -y_src_r)/(y_src_l -y_src_r)
78  real(kind=r8_kind), dimension(:), allocatable :: lon_in !< the coordinates of the source grid
79  real(kind=r8_kind), dimension(:), allocatable :: lat_in !< the coordinates of the source grid
80  real(kind=r8_kind), dimension(:), allocatable :: area_frac_dst !< area fraction in destination grid.
81  real(kind=r8_kind), dimension(:,:), allocatable :: mask_in
82  real(kind=r8_kind) :: max_src_dist
83  logical :: is_allocated = .false. !< set to true upon field allocation
84 
86 
87 !> holds real(4) pointers for use in horiz_interp_type
89  real(kind=r4_kind), dimension(:,:), allocatable :: faci !< weights for conservative scheme
90  real(kind=r4_kind), dimension(:,:), allocatable :: facj !< weights for conservative scheme
91  real(kind=r4_kind), dimension(:,:), allocatable :: area_src !< area of the source grid
92  real(kind=r4_kind), dimension(:,:), allocatable :: area_dst !< area of the destination grid
93  real(kind=r4_kind), dimension(:,:,:), allocatable :: wti !< weights for bilinear interpolation
94  !! wti ist used for derivative "weights" in bicubic
95  real(kind=r4_kind), dimension(:,:,:), allocatable :: wtj !< weights for bilinear interpolation
96  !! wti ist used for derivative "weights" in bicubic
97  real(kind=r4_kind), dimension(:,:,:), allocatable :: src_dist !< distance between destination grid and
98  !! neighbor source grid.
99  real(kind=r4_kind), dimension(:,:), allocatable :: rat_x !< the ratio of coordinates of the dest grid
100  !! (x_dest -x_src_r)/(x_src_l -x_src_r)
101  !! and (y_dest -y_src_r)/(y_src_l -y_src_r)
102  real(kind=r4_kind), dimension(:,:), allocatable :: rat_y !< the ratio of coordinates of the dest grid
103  !! (x_dest -x_src_r)/(x_src_l -x_src_r)
104  !! and (y_dest -y_src_r)/(y_src_l -y_src_r)
105  real(kind=r4_kind), dimension(:), allocatable :: lon_in !< the coordinates of the source grid
106  real(kind=r4_kind), dimension(:), allocatable :: lat_in !< the coordinates of the source grid
107  real(kind=r4_kind), dimension(:), allocatable :: area_frac_dst !< area fraction in destination grid.
108  real(kind=r4_kind), dimension(:,:), allocatable :: mask_in
109  real(kind=r4_kind) :: max_src_dist
110  logical :: is_allocated = .false. !< set to true upon field allocation
111 
112 end type horizinterpreals4_type
113 
114 !> Holds data pointers and metadata for horizontal interpolations, passed between the horiz_interp modules
115 !> @ingroup horiz_interp_type_mod
117  integer, dimension(:,:), allocatable :: ilon !< indices for conservative scheme
118  integer, dimension(:,:), allocatable :: jlat !< indices for conservative scheme
119  !! wti ist used for derivative "weights" in bicubic
120  integer, dimension(:,:,:), allocatable :: i_lon !< indices for bilinear interpolation
121  !! and spherical regrid
122  integer, dimension(:,:,:), allocatable :: j_lat !< indices for bilinear interpolation
123  !! and spherical regrid
124  logical, dimension(:,:), allocatable :: found_neighbors !< indicate whether destination grid
125  !! has some source grid around it.
126  integer, dimension(:,:), allocatable :: num_found
127  integer :: nlon_src !< size of source grid
128  integer :: nlat_src !< size of source grid
129  integer :: nlon_dst !< size of destination grid
130  integer :: nlat_dst !< size of destination grid
131  integer :: interp_method !< interpolation method.
132  !! =1, conservative scheme
133  !! =2, bilinear interpolation
134  !! =3, spherical regrid
135  !! =4, bicubic regrid
136  logical :: i_am_initialized=.false.
137  integer :: version !< indicate conservative
138  !! interpolation version with value 1 or 2
139  !--- The following are for conservative interpolation scheme version 2 ( through xgrid)
140  integer :: nxgrid !< number of exchange grid
141  !! between src and dst grid.
142  integer, dimension(:), allocatable :: i_src !< indices in source grid.
143  integer, dimension(:), allocatable :: j_src !< indices in source grid.
144  integer, dimension(:), allocatable :: i_dst !< indices in destination grid.
145  integer, dimension(:), allocatable :: j_dst !< indices in destination grid.
146  type(horizinterpreals8_type) :: horizinterpreals8_type !< derived type holding kind 8 real data pointers
147  !! if compiled with r8_kind
148  type(horizinterpreals4_type) :: horizinterpreals4_type !< derived type holding kind 4 real data pointers
149  !! if compiled with r8_kind
150  end type
151 
152 !> @addtogroup horiz_interp_type_mod
153 !> @{
154 contains
155 
156 !######################################################################################################################
157 !> @brief horiz_interp_type_eq creates a copy of the horiz_interp_type object
158  subroutine horiz_interp_type_eq(horiz_interp_out, horiz_interp_in)
159  type(horiz_interp_type), intent(inout) :: horiz_interp_out !< Output object being set
160  type(horiz_interp_type), intent(in) :: horiz_interp_in !< Input object being copied
161 
162  if(.not.horiz_interp_in%I_am_initialized) then
163  call mpp_error(fatal,'horiz_interp_type_eq: horiz_interp_type variable on right hand side is unassigned')
164  endif
165 
166  if( allocated(horiz_interp_in%ilon )) &
167  horiz_interp_out%ilon = horiz_interp_in%ilon
168 
169  if( allocated(horiz_interp_in%jlat )) &
170  horiz_interp_out%jlat = horiz_interp_in%jlat
171 
172  if( allocated(horiz_interp_in%i_lon )) &
173  horiz_interp_out%i_lon = horiz_interp_in%i_lon
174 
175  if( allocated(horiz_interp_in%j_lat )) &
176  horiz_interp_out%j_lat = horiz_interp_in%j_lat
177 
178  if( allocated(horiz_interp_in%found_neighbors )) &
179  horiz_interp_out%found_neighbors = horiz_interp_in%found_neighbors
180 
181  if( allocated(horiz_interp_in%num_found )) &
182  horiz_interp_out%num_found = horiz_interp_in%num_found
183 
184  if( allocated(horiz_interp_in%i_src )) &
185  horiz_interp_out%i_src = horiz_interp_in%i_src
186 
187  if( allocated(horiz_interp_in%j_src )) &
188  horiz_interp_out%j_src = horiz_interp_in%j_src
189 
190  if( allocated(horiz_interp_in%i_dst )) &
191  horiz_interp_out%i_dst = horiz_interp_in%i_dst
192 
193  if( allocated(horiz_interp_in%j_dst )) &
194  horiz_interp_out%j_dst = horiz_interp_in%j_dst
195 
196  horiz_interp_out%nlon_src = horiz_interp_in%nlon_src
197  horiz_interp_out%nlat_src = horiz_interp_in%nlat_src
198  horiz_interp_out%nlon_dst = horiz_interp_in%nlon_dst
199  horiz_interp_out%nlat_dst = horiz_interp_in%nlat_dst
200  horiz_interp_out%interp_method = horiz_interp_in%interp_method
201  horiz_interp_out%I_am_initialized = .true.
202 
203  if(horiz_interp_in%horizInterpReals8_type%is_allocated) then
204 
205  if( allocated(horiz_interp_in%horizInterpReals8_type%faci)) &
206  horiz_interp_out%horizInterpReals8_type%faci = horiz_interp_in%horizInterpReals8_type%faci
207 
208  if( allocated( horiz_interp_in%horizInterpReals8_type%facj)) &
209  horiz_interp_out%horizInterpReals8_type%facj = horiz_interp_in%horizInterpReals8_type%facj
210 
211  if( allocated( horiz_interp_in%horizInterpReals8_type%area_src)) &
212  horiz_interp_out%horizInterpReals8_type%area_src = horiz_interp_in%horizInterpReals8_type%area_src
213 
214  if( allocated( horiz_interp_in%horizInterpReals8_type%area_dst)) &
215  horiz_interp_out%horizInterpReals8_type%area_dst = horiz_interp_in%horizInterpReals8_type%area_dst
216 
217  if( allocated( horiz_interp_in%horizInterpReals8_type%wti)) &
218  horiz_interp_out%horizInterpReals8_type%wti = horiz_interp_in%horizInterpReals8_type%wti
219 
220  if( allocated( horiz_interp_in%horizInterpReals8_type%wtj)) &
221  horiz_interp_out%horizInterpReals8_type%wtj = horiz_interp_in%horizInterpReals8_type%wtj
222 
223  if( allocated( horiz_interp_in%horizInterpReals8_type%src_dist)) &
224  horiz_interp_out%horizInterpReals8_type%src_dist = horiz_interp_in%horizInterpReals8_type%src_dist
225 
226  if( allocated( horiz_interp_in%horizInterpReals8_type%rat_x)) &
227  horiz_interp_out%horizInterpReals8_type%rat_x = horiz_interp_in%horizInterpReals8_type%rat_x
228 
229  if( allocated( horiz_interp_in%horizInterpReals8_type%rat_y)) &
230  horiz_interp_out%horizInterpReals8_type%rat_y = horiz_interp_in%horizInterpReals8_type%rat_y
231 
232  if( allocated( horiz_interp_in%horizInterpReals8_type%lon_in)) &
233  horiz_interp_out%horizInterpReals8_type%lon_in = horiz_interp_in%horizInterpReals8_type%lon_in
234 
235  if( allocated( horiz_interp_in%horizInterpReals8_type%lat_in)) &
236  horiz_interp_out%horizInterpReals8_type%lat_in = horiz_interp_in%horizInterpReals8_type%lat_in
237 
238  if( allocated( horiz_interp_in%horizInterpReals8_type%area_frac_dst)) &
239  horiz_interp_out%horizInterpReals8_type%area_frac_dst = horiz_interp_in%horizInterpReals8_type%area_frac_dst
240 
241  horiz_interp_out%horizInterpReals8_type%max_src_dist = horiz_interp_in%horizInterpReals8_type%max_src_dist
242 
243  horiz_interp_out%horizInterpReals8_type%is_allocated = .true.
244  ! this was left out previous to mixed mode
245  if( allocated(horiz_interp_in%horizInterpReals8_type%mask_in)) &
246  horiz_interp_out%horizInterpReals8_type%mask_in = horiz_interp_in%horizInterpReals8_type%mask_in
247 
248  else if (horiz_interp_in%horizInterpReals4_type%is_allocated) then
249  if( allocated(horiz_interp_in%horizInterpReals4_type%faci)) &
250  horiz_interp_out%horizInterpReals4_type%faci = horiz_interp_in%horizInterpReals4_type%faci
251 
252  if( allocated( horiz_interp_in%horizInterpReals4_type%facj)) &
253  horiz_interp_out%horizInterpReals4_type%facj = horiz_interp_in%horizInterpReals4_type%facj
254 
255  if( allocated( horiz_interp_in%horizInterpReals4_type%area_src)) &
256  horiz_interp_out%horizInterpReals4_type%area_src = horiz_interp_in%horizInterpReals4_type%area_src
257 
258  if( allocated( horiz_interp_in%horizInterpReals4_type%area_dst)) &
259  horiz_interp_out%horizInterpReals4_type%area_dst = horiz_interp_in%horizInterpReals4_type%area_dst
260 
261  if( allocated( horiz_interp_in%horizInterpReals4_type%wti)) &
262  horiz_interp_out%horizInterpReals4_type%wti = horiz_interp_in%horizInterpReals4_type%wti
263 
264  if( allocated( horiz_interp_in%horizInterpReals4_type%wtj)) &
265  horiz_interp_out%horizInterpReals4_type%wtj = horiz_interp_in%horizInterpReals4_type%wtj
266 
267  if( allocated( horiz_interp_in%horizInterpReals4_type%src_dist)) &
268  horiz_interp_out%horizInterpReals4_type%src_dist = horiz_interp_in%horizInterpReals4_type%src_dist
269 
270  if( allocated( horiz_interp_in%horizInterpReals4_type%rat_x)) &
271  horiz_interp_out%horizInterpReals4_type%rat_x = horiz_interp_in%horizInterpReals4_type%rat_x
272 
273  if( allocated( horiz_interp_in%horizInterpReals4_type%rat_y)) &
274  horiz_interp_out%horizInterpReals4_type%rat_y = horiz_interp_in%horizInterpReals4_type%rat_y
275 
276  if( allocated( horiz_interp_in%horizInterpReals4_type%lon_in)) &
277  horiz_interp_out%horizInterpReals4_type%lon_in = horiz_interp_in%horizInterpReals4_type%lon_in
278 
279  if( allocated( horiz_interp_in%horizInterpReals4_type%lat_in)) &
280  horiz_interp_out%horizInterpReals4_type%lat_in = horiz_interp_in%horizInterpReals4_type%lat_in
281 
282  if( allocated( horiz_interp_in%horizInterpReals4_type%area_frac_dst)) &
283  horiz_interp_out%horizInterpReals4_type%area_frac_dst = horiz_interp_in%horizInterpReals4_type%area_frac_dst
284 
285  horiz_interp_out%horizInterpReals4_type%max_src_dist = horiz_interp_in%horizInterpReals4_type%max_src_dist
286 
287  horiz_interp_out%horizInterpReals4_type%is_allocated = .true.
288  ! this was left out previous to mixed mode
289  if( allocated(horiz_interp_in%horizInterpReals4_type%mask_in)) &
290  horiz_interp_out%horizInterpReals4_type%mask_in = horiz_interp_in%horizInterpReals4_type%mask_in
291 
292  else
293  call mpp_error(fatal, "horiz_interp_type_eq: cannot assign unallocated real values from horiz_interp_in")
294  endif
295 
296  if(horiz_interp_in%interp_method == conserve) then
297  horiz_interp_out%version = horiz_interp_in%version
298  if(horiz_interp_in%version==2) horiz_interp_out%nxgrid = horiz_interp_in%nxgrid
299  end if
300 
301  end subroutine horiz_interp_type_eq
302 !######################################################################################################################
303 
304 #include "horiz_interp_type_r4.fh"
305 #include "horiz_interp_type_r8.fh"
306 
307 end module horiz_interp_type_mod
308 !> @}
309 ! close documentation grouping
subroutine horiz_interp_type_eq(horiz_interp_out, horiz_interp_in)
horiz_interp_type_eq creates a copy of the horiz_interp_type object
Holds data pointers and metadata for horizontal interpolations, passed between the horiz_interp modul...
subroutine mpp_sync_self(pelist, check, request, msg_size, msg_type)
This is to check if current PE's outstanding puts are complete but we can't use shmem_fence because w...
integer function mpp_npes()
Returns processor count for current pelist.
Definition: mpp_util.inc:420
integer function mpp_pe()
Returns processor ID.
Definition: mpp_util.inc:406
Error handler.
Definition: mpp.F90:381
Receive data from another PE.
Definition: mpp.F90:950
Send data to a receiving PE.
Definition: mpp.F90:1017
holds real(4) pointers for use in horiz_interp_type
real(8) pointers for use in horiz_interp_type