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_type_mod horiz_interp_type_mod
20 !> @ingroup horiz_interp
21 !> @brief define derived data type that contains indices and weights used for subsequent
27 module horiz_interp_type_mod
29 use mpp_mod, only : mpp_send, mpp_recv, mpp_sync_self, mpp_error, FATAL
30 use mpp_mod, only : mpp_pe, mpp_root_pe, mpp_npes
31 use mpp_mod, only : COMM_TAG_1, COMM_TAG_2
32 use platform_mod, only: r4_kind, r8_kind
38 ! parameter to determine interpolation method
39 integer, parameter :: CONSERVE = 1
40 integer, parameter :: BILINEAR = 2
41 integer, parameter :: SPHERICA = 3
42 integer, parameter :: BICUBIC = 4
44 public :: CONSERVE, BILINEAR, SPHERICA, BICUBIC
45 public :: horiz_interp_type, stats, assignment(=)
49 !> @ingroup horiz_interp_type_mod
50 interface assignment(=)
51 module procedure horiz_interp_type_eq
54 !> @ingroup horiz_interp_type_mod
56 module procedure stats_r4
57 module procedure stats_r8
61 !> real(8) pointers for use in horiz_interp_type
62 type horizInterpReals8_type
63 real(kind=r8_kind), dimension(:,:), allocatable :: faci !< weights for conservative scheme
64 real(kind=r8_kind), dimension(:,:), allocatable :: facj !< weights for conservative scheme
65 real(kind=r8_kind), dimension(:,:), allocatable :: area_src !< area of the source grid
66 real(kind=r8_kind), dimension(:,:), allocatable :: area_dst !< area of the destination grid
67 real(kind=r8_kind), dimension(:,:,:), allocatable :: wti !< weights for bilinear interpolation
68 !! wti ist used for derivative "weights" in bicubic
69 real(kind=r8_kind), dimension(:,:,:), allocatable :: wtj !< weights for bilinear interpolation
70 !! wti ist used for derivative "weights" in bicubic
71 real(kind=r8_kind), dimension(:,:,:), allocatable :: src_dist !< distance between destination grid and
72 !! neighbor source grid.
73 real(kind=r8_kind), dimension(:,:), allocatable :: rat_x !< the ratio of coordinates of the dest grid
74 !! (x_dest -x_src_r)/(x_src_l -x_src_r)
75 !! and (y_dest -y_src_r)/(y_src_l -y_src_r)
76 real(kind=r8_kind), dimension(:,:), allocatable :: rat_y !< the ratio of coordinates of the dest grid
77 !! (x_dest -x_src_r)/(x_src_l -x_src_r)
78 !! and (y_dest -y_src_r)/(y_src_l -y_src_r)
79 real(kind=r8_kind), dimension(:), allocatable :: lon_in !< the coordinates of the source grid
80 real(kind=r8_kind), dimension(:), allocatable :: lat_in !< the coordinates of the source grid
81 real(kind=r8_kind), dimension(:), allocatable :: area_frac_dst !< area fraction in destination grid.
82 real(kind=r8_kind), dimension(:,:), allocatable :: mask_in
83 real(kind=r8_kind) :: max_src_dist
84 logical :: is_allocated = .false. !< set to true upon field allocation
86 end type horizInterpReals8_type
88 !> holds real(4) pointers for use in horiz_interp_type
89 type horizInterpReals4_type
90 real(kind=r4_kind), dimension(:,:), allocatable :: faci !< weights for conservative scheme
91 real(kind=r4_kind), dimension(:,:), allocatable :: facj !< weights for conservative scheme
92 real(kind=r4_kind), dimension(:,:), allocatable :: area_src !< area of the source grid
93 real(kind=r4_kind), dimension(:,:), allocatable :: area_dst !< area of the destination grid
94 real(kind=r4_kind), dimension(:,:,:), allocatable :: wti !< weights for bilinear interpolation
95 !! wti ist used for derivative "weights" in bicubic
96 real(kind=r4_kind), dimension(:,:,:), allocatable :: wtj !< weights for bilinear interpolation
97 !! wti ist used for derivative "weights" in bicubic
98 real(kind=r4_kind), dimension(:,:,:), allocatable :: src_dist !< distance between destination grid and
99 !! neighbor source grid.
100 real(kind=r4_kind), dimension(:,:), allocatable :: rat_x !< the ratio of coordinates of the dest grid
101 !! (x_dest -x_src_r)/(x_src_l -x_src_r)
102 !! and (y_dest -y_src_r)/(y_src_l -y_src_r)
103 real(kind=r4_kind), dimension(:,:), allocatable :: rat_y !< the ratio of coordinates of the dest grid
104 !! (x_dest -x_src_r)/(x_src_l -x_src_r)
105 !! and (y_dest -y_src_r)/(y_src_l -y_src_r)
106 real(kind=r4_kind), dimension(:), allocatable :: lon_in !< the coordinates of the source grid
107 real(kind=r4_kind), dimension(:), allocatable :: lat_in !< the coordinates of the source grid
108 real(kind=r4_kind), dimension(:), allocatable :: area_frac_dst !< area fraction in destination grid.
109 real(kind=r4_kind), dimension(:,:), allocatable :: mask_in
110 real(kind=r4_kind) :: max_src_dist
111 logical :: is_allocated = .false. !< set to true upon field allocation
113 end type horizInterpReals4_type
115 !> Holds data pointers and metadata for horizontal interpolations, passed between the horiz_interp modules
116 !> @ingroup horiz_interp_type_mod
117 type horiz_interp_type
118 integer, dimension(:,:), allocatable :: ilon !< indices for conservative scheme
119 integer, dimension(:,:), allocatable :: jlat !< indices for conservative scheme
120 !! wti ist used for derivative "weights" in bicubic
121 integer, dimension(:,:,:), allocatable :: i_lon !< indices for bilinear interpolation
122 !! and spherical regrid
123 integer, dimension(:,:,:), allocatable :: j_lat !< indices for bilinear interpolation
124 !! and spherical regrid
125 logical, dimension(:,:), allocatable :: found_neighbors !< indicate whether destination grid
126 !! has some source grid around it.
127 integer, dimension(:,:), allocatable :: num_found
128 integer :: nlon_src !< size of source grid
129 integer :: nlat_src !< size of source grid
130 integer :: nlon_dst !< size of destination grid
131 integer :: nlat_dst !< size of destination grid
132 integer :: interp_method !< interpolation method.
133 !! =1, conservative scheme
134 !! =2, bilinear interpolation
135 !! =3, spherical regrid
136 !! =4, bicubic regrid
137 logical :: I_am_initialized=.false.
138 integer :: version !< indicate conservative
139 !! interpolation version with value 1 or 2
140 !--- The following are for conservative interpolation scheme version 2 ( through xgrid)
141 integer :: nxgrid !< number of exchange grid
142 !! between src and dst grid.
143 integer, dimension(:), allocatable :: i_src !< indices in source grid.
144 integer, dimension(:), allocatable :: j_src !< indices in source grid.
145 integer, dimension(:), allocatable :: i_dst !< indices in destination grid.
146 integer, dimension(:), allocatable :: j_dst !< indices in destination grid.
147 type(horizInterpReals8_type) :: horizInterpReals8_type !< derived type holding kind 8 real data pointers
148 !! if compiled with r8_kind
149 type(horizInterpReals4_type) :: horizInterpReals4_type !< derived type holding kind 4 real data pointers
150 !! if compiled with r8_kind
153 !> @addtogroup horiz_interp_type_mod
157 !######################################################################################################################
158 !> @brief horiz_interp_type_eq creates a copy of the horiz_interp_type object
159 subroutine horiz_interp_type_eq(horiz_interp_out, horiz_interp_in)
160 type(horiz_interp_type), intent(inout) :: horiz_interp_out !< Output object being set
161 type(horiz_interp_type), intent(in) :: horiz_interp_in !< Input object being copied
163 if(.not.horiz_interp_in%I_am_initialized) then
164 call mpp_error(FATAL,'horiz_interp_type_eq: horiz_interp_type variable on right hand side is unassigned')
167 horiz_interp_out%ilon = horiz_interp_in%ilon
168 horiz_interp_out%jlat = horiz_interp_in%jlat
169 horiz_interp_out%i_lon = horiz_interp_in%i_lon
170 horiz_interp_out%j_lat = horiz_interp_in%j_lat
171 horiz_interp_out%found_neighbors = horiz_interp_in%found_neighbors
172 horiz_interp_out%num_found = horiz_interp_in%num_found
173 horiz_interp_out%nlon_src = horiz_interp_in%nlon_src
174 horiz_interp_out%nlat_src = horiz_interp_in%nlat_src
175 horiz_interp_out%nlon_dst = horiz_interp_in%nlon_dst
176 horiz_interp_out%nlat_dst = horiz_interp_in%nlat_dst
177 horiz_interp_out%interp_method = horiz_interp_in%interp_method
178 horiz_interp_out%I_am_initialized = .true.
179 horiz_interp_out%i_src = horiz_interp_in%i_src
180 horiz_interp_out%j_src = horiz_interp_in%j_src
181 horiz_interp_out%i_dst = horiz_interp_in%i_dst
182 horiz_interp_out%j_dst = horiz_interp_in%j_dst
184 if(horiz_interp_in%horizInterpReals8_type%is_allocated) then
185 horiz_interp_out%horizInterpReals8_type%faci = horiz_interp_in%horizInterpReals8_type%faci
186 horiz_interp_out%horizInterpReals8_type%facj = horiz_interp_in%horizInterpReals8_type%facj
187 horiz_interp_out%horizInterpReals8_type%area_src = horiz_interp_in%horizInterpReals8_type%area_src
188 horiz_interp_out%horizInterpReals8_type%area_dst = horiz_interp_in%horizInterpReals8_type%area_dst
189 horiz_interp_out%horizInterpReals8_type%wti = horiz_interp_in%horizInterpReals8_type%wti
190 horiz_interp_out%horizInterpReals8_type%wtj = horiz_interp_in%horizInterpReals8_type%wtj
191 horiz_interp_out%horizInterpReals8_type%src_dist = horiz_interp_in%horizInterpReals8_type%src_dist
192 horiz_interp_out%horizInterpReals8_type%rat_x = horiz_interp_in%horizInterpReals8_type%rat_x
193 horiz_interp_out%horizInterpReals8_type%rat_y = horiz_interp_in%horizInterpReals8_type%rat_y
194 horiz_interp_out%horizInterpReals8_type%lon_in = horiz_interp_in%horizInterpReals8_type%lon_in
195 horiz_interp_out%horizInterpReals8_type%lat_in = horiz_interp_in%horizInterpReals8_type%lat_in
196 horiz_interp_out%horizInterpReals8_type%area_frac_dst = horiz_interp_in%horizInterpReals8_type%area_frac_dst
197 horiz_interp_out%horizInterpReals8_type%max_src_dist = horiz_interp_in%horizInterpReals8_type%max_src_dist
198 horiz_interp_out%horizInterpReals8_type%is_allocated = .true.
199 ! this was left out previous to mixed mode
200 horiz_interp_out%horizInterpReals8_type%mask_in = horiz_interp_in%horizInterpReals8_type%mask_in
202 else if (horiz_interp_in%horizInterpReals4_type%is_allocated) then
203 horiz_interp_out%horizInterpReals4_type%faci = horiz_interp_in%horizInterpReals4_type%faci
204 horiz_interp_out%horizInterpReals4_type%facj = horiz_interp_in%horizInterpReals4_type%facj
205 horiz_interp_out%horizInterpReals4_type%area_src = horiz_interp_in%horizInterpReals4_type%area_src
206 horiz_interp_out%horizInterpReals4_type%area_dst = horiz_interp_in%horizInterpReals4_type%area_dst
207 horiz_interp_out%horizInterpReals4_type%wti = horiz_interp_in%horizInterpReals4_type%wti
208 horiz_interp_out%horizInterpReals4_type%wtj = horiz_interp_in%horizInterpReals4_type%wtj
209 horiz_interp_out%horizInterpReals4_type%src_dist = horiz_interp_in%horizInterpReals4_type%src_dist
210 horiz_interp_out%horizInterpReals4_type%rat_x = horiz_interp_in%horizInterpReals4_type%rat_x
211 horiz_interp_out%horizInterpReals4_type%rat_y = horiz_interp_in%horizInterpReals4_type%rat_y
212 horiz_interp_out%horizInterpReals4_type%lon_in = horiz_interp_in%horizInterpReals4_type%lon_in
213 horiz_interp_out%horizInterpReals4_type%lat_in = horiz_interp_in%horizInterpReals4_type%lat_in
214 horiz_interp_out%horizInterpReals4_type%area_frac_dst = horiz_interp_in%horizInterpReals4_type%area_frac_dst
215 horiz_interp_out%horizInterpReals4_type%max_src_dist = horiz_interp_in%horizInterpReals4_type%max_src_dist
216 horiz_interp_out%horizInterpReals4_type%is_allocated = .true.
217 ! this was left out previous to mixed mode
218 horiz_interp_out%horizInterpReals4_type%mask_in = horiz_interp_in%horizInterpReals4_type%mask_in
221 call mpp_error(FATAL, "horiz_interp_type_eq: cannot assign unallocated real values from horiz_interp_in")
224 if(horiz_interp_in%interp_method == CONSERVE) then
225 horiz_interp_out%version = horiz_interp_in%version
226 if(horiz_interp_in%version==2) horiz_interp_out%nxgrid = horiz_interp_in%nxgrid
229 end subroutine horiz_interp_type_eq
230 !######################################################################################################################
232 #include "horiz_interp_type_r4.fh"
233 #include "horiz_interp_type_r8.fh"
235 end module horiz_interp_type_mod
237 ! close documentation grouping