MODULE transfert_mpi_collectives_mod IMPLICIT NONE INTERFACE bcast_mpi MODULE PROCEDURE bcast_mpi_c MODULE PROCEDURE bcast_mpi_i, bcast_mpi_i1, bcast_mpi_i2, bcast_mpi_i3, bcast_mpi_i4 MODULE PROCEDURE bcast_mpi_r, bcast_mpi_r1, bcast_mpi_r2, bcast_mpi_r3, bcast_mpi_r4 MODULE PROCEDURE bcast_mpi_l, bcast_mpi_l1, bcast_mpi_l2, bcast_mpi_l3, bcast_mpi_l4 END INTERFACE CONTAINS SUBROUTINE gather_field(field_loc,field_glo) USE field_mod USE domain_mod USE mpi_mod USE mpipara TYPE(t_field),POINTER :: field_loc(:) TYPE(t_field),POINTER :: field_glo(:) INTEGER, ALLOCATABLE :: mpi_req(:) INTEGER, ALLOCATABLE :: status(:,:) INTEGER :: ireq,nreq INTEGER :: ind_glo,ind_loc IF (.NOT. using_mpi) THEN DO ind_loc=1,ndomain field_glo(ind_loc)%rval4d=field_loc(ind_loc)%rval4d ENDDO ELSE nreq=ndomain IF (mpi_rank==0) nreq=nreq+ndomain_glo ALLOCATE(mpi_req(nreq)) ALLOCATE(status(MPI_STATUS_SIZE,nreq)) ireq=0 IF (mpi_rank==0) THEN DO ind_glo=1,ndomain_glo ireq=ireq+1 CALL MPI_IRECV(field_glo(ind_glo)%rval4d,size(field_glo(ind_glo)%rval4d) , MPI_REAL8 , & domglo_rank(ind_glo),domglo_loc_ind(ind_glo), comm_icosa, mpi_req(ireq), ierr) ENDDO ENDIF DO ind_loc=1,ndomain ireq=ireq+1 CALL MPI_ISEND(field_loc(ind_loc)%rval4d,size(field_loc(ind_loc)%rval4d) , MPI_REAL8 , & 0, ind_loc, comm_icosa, mpi_req(ireq), ierr) ENDDO CALL MPI_WAITALL(nreq,mpi_req,status,ierr) ENDIF END SUBROUTINE gather_field SUBROUTINE bcast_field(field_glo) USE field_mod USE domain_mod USE mpi_mod USE mpipara TYPE(t_field),POINTER :: field_glo(:) INTEGER :: ind_glo IF (using_mpi) THEN DO ind_glo=1,ndomain_glo CALL MPI_BCAST(field_glo(ind_glo)%rval4d, size(field_glo(ind_glo)%rval4d) , MPI_REAL8, 0, comm_icosa, ierr) ENDDO ENDIF END SUBROUTINE bcast_field SUBROUTINE scatter_field(field_glo,field_loc) USE field_mod USE domain_mod USE mpi_mod USE mpipara TYPE(t_field),POINTER :: field_glo(:) TYPE(t_field),POINTER :: field_loc(:) INTEGER, ALLOCATABLE :: mpi_req(:) INTEGER, ALLOCATABLE :: status(:,:) INTEGER :: ireq,nreq INTEGER :: ind_glo,ind_loc IF (.NOT. using_mpi) THEN DO ind_loc=1,ndomain field_loc(ind_loc)%rval4d=field_glo(ind_loc)%rval4d ENDDO ELSE nreq=ndomain IF (mpi_rank==0) nreq=nreq+ndomain_glo ALLOCATE(mpi_req(nreq)) ALLOCATE(status(MPI_STATUS_SIZE,nreq)) ireq=0 IF (mpi_rank==0) THEN DO ind_glo=1,ndomain_glo ireq=ireq+1 CALL MPI_ISEND(field_glo(ind_glo)%rval4d,size(field_glo(ind_glo)%rval4d) , MPI_REAL8 , & domglo_rank(ind_glo),domglo_loc_ind(ind_glo), comm_icosa, mpi_req(ireq), ierr) ENDDO ENDIF DO ind_loc=1,ndomain ireq=ireq+1 CALL MPI_IRECV(field_loc(ind_loc)%rval4d,size(field_loc(ind_loc)%rval4d) , MPI_REAL8 , & 0, ind_loc, comm_icosa, mpi_req(ireq), ierr) ENDDO CALL MPI_WAITALL(nreq,mpi_req,status,ierr) ENDIF END SUBROUTINE scatter_field !=================== Broadcast routines for strings ==================! SUBROUTINE bcast_mpi_cgen(var,nb) USE mpi_mod USE mpipara CHARACTER(LEN=*),INTENT(INOUT) :: Var INTEGER,INTENT(IN) :: nb IF (.NOT. using_mpi) RETURN CALL MPI_BCAST(Var,nb,MPI_CHARACTER,mpi_master,comm_icosa,ierr) END SUBROUTINE bcast_mpi_cgen SUBROUTINE bcast_mpi_c(var1) CHARACTER(LEN=*),INTENT(INOUT) :: Var1 CALL bcast_mpi_cgen(Var1,len(Var1)) END SUBROUTINE bcast_mpi_c !============= Broadcast routines for INTEGER scalars and arrays ============! SUBROUTINE bcast_mpi_igen(var,nb) USE mpi_mod USE mpipara INTEGER, INTENT(IN) :: nb INTEGER, DIMENSION(nb), INTENT(INOUT) :: var IF (using_mpi) CALL MPI_BCAST(Var, nb, MPI_INTEGER, mpi_master, comm_icosa, ierr) END SUBROUTINE bcast_mpi_igen SUBROUTINE bcast_mpi_i(var) USE mpipara INTEGER, INTENT(INOUT) :: var INTEGER :: var_tmp(1) IF (is_mpi_master) var_tmp(1)=var CALL bcast_mpi_igen(var_tmp,1) var=var_tmp(1) END SUBROUTINE bcast_mpi_i SUBROUTINE bcast_mpi_i1(var) INTEGER, INTENT(INOUT) :: var(:) CALL bcast_mpi_igen(var,size(var)) END SUBROUTINE bcast_mpi_i1 SUBROUTINE bcast_mpi_i2(var) INTEGER, INTENT(INOUT) :: var(:,:) CALL bcast_mpi_igen(var,size(var)) END SUBROUTINE bcast_mpi_i2 SUBROUTINE bcast_mpi_i3(var) INTEGER, INTENT(INOUT) :: var(:,:,:) CALL bcast_mpi_igen(var,size(var)) END SUBROUTINE bcast_mpi_i3 SUBROUTINE bcast_mpi_i4(var) INTEGER, INTENT(INOUT) :: var(:,:,:,:) CALL bcast_mpi_igen(var,size(var)) END SUBROUTINE bcast_mpi_i4 !============= Broadcast routines for REAL scalars and arrays ============! SUBROUTINE bcast_mpi_rgen(var,nb) USE mpi_mod USE mpipara INTEGER, INTENT(IN) :: nb REAL, DIMENSION(nb), INTENT(INOUT) :: var IF (using_mpi) CALL MPI_BCAST(Var, nb, MPI_REAL8, mpi_master, comm_icosa, ierr) END SUBROUTINE bcast_mpi_rgen SUBROUTINE bcast_mpi_r(var) USE mpipara REAL, INTENT(INOUT) :: var REAL :: var_tmp(1) IF (is_mpi_master) var_tmp(1)=var CALL bcast_mpi_rgen(var_tmp,1) var=var_tmp(1) END SUBROUTINE bcast_mpi_r SUBROUTINE bcast_mpi_r1(var) REAL, INTENT(INOUT) :: var(:) CALL bcast_mpi_rgen(var,size(var)) END SUBROUTINE bcast_mpi_r1 SUBROUTINE bcast_mpi_r2(var) REAL, INTENT(INOUT) :: var(:,:) CALL bcast_mpi_rgen(var,size(var)) END SUBROUTINE bcast_mpi_r2 SUBROUTINE bcast_mpi_r3(var) REAL, INTENT(INOUT) :: var(:,:,:) CALL bcast_mpi_rgen(var,size(var)) END SUBROUTINE bcast_mpi_r3 SUBROUTINE bcast_mpi_r4(var) REAL, INTENT(INOUT) :: var(:,:,:,:) CALL bcast_mpi_rgen(var,size(var)) END SUBROUTINE bcast_mpi_r4 !============= Broadcast routines for LOGICAL scalars and arrays ============! SUBROUTINE bcast_mpi_lgen(var,nb) USE mpi_mod USE mpipara INTEGER, INTENT(IN) :: nb LOGICAL, DIMENSION(nb), INTENT(INOUT) :: var IF (using_mpi) CALL MPI_BCAST(Var, nb, MPI_LOGICAL, mpi_master, comm_icosa, ierr) END SUBROUTINE bcast_mpi_lgen SUBROUTINE bcast_mpi_l(var) USE mpipara LOGICAL, INTENT(INOUT) :: var LOGICAL :: var_tmp(1) IF (is_mpi_master) var_tmp(1)=var CALL bcast_mpi_lgen(var_tmp,1) var=var_tmp(1) END SUBROUTINE bcast_mpi_l SUBROUTINE bcast_mpi_l1(var) LOGICAL, INTENT(INOUT) :: var(:) CALL bcast_mpi_lgen(var,size(var)) END SUBROUTINE bcast_mpi_l1 SUBROUTINE bcast_mpi_l2(var) LOGICAL, INTENT(INOUT) :: var(:,:) CALL bcast_mpi_lgen(var,size(var)) END SUBROUTINE bcast_mpi_l2 SUBROUTINE bcast_mpi_l3(var) LOGICAL, INTENT(INOUT) :: var(:,:,:) CALL bcast_mpi_lgen(var,size(var)) END SUBROUTINE bcast_mpi_l3 SUBROUTINE bcast_mpi_l4(var) LOGICAL, INTENT(INOUT) :: var(:,:,:,:) CALL bcast_mpi_lgen(var,size(var)) END SUBROUTINE bcast_mpi_l4 END MODULE transfert_mpi_collectives_mod