SABER
tools_samp.F90
Go to the documentation of this file.
1 # 1 "/Users/miesch/JEDI/code/working_copy/public/fv3-bundle/saber/src/saber/bump/tools_samp.fypp"
2 # 1 "/Users/miesch/JEDI/code/working_copy/public/fv3-bundle/saber/src/saber/bump/../instrumentation.fypp" 1
3 # 1 "/Users/miesch/JEDI/code/working_copy/public/fv3-bundle/saber/src/saber/bump/../subr_list.fypp" 1
4 !----------------------------------------------------------------------
5 ! Header: subr_list
6 !> Subroutines/functions list
7 ! Author: Benjamin Menetrier
8 ! Licensing: this code is distributed under the CeCILL-C license
9 ! Copyright 2015-... UCAR, CERFACS, METEO-FRANCE and IRIT
10 !----------------------------------------------------------------------
11 
12 # 726 "/Users/miesch/JEDI/code/working_copy/public/fv3-bundle/saber/src/saber/bump/../subr_list.fypp"
13 # 2 "/Users/miesch/JEDI/code/working_copy/public/fv3-bundle/saber/src/saber/bump/../instrumentation.fypp" 2
14 !----------------------------------------------------------------------
15 ! Header: instrumentation
16 !> Instrumentation functions
17 ! Author: Benjamin Menetrier
18 ! Licensing: this code is distributed under the CeCILL-C license
19 ! Copyright 2015-... UCAR, CERFACS, METEO-FRANCE and IRIT
20 !----------------------------------------------------------------------
21 
22 # 112 "/Users/miesch/JEDI/code/working_copy/public/fv3-bundle/saber/src/saber/bump/../instrumentation.fypp"
23 # 2 "/Users/miesch/JEDI/code/working_copy/public/fv3-bundle/saber/src/saber/bump/tools_samp.fypp" 2
24 !----------------------------------------------------------------------
25 ! Module: tools_samp
26 !> Sampling functions
27 ! Author: Benjamin Menetrier
28 ! Licensing: this code is distributed under the CeCILL-C license
29 ! Copyright 2015-... UCAR, CERFACS, METEO-FRANCE and IRIT
30 !----------------------------------------------------------------------
31 module tools_samp
32 
33 use atlas_module, only: atlas_structuredgrid
34 use fckit_mpi_module, only: fckit_mpi_sum,fckit_mpi_status
38 use tools_qsort, only: qsort
39 use tools_repro, only: repro,inf,sup,eq
40 use type_mesh, only: mesh_type
41 use type_mpl, only: mpl_type
42 
43 use type_rng, only: rng_type
44 use type_tree, only: tree_type
45 
46 implicit none
47 
49  module procedure samp_initialize_sampling
50 end interface
52  module procedure samp_initialize_sampling_local
53 end interface
55  module procedure samp_initialize_sampling_global
56 end interface
57 
58 private
60 
61 contains
62 
63 !----------------------------------------------------------------------
64 ! Subroutine: samp_initialize_sampling
65 !> Intialize sampling
66 !----------------------------------------------------------------------
67 subroutine samp_initialize_sampling(mpl,rng,area,n_loc,lon_loc,lat_loc,mask_loc,rh_loc,loc_to_glb,ntry,nrep,ns2_glb,sam2_glb, &
68  & fast,verbosity,n_uni,uni_to_loc,mesh_uni,tree_uni,inside_type)
69 
70 implicit none
71 
72 ! Passed variables
73 type(mpl_type),intent(inout) :: mpl !< MPI data
74 type(rng_type),intent(inout) :: rng !< Random number generator
75 real(kind_real),intent(in) :: area !< Global domain area
76 integer,intent(in) :: n_loc !< Number of points (local)
77 real(kind_real),intent(in) :: lon_loc(n_loc) !< Longitudes (local)
78 real(kind_real),intent(in) :: lat_loc(n_loc) !< Latitudes (local)
79 logical,intent(in) :: mask_loc(n_loc) !< Mask (local)
80 real(kind_real),intent(in) :: rh_loc(n_loc) !< Horizontal support radius (local)
81 integer,intent(in) :: loc_to_glb(n_loc) !< Local to global index
82 integer,intent(in) :: ntry !< Number of tries
83 integer,intent(in) :: nrep !< Number of replacements
84 integer,intent(in) :: ns2_glb !< Number of samplings points (global)
85 integer,intent(out) :: sam2_glb(ns2_glb) !< Horizontal sampling index (global)
86 logical,intent(in),optional :: fast !< Fast sampling flag
87 logical,intent(in),optional :: verbosity !< Verbosity flag
88 integer,intent(in),optional :: n_uni !< Universe size
89 integer,intent(in),optional :: uni_to_loc(:) !< Universe to local index
90 type(mesh_type),intent(in),optional :: mesh_uni !< Universe mesh
91 type(tree_type),intent(in),optional :: tree_uni !< Universe KD-tree
92 character(len=*),intent(in),optional :: inside_type !< Inside operator ('mesh_based' or 'tree_based')
93 
94 ! Local variables
95 integer :: n_glb,n_loc_eff,n_glb_eff,i_glb,i_loc,is2_glb,ns1_glb_eff
96 integer,allocatable :: sam1_glb_eff(:),order(:)
97 real(kind_real),allocatable :: hash_glb(:),hash_loc(:)
98 real(kind_real),allocatable :: lon1_glb_eff(:),lat1_glb_eff(:),rh1_glb_eff(:)
99 real(kind_real),allocatable :: list(:)
100 logical :: lfast,lverbosity
101 logical,allocatable :: mask_glb(:)
102 
103 ! Set name
104 
105 
106 ! Probe in
107 
108 
109 ! Local flags
110 lfast = .false.
111 lverbosity = .true.
112 if (present(fast)) lfast = fast
113 if (present(verbosity)) lverbosity = verbosity
114 
115 ! Global size
116 call mpl%f_comm%allreduce(n_loc,n_glb,fckit_mpi_sum())
117 
118 ! Number of effective points
119 n_loc_eff = count(mask_loc)
120 call mpl%f_comm%allreduce(n_loc_eff,n_glb_eff,fckit_mpi_sum())
121 
122 ! Check mask size
123 if (n_glb_eff==0) then
124  call mpl%abort('samp_initialize_sampling','empty mask in initialize sampling')
125 elseif (n_glb_eff<ns2_glb) then
126  call mpl%abort('samp_initialize_sampling','ns2_glb greater than n_glb_eff in initialize_sampling')
127 elseif (n_glb_eff==ns2_glb) then
128  write(mpl%info,'(a)') ' all points are used'
129  if (lverbosity) call mpl%flush
130 
131  ! Allocation
132  allocate(hash_loc(n_loc))
133  if (mpl%main) then
134  allocate(hash_glb(n_glb))
135  allocate(mask_glb(n_glb))
136  allocate(list(ns2_glb))
137  allocate(order(ns2_glb))
138  else
139  allocate(hash_glb(0))
140  allocate(mask_glb(0))
141  end if
142 
143  ! Compute hash
144  do i_loc=1,n_loc
145  hash_loc(i_loc) = lonlathash(lon_loc(i_loc),lat_loc(i_loc))
146  end do
147 
148  ! Communication
149  call mpl%loc_to_glb(n_loc,n_glb,loc_to_glb,hash_loc,hash_glb)
150  call mpl%loc_to_glb(n_loc,n_glb,loc_to_glb,mask_loc,mask_glb)
151 
152  if (mpl%main) then
153  ! Use all valid points
154  is2_glb = 0
155  do i_glb=1,n_glb
156  if (mask_glb(i_glb)) then
157  is2_glb = is2_glb+1
158  sam2_glb(is2_glb) = i_glb
159  list(is2_glb) = hash_glb(i_glb)
160  end if
161  end do
162 
163  ! Define points order
164  call qsort(ns2_glb,list,order)
165 
166  ! Reorder sampling
167  sam2_glb = sam2_glb(order)
168  end if
169 
170  ! Release memory
171  deallocate(hash_loc)
172  deallocate(hash_glb)
173  deallocate(mask_glb)
174  if (mpl%main) then
175  deallocate(list)
176  deallocate(order)
177  end if
178 else
179  ! First subsampling, local
180  if (present(n_uni).and.present(uni_to_loc).and.present(mesh_uni).and.present(tree_uni).and.present(inside_type)) then
181  call initialize_sampling_local(mpl,area,n_loc,n_loc_eff,n_glb_eff,lon_loc,lat_loc,mask_loc,rh_loc,loc_to_glb,ns2_glb, &
182  & ns1_glb_eff,lon1_glb_eff,lat1_glb_eff,rh1_glb_eff,sam1_glb_eff,lverbosity,n_uni,uni_to_loc,mesh_uni,tree_uni,inside_type)
183  else
184  call initialize_sampling_local(mpl,area,n_loc,n_loc_eff,n_glb_eff,lon_loc,lat_loc,mask_loc,rh_loc,loc_to_glb,ns2_glb, &
185  & ns1_glb_eff,lon1_glb_eff,lat1_glb_eff,rh1_glb_eff,sam1_glb_eff,lverbosity)
186  end if
187 
188  if (mpl%main) then
189  ! Second subsampling, global
190  call initialize_sampling_global(mpl,rng,ns1_glb_eff,lon1_glb_eff,lat1_glb_eff,rh1_glb_eff,sam1_glb_eff, &
191  & ntry,nrep,ns2_glb,sam2_glb,lfast,lverbosity)
192 
193  ! Release memory
194  deallocate(lon1_glb_eff)
195  deallocate(lat1_glb_eff)
196  deallocate(rh1_glb_eff)
197  deallocate(sam1_glb_eff)
198  end if
199 end if
200 
201 ! Broadcast
202 call mpl%f_comm%broadcast(sam2_glb,mpl%rootproc-1)
203 
204 ! Probe out
205 
206 
207 end subroutine samp_initialize_sampling
208 
209 !----------------------------------------------------------------------
210 ! Subroutine: samp_initialize_sampling_local
211 !> Intialize sampling, local, based on ATLAS octahedral grid
212 !----------------------------------------------------------------------
213 subroutine samp_initialize_sampling_local(mpl,area,n_loc,n_loc_eff,n_glb_eff,lon_loc,lat_loc,mask_loc,rh_loc,loc_to_glb,ns2_glb, &
214  & ns1_glb_eff,lon1_glb_eff,lat1_glb_eff,rh1_glb_eff,sam1_glb_eff,verbosity,n_uni,uni_to_loc,mesh_uni,tree_uni,inside_type)
215 
216 implicit none
217 
218 ! Passed variables
219 type(mpl_type),intent(inout) :: mpl !< MPI data
220 real(kind_real),intent(in) :: area !< Global domain area
221 integer,intent(in) :: n_loc !< Number of points (local)
222 integer,intent(in) :: n_loc_eff !< Number of points (local, effective)
223 integer,intent(in) :: n_glb_eff !< Number of points (global, effective)
224 real(kind_real),intent(in) :: lon_loc(n_loc) !< Longitudes (local)
225 real(kind_real),intent(in) :: lat_loc(n_loc) !< Latitudes (local)
226 logical,intent(in) :: mask_loc(n_loc) !< Mask (local)
227 real(kind_real),intent(in) :: rh_loc(n_loc) !< Horizontal support radius (local)
228 integer,intent(in) :: loc_to_glb(n_loc) !< Local to global index
229 integer,intent(in) :: ns2_glb !< Number of samplings points (global)
230 integer,intent(out) :: ns1_glb_eff !< Number of first sampling points (global, effective)
231 real(kind_real),allocatable,intent(out) :: lon1_glb_eff(:) !< First sampling longitudes (global, effective)
232 real(kind_real),allocatable,intent(out) :: lat1_glb_eff(:) !< First sampling latitudes (global, effective)
233 real(kind_real),allocatable,intent(out) :: rh1_glb_eff(:) !< First sampling horizontal support radius (global, effective)
234 integer,allocatable,intent(out) :: sam1_glb_eff(:) !< First sampling indices (global, effective)
235 logical,intent(in),optional :: verbosity !< Verbosity flag
236 integer,intent(in),optional :: n_uni !< Universe size
237 integer,intent(in),optional :: uni_to_loc(:) !< Universe to local index
238 type(mesh_type),intent(in),optional :: mesh_uni !< Universe mesh
239 type(tree_type),intent(in),optional :: tree_uni !< Universe KD-tree
240 character(len=*),intent(in),optional :: inside_type !< Inside operator ('mesh_based' or 'tree_based')
241 
242 ! Local variables
243 integer :: i_loc,i_loc_eff,n,ix,iy,iproc,is1_glb,ns1_loc,is1_loc,nfac,ns1_loc_tmp,nn_index(1),ns1_glb,is1_glb_eff,offset
244 integer,allocatable :: s1_loc_to_glb(:),s1_loc_to_glb_tmp(:),sam1_loc(:),sam1_loc_tmp(:),sam1_glb(:)
245 real(kind_real) :: lonlat(2),nn_dist(1)
246 real(kind_real),allocatable :: lon1_loc(:),lat1_loc(:),dist1_loc(:),rh1_loc(:)
247 real(kind_real),allocatable :: lon1_loc_tmp(:),lat1_loc_tmp(:),dist1_loc_tmp(:),rh1_loc_tmp(:)
248 real(kind_real),allocatable :: lon1_glb(:),lat1_glb(:),dist1_glb(:),rh1_glb(:)
249 logical :: lfull_grid,lverbosity,update,retry,valid
250 logical,allocatable :: lmask(:)
251 character(len=6) :: gridid
252 type(fckit_mpi_status) :: status
253 type(atlas_structuredgrid) :: agrid
254 
255 ! Set name
256 
257 
258 ! Probe in
259 
260 
261 ! Local flags
262 lverbosity = .true.
263 if (present(verbosity)) lverbosity = verbosity
264 lfull_grid = present(n_uni).and.present(uni_to_loc).and.present(mesh_uni).and.present(tree_uni).and.present(inside_type)
265 
266 ! Initialization
267 nfac = 1
268 retry = .true.
269 
270 do while (retry)
271  ! Update nfac
272  nfac = 2*nfac
273 
274  if (lfull_grid) then
275  ! Number of required points
276  ns1_glb = nfac*int(real(ns2_glb,kind_real)*four*pi/area)
277 
278  ! Octahedral grid
279  n = int(-4.5_kind_real+sqrt(20.25_kind_real+quarter*real(ns1_glb,kind_real)))+1
280  write(gridid,'(a,i5.5)') 'O',n
281  agrid = atlas_structuredgrid(gridid)
282  ns1_glb = 4*n**2+36*n
283 
284  ! Allocation
285  allocate(lon1_loc(ns1_glb))
286  allocate(lat1_loc(ns1_glb))
287  allocate(rh1_loc(ns1_glb))
288  allocate(sam1_loc(ns1_glb))
289  allocate(dist1_loc(ns1_glb))
290  allocate(s1_loc_to_glb(ns1_glb))
291  allocate(lmask(n_loc))
292 
293  ! Initialization
294  lmask = mask_loc
295  if (lverbosity) call mpl%prog_init(ns1_glb)
296 
297  ! Loop over octahedral grid points
298  is1_loc = 0
299  is1_glb = 0
300  do iy=1,int(agrid%ny(),kind_int)
301  do ix=1,int(agrid%nx(iy),kind_int)
302  ! Global index
303  is1_glb = is1_glb+1
304 
305  ! Get longitude/latitude
306  lonlat = agrid%lonlat(ix,iy)*deg2rad
307  call lonlatmod(lonlat(1),lonlat(2))
308 
309  ! Check if point is inside the hull
310  if (trim(inside_type)=='mesh_based') then
311  call mesh_uni%inside(mpl,lonlat(1),lonlat(2),valid)
312  elseif (trim(inside_type)=='tree_based') then
313  call tree_uni%inside(lonlat(1),lonlat(2),valid)
314  end if
315 
316  if (valid) then
317  ! Find nearest neighbor in universe
318  call tree_uni%find_nearest_neighbors(lonlat(1),lonlat(2),1,nn_index,nn_dist)
319  i_loc = uni_to_loc(nn_index(1))
320 
321  if (mpl%msv%isnot(i_loc)) then
322  ! Keep valid points
323  if (lmask(i_loc)) then
324  is1_loc = is1_loc+1
325  lon1_loc(is1_loc) = lon_loc(i_loc)
326  lat1_loc(is1_loc) = lat_loc(i_loc)
327  rh1_loc(is1_loc) = rh_loc(i_loc)
328  sam1_loc(is1_loc) = loc_to_glb(i_loc)
329  dist1_loc(is1_loc) = nn_dist(1)
330  s1_loc_to_glb(is1_loc) = is1_glb
331  lmask(i_loc) = .false.
332  end if
333  end if
334  end if
335 
336  ! Update
337  if (lverbosity) call mpl%prog_print(is1_glb)
338  end do
339  end do
340  if (lverbosity) call mpl%prog_final(.false.)
341 
342  ! Number of local valid points
343  ns1_loc = is1_loc
344 
345  ! Release memory
346  deallocate(lmask)
347  else
348  ! All point are used
349  ns1_glb = n_glb_eff
350  ns1_loc = n_loc_eff
351 
352  ! Allocation
353  allocate(lon1_loc(ns1_loc))
354  allocate(lat1_loc(ns1_loc))
355  allocate(rh1_loc(ns1_loc))
356  allocate(sam1_loc(ns1_loc))
357 
358  ! Copy
359  i_loc_eff = 0
360  do i_loc=1,n_loc
361  if (mask_loc(i_loc)) then
362  i_loc_eff = i_loc_eff+1
363  lon1_loc(i_loc_eff) = lon_loc(i_loc)
364  lat1_loc(i_loc_eff) = lat_loc(i_loc)
365  rh1_loc(i_loc_eff) = rh_loc(i_loc)
366  sam1_loc(i_loc_eff) = loc_to_glb(i_loc)
367  end if
368  end do
369  end if
370 
371  if (mpl%main) then
372  if (lfull_grid) then
373  ! Continue printing
374  write(mpl%info,'(a)') ' => '
375  if (lverbosity) call mpl%flush(.false.)
376  end if
377 
378  ! Allocation
379  allocate(lon1_glb(ns1_glb))
380  allocate(lat1_glb(ns1_glb))
381  allocate(rh1_glb(ns1_glb))
382  allocate(sam1_glb(ns1_glb))
383  if (lfull_grid) allocate(dist1_glb(ns1_glb))
384 
385  ! Initialization
386  lon1_glb = mpl%msv%valr
387  lat1_glb = mpl%msv%valr
388  rh1_glb = mpl%msv%valr
389  sam1_glb = mpl%msv%vali
390  if (lfull_grid) then
391  dist1_glb = mpl%msv%valr
392  else
393  offset = 0
394  end if
395 
396  ! Receive and copy data on rootproc
397  do iproc=1,mpl%nproc
398  if (iproc==mpl%rootproc) then
399  ! Copy dimension
400  ns1_loc_tmp = ns1_loc
401  else
402  ! Receive dimension
403  call mpl%f_comm%receive(ns1_loc_tmp,iproc-1,mpl%tag,status)
404  end if
405 
406  ! Allocation
407  allocate(lon1_loc_tmp(ns1_loc_tmp))
408  allocate(lat1_loc_tmp(ns1_loc_tmp))
409  allocate(rh1_loc_tmp(ns1_loc_tmp))
410  allocate(sam1_loc_tmp(ns1_loc_tmp))
411  if (lfull_grid) then
412  allocate(dist1_loc_tmp(ns1_loc_tmp))
413  allocate(s1_loc_to_glb_tmp(ns1_loc_tmp))
414  end if
415 
416  if (iproc==mpl%rootproc) then
417  ! Copy data
418  lon1_loc_tmp = lon1_loc(1:ns1_loc)
419  lat1_loc_tmp = lat1_loc(1:ns1_loc)
420  rh1_loc_tmp = rh1_loc(1:ns1_loc)
421  sam1_loc_tmp = sam1_loc(1:ns1_loc)
422  if (lfull_grid) then
423  dist1_loc_tmp = dist1_loc(1:ns1_loc)
424  s1_loc_to_glb_tmp = s1_loc_to_glb(1:ns1_loc)
425  end if
426  else
427  ! Receive data
428  call mpl%f_comm%receive(lon1_loc_tmp,iproc-1,mpl%tag+1,status)
429  call mpl%f_comm%receive(lat1_loc_tmp,iproc-1,mpl%tag+2,status)
430  call mpl%f_comm%receive(rh1_loc_tmp,iproc-1,mpl%tag+3,status)
431  call mpl%f_comm%receive(sam1_loc_tmp,iproc-1,mpl%tag+4,status)
432  if (lfull_grid) then
433  call mpl%f_comm%receive(dist1_loc_tmp,iproc-1,mpl%tag+5,status)
434  call mpl%f_comm%receive(s1_loc_to_glb_tmp,iproc-1,mpl%tag+6,status)
435  end if
436  end if
437 
438  ! Update global array
439  if (lfull_grid) then
440  do is1_loc=1,ns1_loc_tmp
441  ! Global index
442  is1_glb = s1_loc_to_glb_tmp(is1_loc)
443 
444  ! Check the global array status
445  update = mpl%msv%is(dist1_glb(is1_glb))
446  if (.not.update) update = (dist1_loc_tmp(is1_loc)<dist1_glb(is1_glb))
447  if (update) then
448  lon1_glb(is1_glb) = lon1_loc_tmp(is1_loc)
449  lat1_glb(is1_glb) = lat1_loc_tmp(is1_loc)
450  rh1_glb(is1_glb) = rh1_loc_tmp(is1_loc)
451  sam1_glb(is1_glb) = sam1_loc_tmp(is1_loc)
452  dist1_glb(is1_glb) = dist1_loc_tmp(is1_loc)
453  end if
454  end do
455  else
456  ! Copy with offset
457  lon1_glb(offset+1:offset+ns1_loc_tmp) = lon1_loc_tmp
458  lat1_glb(offset+1:offset+ns1_loc_tmp) = lat1_loc_tmp
459  rh1_glb(offset+1:offset+ns1_loc_tmp) = rh1_loc_tmp
460  sam1_glb(offset+1:offset+ns1_loc_tmp) = sam1_loc_tmp
461  offset = offset+ns1_loc_tmp
462  end if
463 
464  ! Release memory
465  deallocate(lon1_loc_tmp)
466  deallocate(lat1_loc_tmp)
467  deallocate(rh1_loc_tmp)
468  deallocate(sam1_loc_tmp)
469  if (lfull_grid) then
470  deallocate(dist1_loc_tmp)
471  deallocate(s1_loc_to_glb_tmp)
472  end if
473  end do
474 
475  ! Number of global valid points
476  if (lfull_grid) then
477  ns1_glb_eff = count(mpl%msv%isnot(dist1_glb))
478  else
479  ns1_glb_eff = ns1_glb
480  end if
481  retry = (ns1_glb_eff<ns2_glb)
482  else
483  ! Send data to rootproc
484  call mpl%f_comm%send(ns1_loc,mpl%rootproc-1,mpl%tag)
485  call mpl%f_comm%send(lon1_loc(1:ns1_loc),mpl%rootproc-1,mpl%tag+1)
486  call mpl%f_comm%send(lat1_loc(1:ns1_loc),mpl%rootproc-1,mpl%tag+2)
487  call mpl%f_comm%send(rh1_loc(1:ns1_loc),mpl%rootproc-1,mpl%tag+3)
488  call mpl%f_comm%send(sam1_loc(1:ns1_loc),mpl%rootproc-1,mpl%tag+4)
489  if (lfull_grid) then
490  call mpl%f_comm%send(dist1_loc(1:ns1_loc),mpl%rootproc-1,mpl%tag+5)
491  call mpl%f_comm%send(s1_loc_to_glb(1:ns1_loc),mpl%rootproc-1,mpl%tag+6)
492  end if
493 
494  ! Number of global valid points
495  ns1_glb_eff = 0
496 
497  ! Stop printing
498  write(mpl%info,'(a)') ''
499  if (lverbosity) call mpl%flush
500  end if
501  if (lfull_grid) then
502  call mpl%update_tag(7)
503  else
504  call mpl%update_tag(5)
505  end if
506 
507  ! Broadcast
508  call mpl%f_comm%broadcast(retry,mpl%rootproc-1)
509 
510  ! Check situation
511  if ((.not.lfull_grid).and.retry) call mpl%abort('samp_initialize_sampling_local','retry activated for limited grid')
512  if ((nfac>=8).and.retry) call mpl%abort('samp_initialize_sampling_local','retry activated for nfac>=8')
513 
514  ! Release memory
515  deallocate(lon1_loc)
516  deallocate(lat1_loc)
517  deallocate(rh1_loc)
518  deallocate(sam1_loc)
519  if (lfull_grid) then
520  deallocate(dist1_loc)
521  deallocate(s1_loc_to_glb)
522  end if
523  if (mpl%main.and.retry) then
524  deallocate(lon1_glb)
525  deallocate(lat1_glb)
526  deallocate(rh1_glb)
527  deallocate(sam1_glb)
528  if (lfull_grid) deallocate(dist1_glb)
529  end if
530 end do
531 
532 if (mpl%main) then
533  ! Allocation
534  allocate(lon1_glb_eff(ns1_glb_eff))
535  allocate(lat1_glb_eff(ns1_glb_eff))
536  allocate(rh1_glb_eff(ns1_glb_eff))
537  allocate(sam1_glb_eff(ns1_glb_eff))
538 
539  ! Copy valid data
540  if (lfull_grid) then
541  is1_glb_eff = 0
542  do is1_glb=1,ns1_glb
543  if (mpl%msv%isnot(dist1_glb(is1_glb))) then
544  is1_glb_eff = is1_glb_eff+1
545  lon1_glb_eff(is1_glb_eff) = lon1_glb(is1_glb)
546  lat1_glb_eff(is1_glb_eff) = lat1_glb(is1_glb)
547  rh1_glb_eff(is1_glb_eff) = rh1_glb(is1_glb)
548  sam1_glb_eff(is1_glb_eff) = sam1_glb(is1_glb)
549  end if
550  end do
551  else
552  lon1_glb_eff = lon1_glb
553  lat1_glb_eff = lat1_glb
554  rh1_glb_eff = rh1_glb
555  sam1_glb_eff = sam1_glb
556  end if
557 
558  ! Release memory
559  deallocate(lon1_glb)
560  deallocate(lat1_glb)
561  deallocate(rh1_glb)
562  deallocate(sam1_glb)
563  if (lfull_grid) deallocate(dist1_glb)
564 end if
565 
566 ! Probe out
567 
568 
569 end subroutine samp_initialize_sampling_local
570 
571 !----------------------------------------------------------------------
572 ! Subroutine: samp_initialize_sampling_global
573 !> Intialize sampling, global
574 !----------------------------------------------------------------------
575 subroutine samp_initialize_sampling_global(mpl,rng,ns1_glb_eff,lon1_glb_eff,lat1_glb_eff,rh1_glb_eff,sam1_glb_eff,ntry,nrep, &
576  & ns2_glb,sam2_glb,fast,verbosity)
577 
578 implicit none
579 
580 ! Passed variables
581 type(mpl_type),intent(inout) :: mpl !< MPI data
582 type(rng_type),intent(inout) :: rng !< Random number generator
583 integer,intent(in) :: ns1_glb_eff
584 real(kind_real),intent(in) :: lon1_glb_eff(ns1_glb_eff)
585 real(kind_real),intent(in) :: lat1_glb_eff(ns1_glb_eff)
586 real(kind_real),intent(in) :: rh1_glb_eff(ns1_glb_eff)
587 integer,intent(in) :: sam1_glb_eff(ns1_glb_eff)
588 integer,intent(in) :: ntry !< Number of tries
589 integer,intent(in) :: nrep !< Number of replacements
590 integer,intent(in) :: ns2_glb !< Number of samplings points (global)
591 integer,intent(out) :: sam2_glb(ns2_glb) !< Horizontal sampling index (global)
592 logical,intent(in),optional :: fast !< Fast sampling flag
593 logical,intent(in),optional :: verbosity !< Verbosity flag
594 
595 ! Local variables
596 integer :: is2_glb,js,irep,irmax,itry,ir,irval,irvalmin,irvalmax,is2_glb_min,nrep_eff,nn_index(2),is1_glb_eff,ns1_glb_val
597 integer :: order(ns1_glb_eff),sam1_glb_tmp(ns1_glb_eff)
598 integer,allocatable :: to_valid(:),sam2_glb_tmp(:)
599 real(kind_real) :: d,distmax,distmin,nn_dist(2),cdf_norm,rr
600 real(kind_real) :: list(ns1_glb_eff),lon1_glb_tmp(ns1_glb_eff),lat1_glb_tmp(ns1_glb_eff),rh1_glb_tmp(ns1_glb_eff)
601 real(kind_real),allocatable :: cdf(:),lon_rep(:),lat_rep(:),dist(:)
602 logical :: lfast,lverbosity
603 logical,allocatable :: lmask(:),smask(:),rmask(:)
604 type(tree_type) :: tree
605 
606 ! Set name
607 
608 
609 ! Probe in
610 
611 
612 ! Local flags
613 lfast = .false.
614 lverbosity = .true.
615 if (present(fast)) lfast = fast
616 if (present(verbosity)) lverbosity = verbosity
617 
618 ! Define points order
619 do is1_glb_eff=1,ns1_glb_eff
620  list(is1_glb_eff) = lonlathash(lon1_glb_eff(is1_glb_eff),lat1_glb_eff(is1_glb_eff))
621 end do
622 call qsort(ns1_glb_eff,list,order)
623 
624 ! Reorder data
625 lon1_glb_tmp = lon1_glb_eff(order)
626 lat1_glb_tmp = lat1_glb_eff(order)
627 rh1_glb_tmp = rh1_glb_eff(order)
628 sam1_glb_tmp = sam1_glb_eff(order)
629 
630 ! Allocation
631 nrep_eff = min(nrep,ns1_glb_eff-ns2_glb)
632 allocate(sam2_glb_tmp(ns2_glb+nrep_eff))
633 allocate(lmask(ns1_glb_eff))
634 allocate(smask(ns1_glb_eff))
635 allocate(to_valid(ns1_glb_eff))
636 
637 ! Initialization
638 sam2_glb_tmp = mpl%msv%vali
639 lmask = .true.
640 smask = .false.
641 to_valid = mpl%msv%vali
642 do is1_glb_eff=1,ns1_glb_eff
643  to_valid(is1_glb_eff) = is1_glb_eff
644 end do
645 ns1_glb_val = ns1_glb_eff
646 if (lverbosity) call mpl%prog_init(ns2_glb+nrep_eff)
647 
648 if (lfast) then
649  ! Define sampling with a cumulative distribution function
650 
651  ! Allocation
652  allocate(cdf(ns1_glb_eff))
653 
654  ! Initialization
655  cdf(1) = zero
656  do is1_glb_eff=2,ns1_glb_eff
657  if (lmask(is1_glb_eff)) then
658  cdf(is1_glb_eff) = cdf(is1_glb_eff-1)+one/rh1_glb_tmp(is1_glb_eff)**2
659  end if
660  end do
661  cdf_norm = one/cdf(ns1_glb_eff)
662  cdf = cdf*cdf_norm
663 
664  do is2_glb=1,ns2_glb+nrep_eff
665  ! Generate random number
666  call rng%rand(zero,one,rr)
667 
668  ! Dichotomy to find the value
669  irvalmin = 1
670  irvalmax = ns1_glb_val
671  do while (irvalmax-irvalmin>1)
672  irval = (irvalmin+irvalmax)/2
673  if ((sup(cdf(irvalmin),rr).and.sup(cdf(irval),rr)).or.(inf(cdf(irvalmin),rr).and.inf(cdf(irval),rr))) then
674  irvalmin = irval
675  else
676  irvalmax = irval
677  end if
678  end do
679 
680  ! New sampling point
681  ir = to_valid(irval)
682  sam2_glb_tmp(is2_glb) = ir
683 
684  ! Shift valid points array
685  if (irval<ns1_glb_val) then
686  cdf(irval:ns1_glb_val-1) = cdf(irval+1:ns1_glb_val)
687  to_valid(irval:ns1_glb_val-1) = to_valid(irval+1:ns1_glb_val)
688  end if
689  ns1_glb_val = ns1_glb_val-1
690 
691  ! Renormalize cdf
692  cdf_norm = one/cdf(ns1_glb_val)
693  cdf(1:ns1_glb_val) = cdf(1:ns1_glb_val)*cdf_norm
694 
695  ! Update
696  if (lverbosity) call mpl%prog_print(is2_glb)
697  end do
698  if (lverbosity) call mpl%prog_final(.false.)
699 
700  ! Release memory
701  deallocate(cdf)
702 else
703  ! Define sampling with a KD-tree
704  do is2_glb=1,ns2_glb+nrep_eff
705  if (is2_glb>2) then
706  ! Allocation
707  call tree%alloc(mpl,ns1_glb_eff,mask=smask)
708 
709  ! Initialization
710  call tree%init(lon1_glb_tmp,lat1_glb_tmp)
711  end if
712 
713  ! Initialization
714  distmax = zero
715  irmax = 0
716  irvalmax = 0
717  itry = 1
718 
719  ! Find a new point
720  do itry=1,ntry
721  ! Generate a random index among valid points
722  call rng%rand(1,ns1_glb_val,irval)
723  ir = to_valid(irval)
724 
725  ! Check point validity
726  if (is2_glb==1) then
727  ! Accept point
728  irvalmax = irval
729  irmax = ir
730  else
731  if (is2_glb==2) then
732  ! Compute distance
733  nn_index(1) = sam2_glb_tmp(1)
734  call sphere_dist(lon1_glb_tmp(ir),lat1_glb_tmp(ir),lon1_glb_tmp(nn_index(1)),lat1_glb_tmp(nn_index(1)), &
735  & nn_dist(1))
736  else
737  ! Find nearest neighbor distance
738  call tree%find_nearest_neighbors(lon1_glb_tmp(ir),lat1_glb_tmp(ir),1,nn_index(1:1),nn_dist(1:1))
739  end if
740  d = nn_dist(1)**2/(rh1_glb_tmp(ir)**2+rh1_glb_tmp(nn_index(1))**2)
741 
742  ! Check distance
743  if (sup(d,distmax)) then
744  distmax = d
745  irvalmax = irval
746  irmax = ir
747  end if
748  end if
749  end do
750 
751  ! Delete tree
752  if (is2_glb>2) call tree%dealloc
753 
754  ! Add point to sampling
755  if (irmax>0) then
756  ! New sampling point
757  sam2_glb_tmp(is2_glb) = irmax
758  lmask(irmax) = .false.
759  smask(irmax) = .true.
760 
761  ! Shift valid points array
762  if (irvalmax<ns1_glb_val) to_valid(irvalmax:ns1_glb_val-1) = to_valid(irvalmax+1:ns1_glb_val)
763  ns1_glb_val = ns1_glb_val-1
764  end if
765 
766  ! Update
767  if (lverbosity) call mpl%prog_print(is2_glb)
768  end do
769  if (lverbosity) call mpl%prog_final(.false.)
770 end if
771 
772 if (nrep_eff>0) then
773  ! Continue printing
774  write(mpl%info,'(a)') ' => '
775  if (lverbosity) call mpl%flush(.false.)
776 
777  ! Allocation
778  allocate(rmask(ns2_glb+nrep_eff))
779  allocate(lon_rep(ns2_glb+nrep_eff))
780  allocate(lat_rep(ns2_glb+nrep_eff))
781  allocate(dist(ns2_glb+nrep_eff))
782 
783  ! Initialization
784  rmask = .true.
785  do is2_glb=1,ns2_glb+nrep_eff
786  lon_rep(is2_glb) = lon1_glb_tmp(sam2_glb_tmp(is2_glb))
787  lat_rep(is2_glb) = lat1_glb_tmp(sam2_glb_tmp(is2_glb))
788  end do
789  dist = mpl%msv%valr
790  if (lverbosity) call mpl%prog_init(nrep_eff)
791 
792  ! Remove closest points
793  do irep=1,nrep_eff
794  ! Allocation
795  call tree%alloc(mpl,ns2_glb+nrep_eff,mask=rmask)
796 
797  ! Initialization
798  call tree%init(lon_rep,lat_rep)
799 
800  ! Get minimum distance
801  do is2_glb=1,ns2_glb+nrep_eff
802  if (rmask(is2_glb)) then
803  ! Find nearest neighbor distance
804  call tree%find_nearest_neighbors(lon1_glb_tmp(sam2_glb_tmp(is2_glb)),lat1_glb_tmp(sam2_glb_tmp(is2_glb)), &
805  & 2,nn_index,nn_dist)
806  if (nn_index(1)==is2_glb) then
807  dist(is2_glb) = nn_dist(2)
808  elseif (nn_index(2)==is2_glb) then
809  dist(is2_glb) = nn_dist(1)
810  else
811  call mpl%abort('samp_initialize_sampling_global','wrong index in replacement')
812  end if
813  dist(is2_glb) = dist(is2_glb)**2/(rh1_glb_tmp(sam2_glb_tmp(nn_index(1)))**2 &
814  & +rh1_glb_tmp(sam2_glb_tmp(nn_index(2)))**2)
815  end if
816  end do
817 
818  ! Delete tree
819  call tree%dealloc
820 
821  ! Remove worst point
822  distmin = huge_real
823  is2_glb_min = mpl%msv%vali
824  do is2_glb=1,ns2_glb+nrep_eff
825  if (rmask(is2_glb)) then
826  if (inf(dist(is2_glb),distmin)) then
827  is2_glb_min = is2_glb
828  distmin = dist(is2_glb)
829  end if
830  end if
831  end do
832  rmask(is2_glb_min) = .false.
833 
834  ! Update
835  if (lverbosity) call mpl%prog_print(irep)
836  end do
837  if (lverbosity) call mpl%prog_final
838 
839  ! Copy sam2_glb
840  js = 0
841  do is2_glb=1,ns2_glb+nrep_eff
842  if (rmask(is2_glb)) then
843  js = js+1
844  sam2_glb(js) = sam2_glb_tmp(is2_glb)
845  end if
846  end do
847 
848  ! Release memory
849  deallocate(rmask)
850  deallocate(lon_rep)
851  deallocate(lat_rep)
852  deallocate(dist)
853 else
854  ! Stop printing
855  write(mpl%info,'(a)') ''
856  if (lverbosity) call mpl%flush
857 
858  ! Copy sam2_glb
859  sam2_glb = sam2_glb_tmp
860 end if
861 
862 ! Apply first sampling step
863 sam2_glb = sam1_glb_tmp(sam2_glb)
864 
865 ! Release memory
866 deallocate(sam2_glb_tmp)
867 deallocate(lmask)
868 deallocate(smask)
869 deallocate(to_valid)
870 
871 ! Probe out
872 
873 
874 end subroutine samp_initialize_sampling_global
875 
876 end module tools_samp
Subroutines/functions list.
Definition: tools_const.F90:31
real(kind_real), parameter, public one
One.
Definition: tools_const.F90:42
real(kind_real), parameter, public pi
Pi.
Definition: tools_const.F90:51
real(kind_real), parameter, public deg2rad
Degree to radian.
Definition: tools_const.F90:52
real(kind_real), parameter, public zero
Zero.
Definition: tools_const.F90:37
real(kind_real), parameter, public four
Four.
Definition: tools_const.F90:45
real(kind_real), parameter, public quarter
Quarter.
Definition: tools_const.F90:40
Subroutines/functions list.
Definition: tools_func.F90:31
Kinds definition.
Definition: tools_kinds.F90:9
integer, parameter, public kind_int
Integer kind.
Definition: tools_kinds.F90:17
integer, parameter, public kind_real
Real kind alias for the whole code.
Definition: tools_kinds.F90:25
real(kind_real), parameter, public huge_real
Real huge.
Definition: tools_kinds.F90:33
Generic ranks, dimensions and types.
Definition: tools_qsort.F90:46
Generic ranks, dimensions and types.
Definition: tools_repro.F90:42
logical, public repro
Reproducibility flag.
Definition: tools_repro.F90:50
Subroutines/functions list.
Definition: tools_samp.F90:31
subroutine samp_initialize_sampling_local(mpl, area, n_loc, n_loc_eff, n_glb_eff, lon_loc, lat_loc, mask_loc, rh_loc, loc_to_glb, ns2_glb, ns1_glb_eff, lon1_glb_eff, lat1_glb_eff, rh1_glb_eff, sam1_glb_eff, verbosity, n_uni, uni_to_loc, mesh_uni, tree_uni, inside_type)
Intialize sampling, local, based on ATLAS octahedral grid.
Definition: tools_samp.F90:215
subroutine samp_initialize_sampling_global(mpl, rng, ns1_glb_eff, lon1_glb_eff, lat1_glb_eff, rh1_glb_eff, sam1_glb_eff, ntry, nrep, ns2_glb, sam2_glb, fast, verbosity)
Intialize sampling, global.
Definition: tools_samp.F90:577
subroutine samp_initialize_sampling(mpl, rng, area, n_loc, lon_loc, lat_loc, mask_loc, rh_loc, loc_to_glb, ntry, nrep, ns2_glb, sam2_glb, fast, verbosity, n_uni, uni_to_loc, mesh_uni, tree_uni, inside_type)
Intialize sampling.
Definition: tools_samp.F90:69
Subroutines/functions list.
Definition: type_mesh.F90:31
Generic ranks, dimensions and types.
Definition: type_mpl.F90:42
Generic ranks, dimensions and types.
Definition: type_rng.F90:42
Subroutines/functions list.
Definition: type_tree.F90:31