8 use atlas_module,
only: atlas_geometry, atlas_indexkdtree
13 use string_utils,
only: replace_string
14 use netcdf_utils_mod,
only: nccheck
16 use fckit_mpi_module,
only: fckit_mpi_comm, fckit_mpi_sum
29 integer :: ngrid_in_glo
30 type(fckit_mpi_comm) :: comm
31 integer,
allocatable :: displs(:)
32 integer,
allocatable :: rcvcnt(:)
33 real(kind=kind_real),
allocatable :: interp_w(:,:)
34 integer ,
allocatable :: interp_i(:,:)
47 #define LISTED_TYPE unstrc_interp
50 #include "oops/util/linkedList_i.f"
60 #include "oops/util/linkedList_c.f"
64 subroutine create_new( self, comm, nn, wtype, ngrid_in, lats_in, lons_in, &
65 ngrid_out, lats_out, lons_out )
67 type(fckit_mpi_comm),
intent(in) :: comm
68 integer,
intent(in) :: nn
69 character(len=*),
intent(in) :: wtype
70 integer,
intent(in) :: ngrid_in
71 real(kind=kind_real),
intent(in) :: lats_in(ngrid_in)
72 real(kind=kind_real),
intent(in) :: lons_in(ngrid_in)
73 integer,
intent(in) :: ngrid_out
74 real(kind=kind_real),
intent(in) :: lats_out(ngrid_out)
75 real(kind=kind_real),
intent(in) :: lons_out(ngrid_out)
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(:)
96 self%ngrid_in = ngrid_in
97 self%ngrid_out = ngrid_out
99 allocate(self%interp_w(self%nn,self%ngrid_out))
100 allocate(self%interp_i(self%nn,self%ngrid_out))
108 ngrid_in_glo = sum(self%rcvcnt)
109 self%ngrid_in_glo = ngrid_in_glo
115 self%displs(j) = self%displs(j-1) + self%rcvcnt(j-1)
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
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)
131 ageometry = atlas_geometry(
"UnitSphere")
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)
140 allocate(nn_dist(nn,ngrid_out))
145 call kd%closestPoints(lons_out_loc(n), lats_out(n), self%nn, self%interp_i(:,n))
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))
160 self%interp_w = 0.0_kind_real
163 self%interp_w = nn_dist
166 allocate(bw(self%nn))
169 bw(:) = 0.0_kind_real
171 wprod = 1.0_kind_real
172 index1 = self%interp_i(jj,n)
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)
181 bw(jj) = 1.0_kind_real / wprod
185 self%interp_w(:,n) = 0.0_kind_real
186 if (minval(nn_dist(:,n)) < 1e-10)
then
189 jj = minloc(nn_dist(:,n),dim=1)
190 self%interp_w(jj,n) = 1.0_kind_real
197 bsw = bsw + (bw(jj) / nn_dist(jj,n))
201 self%interp_w(jj,n) = ( bw(jj) / nn_dist(jj,n) ) / bsw
210 call abor1_ftn(
"unstrc_interpolation.create: wrong choice of weight output")
217 deallocate(lats_in_glo, lons_in_glo)
225 type(fckit_mpi_comm),
intent(in) :: comm
226 character(len=*),
intent(in) :: filename_in
228 character(len=2048) :: filename
229 integer :: ncid, varid
234 filename = filename_in
238 call nccheck ( nf90_open( trim(filename), nf90_nowrite, ncid ),
"nf90_open"//trim(filename) )
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" )
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" )
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" )
252 allocate(self%interp_w(self%nn,self%ngrid_out))
253 allocate(self%interp_i(self%nn,self%ngrid_out))
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" )
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" )
264 call nccheck ( nf90_close(ncid),
"nf90_close" )
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)
282 subroutine apply(self, field_in, field_out, field_nn_out )
284 real(kind=kind_real),
intent(in) :: field_in(self%ngrid_in)
285 real(kind=kind_real),
intent(out) :: field_out(self%ngrid_out)
286 real(kind=kind_real),
optional,
intent(out) :: field_nn_out(self%nn,self%ngrid_out)
290 real(kind=kind_real),
allocatable :: field_nn(:,:), field_in_glo(:)
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)
298 allocate(field_nn(self%nn,self%ngrid_out))
299 do n = 1, self%ngrid_out
301 field_nn(kk,n) = field_in_glo(self%interp_i(kk,n))
306 deallocate(field_in_glo)
310 do n = 1, self%ngrid_out
311 field_out(n) = 0.0_kind_real
313 field_out(n) = field_out(n) + self%interp_w(kk,n) * field_nn(kk,n)
318 if (
present(field_nn_out)) field_nn_out = field_nn
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)
331 integer :: n, kk, nstart, nfinish
332 real(kind=kind_real),
allocatable :: field_in_glo(:), field_in_glo_all(:)
337 allocate(field_in_glo(self%ngrid_in_glo))
339 do n = 1, self%ngrid_out
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)
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())
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)
364 character(len=*),
intent(in) :: filename_in
366 integer :: ncid, ni_dimid, no_dimid, nn_dimid, vc(2)
367 character(len=2048) :: filename
370 filename = filename_in
374 call nccheck( nf90_create( trim(filename), ior(nf90_netcdf4, nf90_clobber), ncid),
"nf90_create" )
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" )
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" );
386 call nccheck( nf90_enddef(ncid),
"nf90_enddef" )
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" );
393 call nccheck ( nf90_close(ncid),
"nf90_close" )
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(:)
408 allocate(rcvcnt(comm%size()))
409 call comm%allgather(ngrid_in, rcvcnt)
412 allocate(displs(comm%size()))
415 displs(j) = displs(j-1) + rcvcnt(j-1)
423 integer,
intent(in) :: procrank
424 character(len=*),
intent(inout) :: filename
426 character(len=6) :: procrankstr
429 if (procrank>999999)
call abor1_ftn(
"unstructured_interpolation_mod not ready for more than 999999 processors")
432 write(procrankstr,
'(I0.6)') procrank
435 filename = replace_string(filename,
'%PROC',procrankstr)
subroutine write(self, filename_in)
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)