OOPS
unstructured_interpolation_mod.F90
Go to the documentation of this file.
1 ! (C) Copyright 2019-2020 UCAR
2 !
3 ! This software is licensed under the terms of the Apache Licence Version 2.0
4 ! which can be obtained at http://www.apache.org/licenses/LICENSE-2.0.
5 
7 
8 use atlas_module, only: atlas_geometry, atlas_indexkdtree
9 
10 use netcdf
11 
12 use kinds
13 use string_utils, only: replace_string
14 use netcdf_utils_mod, only: nccheck
15 
16 use fckit_mpi_module, only: fckit_mpi_comm, fckit_mpi_sum
17 
18 implicit none
19 private
20 public unstrc_interp
22 
23 !---------------------------------------------------------------------------------------------------
24 
26  integer :: nn ! Number of neighbours
27  integer :: ngrid_in ! Number of input grid points
28  integer :: ngrid_out ! Number of output grid points
29  integer :: ngrid_in_glo ! Number of global input grid points
30  type(fckit_mpi_comm) :: comm ! Communicator
31  integer, allocatable :: displs(:) ! Displacement for global gather
32  integer, allocatable :: rcvcnt(:) ! Receive count for global gather
33  real(kind=kind_real), allocatable :: interp_w(:,:) ! Interpolation weights
34  integer , allocatable :: interp_i(:,:) ! Interpolation index
35  contains
36  generic, public :: create => create_new, create_read
37  procedure, private :: create_new, create_read
38  procedure, public :: delete
39  procedure, public :: apply
40  procedure, public :: apply_ad
41  procedure, public :: write
42 endtype unstrc_interp
43 
44 ! ------------------------------------------------------------------------------
45 !> Registry for unstrc_interpo objects
46 
47 #define LISTED_TYPE unstrc_interp
48 
49 !> Linked list interface - defines registry_t type
50 #include "oops/util/linkedList_i.f"
51 
52 !> Global registry
53 type(registry_t) :: unstrc_interp_registry
54 
55 !---------------------------------------------------------------------------------------------------
56 
57 contains
58 !-------------------------------------------------------------------------------
59 !> Linked list implementation
60 #include "oops/util/linkedList_c.f"
61 
62 !---------------------------------------------------------------------------------------------------
63 
64 subroutine create_new( self, comm, nn, wtype, ngrid_in, lats_in, lons_in, &
65  ngrid_out, lats_out, lons_out )
66 class(unstrc_interp), intent(inout) :: self ! Myself
67 type(fckit_mpi_comm), intent(in) :: comm ! Communicator
68 integer, intent(in) :: nn ! Number of neighbours
69 character(len=*), intent(in) :: wtype ! Weight types to return
70 integer, intent(in) :: ngrid_in ! Number of grid points on input grid
71 real(kind=kind_real), intent(in) :: lats_in(ngrid_in) ! Input grid latitude in degrees
72 real(kind=kind_real), intent(in) :: lons_in(ngrid_in) ! Input grid longitude in degrees
73 integer, intent(in) :: ngrid_out ! Number of grid points on output grid
74 real(kind=kind_real), intent(in) :: lats_out(ngrid_out) ! Output grid latitide in degrees
75 real(kind=kind_real), intent(in) :: lons_out(ngrid_out) ! Output grid longitude in degrees
76 
77 !Locals
78 type(atlas_indexkdtree) :: kd
79 type(atlas_geometry) :: ageometry
80 integer :: j, n, jj, kk, ngrid_in_glo
81 integer :: index1, index2, nindex
82 real(kind=kind_real) :: dist, wprod, bsw
83 real(kind=kind_real) :: lons_in_loc(ngrid_in), lons_out_loc(ngrid_out)
84 real(kind=kind_real), allocatable :: lats_in_glo(:), lons_in_glo(:)
85 real(kind=kind_real), allocatable :: nn_dist(:,:)
86 real(kind=kind_real), allocatable :: bw(:)
87 
88 ! wtype options
89 ! -------------
90 ! 1. (none) No weights needed
91 ! 2. (distance) Distance to nearest neighbours
92 ! 3. (barycent) Barycentric weights
93 
94 ! Allocate self type
95 self%nn = nn
96 self%ngrid_in = ngrid_in
97 self%ngrid_out = ngrid_out
98 self%comm = comm
99 allocate(self%interp_w(self%nn,self%ngrid_out))
100 allocate(self%interp_i(self%nn,self%ngrid_out))
101 
102 ! Get input grid proc counts and displacement
103 call input_grid_share(self%comm, self%ngrid_in, self%rcvcnt, self%displs)
104 
105 ! Global number of grid points
106 
107 ! ----------------------------
108 ngrid_in_glo = sum(self%rcvcnt)
109 self%ngrid_in_glo = ngrid_in_glo
110 
111 ! Diplacement for each processor
112 ! ------------------------------
113 self%displs(1) = 0
114 do j = 2,comm%size()
115  self%displs(j) = self%displs(j-1) + self%rcvcnt(j-1)
116 enddo
117 
118 ! Rotate longitudes less than 0
119 lons_in_loc = lons_in
120 lons_out_loc = lons_out
121 where (lons_in_loc < 0.0) lons_in_loc = lons_in_loc + 360.0_kind_real
122 where (lons_out_loc < 0.0) lons_out_loc = lons_out_loc + 360.0_kind_real
123 
124 ! Gather the lat/lons to all
125 allocate(lats_in_glo(ngrid_in_glo))
126 allocate(lons_in_glo(ngrid_in_glo))
127 call comm%allgather(lats_in ,lats_in_glo,ngrid_in,self%rcvcnt,self%displs)
128 call comm%allgather(lons_in_loc,lons_in_glo,ngrid_in,self%rcvcnt,self%displs)
129 
130 ! Create UnitSphere geometry
131 ageometry = atlas_geometry("UnitSphere")
132 
133 ! Create KDTree
134 kd = atlas_indexkdtree(ageometry)
135 call kd%reserve(ngrid_in_glo)
136 call kd%build(ngrid_in_glo,lons_in_glo,lats_in_glo)
137 
138 ! Loop over observations
139 ! ----------------------
140 allocate(nn_dist(nn,ngrid_out))
141 
142 do n = 1,ngrid_out
143 
144  ! Get nearest neighbours
145  call kd%closestPoints(lons_out_loc(n), lats_out(n), self%nn, self%interp_i(:,n))
146 
147  ! Compute distances
148  do kk = 1, nn
149  nindex = self%interp_i(kk,n)
150  nn_dist(kk,n) = ageometry%distance(lons_out_loc(n), lats_out(n), lons_in_glo(nindex), lats_in_glo(nindex))
151  enddo
152 
153 enddo
154 
155 
156 ! Set weights based on user choice
157 ! --------------------------------
158 select case (wtype)
159  case ('none')
160  self%interp_w = 0.0_kind_real
161 
162  case ('distance')
163  self%interp_w = nn_dist
164 
165  case ('barycent')
166  allocate(bw(self%nn))
167  do n = 1,ngrid_out
168  !Barycentric weights formula
169  bw(:) = 0.0_kind_real
170  do jj = 1, nn
171  wprod = 1.0_kind_real
172  index1 = self%interp_i(jj,n)
173  do kk = 1, nn
174  if (jj.ne.kk) then
175  index2 = self%interp_i(kk,n)
176  dist = ageometry%distance(lons_in_glo(index1),lats_in_glo(index1),&
177  lons_in_glo(index2),lats_in_glo(index2))
178  wprod = wprod * max(dist, 1e-10)
179  endif
180  enddo
181  bw(jj) = 1.0_kind_real / wprod
182  enddo
183 
184  !Barycentric weights
185  self%interp_w(:,n) = 0.0_kind_real
186  if (minval(nn_dist(:,n)) < 1e-10) then
187 
188  ! special case if very close to one grid point
189  jj = minloc(nn_dist(:,n),dim=1)
190  self%interp_w(jj,n) = 1.0_kind_real
191 
192  else
193 
194  !otherwise continue with the normal algorithm
195  bsw = 0.0_kind_real
196  do jj = 1,nn
197  bsw = bsw + (bw(jj) / nn_dist(jj,n))
198  enddo
199 
200  do jj = 1,nn
201  self%interp_w(jj,n) = ( bw(jj) / nn_dist(jj,n) ) / bsw
202  enddo
203  end if
204 
205  enddo
206 
207  deallocate(bw)
208 
209  case default
210  call abor1_ftn("unstrc_interpolation.create: wrong choice of weight output")
211 
212 end select
213 
214 !Deallocate
215 call kd%final()
216 deallocate(nn_dist)
217 deallocate(lats_in_glo, lons_in_glo)
218 
219 end subroutine create_new
220 
221 !---------------------------------------------------------------------------------------------------
222 
223 subroutine create_read(self, comm, filename_in)
224 class(unstrc_interp), intent(inout) :: self
225 type(fckit_mpi_comm), intent(in) :: comm
226 character(len=*), intent(in) :: filename_in
227 
228 character(len=2048) :: filename
229 integer :: ncid, varid
230 
231 self%comm = comm
232 
233 ! One file per processor
234 filename = filename_in
235 call getfilename(self%comm%rank(),filename)
236 
237 ! Open file
238 call nccheck ( nf90_open( trim(filename), nf90_nowrite, ncid ), "nf90_open"//trim(filename) )
239 
240 ! Get the dimensions
241 ! ------------------
242 call nccheck ( nf90_inq_dimid(ncid, "ngrid_in", varid), "nf90_inq_dimid ngrid_in" )
243 call nccheck ( nf90_inquire_dimension(ncid, varid, len = self%ngrid_in), "nf90_inquire_dimension ngrid_in" )
244 
245 call nccheck ( nf90_inq_dimid(ncid, "ngrid_out", varid), "nf90_inq_dimid ngrid_out" )
246 call nccheck ( nf90_inquire_dimension(ncid, varid, len = self%ngrid_out), "nf90_inquire_dimension ngrid_out" )
247 
248 call nccheck ( nf90_inq_dimid(ncid, "nn", varid), "nf90_inq_dimid nn" )
249 call nccheck ( nf90_inquire_dimension(ncid, varid, len = self%nn), "nf90_inquire_dimension nn" )
250 
251 ! Allocate arrays
252 allocate(self%interp_w(self%nn,self%ngrid_out))
253 allocate(self%interp_i(self%nn,self%ngrid_out))
254 
255 ! Get the interpolation weights and indices
256 ! ----------------------------------------
257 call nccheck ( nf90_inq_varid(ncid, "interp_weights", varid), "nf90_inq_varid interp_weights" )
258 call nccheck ( nf90_get_var(ncid, varid, self%interp_w), "nf90_get_var interp_weights" )
259 
260 call nccheck ( nf90_inq_varid(ncid, "interp_indices", varid), "nf90_inq_varid interp_indices" )
261 call nccheck ( nf90_get_var(ncid, varid, self%interp_i), "nf90_get_var interp_indices" )
262 
263 ! Close the file
264 call nccheck ( nf90_close(ncid), "nf90_close" )
265 
266 end subroutine create_read
267 
268 !---------------------------------------------------------------------------------------------------
269 
270 subroutine delete( self )
271 class(unstrc_interp), intent(inout) :: self
272 
273 if (allocated(self%interp_w)) deallocate(self%interp_w)
274 if (allocated(self%interp_i)) deallocate(self%interp_i)
275 if (allocated(self%rcvcnt)) deallocate(self%rcvcnt)
276 if (allocated(self%displs)) deallocate(self%displs)
277 
278 end subroutine delete
279 
280 !---------------------------------------------------------------------------------------------------
281 
282 subroutine apply(self, field_in, field_out, field_nn_out )
283 class(unstrc_interp), intent(in) :: self ! Myself
284 real(kind=kind_real), intent(in) :: field_in(self%ngrid_in) ! Input field
285 real(kind=kind_real), intent(out) :: field_out(self%ngrid_out) ! Result of interpolation
286 real(kind=kind_real), optional, intent(out) :: field_nn_out(self%nn,self%ngrid_out) ! Neighbours
287 
288 !Locals
289 integer :: n, kk
290 real(kind=kind_real), allocatable :: field_nn(:,:), field_in_glo(:)
291 
292 
293 ! Gather the field
294 allocate(field_in_glo(sum(self%rcvcnt)))
295 call self%comm%allgather(field_in,field_in_glo,self%ngrid_in,self%rcvcnt,self%displs)
296 
297 ! Get output neighbours
298 allocate(field_nn(self%nn,self%ngrid_out))
299 do n = 1, self%ngrid_out
300  do kk = 1, self%nn
301  field_nn(kk,n) = field_in_glo(self%interp_i(kk,n))
302  enddo
303 enddo
304 
305 ! Deallocate global field
306 deallocate(field_in_glo)
307 
308 ! Apply weights
309 ! -------------
310 do n = 1, self%ngrid_out
311  field_out(n) = 0.0_kind_real
312  do kk = 1, self%nn
313  field_out(n) = field_out(n) + self%interp_w(kk,n) * field_nn(kk,n)
314  enddo
315 enddo
316 
317 ! Return neighbours if requested
318 if (present(field_nn_out)) field_nn_out = field_nn
319 
320 deallocate(field_nn)
321 
322 end subroutine apply
323 
324 !---------------------------------------------------------------------------------------------------
325 
326 subroutine apply_ad(self, field_in, field_out )
327 class(unstrc_interp), intent(in) :: self
328 real(kind=kind_real), intent(out) :: field_in(self%ngrid_in)
329 real(kind=kind_real), optional, intent(in) :: field_out(self%ngrid_out)
330 
331 integer :: n, kk, nstart, nfinish
332 real(kind=kind_real), allocatable :: field_in_glo(:), field_in_glo_all(:)
333 
334 
335 ! Apply backward interpolation
336 ! ----------------------------
337 allocate(field_in_glo(self%ngrid_in_glo))
338 field_in_glo = 0.0
339 do n = 1, self%ngrid_out
340  do kk = 1, self%nn
341  field_in_glo(self%interp_i(kk,n)) = field_in_glo(self%interp_i(kk,n)) + self%interp_w(kk,n)*field_out(n)
342  enddo
343 enddo
344 
345 ! Global sum
346 ! ----------
347 allocate(field_in_glo_all(self%ngrid_in_glo))
348 call self%comm%allreduce(field_in_glo, field_in_glo_all, fckit_mpi_sum())
349 
350 ! Return local field
351 ! ------------------
352 nstart = self%displs(self%comm%rank()+1) + 1
353 nfinish = self%displs(self%comm%rank()+1) + self%rcvcnt(self%comm%rank()+1)
354 field_in(1:self%ngrid_in) = field_in_glo_all(nstart:nfinish)
355 deallocate(field_in_glo)
356 deallocate(field_in_glo_all)
357 
358 end subroutine apply_ad
359 
360 !---------------------------------------------------------------------------------------------------
361 
362 subroutine write(self,filename_in)
363 class(unstrc_interp), intent(in) :: self
364 character(len=*), intent(in) :: filename_in
365 
366 integer :: ncid, ni_dimid, no_dimid, nn_dimid, vc(2)
367 character(len=2048) :: filename
368 
369 ! One file per processor
370 filename = filename_in
371 call getfilename(self%comm%rank(),filename)
372 
373 ! Create file
374 call nccheck( nf90_create( trim(filename), ior(nf90_netcdf4, nf90_clobber), ncid), "nf90_create" )
375 
376 ! Define dimensions
377 call nccheck( nf90_def_dim(ncid, "ngrid_in", self%ngrid_in, ni_dimid), "nf90_def_dim ng" )
378 call nccheck( nf90_def_dim(ncid, "ngrid_out", self%ngrid_out, no_dimid), "nf90_def_dim ng" )
379 call nccheck( nf90_def_dim(ncid, "nn", self%nn, nn_dimid), "nf90_def_dim nn" )
380 
381 ! Define variables
382 call nccheck( nf90_def_var(ncid, "interp_weights", nf90_double, (/ nn_dimid, no_dimid /), vc(1)), "nf90_def_var interp_weights" );
383 call nccheck( nf90_def_var(ncid, "interp_indices", nf90_int, (/ nn_dimid, no_dimid /), vc(2)), "nf90_def_var interp_indices" );
384 
385 ! End define
386 call nccheck( nf90_enddef(ncid), "nf90_enddef" )
387 
388 ! Write variables
389 call nccheck( nf90_put_var( ncid, vc(1), self%interp_w ), "nf90_put_var interp_weights" );
390 call nccheck( nf90_put_var( ncid, vc(2), self%interp_I ), "nf90_put_var interp_indices" );
391 
392 ! Close file
393 call nccheck ( nf90_close(ncid), "nf90_close" )
394 
395 end subroutine write
396 
397 !---------------------------------------------------------------------------------------------------
398 
399 subroutine input_grid_share(comm, ngrid_in, rcvcnt, displs)
400 type(fckit_mpi_comm), intent(in) :: comm
401 integer, intent(in) :: ngrid_in
402 integer, allocatable, intent(inout) :: rcvcnt(:)
403 integer, allocatable, intent(inout) :: displs(:)
404 
405 integer :: j
406 
407 ! Gather receive count from each processor
408 allocate(rcvcnt(comm%size()))
409 call comm%allgather(ngrid_in, rcvcnt)
410 
411 ! Displacement for each processor
412 allocate(displs(comm%size()))
413 displs(1) = 0
414 do j = 2,comm%size()
415  displs(j) = displs(j-1) + rcvcnt(j-1)
416 enddo
417 
418 end subroutine input_grid_share
419 
420 !---------------------------------------------------------------------------------------------------
421 
422 subroutine getfilename(procrank, filename)
423 integer, intent(in) :: procrank
424 character(len=*), intent(inout) :: filename
425 
426 character(len=6) :: procrankstr
427 
428 ! Error check
429 if (procrank>999999) call abor1_ftn("unstructured_interpolation_mod not ready for more than 999999 processors")
430 
431 ! String processor number
432 write(procrankstr,'(I0.6)') procrank
433 
434 ! Get new filename
435 filename = replace_string(filename,'%PROC',procrankstr)
436 
437 end subroutine getfilename
438 
439 ! --------------------------------------------------------------------------------------------------
440 
subroutine getfilename(procrank, filename)
subroutine input_grid_share(comm, ngrid_in, rcvcnt, displs)
subroutine apply(self, field_in, field_out, field_nn_out)
type(registry_t), public unstrc_interp_registry
Registry for unstrc_interpo objects.
subroutine create_new(self, comm, nn, wtype, ngrid_in, lats_in, lons_in, ngrid_out, lats_out, lons_out)
Linked list implementation.
subroutine apply_ad(self, field_in, field_out)
subroutine create_read(self, comm, filename_in)