27 module horiz_interp_type_mod
31 use mpp_mod,
only : comm_tag_1, comm_tag_2
32 use platform_mod,
only: r4_kind, r8_kind
39 integer,
parameter :: CONSERVE = 1
40 integer,
parameter :: BILINEAR = 2
41 integer,
parameter :: SPHERICAL = 3
42 integer,
parameter :: BICUBIC = 4
44 public :: conserve, bilinear, spherical, bicubic
50 interface assignment(=)
56 module procedure stats_r4
57 module procedure stats_r8
63 real(kind=r8_kind),
dimension(:,:),
allocatable :: faci
64 real(kind=r8_kind),
dimension(:,:),
allocatable :: facj
65 real(kind=r8_kind),
dimension(:,:),
allocatable :: area_src
66 real(kind=r8_kind),
dimension(:,:),
allocatable :: area_dst
67 real(kind=r8_kind),
dimension(:,:,:),
allocatable :: wti
69 real(kind=r8_kind),
dimension(:,:,:),
allocatable :: wtj
71 real(kind=r8_kind),
dimension(:,:,:),
allocatable :: src_dist
73 real(kind=r8_kind),
dimension(:,:),
allocatable :: rat_x
76 real(kind=r8_kind),
dimension(:,:),
allocatable :: rat_y
79 real(kind=r8_kind),
dimension(:),
allocatable :: lon_in
80 real(kind=r8_kind),
dimension(:),
allocatable :: lat_in
81 real(kind=r8_kind),
dimension(:),
allocatable :: area_frac_dst
82 real(kind=r8_kind),
dimension(:,:),
allocatable :: mask_in
83 real(kind=r8_kind) :: max_src_dist
84 logical :: is_allocated = .false.
90 real(kind=r4_kind),
dimension(:,:),
allocatable :: faci
91 real(kind=r4_kind),
dimension(:,:),
allocatable :: facj
92 real(kind=r4_kind),
dimension(:,:),
allocatable :: area_src
93 real(kind=r4_kind),
dimension(:,:),
allocatable :: area_dst
94 real(kind=r4_kind),
dimension(:,:,:),
allocatable :: wti
96 real(kind=r4_kind),
dimension(:,:,:),
allocatable :: wtj
98 real(kind=r4_kind),
dimension(:,:,:),
allocatable :: src_dist
100 real(kind=r4_kind),
dimension(:,:),
allocatable :: rat_x
103 real(kind=r4_kind),
dimension(:,:),
allocatable :: rat_y
106 real(kind=r4_kind),
dimension(:),
allocatable :: lon_in
107 real(kind=r4_kind),
dimension(:),
allocatable :: lat_in
108 real(kind=r4_kind),
dimension(:),
allocatable :: area_frac_dst
109 real(kind=r4_kind),
dimension(:,:),
allocatable :: mask_in
110 real(kind=r4_kind) :: max_src_dist
111 logical :: is_allocated = .false.
118 integer,
dimension(:,:),
allocatable :: ilon
119 integer,
dimension(:,:),
allocatable :: jlat
121 integer,
dimension(:,:,:),
allocatable :: i_lon
123 integer,
dimension(:,:,:),
allocatable :: j_lat
125 logical,
dimension(:,:),
allocatable :: found_neighbors
127 integer,
dimension(:,:),
allocatable :: num_found
132 integer :: interp_method
137 logical :: i_am_initialized=.false.
143 integer,
dimension(:),
allocatable :: i_src
144 integer,
dimension(:),
allocatable :: j_src
145 integer,
dimension(:),
allocatable :: i_dst
146 integer,
dimension(:),
allocatable :: j_dst
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.
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.
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
232 #include "horiz_interp_type_r4.fh"
233 #include "horiz_interp_type_r8.fh"
235 end module horiz_interp_type_mod
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.
integer function mpp_pe()
Returns processor ID.
Recieve data from another PE.
Send data to a receiving PE.
holds real(4) pointers for use in horiz_interp_type
real(8) pointers for use in horiz_interp_type