16 use fckit_mpi_module,
only: fckit_mpi_status
23 integer,
allocatable :: ext_to_proc(:)
24 integer,
allocatable :: ext_to_red(:)
27 character(len=1024) :: prefix
30 integer,
allocatable :: red_to_ext(:)
33 integer,
allocatable :: jhalocounts(:)
34 integer,
allocatable :: jexclcounts(:)
35 integer,
allocatable :: jhalodispl(:)
36 integer,
allocatable :: jexcldispl(:)
37 integer,
allocatable :: halo(:)
38 integer,
allocatable :: excl(:)
66 class(com_type),
intent(inout) :: com
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)
90 class(com_type),
intent(in) :: com
91 type(mpl_type),
intent(in) :: mpl
92 real(kind_real),
intent(in) :: vec_red(com%nred)
93 real(kind_real),
intent(out) :: vec_ext(com%next)
96 integer :: iexcl,ired,ihalo
97 real(kind_real) :: sbuf(com%nexcl),rbuf(com%nhalo)
102 sbuf(iexcl) = vec_red(com%excl(iexcl))
107 call mpl%f_comm%alltoallv(sbuf,com%jexclcounts,com%jexcldispl,rbuf,com%jhalocounts,com%jhalodispl)
112 vec_ext(com%red_to_ext(ired)) = vec_red(ired)
119 vec_ext(com%halo(ihalo)) = rbuf(ihalo)
129 subroutine com_ext_2d(com,mpl,nl,vec_red,vec_ext)
134 class(com_type),
intent(in) :: com
135 type(mpl_type),
intent(in) :: mpl
136 integer,
intent(in) :: nl
137 real(kind_real),
intent(in) :: vec_red(com%nred,nl)
138 real(kind_real),
intent(out) :: vec_ext(com%next,nl)
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)
149 sbuf((iexcl-1)*nl+il) = vec_red(com%excl(iexcl),il)
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)
165 vec_ext(com%red_to_ext(ired),il) = vec_red(ired,il)
174 vec_ext(com%halo(ihalo),il) = rbuf((ihalo-1)*nl+il)
185 subroutine com_red_1d(com,mpl,vec_ext,vec_red)
190 class(com_type),
intent(in) :: com
191 type(mpl_type),
intent(in) :: mpl
192 real(kind_real),
intent(in) :: vec_ext(com%next)
193 real(kind_real),
intent(out) :: vec_red(com%nred)
196 integer :: ihalo,ired,iexcl,ithread
197 real(kind_real) :: sbuf(com%nhalo),rbuf(com%nexcl),vec_red_arr(com%nred,mpl%nthread)
202 sbuf(ihalo) = vec_ext(com%halo(ihalo))
207 call mpl%f_comm%alltoallv(sbuf,com%jhalocounts,com%jhalodispl,rbuf,com%jexclcounts,com%jexcldispl)
212 vec_red(ired) = vec_ext(com%red_to_ext(ired))
222 vec_red_arr(com%excl(iexcl),ithread) = vec_red_arr(com%excl(iexcl),ithread)+rbuf(iexcl)
227 do ithread=1,mpl%nthread
229 vec_red(ired) = vec_red(ired)+vec_red_arr(ired,ithread)
239 subroutine com_red_2d(com,mpl,nl,vec_ext,vec_red)
244 class(com_type),
intent(in) :: com
245 type(mpl_type),
intent(in) :: mpl
246 integer,
intent(in) :: nl
247 real(kind_real),
intent(in) :: vec_ext(com%next,nl)
248 real(kind_real),
intent(out) :: vec_red(com%nred,nl)
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)
258 sbuf((ihalo-1)*nl+il) = vec_ext(com%halo(ihalo),il)
264 call mpl%f_comm%alltoallv(sbuf,com%jhalocounts*nl,com%jhalodispl*nl,rbuf,com%jexclcounts*nl,com%jexcldispl*nl)
270 vec_red(ired,il) = vec_ext(com%red_to_ext(ired),il)
282 vec_red_arr(com%excl(iexcl),il,ithread) = vec_red_arr(com%excl(iexcl),il,ithread)+rbuf((iexcl-1)*nl+il)
288 do ithread=1,mpl%nthread
291 vec_red(ired,il) = vec_red(ired,il)+vec_red_arr(ired,il,ithread)
302 subroutine com_read(com,mpl,ncid,prefix)
307 class(com_type),
intent(inout) :: com
308 type(mpl_type),
intent(in) :: mpl
309 integer,
intent(in) :: ncid
310 character(len=*),
intent(in) :: prefix
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' 319 com%prefix = trim(prefix)
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))
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))
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))
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))
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))
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))
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))
385 class(com_type),
intent(in) :: com
386 type(mpl_type),
intent(in) :: mpl
387 integer,
intent(in) :: ncid
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' 396 call mpl%ncerr(subr,nf90_redef(ncid))
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))
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))
416 call mpl%ncerr(subr,nf90_enddef(ncid))
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))
433 subroutine com_setup(com_out,mpl,prefix,nglb,nred,next,ext_to_glb,red_to_ext,glb_to_proc,glb_to_red)
438 class(com_type),
intent(inout) :: com_out
439 type(mpl_type),
intent(inout) :: mpl
440 character(len=*),
intent(in) :: prefix
441 integer,
intent(in) :: nglb
442 integer,
intent(in) :: nred
443 integer,
intent(in) :: next
444 integer,
intent(in) :: ext_to_glb(next)
445 integer,
intent(in) :: red_to_ext(nred)
446 integer,
intent(in) :: glb_to_proc(nglb)
447 integer,
intent(in) :: glb_to_red(nglb)
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
458 if (iproc==mpl%ioproc)
then 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)
469 allocate(ext_to_glb_tmp(next_tmp))
470 allocate(red_to_ext_tmp(nred_tmp))
473 if (iproc==mpl%ioproc)
then 475 ext_to_glb_tmp = ext_to_glb
476 red_to_ext_tmp = red_to_ext
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)
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))
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
497 com_in(iproc)%red_to_ext = red_to_ext_tmp
500 deallocate(ext_to_glb_tmp)
501 deallocate(red_to_ext_tmp)
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)
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)
512 call mpl%update_tag(4)
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))
525 com_in(iproc)%jhalocounts = 0
526 com_in(iproc)%jexclcounts = 0
527 com_in(iproc)%jhalodispl = 0
528 com_in(iproc)%jexcldispl = 0
533 do iext=1,com_in(iproc)%next
534 jproc = com_in(iproc)%ext_to_proc(iext)
535 if (jproc/=iproc)
then 537 com_in(iproc)%jhalocounts(jproc) = com_in(iproc)%jhalocounts(jproc)+1
540 com_in(jproc)%jexclcounts(iproc) = com_in(jproc)%jexclcounts(iproc)+1
547 com_in(iproc)%jhalodispl(1) = 0
548 com_in(iproc)%jexcldispl(1) = 0
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)
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))
565 com_in(iproc)%jhalocounts = 0
568 do iext=1,com_in(iproc)%next
570 jproc = com_in(iproc)%ext_to_proc(iext)
571 if (jproc/=iproc)
then 573 com_in(iproc)%jhalocounts(jproc) = com_in(iproc)%jhalocounts(jproc)+1
576 com_in(iproc)%halo(com_in(iproc)%jhalodispl(jproc)+com_in(iproc)%jhalocounts(jproc)) = iext
585 do icount=1,com_in(iproc)%jhalocounts(jproc)
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))
595 if (iproc==mpl%ioproc)
then 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
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)
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)
616 call mpl%update_tag(4)
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))
630 if (iproc==mpl%ioproc)
then 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
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)
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)
660 call mpl%update_tag(7)
663 com_out%prefix = trim(prefix)
subroutine com_read(com, mpl, ncid, prefix)
subroutine com_red_1d(com, mpl, vec_ext, vec_red)
subroutine com_dealloc(com)
subroutine com_ext_2d(com, mpl, nl, vec_red, vec_ext)
subroutine com_red_2d(com, mpl, nl, vec_ext, vec_red)
subroutine com_setup(com_out, mpl, prefix, nglb, nred, next, ext_to_glb, red_to_ext, glb_to_proc, glb_to_red)
subroutine com_write(com, mpl, ncid)
integer, parameter, public kind_real
subroutine com_ext_1d(com, mpl, vec_red, vec_ext)