MODULE transfert_mpi_collectives_mod IMPLICIT NONE INTERFACE bcast_mpi MODULE PROCEDURE bcast_mpi_c, & bcast_mpi_i,bcast_mpi_i1,bcast_mpi_i2,bcast_mpi_i3,bcast_mpi_i4, & bcast_mpi_r,bcast_mpi_r1,bcast_mpi_r2,bcast_mpi_r3,bcast_mpi_r4, & 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 IMPLICIT NONE 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 IF (field_loc(ind_loc)%ndim==2) field_glo(ind_loc)%rval2d=field_loc(ind_loc)%rval2d IF (field_loc(ind_loc)%ndim==3) field_glo(ind_loc)%rval3d=field_loc(ind_loc)%rval3d IF (field_loc(ind_loc)%ndim==4) 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 IF (field_glo(ind_glo)%ndim==2) THEN CALL MPI_IRECV(field_glo(ind_glo)%rval2d,size(field_glo(ind_glo)%rval2d) , MPI_REAL8 , & domglo_rank(ind_glo),domglo_loc_ind(ind_glo), comm_icosa, mpi_req(ireq), ierr) ELSE IF (field_glo(ind_glo)%ndim==3) THEN CALL MPI_IRECV(field_glo(ind_glo)%rval3d,size(field_glo(ind_glo)%rval3d) , MPI_REAL8 , & domglo_rank(ind_glo),domglo_loc_ind(ind_glo), comm_icosa, mpi_req(ireq), ierr) ELSE IF (field_glo(ind_glo)%ndim==4) THEN 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) ENDIF ENDDO ENDIF DO ind_loc=1,ndomain ireq=ireq+1 IF (field_loc(ind_loc)%ndim==2) THEN CALL MPI_ISEND(field_loc(ind_loc)%rval2d,size(field_loc(ind_loc)%rval2d) , MPI_REAL8 , & 0, ind_loc, comm_icosa, mpi_req(ireq), ierr) ELSE IF (field_loc(ind_loc)%ndim==3) THEN CALL MPI_ISEND(field_loc(ind_loc)%rval3d,size(field_loc(ind_loc)%rval3d) , MPI_REAL8 , & 0, ind_loc, comm_icosa, mpi_req(ireq), ierr) ELSE IF (field_loc(ind_loc)%ndim==4) THEN 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) ENDIF 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 IMPLICIT NONE TYPE(t_field),POINTER :: field_glo(:) INTEGER :: ind_glo IF (.NOT. using_mpi) THEN ! nothing to do ELSE DO ind_glo=1,ndomain_glo IF (field_glo(ind_glo)%ndim==2) THEN CALL MPI_BCAST(field_glo(ind_glo)%rval2d, size(field_glo(ind_glo)%rval2d) , MPI_REAL8, 0, comm_icosa, ierr) ELSE IF (field_glo(ind_glo)%ndim==3) THEN CALL MPI_BCAST(field_glo(ind_glo)%rval3d, size(field_glo(ind_glo)%rval3d) , MPI_REAL8, 0, comm_icosa, ierr) ELSE IF (field_glo(ind_glo)%ndim==4) THEN CALL MPI_BCAST(field_glo(ind_glo)%rval4d, size(field_glo(ind_glo)%rval4d) , MPI_REAL8, 0, comm_icosa, ierr) ENDIF ENDDO ENDIF END SUBROUTINE bcast_field SUBROUTINE scatter_field(field_glo,field_loc) USE field_mod USE domain_mod USE mpi_mod USE mpipara IMPLICIT NONE 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 IF (field_loc(ind_loc)%ndim==2) field_loc(ind_loc)%rval2d=field_glo(ind_loc)%rval2d IF (field_loc(ind_loc)%ndim==3) field_loc(ind_loc)%rval3d=field_glo(ind_loc)%rval3d IF (field_loc(ind_loc)%ndim==4) 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 IF (field_glo(ind_glo)%ndim==2) THEN CALL MPI_ISEND(field_glo(ind_glo)%rval2d,size(field_glo(ind_glo)%rval2d) , MPI_REAL8 , & domglo_rank(ind_glo),domglo_loc_ind(ind_glo), comm_icosa, mpi_req(ireq), ierr) ELSE IF (field_glo(ind_glo)%ndim==3) THEN CALL MPI_ISEND(field_glo(ind_glo)%rval3d,size(field_glo(ind_glo)%rval3d) , MPI_REAL8 , & domglo_rank(ind_glo),domglo_loc_ind(ind_glo), comm_icosa, mpi_req(ireq), ierr) ELSE IF (field_glo(ind_glo)%ndim==4) THEN 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) ENDIF ENDDO ENDIF DO ind_loc=1,ndomain ireq=ireq+1 IF (field_loc(ind_loc)%ndim==2) THEN CALL MPI_IRECV(field_loc(ind_loc)%rval2d,size(field_loc(ind_loc)%rval2d) , MPI_REAL8 , & 0, ind_loc, comm_icosa, mpi_req(ireq), ierr) ELSE IF (field_loc(ind_loc)%ndim==3) THEN CALL MPI_IRECV(field_loc(ind_loc)%rval3d,size(field_loc(ind_loc)%rval3d) , MPI_REAL8 , & 0, ind_loc, comm_icosa, mpi_req(ireq), ierr) ELSE IF (field_loc(ind_loc)%ndim==4) THEN 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) ENDIF ENDDO CALL MPI_WAITALL(nreq,mpi_req,status,ierr) ENDIF END SUBROUTINE scatter_field !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! !! Definition des Broadcast --> 4D !! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! !! -- Les chaine de charact�re -- !! SUBROUTINE bcast_mpi_c(var1) IMPLICIT NONE CHARACTER(LEN=*),INTENT(INOUT) :: Var1 CALL bcast_mpi_cgen(Var1,len(Var1)) END SUBROUTINE bcast_mpi_c !! -- Les entiers -- !! SUBROUTINE bcast_mpi_i(var) USE mpipara IMPLICIT NONE 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) IMPLICIT NONE INTEGER,INTENT(INOUT) :: Var(:) CALL bcast_mpi_igen(Var,size(Var)) END SUBROUTINE bcast_mpi_i1 SUBROUTINE bcast_mpi_i2(var) IMPLICIT NONE INTEGER,INTENT(INOUT) :: Var(:,:) CALL bcast_mpi_igen(Var,size(Var)) END SUBROUTINE bcast_mpi_i2 SUBROUTINE bcast_mpi_i3(var) IMPLICIT NONE INTEGER,INTENT(INOUT) :: Var(:,:,:) CALL bcast_mpi_igen(Var,size(Var)) END SUBROUTINE bcast_mpi_i3 SUBROUTINE bcast_mpi_i4(var) IMPLICIT NONE INTEGER,INTENT(INOUT) :: Var(:,:,:,:) CALL bcast_mpi_igen(Var,size(Var)) END SUBROUTINE bcast_mpi_i4 !! -- Les reels -- !! SUBROUTINE bcast_mpi_r(var) USE mpipara IMPLICIT NONE 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) IMPLICIT NONE REAL,INTENT(INOUT) :: Var(:) CALL bcast_mpi_rgen(Var,size(Var)) END SUBROUTINE bcast_mpi_r1 SUBROUTINE bcast_mpi_r2(var) IMPLICIT NONE REAL,INTENT(INOUT) :: Var(:,:) CALL bcast_mpi_rgen(Var,size(Var)) END SUBROUTINE bcast_mpi_r2 SUBROUTINE bcast_mpi_r3(var) IMPLICIT NONE REAL,INTENT(INOUT) :: Var(:,:,:) CALL bcast_mpi_rgen(Var,size(Var)) END SUBROUTINE bcast_mpi_r3 SUBROUTINE bcast_mpi_r4(var) IMPLICIT NONE REAL,INTENT(INOUT) :: Var(:,:,:,:) CALL bcast_mpi_rgen(Var,size(Var)) END SUBROUTINE bcast_mpi_r4 !! -- Les booleans -- !! SUBROUTINE bcast_mpi_l(var) USE mpipara IMPLICIT NONE 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) IMPLICIT NONE LOGICAL,INTENT(INOUT) :: Var(:) CALL bcast_mpi_lgen(Var,size(Var)) END SUBROUTINE bcast_mpi_l1 SUBROUTINE bcast_mpi_l2(var) IMPLICIT NONE LOGICAL,INTENT(INOUT) :: Var(:,:) CALL bcast_mpi_lgen(Var,size(Var)) END SUBROUTINE bcast_mpi_l2 SUBROUTINE bcast_mpi_l3(var) IMPLICIT NONE LOGICAL,INTENT(INOUT) :: Var(:,:,:) CALL bcast_mpi_lgen(Var,size(Var)) END SUBROUTINE bcast_mpi_l3 SUBROUTINE bcast_mpi_l4(var) IMPLICIT NONE LOGICAL,INTENT(INOUT) :: Var(:,:,:,:) CALL bcast_mpi_lgen(Var,size(Var)) END SUBROUTINE bcast_mpi_l4 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! !! DEFINITION DES FONCTIONS DE TRANSFERT GENERIQUES ! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! SUBROUTINE bcast_mpi_cgen(var,nb) USE mpi_mod USE mpipara IMPLICIT NONE 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_igen(var,nb) USE mpi_mod USE mpipara IMPLICIT NONE INTEGER,INTENT(IN) :: nb INTEGER,DIMENSION(nb),INTENT(INOUT) :: Var IF (.NOT. using_mpi) RETURN CALL MPI_BCAST(Var,nb,MPI_INTEGER,mpi_master,comm_icosa,ierr) END SUBROUTINE bcast_mpi_igen SUBROUTINE bcast_mpi_rgen(var,nb) USE mpi_mod USE mpipara IMPLICIT NONE INTEGER,INTENT(IN) :: nb REAL,DIMENSION(nb),INTENT(INOUT) :: Var IF (.NOT. using_mpi) RETURN CALL MPI_BCAST(Var,nb,MPI_REAL8,mpi_master,comm_icosa,ierr) END SUBROUTINE bcast_mpi_rgen SUBROUTINE bcast_mpi_lgen(var,nb) USE mpi_mod USE mpipara IMPLICIT NONE INTEGER,INTENT(IN) :: nb LOGICAL,DIMENSION(nb),INTENT(INOUT) :: Var IF (.NOT. using_mpi) RETURN CALL MPI_BCAST(Var,nb,MPI_LOGICAL,mpi_master,comm_icosa,ierr) END SUBROUTINE bcast_mpi_lgen END MODULE transfert_mpi_collectives_mod