FV3 Bundle
type_com.F90
Go to the documentation of this file.
1 !----------------------------------------------------------------------
2 ! Module: type_com
3 ! Purpose: communications derived type
4 ! Author: Benjamin Menetrier
5 ! Licensing: this code is distributed under the CeCILL-C license
6 ! Copyright © 2015-... UCAR, CERFACS, METEO-FRANCE and IRIT
7 !----------------------------------------------------------------------
8 module type_com
9 
10 use netcdf
11 !$ use omp_lib
12 use tools_kinds, only: kind_real
13 use tools_missing, only: msi,msr
14 use tools_nc, only: ncfloat
15 use type_mpl, only: mpl_type
16 use fckit_mpi_module, only: fckit_mpi_status
17 
18 implicit none
19 
20 ! Communication derived type
22  ! Setup data
23  integer,allocatable :: ext_to_proc(:) ! Extended index to processor
24  integer,allocatable :: ext_to_red(:) ! Extended index to reduced index
25 
26  ! Communication data
27  character(len=1024) :: prefix ! Communication prefix
28  integer :: nred ! Reduction size
29  integer :: next ! Extension size
30  integer,allocatable :: red_to_ext(:) ! Indices conversion
31  integer :: nhalo ! Halo buffer size
32  integer :: nexcl ! Exclusive interior buffer size
33  integer,allocatable :: jhalocounts(:) ! Halo counts
34  integer,allocatable :: jexclcounts(:) ! Exclusive interior counts
35  integer,allocatable :: jhalodispl(:) ! Halo displacement
36  integer,allocatable :: jexcldispl(:) ! Exclusive interior displacement
37  integer,allocatable :: halo(:) ! Halo buffer
38  integer,allocatable :: excl(:) ! Exclusive interior buffer
39 contains
40  procedure :: dealloc => com_dealloc
41  procedure :: com_ext_1d
42  procedure :: com_ext_2d
43  generic :: ext => com_ext_1d,com_ext_2d
44  procedure :: com_red_1d
45  procedure :: com_red_2d
46  generic :: red => com_red_1d,com_red_2d
47  procedure :: read => com_read
48  procedure :: write => com_write
49  procedure :: setup => com_setup
50 end type com_type
51 
52 private
53 public :: com_type
54 
55 contains
56 
57 !----------------------------------------------------------------------
58 ! Subroutine: com_dealloc
59 ! Purpose: communications data deallocation
60 !----------------------------------------------------------------------
61 subroutine com_dealloc(com)
62 
63 implicit none
64 
65 ! Passed variables
66 class(com_type),intent(inout) :: com ! Communication data
67 
68 ! Release memory
69 if (allocated(com%ext_to_proc)) deallocate(com%ext_to_proc)
70 if (allocated(com%ext_to_red)) deallocate(com%ext_to_red)
71 if (allocated(com%red_to_ext)) deallocate(com%red_to_ext)
72 if (allocated(com%jhalocounts)) deallocate(com%jhalocounts)
73 if (allocated(com%jexclcounts)) deallocate(com%jexclcounts)
74 if (allocated(com%jhalodispl)) deallocate(com%jhalodispl)
75 if (allocated(com%jexcldispl)) deallocate(com%jexcldispl)
76 if (allocated(com%halo)) deallocate(com%halo)
77 if (allocated(com%excl)) deallocate(com%excl)
78 
79 end subroutine com_dealloc
80 
81 !----------------------------------------------------------------------
82 ! Subroutine: com_ext_1d
83 ! Purpose: communicate field to halo (extension), 1d
84 !----------------------------------------------------------------------
85 subroutine com_ext_1d(com,mpl,vec_red,vec_ext)
86 
87 implicit none
88 
89 ! Passed variables
90 class(com_type),intent(in) :: com ! Communication data
91 type(mpl_type),intent(in) :: mpl ! MPI data
92 real(kind_real),intent(in) :: vec_red(com%nred) ! Reduced vector
93 real(kind_real),intent(out) :: vec_ext(com%next) ! Extended vector
94 
95 ! Local variables
96 integer :: iexcl,ired,ihalo
97 real(kind_real) :: sbuf(com%nexcl),rbuf(com%nhalo)
98 
99 ! Prepare buffers to send
100 !$omp parallel do schedule(static) private(iexcl)
101 do iexcl=1,com%nexcl
102  sbuf(iexcl) = vec_red(com%excl(iexcl))
103 end do
104 !$omp end parallel do
105 
106 ! Communication
107 call mpl%f_comm%alltoallv(sbuf,com%jexclcounts,com%jexcldispl,rbuf,com%jhalocounts,com%jhalodispl)
108 
109 ! Copy interior
110 !$omp parallel do schedule(static) private(ired)
111 do ired=1,com%nred
112  vec_ext(com%red_to_ext(ired)) = vec_red(ired)
113 end do
114 !$omp end parallel do
115 
116 ! Copy halo
117 !$omp parallel do schedule(static) private(ihalo)
118 do ihalo=1,com%nhalo
119  vec_ext(com%halo(ihalo)) = rbuf(ihalo)
120 end do
121 !$omp end parallel do
122 
123 end subroutine com_ext_1d
124 
125 !----------------------------------------------------------------------
126 ! Subroutine: com_ext_2d
127 ! Purpose: communicate field to halo (extension), 2d
128 !----------------------------------------------------------------------
129 subroutine com_ext_2d(com,mpl,nl,vec_red,vec_ext)
131 implicit none
132 
133 ! Passed variables
134 class(com_type),intent(in) :: com ! Communication data
135 type(mpl_type),intent(in) :: mpl ! MPI data
136 integer,intent(in) :: nl ! Number of levels
137 real(kind_real),intent(in) :: vec_red(com%nred,nl) ! Reduced vector
138 real(kind_real),intent(out) :: vec_ext(com%next,nl) ! Extended vector
139 
140 ! Local variables
141 integer :: il,iexcl,ired,ihalo
142 integer :: jexclcounts(mpl%nproc),jexcldispl(mpl%nproc),jhalocounts(mpl%nproc),jhalodispl(mpl%nproc)
143 real(kind_real) :: sbuf(com%nexcl*nl),rbuf(com%nhalo*nl)
144 
145 ! Prepare buffers to send
146 !$omp parallel do schedule(static) private(il,iexcl)
147 do il=1,nl
148  do iexcl=1,com%nexcl
149  sbuf((iexcl-1)*nl+il) = vec_red(com%excl(iexcl),il)
150  end do
151 end do
152 !$omp end parallel do
153 
154 ! Communication
155 jexclcounts = com%jexclcounts*nl
156 jexcldispl = com%jexcldispl*nl
157 jhalocounts = com%jhalocounts*nl
158 jhalodispl = com%jhalodispl*nl
159 call mpl%f_comm%alltoallv(sbuf,jexclcounts,jexcldispl,rbuf,jhalocounts,jhalodispl)
160 
161 ! Copy interior
162 !$omp parallel do schedule(static) private(il,ired)
163 do il=1,nl
164  do ired=1,com%nred
165  vec_ext(com%red_to_ext(ired),il) = vec_red(ired,il)
166  end do
167 end do
168 !$omp end parallel do
169 
170 ! Copy halo
171 !$omp parallel do schedule(static) private(il,ihalo)
172 do il=1,nl
173  do ihalo=1,com%nhalo
174  vec_ext(com%halo(ihalo),il) = rbuf((ihalo-1)*nl+il)
175  end do
176 end do
177 !$omp end parallel do
178 
179 end subroutine com_ext_2d
180 
181 !----------------------------------------------------------------------
182 ! Subroutine: com_red_1d
183 ! Purpose: communicate vector from halo (reduction)
184 !----------------------------------------------------------------------
185 subroutine com_red_1d(com,mpl,vec_ext,vec_red)
187 implicit none
188 
189 ! Passed variables
190 class(com_type),intent(in) :: com ! Communication data
191 type(mpl_type),intent(in) :: mpl ! MPI data
192 real(kind_real),intent(in) :: vec_ext(com%next) ! Extended vector
193 real(kind_real),intent(out) :: vec_red(com%nred) ! Reduced vector
194 
195 ! Local variables
196 integer :: ihalo,ired,iexcl,ithread
197 real(kind_real) :: sbuf(com%nhalo),rbuf(com%nexcl),vec_red_arr(com%nred,mpl%nthread)
198 
199 ! Prepare buffers to send
200 !$omp parallel do schedule(static) private(ihalo)
201 do ihalo=1,com%nhalo
202  sbuf(ihalo) = vec_ext(com%halo(ihalo))
203 end do
204 !$omp end parallel do
205 
206 ! Communication
207 call mpl%f_comm%alltoallv(sbuf,com%jhalocounts,com%jhalodispl,rbuf,com%jexclcounts,com%jexcldispl)
208 
209 ! Copy interior
210 !$omp parallel do schedule(static) private(ired)
211 do ired=1,com%nred
212  vec_red(ired) = vec_ext(com%red_to_ext(ired))
213 end do
214 !$omp end parallel do
215 
216 ! Copy halo
217 vec_red_arr = 0.0
218 !$omp parallel do schedule(static) private(iexcl,ithread)
219 do iexcl=1,com%nexcl
220  ithread = 1
221 !$ ithread = omp_get_thread_num()+1
222  vec_red_arr(com%excl(iexcl),ithread) = vec_red_arr(com%excl(iexcl),ithread)+rbuf(iexcl)
223 end do
224 !$omp end parallel do
225 
226 ! Sum over threads
227 do ithread=1,mpl%nthread
228  do ired=1,com%nred
229  vec_red(ired) = vec_red(ired)+vec_red_arr(ired,ithread)
230  end do
231 end do
232 
233 end subroutine com_red_1d
234 
235 !----------------------------------------------------------------------
236 ! Subroutine: com_red_2d
237 ! Purpose: communicate vector from halo (reduction)
238 !----------------------------------------------------------------------
239 subroutine com_red_2d(com,mpl,nl,vec_ext,vec_red)
241 implicit none
242 
243 ! Passed variables
244 class(com_type),intent(in) :: com ! Communication data
245 type(mpl_type),intent(in) :: mpl ! MPI data
246 integer,intent(in) :: nl ! Number of levels
247 real(kind_real),intent(in) :: vec_ext(com%next,nl) ! Extended vector
248 real(kind_real),intent(out) :: vec_red(com%nred,nl) ! Reduced vector
249 
250 ! Local variables
251 integer :: il,ihalo,ired,iexcl,ithread
252 real(kind_real) :: sbuf(com%nhalo*nl),rbuf(com%nexcl*nl),vec_red_arr(com%nred,nl,mpl%nthread)
253 
254 ! Prepare buffers to send
255 !$omp parallel do schedule(static) private(il,ihalo)
256 do il=1,nl
257  do ihalo=1,com%nhalo
258  sbuf((ihalo-1)*nl+il) = vec_ext(com%halo(ihalo),il)
259  end do
260 end do
261 !$omp end parallel do
262 
263 ! Communication
264 call mpl%f_comm%alltoallv(sbuf,com%jhalocounts*nl,com%jhalodispl*nl,rbuf,com%jexclcounts*nl,com%jexcldispl*nl)
265 
266 ! Copy interior
267 !$omp parallel do schedule(static) private(il,ired)
268 do il=1,nl
269  do ired=1,com%nred
270  vec_red(ired,il) = vec_ext(com%red_to_ext(ired),il)
271  end do
272 end do
273 !$omp end parallel do
274 
275 ! Copy halo
276 vec_red_arr = 0.0
277 !$omp parallel do schedule(static) private(il,iexcl,ithread)
278 do il=1,nl
279  do iexcl=1,com%nexcl
280  ithread = 1
281 !$ ithread = omp_get_thread_num()+1
282  vec_red_arr(com%excl(iexcl),il,ithread) = vec_red_arr(com%excl(iexcl),il,ithread)+rbuf((iexcl-1)*nl+il)
283  end do
284 end do
285 !$omp end parallel do
286 
287 ! Sum over threads
288 do ithread=1,mpl%nthread
289  do il=1,nl
290  do ired=1,com%nred
291  vec_red(ired,il) = vec_red(ired,il)+vec_red_arr(ired,il,ithread)
292  end do
293  end do
294 end do
295 
296 end subroutine com_red_2d
297 
298 !----------------------------------------------------------------------
299 ! Subroutine: com_read
300 ! Purpose: read communications from a NetCDF file
301 !----------------------------------------------------------------------
302 subroutine com_read(com,mpl,ncid,prefix)
304 implicit none
305 
306 ! Passed variables
307 class(com_type),intent(inout) :: com ! Communication data
308 type(mpl_type),intent(in) :: mpl ! MPI data
309 integer,intent(in) :: ncid ! NetCDF file id
310 character(len=*),intent(in) :: prefix ! Communication prefix
311 
312 ! Local variables
313 integer :: info
314 integer :: nred_id,next_id,red_to_ext_id,nhalo_id,nexcl_id
315 integer :: jhalocounts_id,jexclcounts_id,jhalodispl_id,jexcldispl_id,halo_id,excl_id
316 character(len=1024) :: subr = 'com_read'
317 
318 ! Copy prefix
319 com%prefix = trim(prefix)
320 
321 ! Get dimensions
322 info = nf90_inq_dimid(ncid,trim(prefix)//'_nred',nred_id)
323 if (info==nf90_noerr) then
324  call mpl%ncerr(subr,nf90_inquire_dimension(ncid,nred_id,len=com%nred))
325 else
326  com%nred = 0
327 end if
328 info = nf90_inq_dimid(ncid,trim(prefix)//'_next',next_id)
329 if (info==nf90_noerr) then
330  call mpl%ncerr(subr,nf90_inquire_dimension(ncid,next_id,len=com%next))
331 else
332  com%next = 0
333 end if
334 info = nf90_inq_dimid(ncid,trim(prefix)//'_nhalo',nhalo_id)
335 if (info==nf90_noerr) then
336  call mpl%ncerr(subr,nf90_inquire_dimension(ncid,nhalo_id,len=com%nhalo))
337 else
338  com%nhalo = 0
339 end if
340 info = nf90_inq_dimid(ncid,trim(prefix)//'_nexcl',nexcl_id)
341 if (info==nf90_noerr) then
342  call mpl%ncerr(subr,nf90_inquire_dimension(ncid,nexcl_id,len=com%nexcl))
343 else
344  com%nexcl = 0
345 end if
346 
347 ! Allocation
348 if (com%nred>0) allocate(com%red_to_ext(com%nred))
349 allocate(com%jhalocounts(mpl%nproc))
350 allocate(com%jexclcounts(mpl%nproc))
351 allocate(com%jhalodispl(mpl%nproc))
352 allocate(com%jexcldispl(mpl%nproc))
353 if (com%nhalo>0) allocate(com%halo(com%nhalo))
354 if (com%nexcl>0) allocate(com%excl(com%nexcl))
355 
356 ! Get variables id
357 if (com%nred>0) call mpl%ncerr(subr,nf90_inq_varid(ncid,trim(prefix)//'_red_to_ext',red_to_ext_id))
358 call mpl%ncerr(subr,nf90_inq_varid(ncid,trim(prefix)//'_jhalocounts',jhalocounts_id))
359 call mpl%ncerr(subr,nf90_inq_varid(ncid,trim(prefix)//'_jexclcounts',jexclcounts_id))
360 call mpl%ncerr(subr,nf90_inq_varid(ncid,trim(prefix)//'_jhalodispl',jhalodispl_id))
361 call mpl%ncerr(subr,nf90_inq_varid(ncid,trim(prefix)//'_jexcldispl',jexcldispl_id))
362 if (com%nhalo>0) call mpl%ncerr(subr,nf90_inq_varid(ncid,trim(prefix)//'_halo',halo_id))
363 if (com%nexcl>0) call mpl%ncerr(subr,nf90_inq_varid(ncid,trim(prefix)//'_excl',excl_id))
364 
365 ! Get variable
366 if (com%nred>0) call mpl%ncerr(subr,nf90_get_var(ncid,red_to_ext_id,com%red_to_ext))
367 call mpl%ncerr(subr,nf90_get_var(ncid,jhalocounts_id,com%jhalocounts))
368 call mpl%ncerr(subr,nf90_get_var(ncid,jexclcounts_id,com%jexclcounts))
369 call mpl%ncerr(subr,nf90_get_var(ncid,jhalodispl_id,com%jhalodispl))
370 call mpl%ncerr(subr,nf90_get_var(ncid,jexcldispl_id,com%jexcldispl))
371 if (com%nhalo>0) call mpl%ncerr(subr,nf90_get_var(ncid,halo_id,com%halo))
372 if (com%nexcl>0) call mpl%ncerr(subr,nf90_get_var(ncid,excl_id,com%excl))
373 
374 end subroutine com_read
375 
376 !----------------------------------------------------------------------
377 ! Subroutine: com_write
378 ! Purpose: write communications to a NetCDF file
379 !----------------------------------------------------------------------
380 subroutine com_write(com,mpl,ncid)
382 implicit none
383 
384 ! Passed variables
385 class(com_type),intent(in) :: com ! Communication data
386 type(mpl_type),intent(in) :: mpl ! MPI data
387 integer,intent(in) :: ncid ! NetCDF file id
388 
389 ! Local variables
390 integer :: info
391 integer :: nproc_id,nred_id,next_id,red_to_ext_id,nhalo_id,nexcl_id
392 integer :: jhalocounts_id,jexclcounts_id,jhalodispl_id,jexcldispl_id,halo_id,excl_id
393 character(len=1024) :: subr = 'com_write'
394 
395 ! Start definition mode
396 call mpl%ncerr(subr,nf90_redef(ncid))
397 
398 ! Define dimensions
399 info = nf90_inq_dimid(ncid,'nproc',nproc_id)
400 if (info/=nf90_noerr) call mpl%ncerr(subr,nf90_def_dim(ncid,'nproc',mpl%nproc,nproc_id))
401 if (com%nred>0) call mpl%ncerr(subr,nf90_def_dim(ncid,trim(com%prefix)//'_nred',com%nred,nred_id))
402 if (com%next>0) call mpl%ncerr(subr,nf90_def_dim(ncid,trim(com%prefix)//'_next',com%next,next_id))
403 if (com%nhalo>0) call mpl%ncerr(subr,nf90_def_dim(ncid,trim(com%prefix)//'_nhalo',com%nhalo,nhalo_id))
404 if (com%nexcl>0) call mpl%ncerr(subr,nf90_def_dim(ncid,trim(com%prefix)//'_nexcl',com%nexcl,nexcl_id))
405 
406 ! Define variables
407 if (com%nred>0) call mpl%ncerr(subr,nf90_def_var(ncid,trim(com%prefix)//'_red_to_ext',nf90_int,(/nred_id/),red_to_ext_id))
408 call mpl%ncerr(subr,nf90_def_var(ncid,trim(com%prefix)//'_jhalocounts',nf90_int,(/nproc_id/),jhalocounts_id))
409 call mpl%ncerr(subr,nf90_def_var(ncid,trim(com%prefix)//'_jexclcounts',nf90_int,(/nproc_id/),jexclcounts_id))
410 call mpl%ncerr(subr,nf90_def_var(ncid,trim(com%prefix)//'_jhalodispl',nf90_int,(/nproc_id/),jhalodispl_id))
411 call mpl%ncerr(subr,nf90_def_var(ncid,trim(com%prefix)//'_jexcldispl',nf90_int,(/nproc_id/),jexcldispl_id))
412 if (com%nhalo>0) call mpl%ncerr(subr,nf90_def_var(ncid,trim(com%prefix)//'_halo',nf90_int,(/nhalo_id/),halo_id))
413 if (com%nexcl>0) call mpl%ncerr(subr,nf90_def_var(ncid,trim(com%prefix)//'_excl',nf90_int,(/nexcl_id/),excl_id))
414 
415 ! End definition mode
416 call mpl%ncerr(subr,nf90_enddef(ncid))
417 
418 ! Put variables
419 if (com%nred>0) call mpl%ncerr(subr,nf90_put_var(ncid,red_to_ext_id,com%red_to_ext))
420 call mpl%ncerr(subr,nf90_put_var(ncid,jhalocounts_id,com%jhalocounts))
421 call mpl%ncerr(subr,nf90_put_var(ncid,jexclcounts_id,com%jexclcounts))
422 call mpl%ncerr(subr,nf90_put_var(ncid,jhalodispl_id,com%jhalodispl))
423 call mpl%ncerr(subr,nf90_put_var(ncid,jexcldispl_id,com%jexcldispl))
424 if (com%nhalo>0) call mpl%ncerr(subr,nf90_put_var(ncid,halo_id,com%halo))
425 if (com%nexcl>0) call mpl%ncerr(subr,nf90_put_var(ncid,excl_id,com%excl))
426 
427 end subroutine com_write
428 
429 !----------------------------------------------------------------------
430 ! Subroutine: com_setup
431 ! Purpose: setup communications
432 !----------------------------------------------------------------------
433 subroutine com_setup(com_out,mpl,prefix,nglb,nred,next,ext_to_glb,red_to_ext,glb_to_proc,glb_to_red)
435 implicit none
436 
437 ! Passed variables
438 class(com_type),intent(inout) :: com_out ! Communication data
439 type(mpl_type),intent(inout) :: mpl ! MPI data
440 character(len=*),intent(in) :: prefix ! Prefix
441 integer,intent(in) :: nglb ! Global size
442 integer,intent(in) :: nred ! Reduced halo size
443 integer,intent(in) :: next ! Extended halo size
444 integer,intent(in) :: ext_to_glb(next) ! Extended halo to global
445 integer,intent(in) :: red_to_ext(nred) ! Reduced halo to extended halo
446 integer,intent(in) :: glb_to_proc(nglb) ! Global to processor
447 integer,intent(in) :: glb_to_red(nglb) ! Global to reduced halo
448 
449 ! Local variables
450 integer :: iproc,jproc,iext,iglb,ired,icount,nred_tmp,next_tmp
451 integer,allocatable :: ext_to_glb_tmp(:),red_to_ext_tmp(:)
452 type(com_type) :: com_in(mpl%nproc)
453 type(fckit_mpi_status) :: status
454 
455 if (mpl%main) then
456  do iproc=1,mpl%nproc
457  ! Communicate dimensions
458  if (iproc==mpl%ioproc) then
459  ! Copy dimensions
460  nred_tmp = nred
461  next_tmp = next
462  else
463  ! Receive dimensions on ioproc
464  call mpl%f_comm%receive(nred_tmp,iproc-1,mpl%tag,status)
465  call mpl%f_comm%receive(next_tmp,iproc-1,mpl%tag+1,status)
466  end if
467 
468  ! Allocation
469  allocate(ext_to_glb_tmp(next_tmp))
470  allocate(red_to_ext_tmp(nred_tmp))
471 
472  ! Communicate data
473  if (iproc==mpl%ioproc) then
474  ! Copy data
475  ext_to_glb_tmp = ext_to_glb
476  red_to_ext_tmp = red_to_ext
477  else
478  ! Receive data on ioproc
479  call mpl%f_comm%receive(ext_to_glb_tmp,iproc-1,mpl%tag+2,status)
480  call mpl%f_comm%receive(red_to_ext_tmp,iproc-1,mpl%tag+3,status)
481  end if
482 
483  ! Allocation
484  com_in(iproc)%nred = nred_tmp
485  com_in(iproc)%next = next_tmp
486  allocate(com_in(iproc)%ext_to_proc(com_in(iproc)%next))
487  allocate(com_in(iproc)%ext_to_red(com_in(iproc)%next))
488  allocate(com_in(iproc)%red_to_ext(com_in(iproc)%nred))
489 
490  ! Communication parameters
491  do iext=1,next_tmp
492  iglb = ext_to_glb_tmp(iext)
493  com_in(iproc)%ext_to_proc(iext) = glb_to_proc(iglb)
494  ired = glb_to_red(iglb)
495  com_in(iproc)%ext_to_red(iext) = ired
496  end do
497  com_in(iproc)%red_to_ext = red_to_ext_tmp
498 
499  ! Release memory
500  deallocate(ext_to_glb_tmp)
501  deallocate(red_to_ext_tmp)
502  end do
503 else
504  ! Send dimensions to ioproc
505  call mpl%f_comm%send(nred,mpl%ioproc-1,mpl%tag)
506  call mpl%f_comm%send(next,mpl%ioproc-1,mpl%tag+1)
507 
508  ! Send data to ioproc
509  call mpl%f_comm%send(ext_to_glb,mpl%ioproc-1,mpl%tag+2)
510  call mpl%f_comm%send(red_to_ext,mpl%ioproc-1,mpl%tag+3)
511 end if
512 call mpl%update_tag(4)
513 
514 if (mpl%main) then
515  ! Allocation
516  do iproc=1,mpl%nproc
517  allocate(com_in(iproc)%jhalocounts(mpl%nproc))
518  allocate(com_in(iproc)%jexclcounts(mpl%nproc))
519  allocate(com_in(iproc)%jhalodispl(mpl%nproc))
520  allocate(com_in(iproc)%jexcldispl(mpl%nproc))
521  end do
522 
523  ! Initialization
524  do iproc=1,mpl%nproc
525  com_in(iproc)%jhalocounts = 0
526  com_in(iproc)%jexclcounts = 0
527  com_in(iproc)%jhalodispl = 0
528  com_in(iproc)%jexcldispl = 0
529  end do
530 
531  ! Compute counts
532  do iproc=1,mpl%nproc
533  do iext=1,com_in(iproc)%next
534  jproc = com_in(iproc)%ext_to_proc(iext)
535  if (jproc/=iproc) then
536  ! Count of points sent from IPROC to JPROC
537  com_in(iproc)%jhalocounts(jproc) = com_in(iproc)%jhalocounts(jproc)+1
538 
539  ! Count of points received on JPROC from IPROC
540  com_in(jproc)%jexclcounts(iproc) = com_in(jproc)%jexclcounts(iproc)+1
541  end if
542  end do
543  end do
544 
545  ! Compute displacement
546  do iproc=1,mpl%nproc
547  com_in(iproc)%jhalodispl(1) = 0
548  com_in(iproc)%jexcldispl(1) = 0
549  do jproc=2,mpl%nproc
550  com_in(iproc)%jhalodispl(jproc) = com_in(iproc)%jhalodispl(jproc-1)+com_in(iproc)%jhalocounts(jproc-1)
551  com_in(iproc)%jexcldispl(jproc) = com_in(iproc)%jexcldispl(jproc-1)+com_in(iproc)%jexclcounts(jproc-1)
552  end do
553  end do
554 
555  ! Allocation
556  do iproc=1,mpl%nproc
557  com_in(iproc)%nhalo = sum(com_in(iproc)%jhalocounts)
558  com_in(iproc)%nexcl = sum(com_in(iproc)%jexclcounts)
559  allocate(com_in(iproc)%halo(com_in(iproc)%nhalo))
560  allocate(com_in(iproc)%excl(com_in(iproc)%nexcl))
561  end do
562 
563  ! Fill halo array
564  do iproc=1,mpl%nproc
565  com_in(iproc)%jhalocounts = 0
566  end do
567  do iproc=1,mpl%nproc
568  do iext=1,com_in(iproc)%next
569  ! Check for halo points
570  jproc = com_in(iproc)%ext_to_proc(iext)
571  if (jproc/=iproc) then
572  ! Count of points sent from IPROC to JPROC
573  com_in(iproc)%jhalocounts(jproc) = com_in(iproc)%jhalocounts(jproc)+1
574 
575  ! Local index of points sent from IPROC to JPROC
576  com_in(iproc)%halo(com_in(iproc)%jhalodispl(jproc)+com_in(iproc)%jhalocounts(jproc)) = iext
577  end if
578  end do
579  end do
580 
581  ! Fill excl array
582  do jproc=1,mpl%nproc
583  ! Loop over processors sending data to JPROC
584  do iproc=1,mpl%nproc
585  do icount=1,com_in(iproc)%jhalocounts(jproc)
586  ! Local index of points received on JPROC from IPROC
587  com_in(jproc)%excl(com_in(jproc)%jexcldispl(iproc)+icount) = &
588  & com_in(iproc)%ext_to_red(com_in(iproc)%halo(com_in(iproc)%jhalodispl(jproc)+icount))
589  end do
590  end do
591  end do
592 
593  ! Communicate dimensions
594  do iproc=1,mpl%nproc
595  if (iproc==mpl%ioproc) then
596  ! Copy dimensions
597  com_out%nred = com_in(iproc)%nred
598  com_out%next = com_in(iproc)%next
599  com_out%nhalo = com_in(iproc)%nhalo
600  com_out%nexcl = com_in(iproc)%nexcl
601  else
602  ! Send dimensions to iproc
603  call mpl%f_comm%send(com_in(iproc)%nred,iproc-1,mpl%tag)
604  call mpl%f_comm%send(com_in(iproc)%next,iproc-1,mpl%tag+1)
605  call mpl%f_comm%send(com_in(iproc)%nhalo,iproc-1,mpl%tag+2)
606  call mpl%f_comm%send(com_in(iproc)%nexcl,iproc-1,mpl%tag+3)
607  end if
608  end do
609 else
610  ! Receive dimensions from ioproc
611  call mpl%f_comm%receive(com_out%nred,mpl%ioproc-1,mpl%tag,status)
612  call mpl%f_comm%receive(com_out%next,mpl%ioproc-1,mpl%tag+1,status)
613  call mpl%f_comm%receive(com_out%nhalo,mpl%ioproc-1,mpl%tag+2,status)
614  call mpl%f_comm%receive(com_out%nexcl,mpl%ioproc-1,mpl%tag+3,status)
615 end if
616 call mpl%update_tag(4)
617 
618 ! Allocation
619 allocate(com_out%red_to_ext(com_out%nred))
620 allocate(com_out%jhalocounts(mpl%nproc))
621 allocate(com_out%jexclcounts(mpl%nproc))
622 allocate(com_out%jhalodispl(mpl%nproc))
623 allocate(com_out%jexcldispl(mpl%nproc))
624 allocate(com_out%halo(com_out%nhalo))
625 allocate(com_out%excl(com_out%nexcl))
626 
627 ! Communicate data
628 if (mpl%main) then
629  do iproc=1,mpl%nproc
630  if (iproc==mpl%ioproc) then
631  ! Copy dimensions
632  com_out%red_to_ext = com_in(iproc)%red_to_ext
633  com_out%jhalocounts = com_in(iproc)%jhalocounts
634  com_out%jexclcounts = com_in(iproc)%jexclcounts
635  com_out%jhalodispl = com_in(iproc)%jhalodispl
636  com_out%jexcldispl = com_in(iproc)%jexcldispl
637  com_out%halo = com_in(iproc)%halo
638  com_out%excl = com_in(iproc)%excl
639  else
640  ! Send dimensions to iproc
641  call mpl%f_comm%send(com_in(iproc)%red_to_ext,iproc-1,mpl%tag)
642  call mpl%f_comm%send(com_in(iproc)%jhalocounts,iproc-1,mpl%tag+1)
643  call mpl%f_comm%send(com_in(iproc)%jexclcounts,iproc-1,mpl%tag+2)
644  call mpl%f_comm%send(com_in(iproc)%jhalodispl,iproc-1,mpl%tag+3)
645  call mpl%f_comm%send(com_in(iproc)%jexcldispl,iproc-1,mpl%tag+4)
646  if (com_in(iproc)%nhalo>0) call mpl%f_comm%send(com_in(iproc)%halo,iproc-1,mpl%tag+5)
647  if (com_in(iproc)%nexcl>0) call mpl%f_comm%send(com_in(iproc)%excl,iproc-1,mpl%tag+6)
648  end if
649  end do
650 else
651  ! Receive dimensions from ioproc
652  call mpl%f_comm%receive(com_out%red_to_ext,mpl%ioproc-1,mpl%tag,status)
653  call mpl%f_comm%receive(com_out%jhalocounts,mpl%ioproc-1,mpl%tag+1,status)
654  call mpl%f_comm%receive(com_out%jexclcounts,mpl%ioproc-1,mpl%tag+2,status)
655  call mpl%f_comm%receive(com_out%jhalodispl,mpl%ioproc-1,mpl%tag+3,status)
656  call mpl%f_comm%receive(com_out%jexcldispl,mpl%ioproc-1,mpl%tag+4,status)
657  if (com_out%nhalo>0) call mpl%f_comm%receive(com_out%halo,mpl%ioproc-1,mpl%tag+5,status)
658  if (com_out%nexcl>0) call mpl%f_comm%receive(com_out%excl,mpl%ioproc-1,mpl%tag+6,status)
659 end if
660 call mpl%update_tag(7)
661 
662 ! Set prefix
663 com_out%prefix = trim(prefix)
664 
665 end subroutine com_setup
666 
667 end module type_com
subroutine com_read(com, mpl, ncid, prefix)
Definition: type_com.F90:303
subroutine com_red_1d(com, mpl, vec_ext, vec_red)
Definition: type_com.F90:186
subroutine com_dealloc(com)
Definition: type_com.F90:62
subroutine com_ext_2d(com, mpl, nl, vec_red, vec_ext)
Definition: type_com.F90:130
subroutine com_red_2d(com, mpl, nl, vec_ext, vec_red)
Definition: type_com.F90:240
subroutine com_setup(com_out, mpl, prefix, nglb, nred, next, ext_to_glb, red_to_ext, glb_to_proc, glb_to_red)
Definition: type_com.F90:434
subroutine com_write(com, mpl, ncid)
Definition: type_com.F90:381
integer, parameter, public kind_real
integer, parameter, public ncfloat
Definition: tools_nc.F90:27
subroutine com_ext_1d(com, mpl, vec_red, vec_ext)
Definition: type_com.F90:86