Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_collective.F90
Go to the documentation of this file.
1 !**************************************************************
2 ! Copyright INRIA
3 ! Authors :
4 ! CALVI project team
5 !
6 ! This code SeLaLib (for Semi-Lagrangian-Library)
7 ! is a parallel library for simulating the plasma turbulence
8 ! in a tokamak.
9 !
10 ! This software is governed by the CeCILL-B license
11 ! under French law and abiding by the rules of distribution
12 ! of free software. You can use, modify and redistribute
13 ! the software under the terms of the CeCILL-B license as
14 ! circulated by CEA, CNRS and INRIA at the following URL
15 ! "http://www.cecill.info".
16 !**************************************************************
17 
18 !------------------------------------------------------------------------------
19 ! Selalib
20 !------------------------------------------------------------------------------
21 ! MODULE: sll_m_collective
22 !
23 ! DESCRIPTION:
107 ! REVISION HISTORY:
108 ! DD Mmm YYYY - Initial Version
109 ! TODO_dd_mmm_yyyy - TODO_describe_appropriate_changes - TODO_name
110 !------------------------------------------------------------------------------
111 
112 ! ***************************************************************************
113 ! sll_m_collective is a module that encapsulates our calls to the MPI library.
114 ! Any module that wants to access distributed multiprocessing capabilities
115 ! must use sll_m_collective instead of MPI directly. If one were to desire
116 ! some MPI capability that is not included here, such capability should not
117 ! be plugged in directly but rather this module should be expanded.
118 !
119 ! By centralizing the interaction with the MPI library we achieve several
120 ! things:
121 ! - obtain an explicitely reduced 'parallel' vocabulary to build our
122 ! application,
123 ! - centralize argument checking and some error handling, and
124 ! - adjust the interface to our wishes.
125 !
126 ! ***************************************************************************
127 
129 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
130 #include "sll_assert.h"
131 #include "sll_memory.h"
132 #include "sll_working_precision.h"
133 
134  use mpi, only: &
135  mpi_allgather, &
136  mpi_allgatherv, &
137  mpi_allreduce, &
138  mpi_alltoall, &
139  mpi_alltoallv, &
140  mpi_barrier, &
141  mpi_bcast, &
142  mpi_comm_rank, &
143  mpi_comm_size, &
144  mpi_comm_split, &
145  mpi_comm_world, &
146  mpi_complex, &
147  mpi_double_complex, &
148  mpi_double_precision, &
149  mpi_finalize, &
150  mpi_gather, &
151  mpi_gatherv, &
152  mpi_init_thread, &
153  mpi_integer, &
154  mpi_integer8, &
155  mpi_logical, &
156  mpi_real, &
157  mpi_real8, &
158  mpi_reduce, &
159  mpi_scatter, &
160  mpi_scatterv, &
161  mpi_success, &
162  mpi_sum, &
163  mpi_thread_funneled, &
164  mpi_byte
165 
166  implicit none
167 
168  public :: &
212 
213  private
214 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
215  ! This is the only place in the prototype that should have to include
216  ! the mpi header file.
217  !include 'mpif.h'
218  !implicit none
219 
220  ! **********************************************************************
221  ! The sll_t_collective_t is essentially a wrapper around an MPI
222  ! communicator. We store some of the most frequently accessed values
223  ! to save a little on the overhead of the corresponding function calls.
224  !
225  ! The original idea was to have this as an opaque type. But how to
226  ! do this in Fortran without forward declarations? Any help would be
227  ! appreciated... In any case, the only means allowed to interact with
228  ! this type, outside of this module, are the functions/subroutines
229  ! defined herein.
230  !
231  !***********************************************************************
232 
235  type(sll_t_collective_t), pointer :: parent => null()
236  sll_int32 :: comm
239  sll_int32 :: color
240  sll_int32 :: key
241  sll_int32 :: rank
242  sll_int32 :: size
243  sll_int32 :: thread_level_required
244  sll_int32 :: thread_level_provided
245  end type sll_t_collective_t
246 
247  ! **********************************************************************
248  ! A few rare global variables.
249  ! The sll_v_world_collective just wraps around MPI_WORLD_COMM. This gets
250  ! initialized in the call to sll_s_boot_collective().
251  !
252  ! **********************************************************************
253 
256 
260  module procedure sll_collective_bcast_int32
261  module procedure sll_collective_bcast_int64
262  module procedure sll_collective_bcast_real32
263  module procedure sll_collective_bcast_real64
264  module procedure sll_s_collective_bcast_real64 ! TODO: delete (use subr above)
265  module procedure sll_collective_bcast_comp32 ! TODO: remove 'size' argument
266  module procedure sll_collective_bcast_comp64 ! TODO: remove 'size' argument
267  module procedure sll_collective_bcast_logical
268  end interface
269 
272  module procedure sll_collective_gather_real32
273  module procedure sll_collective_gather_real64
274  module procedure sll_collective_gather_logical
275  end interface
276 
280  module procedure sll_collective_allgather_int, &
282  end interface
283 
287  module procedure sll_collective_allgatherv_real32, &
289  end interface
290 
293  module procedure sll_s_collective_gatherv_real, &
295  end interface
296 
300  module procedure sll_collective_scatter_real
301  module procedure sll_collective_scatter_real64
302  module procedure sll_collective_scatter_int32
303  end interface
304 
307  module procedure sll_s_collective_scatterv_real
308  end interface
309 
313  module procedure sll_s_collective_allreduce_real32, &
320  end interface
321 
324  module procedure sll_s_collective_reduce_real32, &
331  end interface
332 
335  module procedure sll_s_collective_alltoall_int, &
338  end interface
339 
344  module procedure sll_s_collective_alltoallv_int, &
348  end interface
349 
352  module procedure sll_s_collective_globalsum_array_real64, &
362  end interface
363 
364 contains !************************** Operations **************************
365 
366  ! First a couple of utilities for this module
369  subroutine sll_check_collective_ptr(ptr)
370  type(sll_t_collective_t), pointer :: ptr
371  if (.not. associated(ptr)) then
372  write (*, '(a)') 'sll_check_collective_ptr: non-associated pointer'
373  stop 'SLL_CHECK_COLLECTIVE_PTR'
374  end if
375  end subroutine sll_check_collective_ptr
376 
384  subroutine sll_s_test_mpi_error(ierr, descriptor)
385  sll_int32, intent(in) :: ierr
386  character(len=*), intent(in) :: descriptor
387 
388  if (ierr .ne. mpi_success) then
389  write (*, '(a, a)') 'MPI error code failure: ', descriptor
390  end if
391  end subroutine sll_s_test_mpi_error
392 
393  ! Following what is exposed to the users of the module.
394 
395  ! The following function is a little tricky. The only test of equality we
396  ! use is the identifier of the underlying communicator. Can the collective
397  ! differ in any other way?
398  function sll_f_collectives_are_same(col1, col2)
399  logical :: sll_f_collectives_are_same
400  type(sll_t_collective_t), pointer :: col1
401  type(sll_t_collective_t), pointer :: col2
402  sll_assert(associated(col1))
403  sll_assert(associated(col2))
404  if (col1%comm == col2%comm) then
406  else
408  end if
409  end function sll_f_collectives_are_same
410 
411  ! sll_s_boot_collective allocates and initializes the global variable
412  ! sll_v_world_collective and boots the MPI environment.
413 
415  subroutine sll_s_boot_collective(required)
416  sll_int32, intent(in), optional :: required
417  sll_int32 :: ierr
418 
420 
421  if (present(required)) then
422  sll_v_world_collective%thread_level_required = required
423  else
424  sll_v_world_collective%thread_level_required = mpi_thread_funneled
425  end if
426 
427  call mpi_init_thread(sll_v_world_collective%thread_level_required, &
428  sll_v_world_collective%thread_level_provided, &
429  ierr)
430  call sll_s_set_communicator_collective(mpi_comm_world)
431  end subroutine sll_s_boot_collective
432 
433  subroutine sll_s_allocate_collective() bind(C,name='sll_s_allocate_collective')
434  sll_int32 :: ierr
435  sll_allocate( sll_v_world_collective, ierr )
436  end subroutine sll_s_allocate_collective
437 
438  subroutine sll_s_set_communicator_collective(communicator_in) bind(C,name='sll_s_set_communicator_collective')
439  sll_int32, intent(in) :: communicator_in
440  sll_int32 :: ierr
441 
442  sll_v_world_collective%comm = communicator_in
443  sll_v_world_collective%color = 0
444  sll_v_world_collective%key = 0
445  call mpi_comm_rank( communicator_in, sll_v_world_collective%rank, ierr )
446  call sll_s_test_mpi_error( ierr, 'sll_s_boot_collective(): MPI_COMM_RANK()' )
447  call mpi_comm_size( communicator_in, sll_v_world_collective%size, ierr )
448  call sll_s_test_mpi_error( ierr, 'sll_s_boot_collective(): MPI_COMM_SIZE()')
449  end subroutine sll_s_set_communicator_collective
450 
452  subroutine sll_s_halt_collective() bind(C,name='sll_s_halt_collective')
453  sll_int32 :: ierr
454  call mpi_barrier(mpi_comm_world, ierr)
455  call sll_s_test_mpi_error(ierr, &
456  'sll_s_halt_collective(): MPI_BARRIER()')
457  call mpi_finalize(ierr)
458  end subroutine sll_s_halt_collective
459 
460  ! sll_f_new_collective follows a somewhat similar calling convention as
461  ! MPI_COMM_SPLIT() and is a wrapper around it. Error checking
462  ! is done internally. The departure from the standard syntax is to permit
463  ! a call like:
464  !
465  ! type(sll_t_collective_t), pointer :: new_col
466  !
467  ! followed by
468  !
469  ! new_col = sll_f_new_collective( parent, color, key )
470  !
471  ! Note that the library's basic types are always pointers. The alternative
472  ! syntax would be:
473  !
474  ! type(sll_t_collective_t), pointer :: new_col
475  !
476  ! call sll_f_new_collective( parent, color, key, new_col )
477  !
478  ! The collective stores the information:
479  ! - a pointer to the parent collective from which it was split
480  ! - the values of color and key used for the splitting
481  ! - the other fields, rank and size are pre-populated for convenience.
482 
488  function sll_f_new_collective(parent, color, key)
489  type(sll_t_collective_t), pointer :: parent
490  sll_int32, intent(in) :: color
491  sll_int32, intent(in) :: key
492  sll_int32 :: ierr
493  type(sll_t_collective_t), pointer :: sll_f_new_collective
494 
495  call sll_check_collective_ptr(parent)
496  sll_allocate(sll_f_new_collective, ierr)
497  sll_f_new_collective%parent => parent
498  sll_f_new_collective%color = color
499  sll_f_new_collective%key = key
500  call mpi_comm_split(parent%comm, color, key, sll_f_new_collective%comm, &
501  ierr)
502  call sll_s_test_mpi_error(ierr, 'sll_f_new_collective(): MPI_COMM_SPLIT()')
503  ! fill out the rest of the fields.
504  call mpi_comm_rank(sll_f_new_collective%comm, sll_f_new_collective%rank, ierr)
505  call sll_s_test_mpi_error(ierr, 'sll_f_new_collective(): MPI_COMM_RANK()')
506  call mpi_comm_size(sll_f_new_collective%comm, sll_f_new_collective%size, ierr)
507  call sll_s_test_mpi_error(ierr, 'sll_f_new_collective(): MPI_COMM_RANK()')
508  end function sll_f_new_collective
509 
512  subroutine sll_delete_collective(col)
513  type(sll_t_collective_t), pointer :: col
514  sll_int32 :: ierr
515  ! Why don't use MPI_COMM_FREE(col%comm,ierr)
516  call sll_check_collective_ptr(col)
517  sll_deallocate(col, ierr)
518  end subroutine sll_delete_collective
519 
523  type(sll_t_collective_t), pointer :: col
524  sll_int32 :: sll_get_collective_comm
525  call sll_check_collective_ptr(col)
526  sll_get_collective_comm = col%comm
527  end function sll_get_collective_comm
528 
533  type(sll_t_collective_t), pointer :: col
534  sll_int32 :: sll_f_get_collective_rank
535  call sll_check_collective_ptr(col)
536  sll_f_get_collective_rank = col%rank
537  end function sll_f_get_collective_rank
538 
543  type(sll_t_collective_t), pointer :: col
544  sll_int32 :: sll_get_collective_color
545  call sll_check_collective_ptr(col)
546  sll_get_collective_color = col%color
547  end function sll_get_collective_color
548 
553  type(sll_t_collective_t), pointer :: col
554  sll_int32 :: sll_f_get_collective_size
555  call sll_check_collective_ptr(col)
556  sll_f_get_collective_size = col%size
557  end function sll_f_get_collective_size
558 
563  type(sll_t_collective_t), pointer :: col
565  call sll_check_collective_ptr(col)
566  sll_get_collective_parent => col%parent
567  end function sll_get_collective_parent
568 
572  subroutine sll_s_collective_barrier(col)
573  type(sll_t_collective_t), pointer :: col
574  sll_int32 :: ierr
575  call sll_check_collective_ptr(col)
576  call mpi_barrier(col%comm, ierr)
577  call sll_s_test_mpi_error(ierr, &
578  'sll_s_collective_barrier(): MPI_BARRIER()')
579  end subroutine sll_s_collective_barrier
580 
586  subroutine sll_collective_bcast_logical(col, buffer, root)
587  type(sll_t_collective_t), pointer :: col
588  logical, intent(inout) :: buffer(:)
589  sll_int32, intent(in) :: root
590  sll_int32 :: ierr
591  call mpi_bcast(buffer, size(buffer), mpi_logical, root, col%comm, ierr)
592  call sll_s_test_mpi_error(ierr, 'sll_collective_bcast_logical(): MPI_BCAST()')
593  end subroutine sll_collective_bcast_logical
594 
600  subroutine sll_collective_bcast_int32(col, buffer, root)
601  type(sll_t_collective_t), pointer :: col
602  sll_int32, intent(inout) :: buffer(:)
603  sll_int32, intent(in) :: root
604 
605  sll_int32 :: ierr
606 
607  call mpi_bcast(buffer, size(buffer), mpi_integer, root, col%comm, ierr)
608  call sll_s_test_mpi_error(ierr, &
609  'sll_collective_bcast_int32(): MPI_BCAST()')
610  end subroutine sll_collective_bcast_int32
611 
617  subroutine sll_collective_bcast_int64(col, buffer, root)
618  type(sll_t_collective_t), pointer :: col
619  sll_int64, intent(inout) :: buffer(:)
620  sll_int32, intent(in) :: root
621 
622  sll_int32 :: ierr
623 
624  call mpi_bcast(buffer, size(buffer), mpi_integer8, root, col%comm, ierr)
625  call sll_s_test_mpi_error(ierr, &
626  'sll_collective_bcast_int64(): MPI_BCAST()')
627  end subroutine sll_collective_bcast_int64
628 
634  subroutine sll_collective_bcast_real64(col, buffer, root)
635  type(sll_t_collective_t), pointer :: col
636  sll_real64, intent(inout) :: buffer(:)
637  sll_int32, intent(in) :: root
638 
639  sll_int32 :: ierr
640 
641  call mpi_bcast(buffer, size(buffer), mpi_real8, root, col%comm, ierr)
642  call sll_s_test_mpi_error(ierr, &
643  'sll_collective_bcast_real64(): MPI_BCAST()')
644  end subroutine sll_collective_bcast_real64
645 
652  subroutine sll_collective_bcast_real32(col, buffer, size, root)
653  type(sll_t_collective_t), pointer :: col
654  sll_real32, intent(inout) :: buffer(:)
655  sll_int32, intent(in) :: size
656  sll_int32, intent(in) :: root
657 
658  sll_int32 :: ierr
659 
660  call mpi_bcast(buffer, size, mpi_real, root, col%comm, ierr)
661  call sll_s_test_mpi_error(ierr, &
662  'sll_collective_bcast_real32(): MPI_BCAST()')
663  end subroutine sll_collective_bcast_real32
664 
671  subroutine sll_s_collective_bcast_real64(col, buffer, size, root)
672  type(sll_t_collective_t), pointer :: col
673  sll_real64, intent(inout) :: buffer(:)
674  sll_int32, intent(in) :: size
675  sll_int32, intent(in) :: root
676 
677  sll_int32 :: ierr
678 
679  call mpi_bcast(buffer, size, mpi_real8, root, col%comm, ierr)
680  call sll_s_test_mpi_error(ierr, &
681  'sll_s_collective_bcast_real64(): MPI_BCAST()')
682  end subroutine sll_s_collective_bcast_real64
683 
690  subroutine sll_collective_bcast_comp64(col, buffer, size, root)
691  type(sll_t_collective_t), pointer :: col
692  sll_comp64, intent(inout) :: buffer(:)
693  sll_int32, intent(in) :: size
694  sll_int32, intent(in) :: root
695 
696  sll_int32 :: ierr
697 
698  call mpi_bcast(buffer, size, mpi_double_complex, root, col%comm, ierr)
699  call sll_s_test_mpi_error(ierr, &
700  'sll_collective_bcast_comp64(): MPI_BCAST()')
701  end subroutine sll_collective_bcast_comp64
702 
709  subroutine sll_collective_bcast_comp32(col, buffer, size, root)
710  type(sll_t_collective_t), pointer :: col
711  sll_comp32, intent(inout) :: buffer(:)
712  sll_int32, intent(in) :: size
713  sll_int32, intent(in) :: root
714 
715  sll_int32 :: ierr
716 
717  call mpi_bcast(buffer, size, mpi_complex, root, col%comm, ierr)
718  call sll_s_test_mpi_error(ierr, &
719  'sll_collective_bcast_comp32(): MPI_BCAST()')
720  end subroutine sll_collective_bcast_comp32
721 
728  subroutine sll_s_collective_bcast_3d_real64(col, buffer, root)
729  type(sll_t_collective_t) :: col
730  sll_real64, dimension(:, :, :), intent(inout) :: buffer
731  sll_int32, intent(in) :: root
732  sll_int32 :: count, ierr
733  count = size(buffer)
734  call mpi_bcast(buffer, count, mpi_double_precision, root, col%comm, ierr)
735  sll_assert(ierr == mpi_success)
736  end subroutine sll_s_collective_bcast_3d_real64
737 
738  ! Consider joining the next 'gather' interfaces into a single one.
745  subroutine sll_collective_gather_real32(col, send_buf, send_sz, root, &
746  rec_buf)
747  type(sll_t_collective_t), pointer :: col
748  sll_real32, dimension(:), intent(in) :: send_buf ! what would change...
749  sll_int32 :: send_sz
750  sll_real32, dimension(:), intent(out) :: rec_buf ! would also change
751  sll_int32, intent(in) :: root
752  sll_int32 :: ierr
753  !sll_int32 :: rec_count ! size of receive buf
754  ! FIXME: add some argument checking here
755  !rec_count = send_sz*col%size
756  !Note that the 5th argument at the root indicates the number of items
757  !it receives from each task. It is not the total number of items received.
758  call mpi_gather(send_buf, send_sz, mpi_real, rec_buf, send_sz, &
759  mpi_real, root, col%comm, ierr)
760  call sll_s_test_mpi_error(ierr, &
761  'sll_collective_gather_real(): MPI_GATHER()')
762  end subroutine sll_collective_gather_real32
763 
764  subroutine sll_collective_gather_real64(col, send_buf, send_sz, root, &
765  rec_buf)
766  type(sll_t_collective_t), pointer :: col
767  sll_real64, dimension(:), intent(in) :: send_buf ! what would change...
768  sll_int32 :: send_sz
769  sll_real64, dimension(:), intent(out) :: rec_buf ! would also change
770  sll_int32, intent(in) :: root
771  sll_int32 :: ierr
772  !sll_int32 :: rec_count ! size of receive buf
773  ! FIXME: add some argument checking here
774  !rec_count = send_sz*col%size
775  !Note that the 5th argument at the root indicates the number of items
776  !it receives from each task. It is not the total number of items received.
777  call mpi_gather(send_buf, send_sz, mpi_double_precision, rec_buf, send_sz, &
778  mpi_double_precision, root, col%comm, ierr)
779  call sll_s_test_mpi_error(ierr, &
780  'sll_collective_gather_real(): MPI_GATHER()')
781  end subroutine sll_collective_gather_real64
782 
783  subroutine sll_collective_gather_logical(col, send_buf, root, rec_buf)
784  type(sll_t_collective_t), pointer :: col
785  logical, intent(in) :: send_buf(:)
786  sll_int32, intent(in) :: root
787  logical, intent(out) :: rec_buf(:)
788 
789  sll_int32 :: send_sz, ierr
790 
791  ! Get message size
792  send_sz = size(send_buf)
793 
794  ! Root: Verify that receive buffer is of correct size
795  if (col%rank == root) then
796  sll_assert(size(rec_buf) == send_sz*col%size)
797  end if
798 
799  !Note that the 5th argument at the root indicates the number of items
800  !it receives from each task. It is not the total number of items received.
801  call mpi_gather(send_buf, send_sz, mpi_logical, rec_buf, send_sz, &
802  mpi_logical, root, col%comm, ierr)
803  call sll_s_test_mpi_error(ierr, &
804  'sll_collective_gather_real(): MPI_GATHER()')
805  end subroutine sll_collective_gather_logical
806 
817  subroutine sll_s_collective_gatherv_real(col, send_buf, send_count, &
818  recvcnts, displs, root, rec_buf)
819  type(sll_t_collective_t), pointer :: col
820  sll_real32, dimension(:), intent(in) :: send_buf ! what would change...
821  sll_int32, intent(in) :: send_count
822  sll_int32, dimension(:), intent(in) :: recvcnts
823  sll_int32, dimension(:), intent(in) :: displs
824  sll_real32, dimension(:), intent(out) :: rec_buf ! would also change
825  sll_int32, intent(in) :: root
826  sll_int32 :: ierr
827  ! FIXME: Argument checking
828  call sll_check_collective_ptr(col)
829  ! displs, rec_buf and recvcnts significant only for root
830  if (col%rank .eq. root) then
831  sll_assert(SIZE(recvcnts) .eq. col%size)
832  sll_assert(SIZE(displs) .eq. col%size)
833  sll_assert(SIZE(rec_buf) .eq. sum(recvcnts))
834  end if
835  call mpi_gatherv( &
836  send_buf, &
837  send_count, &
838  mpi_real, &
839  rec_buf, &
840  recvcnts, &
841  displs, &
842  mpi_real, &
843  root, &
844  col%comm, &
845  ierr)
846  call sll_s_test_mpi_error(ierr, &
847  'sll_s_collective_gatherv_real(): MPI_GATHERV()')
848  end subroutine sll_s_collective_gatherv_real
849 
860  subroutine sll_s_collective_gatherv_real64(col, send_buf, send_count, &
861  recvcnts, displs, root, rec_buf)
862  type(sll_t_collective_t), pointer :: col
863  sll_real64, dimension(:), intent(in) :: send_buf ! what would change...
864  sll_int32, intent(in) :: send_count
865  sll_int32, dimension(:), intent(in) :: recvcnts
866  sll_int32, dimension(:), intent(in) :: displs
867  sll_real64, dimension(:), intent(out) :: rec_buf ! would also change
868  sll_int32, intent(in) :: root
869  sll_int32 :: ierr
870  ! FIXME: Argument checking
871  call sll_check_collective_ptr(col)
872  ! displs, rec_buf and recvcnts significant only for root
873  if (col%rank .eq. root) then
874  sll_assert(SIZE(recvcnts) .eq. col%size)
875  sll_assert(SIZE(displs) .eq. col%size)
876  sll_assert(SIZE(rec_buf) .eq. sum(recvcnts))
877  end if
878  call mpi_gatherv(send_buf, send_count, mpi_real8, rec_buf, recvcnts, &
879  displs, mpi_real8, root, col%comm, ierr)
880  call sll_s_test_mpi_error(ierr, &
881  'sll_s_collective_gatherv_real64(): MPI_GATHERV()')
882  end subroutine sll_s_collective_gatherv_real64
883 
891  subroutine sll_collective_allgather_int(col, send_buf, send_sz, &
892  recv_buf, recv_sz)
893  type(sll_t_collective_t), pointer :: col
894  sll_int32, dimension(:), intent(in) :: send_buf ! what would change...
895  sll_int32, intent(in) :: send_sz
896  sll_int32, dimension(:), intent(inout) :: recv_buf ! would also change
897  sll_int32, intent(in) :: recv_sz
898  sll_int32 :: ierr
899  ! FIXME: Argument checking
900  call sll_check_collective_ptr(col)
901  call mpi_allgather(send_buf(:), send_sz, mpi_integer, &
902  recv_buf(:), recv_sz, mpi_integer, col%comm, ierr)
903  call sll_s_test_mpi_error(ierr, &
904  'sll_collective_allgather_int(): MPI_ALLGATHER()')
905 
906  end subroutine
907 
909  col, &
910  send_buf, &
911  send_sz, &
912  recv_buf, &
913  recv_sz)
914 
915  type(sll_t_collective_t), pointer :: col
916  sll_real64, dimension(:), intent(in) :: send_buf
917  sll_int32, intent(in) :: send_sz
918  sll_real64, dimension(:), intent(out) :: recv_buf ! would change
919  sll_int32, intent(in) :: recv_sz
920  sll_int32 :: ierr
921  ! FIXME: Argument checking
922  call sll_check_collective_ptr(col)
923  call mpi_allgather(send_buf(:), send_sz, mpi_real8, &
924  recv_buf(:), recv_sz, mpi_real8, &
925  col%comm, ierr)
926  call sll_s_test_mpi_error(ierr, &
927  'sll_collective_allgather_int(): MPI_ALLGATHER()')
928 
929  end subroutine sll_collective_allgather_real64
930 
940  subroutine sll_collective_allgatherv_real32(col, send_buf, send_cnt, &
941  rec_cnt, displs, rec_buf)
942  type(sll_t_collective_t), pointer :: col
943  sll_real32, dimension(:), intent(in) :: send_buf ! what would change...
944  sll_int32, intent(in) :: send_cnt
945  sll_int32, dimension(:), intent(in) :: rec_cnt
946  sll_int32, dimension(:), intent(in) :: displs
947  sll_real32, dimension(:), intent(out) :: rec_buf ! would also change
948  sll_int32 :: ierr
949  ! FIXME: argument checking
950  sll_assert(col%size .eq. SIZE(displs))
951  sll_assert(col%size .eq. SIZE(rec_cnt))
952  call mpi_allgatherv(send_buf, send_cnt, mpi_real, rec_buf, rec_cnt, &
953  displs, mpi_real, col%comm, ierr)
954  call sll_s_test_mpi_error(ierr, &
955  'sll_collective_allgatherv_real32(): MPI_ALLGATHERV()')
956  end subroutine sll_collective_allgatherv_real32
957 
959  col, &
960  send_buf, &
961  send_cnt, &
962  rec_cnt, &
963  displs, &
964  rec_buf)
965 
966  type(sll_t_collective_t), pointer :: col
967  sll_real64, dimension(:), intent(in) :: send_buf ! what would change...
968  sll_int32, intent(in) :: send_cnt
969  sll_int32, dimension(:), intent(in) :: rec_cnt
970  sll_int32, dimension(:), intent(in) :: displs
971  sll_real64, dimension(:), intent(out) :: rec_buf ! would also change
972  sll_int32 :: ierr
973  ! FIXME: argument checking
974  sll_assert(col%size .eq. SIZE(displs))
975  sll_assert(col%size .eq. SIZE(rec_cnt))
976  call mpi_allgatherv( &
977  send_buf, &
978  send_cnt, &
979  mpi_double_precision, &
980  rec_buf, &
981  rec_cnt, &
982  displs, &
983  mpi_double_precision, &
984  col%comm, &
985  ierr)
986  call sll_s_test_mpi_error(ierr, &
987  'sll_s_collective_allgatherv_real64(): MPI_ALLGATHERV()')
989 
997  subroutine sll_collective_scatter_real(col, send_buf, send_count, root, &
998  rec_buf)
999  type(sll_t_collective_t), pointer :: col
1000  sll_real32, dimension(:), intent(in) :: send_buf ! what would change...
1001  sll_int32, intent(in) :: send_count
1002  sll_int32, intent(in) :: root
1003  sll_real32, dimension(:), intent(in) :: rec_buf ! would also change
1004  sll_int32 :: ierr
1005  !sll_int32 :: recvcount
1006  ! FIXME: argument checking
1007  ! recvcount = size/col%size
1008  !send_buf and send_count significant only at root
1009  if (col%rank .eq. root) then
1010  sll_assert(SIZE(send_buf) .eq. (send_count*(col%size)))
1011  end if
1012  call mpi_scatter(send_buf, send_count, mpi_real, rec_buf, send_count, &
1013  mpi_real, root, col%comm, ierr)
1014  call sll_s_test_mpi_error(ierr, &
1015  'sll_collective_scatter_real(): MPI_SCATTER()')
1016  end subroutine sll_collective_scatter_real
1017 
1025  subroutine sll_collective_scatter_real64(col, send_buf, send_count, root, &
1026  rec_buf)
1027  type(sll_t_collective_t), pointer :: col
1028  sll_real64, dimension(:), intent(in) :: send_buf ! what would change...
1029  sll_int32, intent(in) :: send_count
1030  sll_int32, intent(in) :: root
1031  sll_real64, dimension(:), intent(in) :: rec_buf ! would also change
1032  sll_int32 :: ierr
1033  !sll_int32 :: recvcount
1034  ! FIXME: argument checking
1035  ! recvcount = size/col%size
1036  !send_buf and send_count significant only at root
1037  if (col%rank .eq. root) then
1038  sll_assert(SIZE(send_buf) .eq. (send_count*(col%size)))
1039  end if
1040  call mpi_scatter(send_buf, send_count, mpi_double_precision, rec_buf, send_count, &
1041  mpi_double_precision, root, col%comm, ierr)
1042  call sll_s_test_mpi_error(ierr, &
1043  'sll_collective_scatter_real64(): MPI_SCATTER()')
1044  end subroutine sll_collective_scatter_real64
1045 
1053  subroutine sll_collective_scatter_int32(col, send_buf, send_count, root, &
1054  rec_buf)
1055  type(sll_t_collective_t), pointer :: col
1056  sll_int32, dimension(:), intent(in) :: send_buf ! what would change...
1057  sll_int32, intent(in) :: send_count
1058  sll_int32, intent(in) :: root
1059  sll_int32, dimension(:), intent(in) :: rec_buf ! would also change
1060  sll_int32 :: ierr
1061  !sll_int32 :: recvcount
1062  ! FIXME: argument checking
1063  ! recvcount = size/col%size
1064  !send_buf and send_count significant only at root
1065  if (col%rank .eq. root) then
1066  sll_assert(SIZE(send_buf) .eq. (send_count*(col%size)))
1067  end if
1068  call mpi_scatter(send_buf, send_count, mpi_double_precision, rec_buf, send_count, &
1069  mpi_integer, root, col%comm, ierr)
1070  call sll_s_test_mpi_error(ierr, &
1071  'sll_collective_scatter_real64(): MPI_SCATTER()')
1072  end subroutine sll_collective_scatter_int32
1073 
1085  subroutine sll_s_collective_scatterv_real(col, send_buf, send_count, &
1086  displs, recv_count, root, rec_buf)
1087  type(sll_t_collective_t), pointer :: col
1088  sll_real32, dimension(:), intent(in) :: send_buf ! what would change...
1089  sll_int32, dimension(:), intent(in) :: send_count
1090  sll_int32, dimension(:), intent(in) :: displs
1091  sll_int32, intent(in) :: recv_count
1092  sll_real32, dimension(:), intent(out) :: rec_buf ! would also change
1093  sll_int32, intent(in) :: root
1094  sll_int32 :: ierr
1095  ! FIXME: ARG CHECKING!
1096  sll_assert(SIZE(send_count) .eq. col%size)
1097  sll_assert(SIZE(displs) .eq. col%size)
1098  call mpi_scatterv(send_buf, send_count, displs, mpi_real, rec_buf, &
1099  recv_count, mpi_real, root, col%comm, ierr)
1100  call sll_s_test_mpi_error(ierr, &
1101  'sll_s_collective_scatterv_real(): MPI_SCATTERV()')
1102  end subroutine sll_s_collective_scatterv_real
1103 
1111  subroutine sll_s_collective_allreduce_real32(col, send_buf, count, op, &
1112  rec_buf)
1113  type(sll_t_collective_t), pointer :: col
1114  sll_real32, dimension(:), intent(in) :: send_buf ! what would change...
1115  sll_int32, intent(in) :: count
1116  sll_int32, intent(in) :: op
1117  sll_real32, dimension(:), intent(out) :: rec_buf ! would also change
1118  sll_int32 :: ierr
1119  ! FIXME: ARG CHECKING!
1120  call sll_check_collective_ptr(col)
1121  call mpi_allreduce( &
1122  send_buf, &
1123  rec_buf, &
1124  count, &
1125  mpi_real, &
1126  op, &
1127  col%comm, &
1128  ierr)
1129  call sll_s_test_mpi_error( &
1130  ierr, &
1131  'sll_s_collective_allreduce_real32(): MPI_ALLREDUCE()')
1132  end subroutine sll_s_collective_allreduce_real32
1133 
1141  subroutine sll_collective_allreduce_real64(col, send_buf, count, op, &
1142  rec_buf)
1143  type(sll_t_collective_t), pointer :: col
1144  sll_real64, dimension(:), intent(in) :: send_buf ! what would change...
1145  sll_int32, intent(in) :: count
1146  sll_int32, intent(in) :: op
1147  sll_real64, dimension(:), intent(out) :: rec_buf ! would also change
1148  sll_int32 :: ierr
1149  ! FIXME: ARG CHECKING!
1150  call sll_check_collective_ptr(col)
1151  call mpi_allreduce( &
1152  send_buf, &
1153  rec_buf, &
1154  count, &
1155  mpi_double_precision, &
1156  op, &
1157  col%comm, &
1158  ierr)
1159  call sll_s_test_mpi_error( &
1160  ierr, &
1161  'sll_collective_allreduce_real64(): MPI_ALLREDUCE()')
1162  end subroutine sll_collective_allreduce_real64
1163 
1171  subroutine sll_collective_allreduce_real64_2darray(col, send_buf, count, op, &
1172  rec_buf)
1173  type(sll_t_collective_t), pointer :: col
1174  sll_real64, dimension(:, :), intent(in) :: send_buf ! what would change...
1175  sll_int32, intent(in) :: count
1176  sll_int32, intent(in) :: op
1177  sll_real64, dimension(:, :), intent(out) :: rec_buf ! would also change
1178  sll_int32 :: ierr
1179  ! FIXME: ARG CHECKING!
1180  call sll_check_collective_ptr(col)
1181  call mpi_barrier(col%comm, ierr)
1182  call mpi_allreduce( &
1183  send_buf, &
1184  rec_buf, &
1185  count, &
1186  mpi_double_precision, &
1187  op, &
1188  col%comm, &
1189  ierr)
1190  call sll_s_test_mpi_error( &
1191  ierr, &
1192  'sll_collective_allreduce_real64(): MPI_ALLREDUCE()')
1194 
1202  subroutine sll_collective_allreduce_comp64(col, send_buf, count, op, &
1203  rec_buf)
1204  type(sll_t_collective_t), pointer :: col
1205  sll_comp64, dimension(:), intent(in) :: send_buf ! what would change...
1206  sll_int32, intent(in) :: count
1207  sll_int32, intent(in) :: op
1208  sll_comp64, dimension(:), intent(out) :: rec_buf ! would also change
1209  sll_int32 :: ierr
1210  ! FIXME: ARG CHECKING!
1211  call sll_check_collective_ptr(col)
1212  call mpi_allreduce( &
1213  send_buf, &
1214  rec_buf, &
1215  count, &
1216  mpi_double_complex, &
1217  op, &
1218  col%comm, &
1219  ierr)
1220  call sll_s_test_mpi_error( &
1221  ierr, &
1222  'sll_collective_allreduce_comp64(): MPI_ALLREDUCE()')
1223  end subroutine sll_collective_allreduce_comp64
1224 
1225  subroutine sll_collective_allreduce_comp32(col, send_buf, count, op, &
1226  rec_buf)
1227  type(sll_t_collective_t), pointer :: col
1228  sll_comp32, dimension(:), intent(in) :: send_buf ! what would change...
1229  sll_int32, intent(in) :: count
1230  sll_int32, intent(in) :: op
1231  sll_comp32, dimension(:), intent(out) :: rec_buf ! would also change
1232  sll_int32 :: ierr
1233  ! FIXME: ARG CHECKING!
1234  call sll_check_collective_ptr(col)
1235  call mpi_allreduce( &
1236  send_buf, &
1237  rec_buf, &
1238  count, &
1239  mpi_complex, &
1240  op, &
1241  col%comm, &
1242  ierr)
1243  call sll_s_test_mpi_error( &
1244  ierr, &
1245  'sll_collective_allreduce_comp32(): MPI_ALLREDUCE()')
1246  end subroutine sll_collective_allreduce_comp32
1247 
1255  subroutine sll_collective_allreduce_int32(col, send_buf, count, op, &
1256  rec_buf)
1257  type(sll_t_collective_t), pointer :: col
1258  sll_int32, dimension(:), intent(in) :: send_buf ! what would change...
1259  sll_int32, intent(in) :: count
1260  sll_int32, intent(in) :: op
1261  sll_int32, dimension(:), intent(out) :: rec_buf ! would also change
1262  sll_int32 :: ierr
1263  ! FIXME: ARG CHECKING!
1264  call sll_check_collective_ptr(col)
1265  call mpi_allreduce( &
1266  send_buf, &
1267  rec_buf, &
1268  count, &
1269  mpi_integer, &
1270  op, &
1271  col%comm, &
1272  ierr)
1273  call sll_s_test_mpi_error( &
1274  ierr, &
1275  'sll_collective_allreduce_int32(): MPI_ALLREDUCE()')
1276  end subroutine sll_collective_allreduce_int32
1277 
1285  subroutine sll_s_collective_allreduce_logical(col, send_buf, count, op, &
1286  rec_buf)
1287  type(sll_t_collective_t), pointer :: col
1288  logical, dimension(:), intent(in) :: send_buf ! what would change...
1289  sll_int32, intent(in) :: count
1290  sll_int32, intent(in) :: op
1291  logical, dimension(:), intent(out) :: rec_buf ! would also change
1292  sll_int32 :: ierr
1293  ! FIXME: MORE ARG CHECKING!
1294  call sll_check_collective_ptr(col)
1295  call mpi_allreduce(send_buf, rec_buf, count, mpi_logical, op, &
1296  col%comm, ierr)
1297  call sll_s_test_mpi_error(ierr, &
1298  'sll_s_collective_allreduce_logical(): MPI_ALLREDUCE()')
1299  end subroutine sll_s_collective_allreduce_logical
1300 
1306  use mpi
1307  type(sll_t_collective_t) :: col
1308  sll_real64, dimension(:, :, :), intent(inout) :: buffer
1309  sll_real64, dimension(:, :, :), allocatable :: temp
1310  sll_int32 :: count, ierr
1311  count = size(buffer)
1312  allocate (temp(size(buffer, dim=1), size(buffer, dim=2), size(buffer, dim=3)))
1313  call mpi_allreduce(buffer, temp, count, mpi_double_precision, mpi_sum, &
1314  col%comm, ierr)
1315  sll_assert(ierr == mpi_success)
1316  buffer = temp
1317  deallocate (temp)
1319 
1327  subroutine sll_s_collective_reduce_int(col, send_buf, size, op, root_rank, &
1328  rec_buf)
1329  type(sll_t_collective_t), pointer :: col
1330  sll_int32, dimension(:), intent(in) :: send_buf ! what would change...
1331  sll_int32, intent(in) :: size
1332  sll_int32, intent(in) :: op
1333  sll_int32, intent(in) :: root_rank
1334  sll_int32, dimension(:), intent(in) :: rec_buf ! would also change
1335 
1336  sll_int32 :: ierr
1337  ! FIXME: ARG CHECKING!
1338  call mpi_reduce(send_buf, rec_buf, size, mpi_integer, op, root_rank, &
1339  col%comm, ierr)
1340  call sll_s_test_mpi_error(ierr, &
1341  'sll_collective_reduce_real(): MPI_REDUCE()')
1342  end subroutine sll_s_collective_reduce_int
1343 
1351  subroutine sll_s_collective_reduce_real32(col, send_buf, size, op, root_rank, &
1352  rec_buf)
1353  type(sll_t_collective_t), pointer :: col
1354  sll_real32, dimension(:), intent(in) :: send_buf ! what would change...
1355  sll_int32, intent(in) :: size
1356  sll_int32, intent(in) :: op
1357  sll_int32, intent(in) :: root_rank
1358  sll_real32, dimension(:), intent(in) :: rec_buf ! would also change
1359 
1360  sll_int32 :: ierr
1361  ! FIXME: ARG CHECKING!
1362  call mpi_reduce(send_buf, rec_buf, size, mpi_real, op, root_rank, &
1363  col%comm, ierr)
1364  call sll_s_test_mpi_error(ierr, &
1365  'sll_collective_reduce_real(): MPI_REDUCE()')
1366  end subroutine sll_s_collective_reduce_real32
1367 
1368  subroutine sll_s_collective_reduce_real64(col, send_buf, size, op, root_rank, &
1369  rec_buf)
1370  type(sll_t_collective_t), pointer :: col
1371  sll_real64, dimension(:), intent(in) :: send_buf
1372  sll_int32, intent(in) :: size
1373  sll_int32, intent(in) :: op
1374  sll_int32, intent(in) :: root_rank
1375  sll_real64, dimension(:), intent(out) :: rec_buf
1376  sll_int32 :: ierr
1377 
1378  ! FIXME: ARG CHECKING!
1379  call mpi_reduce(send_buf, rec_buf, size, mpi_double_precision, op, &
1380  root_rank, col%comm, ierr)
1381  call sll_s_test_mpi_error(ierr, &
1382  'sll_s_collective_reduce_real64(): MPI_REDUCE()')
1383  end subroutine sll_s_collective_reduce_real64
1384 
1385  subroutine sll_collective_reduce_real64_2d(col, send_buf, size, op, root_rank, &
1386  rec_buf)
1387  type(sll_t_collective_t), pointer :: col
1388  sll_real64, dimension(:, :), intent(in) :: send_buf
1389  sll_int32, intent(in) :: size
1390  sll_int32, intent(in) :: op
1391  sll_int32, intent(in) :: root_rank
1392  sll_real64, dimension(:, :), intent(out) :: rec_buf
1393  sll_int32 :: ierr
1394 
1395  ! FIXME: ARG CHECKING!
1396  call mpi_reduce(send_buf, rec_buf, size, mpi_double_precision, op, &
1397  root_rank, col%comm, ierr)
1398  call sll_s_test_mpi_error(ierr, &
1399  'sll_collective_reduce_real64_2d(): MPI_REDUCE()')
1400  end subroutine sll_collective_reduce_real64_2d
1401 
1409  subroutine sll_collective_reduce_comp64(col, send_buf, size, op, root_rank, &
1410  rec_buf)
1411  type(sll_t_collective_t), pointer :: col
1412  sll_comp64, dimension(:), intent(in) :: send_buf
1413  sll_int32, intent(in) :: size
1414  sll_int32, intent(in) :: op
1415  sll_int32, intent(in) :: root_rank
1416  sll_comp64, dimension(:), intent(out) :: rec_buf
1417  sll_int32 :: ierr
1418 
1419  ! FIXME: ARG CHECKING!
1420  call mpi_reduce(send_buf, rec_buf, size, mpi_double_complex, op, &
1421  root_rank, col%comm, ierr)
1422  call sll_s_test_mpi_error(ierr, &
1423  'sll_collective_reduce_comp64(): MPI_REDUCE()')
1424  end subroutine sll_collective_reduce_comp64
1425 
1426  subroutine sll_collective_reduce_comp32(col, send_buf, size, op, root_rank, &
1427  rec_buf)
1428  type(sll_t_collective_t), pointer :: col
1429  sll_comp32, dimension(:), intent(in) :: send_buf
1430  sll_int32, intent(in) :: size
1431  sll_int32, intent(in) :: op
1432  sll_int32, intent(in) :: root_rank
1433  sll_comp32, dimension(:), intent(out) :: rec_buf
1434  sll_int32 :: ierr
1435 
1436  ! FIXME: ARG CHECKING!
1437  call mpi_reduce(send_buf, rec_buf, size, mpi_complex, op, &
1438  root_rank, col%comm, ierr)
1439  call sll_s_test_mpi_error(ierr, &
1440  'sll_collective_reduce_comp32(): MPI_REDUCE()')
1441  end subroutine sll_collective_reduce_comp32
1442 
1450  subroutine sll_s_collective_reduce_logical(col, send_buf, size, op, root_rank, &
1451  rec_buf)
1452  type(sll_t_collective_t), pointer :: col
1453  LOGICAL, DIMENSION(:), intent(in) :: send_buf ! what would change...
1454  sll_int32, intent(in) :: size
1455  sll_int32, intent(in) :: op
1456  sll_int32, intent(in) :: root_rank
1457  LOGICAL, DIMENSION(:), intent(in) :: rec_buf ! would also change
1458 
1459  sll_int32 :: ierr
1460  ! FIXME: ARG CHECKING!
1461  call mpi_reduce(send_buf, rec_buf, size, mpi_logical, op, root_rank, &
1462  col%comm, ierr)
1463  call sll_s_test_mpi_error(ierr, &
1464  'sll_collective_reduce_real(): MPI_REDUCE()')
1465  end subroutine sll_s_collective_reduce_logical
1466 
1473  subroutine sll_s_collective_alltoall_int(send_buf, send_count, &
1474  recv_count, recv_buf, col)
1475  sll_int32, dimension(:), intent(in) :: send_buf
1476  sll_int32, intent(in) :: send_count
1477  sll_int32, intent(in) :: recv_count
1478  sll_int32, dimension(:), intent(out) :: recv_buf
1479  type(sll_t_collective_t), pointer :: col
1480  sll_int32 :: ierr
1481  call sll_check_collective_ptr(col)
1482  ! FIXME: MORE ARG CHECKING
1483  call mpi_alltoall(send_buf, send_count, mpi_integer, &
1484  recv_buf, recv_count, mpi_integer, col%comm, ierr)
1485  call sll_s_test_mpi_error(ierr, &
1486  'sll_s_collective_alltoall_int(): MPI_ALLTOALLV()')
1487 
1488  end subroutine sll_s_collective_alltoall_int
1489 
1490  subroutine sll_collective_alltoall_double(send_buf, send_count, &
1491  recv_count, recv_buf, col)
1492  sll_real64, dimension(:), intent(in) :: send_buf
1493  sll_int32, intent(in) :: send_count
1494  sll_int32, intent(in) :: recv_count
1495  sll_real64, dimension(:), intent(out) :: recv_buf
1496  type(sll_t_collective_t), pointer :: col
1497  sll_int32 :: ierr
1498  call sll_check_collective_ptr(col)
1499  ! FIXME: MORE ARG CHECKING
1500  call mpi_alltoall(send_buf, send_count, mpi_double_precision, &
1501  recv_buf, recv_count, mpi_double_precision, &
1502  col%comm, ierr)
1503  call sll_s_test_mpi_error(ierr, &
1504  'sll_collective_alltoall_double(): MPI_ALLTOALL()')
1505  end subroutine sll_collective_alltoall_double
1506 
1507  subroutine sll_collective_alltoall_complex_double(send_buf, send_count, &
1508  recv_count, recv_buf, col)
1509  sll_comp64, dimension(:), intent(in) :: send_buf
1510  sll_int32, intent(in) :: send_count
1511  sll_int32, intent(in) :: recv_count
1512  sll_comp64, dimension(:), intent(out) :: recv_buf
1513  type(sll_t_collective_t), pointer :: col
1514  sll_int32 :: ierr
1515  call sll_check_collective_ptr(col)
1516  ! FIXME: MORE ARG CHECKING
1517  call mpi_alltoall(send_buf, send_count, mpi_double_complex, &
1518  recv_buf, recv_count, mpi_double_complex, &
1519  col%comm, ierr)
1520  call sll_s_test_mpi_error(ierr, &
1521  'sll_collective_alltoall_complex_double(): MPI_ALLTOALL()')
1522 
1524 
1525  ! Explore making this effectively typeless... need a macro implementation;
1526  ! pointer arguments won't work...
1546  subroutine sll_s_collective_alltoallv_real(send_buf, send_cnts, &
1547  send_displs, &
1548  recv_buf, recv_cnts, &
1549  recv_displs, col)
1550  sll_real32, dimension(:), intent(in) :: send_buf
1551  sll_int32, dimension(:), intent(in) :: send_cnts
1552  sll_int32, dimension(:), intent(in) :: send_displs
1553  sll_real32, dimension(:), intent(out) :: recv_buf
1554  sll_int32, dimension(:), intent(in) :: recv_cnts
1555  sll_int32, dimension(:), intent(in) :: recv_displs
1556  type(sll_t_collective_t), pointer :: col
1557  sll_int32 :: ierr
1558  ! FIXME: ARG CHECKING!
1559  call sll_check_collective_ptr(col)
1560  call mpi_alltoallv(send_buf, send_cnts, send_displs, mpi_real, &
1561  recv_buf, recv_cnts, recv_displs, mpi_real, &
1562  col%comm, ierr)
1563  call sll_s_test_mpi_error(ierr, &
1564  'sll_s_collective_alltoallv_real(): MPI_ALLTOALLV()')
1565 
1566  end subroutine sll_s_collective_alltoallv_real
1567 
1568  subroutine sll_s_collective_alltoallv_double(send_buf, send_cnts, &
1569  send_displs, &
1570  recv_buf, recv_cnts, &
1571  recv_displs, col)
1572  sll_real64, dimension(:), intent(in) :: send_buf
1573  sll_int32, dimension(:), intent(in) :: send_cnts
1574  sll_int32, dimension(:), intent(in) :: send_displs
1575  sll_real64, dimension(:), intent(out) :: recv_buf
1576  sll_int32, dimension(:), intent(in) :: recv_cnts
1577  sll_int32, dimension(:), intent(in) :: recv_displs
1578  type(sll_t_collective_t), pointer :: col
1579  sll_int32 :: ierr
1580  ! FIXME: ARG CHECKING!
1581  call sll_check_collective_ptr(col)
1582  call mpi_alltoallv(send_buf, send_cnts, send_displs, &
1583  mpi_double_precision, &
1584  recv_buf, recv_cnts, recv_displs, &
1585  mpi_double_precision, &
1586  col%comm, ierr)
1587  call sll_s_test_mpi_error(ierr, &
1588  'sll_s_collective_alltoallv_double(): MPI_ALLTOALLV()')
1589 
1590  end subroutine sll_s_collective_alltoallv_double
1591 
1592  subroutine sll_collective_alltoallv_complex_double(send_buf, send_cnts, &
1593  send_displs, &
1594  recv_buf, recv_cnts, &
1595  recv_displs, col)
1596  sll_comp64, dimension(:), intent(in) :: send_buf
1597  sll_int32, dimension(:), intent(in) :: send_cnts
1598  sll_int32, dimension(:), intent(in) :: send_displs
1599  sll_comp64, dimension(:), intent(out) :: recv_buf
1600  sll_int32, dimension(:), intent(in) :: recv_cnts
1601  sll_int32, dimension(:), intent(in) :: recv_displs
1602  type(sll_t_collective_t), pointer :: col
1603  sll_int32 :: ierr
1604  ! FIXME: ARG CHECKING!
1605  call sll_check_collective_ptr(col)
1606  call mpi_alltoallv(send_buf, send_cnts, send_displs, &
1607  mpi_double_complex, &
1608  recv_buf, recv_cnts, recv_displs, &
1609  mpi_double_complex, &
1610  col%comm, ierr)
1611  call sll_s_test_mpi_error(ierr, &
1612  'sll_collective_alltoallV_complex_double(): MPI_ALLTOALLV()')
1613 
1615 
1635  subroutine sll_s_collective_alltoallv_int(send_buf, send_cnts, &
1636  send_displs, &
1637  recv_buf, recv_cnts, &
1638  recv_displs, col)
1639  sll_int32, dimension(:), intent(in) :: send_buf
1640  sll_int32, dimension(:), intent(in) :: send_cnts
1641  sll_int32, dimension(:), intent(in) :: send_displs
1642  sll_int32, dimension(:), intent(out) :: recv_buf
1643  sll_int32, dimension(:), intent(in) :: recv_cnts
1644  sll_int32, dimension(:), intent(in) :: recv_displs
1645  type(sll_t_collective_t), pointer :: col
1646  sll_int32 :: ierr
1647  ! FIXME: ARG CHECKING!
1648  call sll_check_collective_ptr(col)
1649  call mpi_alltoallv(send_buf(:), send_cnts(:), send_displs(:), mpi_integer, &
1650  recv_buf(:), recv_cnts(:), recv_displs(:), mpi_integer, &
1651  col%comm, ierr)
1652  call sll_s_test_mpi_error(ierr, &
1653  'sll_s_collective_alltoallv_int(): MPI_ALLTOALLV()')
1654 
1655  end subroutine sll_s_collective_alltoallv_int
1656 
1657 !
1658 !FOR SOME REASON THIS ARGUMENT OF THIS FUNCTION IS DECLARED ALLOCATBLE AND
1659 !POSING PROBLEM WITH F95 STANDARD. WE SHOULD COME BACK TO THIS
1660 !
1661 !warning sll_s_collective_alltoallv_int_simple is not fixed
1662 
1663  ! This toutine is a simpler version of the sll_s_collective_alltoallv_int subroutine
1664  subroutine sll_s_collective_alltoallv_int_simple(send_buf, send_cnts, &
1665  recv_buf, col)
1666  sll_int32, dimension(:), intent(in) :: send_buf
1667  sll_int32, dimension(:), intent(in) :: send_cnts
1668  !sll_int32, ,allocatable dimension(:) :: recv_buf
1669  sll_int32, pointer, dimension(:) :: recv_buf
1670  sll_int32, allocatable, dimension(:) :: send_displs
1671  sll_int32, allocatable, dimension(:) :: recv_cnts
1672  sll_int32, allocatable, dimension(:) :: recv_displs
1673  type(sll_t_collective_t), pointer :: col
1674  sll_int32 :: ierr, size_comm, i, sendcnts_size
1675  sll_int32 :: dum
1676  sendcnts_size = size(send_cnts)
1677 
1678  ! FIXME: ARG CHECKING!
1679  call sll_check_collective_ptr(col)
1680 
1681  size_comm = sll_f_get_collective_size(col)
1682  sll_assert(sendcnts_size .eq. size_comm)
1683 
1684  ! Define RECV_CNTS
1685  sll_allocate(recv_cnts(size_comm), ierr)
1686  call sll_s_collective_alltoall_int(send_cnts, 1, 1, &
1687  recv_cnts, col)
1688 
1689  ! Define RECV_BUF
1690  dum = sum(recv_cnts)
1691  sll_allocate(recv_buf(dum), ierr)
1692 
1693  ! Define SEND_DISPLS
1694  sll_allocate(send_displs(size_comm), ierr)
1695  send_displs(1) = 0
1696  do i = 2, size_comm
1697  send_displs(i) = send_displs(i - 1) + send_cnts(i - 1)
1698  end do
1699 
1700  ! Define RECV_DISPLS
1701  sll_allocate(recv_displs(size_comm), ierr)
1702  recv_displs(1) = 0
1703  do i = 2, size_comm
1704  recv_displs(i) = recv_displs(i - 1) + recv_cnts(i - 1)
1705  end do
1706 
1707  call mpi_alltoallv(send_buf(:), send_cnts(:), send_displs(:), mpi_integer, &
1708  recv_buf(:), recv_cnts(:), recv_displs(:), mpi_integer, &
1709  col%comm, ierr)
1710  call sll_s_test_mpi_error(ierr, &
1711  'sll_s_collective_alltoallv_int(): MPI_ALLTOALLV()')
1712  sll_deallocate_array(recv_displs, ierr)
1713  sll_deallocate_array(send_displs, ierr)
1714  sll_deallocate_array(recv_cnts, ierr)
1716 
1717  ! Explore if the Irecv calls can be made into collective calls in this module
1718  ! subroutine sll_collective_Irecv( )
1719 
1727  subroutine sll_s_collective_globalsum_array_real64(col, summand, root_rank)
1728  type(sll_t_collective_t), pointer :: col
1729  sll_real64, dimension(:), intent(inout) :: summand
1730  sll_int32, intent(in), optional :: root_rank
1731  sll_real64, dimension(:), allocatable :: recvsum
1732  sll_int32:: summand_size
1733  sll_int32 :: ierr
1734 
1735  summand_size = size(summand)
1736  sll_allocate(recvsum(1:summand_size), ierr)
1737 
1738  if (present(root_rank)) then
1739  !Write the result only to node with root_rank
1740  call sll_s_collective_reduce_real64(col, summand, summand_size, mpi_sum, root_rank, &
1741  recvsum)
1742  else
1743  !Write the result to all nodes
1744  call sll_collective_allreduce_real64(col, summand, summand_size, mpi_sum, &
1745  recvsum)
1746  end if
1747  summand = recvsum
1748  sll_deallocate_array(recvsum, ierr)
1750 
1758  subroutine sll_s_collective_globalsum_array_comp64(col, summand, root_rank)
1759  type(sll_t_collective_t), pointer :: col
1760  sll_comp64, dimension(:), intent(inout) :: summand
1761  sll_int32, intent(in), optional :: root_rank
1762  sll_comp64, dimension(:), allocatable :: recvsum
1763  sll_int32:: summand_size
1764  sll_int32 :: ierr
1765 
1766  summand_size = size(summand)
1767  sll_allocate(recvsum(1:summand_size), ierr)
1768 
1769  if (present(root_rank)) then
1770  !Write the result only to node with root_rank
1771  call sll_collective_reduce_comp64(col, summand, summand_size, mpi_sum, root_rank, &
1772  recvsum)
1773  else
1774  !Write the result to all nodes
1775  call sll_collective_allreduce_comp64(col, summand, summand_size, mpi_sum, &
1776  recvsum)
1777  end if
1778  summand = recvsum
1779  sll_deallocate_array(recvsum, ierr)
1781 
1787  subroutine sll_collective_globalsum_real64(col, summand, root_rank)
1788  type(sll_t_collective_t), pointer :: col
1789  sll_real64, intent(inout) :: summand
1790  sll_int32, optional, intent(in) :: root_rank
1791  sll_real64, dimension(1) :: summand_tmp
1792 
1793  !Promote the scalar to an array
1794  summand_tmp(1) = summand
1795  if (present(root_rank)) then
1796  call sll_s_collective_globalsum_array_real64(col, summand_tmp, root_rank)
1797  else
1798  call sll_s_collective_globalsum_array_real64(col, summand_tmp)
1799  end if
1800 
1801  summand = summand_tmp(1)
1802  end subroutine sll_collective_globalsum_real64
1803 
1811  subroutine sll_collective_globalsum_array_real32(col, summand, root_rank)
1812  type(sll_t_collective_t), pointer :: col
1813  sll_real32, dimension(:), intent(inout) :: summand
1814  sll_int32, intent(in), optional :: root_rank
1815  sll_real32, dimension(:), allocatable :: recvsum
1816  sll_int32:: summand_size
1817  sll_int32 :: ierr
1818 
1819  summand_size = size(summand)
1820  sll_allocate(recvsum(1:summand_size), ierr)
1821 
1822  if (present(root_rank)) then
1823  !Write the result only to node with root_rank
1824  call sll_s_collective_reduce_real32(col, summand, summand_size, mpi_sum, root_rank, &
1825  recvsum)
1826  else
1827  !Write the result to all nodes
1828  call sll_s_collective_allreduce_real32(col, summand, summand_size, mpi_sum, &
1829  recvsum)
1830  end if
1831  summand = recvsum
1832  sll_deallocate_array(recvsum, ierr)
1834 
1842  subroutine sll_collective_globalsum_array_comp32(col, summand, root_rank)
1843  type(sll_t_collective_t), pointer :: col
1844  sll_comp32, dimension(:), intent(inout) :: summand
1845  sll_int32, intent(in), optional :: root_rank
1846  sll_comp32, dimension(:), allocatable :: recvsum
1847  sll_int32:: summand_size
1848  sll_int32 :: ierr
1849 
1850  summand_size = size(summand)
1851  sll_allocate(recvsum(1:summand_size), ierr)
1852 
1853  if (present(root_rank)) then
1854  !Write the result only to node with root_rank
1855  call sll_collective_reduce_comp32(col, summand, summand_size, mpi_sum, root_rank, &
1856  recvsum)
1857  else
1858  !Write the result to all nodes
1859  call sll_collective_allreduce_comp32(col, summand, summand_size, mpi_sum, &
1860  recvsum)
1861  end if
1862  summand = recvsum
1863  sll_deallocate_array(recvsum, ierr)
1865 
1873  subroutine sll_collective_globalsum_array_int32(col, summand, root_rank)
1874  type(sll_t_collective_t), pointer :: col
1875  sll_int32, dimension(:), intent(inout) :: summand
1876  sll_int32, intent(in), optional :: root_rank
1877  sll_int32, dimension(:), allocatable :: recvsum
1878  sll_int32:: summand_size
1879  sll_int32 :: ierr
1880 
1881  summand_size = size(summand)
1882  sll_allocate(recvsum(1:summand_size), ierr)
1883 
1884  if (present(root_rank)) then
1885  !Write the result only to node with root_rank
1886  call sll_s_collective_reduce_int(col, summand, summand_size, mpi_sum, root_rank, &
1887  recvsum)
1888  else
1889  !Write the result to all nodes
1890  call sll_collective_allreduce_int32(col, summand, summand_size, mpi_sum, &
1891  recvsum)
1892  end if
1893  summand = recvsum
1894  sll_deallocate_array(recvsum, ierr)
1896 
1902  subroutine sll_collective_globalsum_comp64(col, summand, root_rank)
1903  type(sll_t_collective_t), pointer :: col
1904  sll_comp64, intent(inout) :: summand
1905  sll_int32, optional, intent(in) :: root_rank
1906  sll_comp64, dimension(1) :: summand_tmp
1907 
1908  !Promote the scalar to an array
1909  summand_tmp(1) = summand
1910  if (present(root_rank)) then
1911  call sll_s_collective_globalsum_array_comp64(col, summand_tmp, root_rank)
1912  else
1913  call sll_s_collective_globalsum_array_comp64(col, summand_tmp)
1914  end if
1915 
1916  summand = summand_tmp(1)
1917  end subroutine sll_collective_globalsum_comp64
1918 
1924  subroutine sll_collective_globalsum_int32(col, summand, root_rank)
1925  type(sll_t_collective_t), pointer :: col
1926  sll_int32, intent(inout) :: summand
1927  sll_int32, optional, intent(in) :: root_rank
1928  sll_int32, dimension(1) :: summand_tmp
1929 
1930  !Promote the scalar to an array
1931  summand_tmp(1) = summand
1932  if (present(root_rank)) then
1933  call sll_collective_globalsum_array_int32(col, summand_tmp, root_rank)
1934  else
1935  call sll_collective_globalsum_array_int32(col, summand_tmp)
1936  end if
1937 
1938  summand = summand_tmp(1)
1939  end subroutine sll_collective_globalsum_int32
1940 
1946  subroutine sll_collective_globalsum_comp32(col, summand, root_rank)
1947  type(sll_t_collective_t), pointer :: col
1948  sll_comp32, intent(inout) :: summand
1949  sll_int32, optional, intent(in) :: root_rank
1950  sll_comp32, dimension(1) :: summand_tmp
1951 
1952  !Promote the scalar to an array
1953  summand_tmp(1) = summand
1954  if (present(root_rank)) then
1955  call sll_collective_globalsum_array_comp32(col, summand_tmp, root_rank)
1956  else
1957  call sll_collective_globalsum_array_comp32(col, summand_tmp)
1958  end if
1959 
1960  summand = summand_tmp(1)
1961  end subroutine sll_collective_globalsum_comp32
1962 
1968  subroutine sll_collective_globalsum_real32(col, summand, root_rank)
1969  type(sll_t_collective_t), pointer :: col
1970  sll_real32, intent(inout) :: summand
1971  sll_int32, optional, intent(in) :: root_rank
1972  sll_real32, dimension(1) :: summand_tmp
1973 
1974  !Promote the scalar to an array
1975  summand_tmp(1) = summand
1976  if (present(root_rank)) then
1977  call sll_collective_globalsum_array_real32(col, summand_tmp, root_rank)
1978  else
1979  call sll_collective_globalsum_array_real32(col, summand_tmp)
1980  end if
1981 
1982  summand = summand_tmp(1)
1983  end subroutine sll_collective_globalsum_real32
1984 
1987  use mpi
1989  sll_int32, intent(in) :: comm
1990 
1991  sll_int32 :: ierr
1992 
1993  sll_allocate(sll_f_create_collective, ierr)
1994  sll_assert(ierr == mpi_success)
1995  sll_f_create_collective%comm = comm
1996  sll_f_create_collective%color = 0
1997  sll_f_create_collective%key = 0
1998  ! plain MPI calls will be fine inside "sll_m_collective.F90"
1999  call mpi_comm_rank(sll_f_create_collective%comm, sll_f_create_collective%rank, ierr)
2000  sll_assert(ierr == mpi_success)
2001  call mpi_comm_size(sll_f_create_collective%comm, sll_f_create_collective%size, ierr)
2002  sll_assert(ierr == mpi_success)
2003  end function sll_f_create_collective
2004 
2005 end module sll_m_collective
Gathers into specified locations from all processes in a group.
Scatters a buffer in parts to all processes in a communicator.
Gathers data from all tasks and distribute the combined data to all tasks.
Gathers data from all tasks and deliver the combined data to all tasks.
Combines values from all processes and distributes the result back to all processes.
Sends data from all to all processes.
Sends data from all to all processes; each process may send a different amount of data and provide di...
Broadcasts a message from the process with rank "root" to all other processes of the communicator.
Gathers together values from a group of processes.
Performs a global sum and writes the result in the optional given node.
Reduces values on all processes to a single value.
Sends data from one process to all other processes in a communicator.
Parallelizing facility.
type(sll_t_collective_t) function, pointer, public sll_f_new_collective(parent, color, key)
Creates new (wrapper around) communicators based on colors and keys.
subroutine sll_collective_bcast_int64(col, buffer, root)
Broadcasts an array of type 'sll_int64' from the process with rank "root" to all other processes of t...
subroutine sll_collective_globalsum_real64(col, summand, root_rank)
Performs a global sum and writes the result in the given node.
subroutine, public sll_s_collective_bcast_3d_real64(col, buffer, root)
Broadcasts a 3d array of type 'sll_real64' from the process with rank "root" to all other processes o...
subroutine sll_collective_globalsum_array_real32(col, summand, root_rank)
Performs a global sum over an array and writes the result in the given node If no node in root_rank i...
subroutine, public sll_s_collective_allreduce_logical(col, send_buf, count, op, rec_buf)
Combines logical values from all processes and distributes the result back to all processes.
subroutine, public sll_s_collective_scatterv_real(col, send_buf, send_count, displs, recv_count, root, rec_buf)
Scatters a buffer in parts to all processes in a communicator.
subroutine sll_collective_bcast_comp64(col, buffer, size, root)
Broadcasts an array of type 'sll_comp64' from the process with rank "root" to all other processes of ...
subroutine sll_collective_allreduce_comp64(col, send_buf, count, op, rec_buf)
Combines complex values from all processes and distributes the result back to all processes.
subroutine, public sll_s_collective_alltoallv_int_simple(send_buf, send_cnts, recv_buf, col)
subroutine, public sll_s_collective_bcast_real64(col, buffer, size, root)
Broadcasts an array of type 'sll_real64' from the process with rank "root" to all other processes of ...
subroutine sll_check_collective_ptr(ptr)
Checks if the pointer ptr is associated to an object.
integer(kind=i32) function, public sll_f_get_collective_rank(col)
Determines the rank of the calling process in the communicator.
integer(kind=i32) function sll_get_collective_comm(col)
Gets the id (integer) of the communicator.
subroutine, public sll_s_collective_reduce_int(col, send_buf, size, op, root_rank, rec_buf)
Reduces integer values on all processes to a single value.
subroutine, public sll_s_collective_alltoallv_real(send_buf, send_cnts, send_displs, recv_buf, recv_cnts, recv_displs, col)
Sends real data from all to all processes; each process may send a different amount of data and provi...
subroutine sll_collective_allreduce_comp32(col, send_buf, count, op, rec_buf)
subroutine sll_collective_allreduce_real64_2darray(col, send_buf, count, op, rec_buf)
Combines real values from all processes and distributes the result back to all processes.
subroutine sll_collective_alltoall_double(send_buf, send_count, recv_count, recv_buf, col)
subroutine, public sll_s_collective_allreduce_real32(col, send_buf, count, op, rec_buf)
Combines real values from all processes and distributes the result back to all processes.
subroutine, public sll_s_halt_collective()
Ends the paralell environment.
subroutine sll_collective_reduce_comp32(col, send_buf, size, op, root_rank, rec_buf)
integer(kind=i32) function sll_get_collective_color(col)
Determines the color of the calling process in the communicator.
subroutine, public sll_s_collective_allreduce_sum_3d_real64(col, buffer)
Sums real values from all processes and distributes the result back to all processes.
subroutine sll_collective_bcast_real64(col, buffer, root)
Broadcasts an array of type 'sll_real64' from the process with rank "root" to all other processes of ...
subroutine sll_collective_scatter_int32(col, send_buf, send_count, root, rec_buf)
Sends data from one process to all other processes in a communicator.
subroutine sll_collective_allgatherv_real32(col, send_buf, send_cnt, rec_cnt, displs, rec_buf)
Gathers real data from all tasks and deliver the combined data to all tasks.
subroutine sll_collective_alltoallv_complex_double(send_buf, send_cnts, send_displs, recv_buf, recv_cnts, recv_displs, col)
type(sll_t_collective_t) function, pointer, public sll_f_create_collective(comm)
Function to derive a "collective" from a plain MPI communicator,.
subroutine, public sll_s_collective_reduce_real64(col, send_buf, size, op, root_rank, rec_buf)
subroutine sll_collective_scatter_real(col, send_buf, send_count, root, rec_buf)
Sends data from one process to all other processes in a communicator.
subroutine sll_collective_allgather_int(col, send_buf, send_sz, recv_buf, recv_sz)
Gathers integer data from all tasks and distribute the combined data to all tasks.
subroutine sll_collective_gather_logical(col, send_buf, root, rec_buf)
subroutine sll_collective_globalsum_comp64(col, summand, root_rank)
Performs a global sum and writes the result in the given node.
subroutine sll_collective_bcast_comp32(col, buffer, size, root)
Broadcasts an array of type 'sll_comp32' from the process with rank "root" to all other processes of ...
subroutine, public sll_s_collective_gatherv_real64(col, send_buf, send_count, recvcnts, displs, root, rec_buf)
Gathers real64 values into specified locations from all processes in a group.
subroutine, public sll_s_boot_collective(required)
Starts the paralell environment.
subroutine sll_collective_gather_real64(col, send_buf, send_sz, root, rec_buf)
subroutine, public sll_s_collective_barrier(col)
Blocks until all processes in the communicator have reached this routine.
subroutine, public sll_s_test_mpi_error(ierr, descriptor)
Checks the good execution of collective instruction.
subroutine sll_collective_allreduce_int32(col, send_buf, count, op, rec_buf)
Combines complex values from all processes and distributes the result back to all processes.
subroutine sll_delete_collective(col)
Marks the communicator object for deallocation.
subroutine, public sll_s_collective_reduce_logical(col, send_buf, size, op, root_rank, rec_buf)
Reduces logical values on all processes to a single value.
subroutine, public sll_s_collective_alltoall_int(send_buf, send_count, recv_count, recv_buf, col)
Sends integer data from all to all processes.
subroutine sll_collective_reduce_real64_2d(col, send_buf, size, op, root_rank, rec_buf)
subroutine, public sll_s_collective_allgatherv_real64(col, send_buf, send_cnt, rec_cnt, displs, rec_buf)
subroutine, public sll_s_collective_alltoallv_int(send_buf, send_cnts, send_displs, recv_buf, recv_cnts, recv_displs, col)
Sends integer data from all to all processes; each process may send a different amount of data and pr...
logical function, public sll_f_collectives_are_same(col1, col2)
subroutine sll_collective_bcast_real32(col, buffer, size, root)
Broadcasts an array of type 'sll_real32' from the process with rank "root" to all other processes of ...
type(sll_t_collective_t), pointer, public sll_v_world_collective
Control of subset assignment. Processes with the same color are in the same new communicator.
subroutine sll_collective_allgather_real64(col, send_buf, send_sz, recv_buf, recv_sz)
subroutine, public sll_s_collective_globalsum_array_real64(col, summand, root_rank)
Performs a global sum over an array and writes the result in the given node If no node in root_rank i...
subroutine sll_collective_globalsum_array_int32(col, summand, root_rank)
Performs a global sum over an array and writes the result in the given node If no node in root_rank i...
subroutine sll_collective_globalsum_array_comp32(col, summand, root_rank)
Performs a global sum over an array and writes the result in the given node If no node in root_rank i...
subroutine sll_collective_globalsum_real32(col, summand, root_rank)
Performs a global sum and writes the result in the given node.
type(sll_t_collective_t) function, pointer sll_get_collective_parent(col)
Gets collective parent.
subroutine sll_collective_scatter_real64(col, send_buf, send_count, root, rec_buf)
Sends data from one process to all other processes in a communicator.
subroutine, public sll_s_collective_alltoallv_double(send_buf, send_cnts, send_displs, recv_buf, recv_cnts, recv_displs, col)
subroutine sll_collective_allreduce_real64(col, send_buf, count, op, rec_buf)
Combines real values from all processes and distributes the result back to all processes.
subroutine sll_collective_reduce_comp64(col, send_buf, size, op, root_rank, rec_buf)
Reduces complex values on all processes to a single value.
subroutine sll_collective_globalsum_int32(col, summand, root_rank)
Performs a global sum and writes the result in the given node.
subroutine sll_collective_globalsum_comp32(col, summand, root_rank)
Performs a global sum and writes the result in the given node.
subroutine, public sll_s_collective_globalsum_array_comp64(col, summand, root_rank)
Performs a global sum over an array and writes the result in the given node If no node in root_rank i...
integer(kind=i32) function, public sll_f_get_collective_size(col)
Determines the size of the group associated with a communicator.
subroutine, public sll_s_allocate_collective()
subroutine, public sll_s_collective_reduce_real32(col, send_buf, size, op, root_rank, rec_buf)
Reduces real values on all processes to a single value.
subroutine sll_collective_bcast_int32(col, buffer, root)
Broadcasts an array of type 'sll_int32' from the process with rank "root" to all other processes of t...
subroutine sll_collective_alltoall_complex_double(send_buf, send_count, recv_count, recv_buf, col)
subroutine sll_collective_bcast_logical(col, buffer, root)
Broadcasts an array of type 'logical' from the process with rank "root" to all other processes of the...
subroutine, public sll_s_collective_gatherv_real(col, send_buf, send_count, recvcnts, displs, root, rec_buf)
Gathers real values into specified locations from all processes in a group.
subroutine, public sll_s_set_communicator_collective(communicator_in)
subroutine sll_collective_gather_real32(col, send_buf, send_sz, root, rec_buf)
Gathers together real values from a group of processes.
Wrapper around the communicator.
    Report Typos and Errors