Skip to content
Snippets Groups Projects
Commit f85f7606 authored by Antoine Cyril David Hoffmann's avatar Antoine Cyril David Hoffmann
Browse files

Update ghosts cells exchange for parallel 2D

parent 639ccdd1
No related branches found
No related tags found
No related merge requests found
!Communicate two extra buffer zones needed for using left/right fd schemes in parallel
!gradients routines
SUBROUTINE update_ghosts_p
USE basic
USE fields
USE grid
USE ppinit
use prec_const
IMPLICIT NONE
INTEGER :: status(MPI_STATUS_SIZE)
complex(dp):: buff_snd_L( 1: 2,ijs_e:ije_e,ikrs:ikre,ikzs:ikze)
complex(dp):: buff_snd_R(-2:-1,ijs_e:ije_e,ikrs:ikre,ikzs:ikze)
complex(dp):: buff_rcv_R( 1: 2,ijs_e:ije_e,ikrs:ikre,ikzs:ikze)
complex(dp):: buff_rcv_L(-2:-1,ijs_e:ije_e,ikrs:ikre,ikzs:ikze)
!! Set up data to send to left neighbor
DO ij=ijs_e,ije_e
DO ikr=ikrs,ikre
DO ikz=ikzs,ikze
! send to left ipe + 1 moments
buff_snd_L(1,ij,ikr,ikz) = moments_e(ips,ij,ikr,ikz,updatetlevel)
! send to left ipe + 2 moments
buff_snd_L(2,ij,ikr,ikz) = moments_e(ips+1,ij,ikr,ikz,updatetlevel)
END DO
END DO
END DO
! Exchange data with left neighbor
CALL mpi_sendrecv(buff_snd_L, count, MPI_DOUBLE_PRECISION, left_neighbor, 0, &
buff_rcv_L, count, MPI_DOUBLE_PRECISION, source, 0, &
commx, status, ierr)
! Write received data
DO ij=ijs_e,ije_e
DO ikr=ikrs,ikre
DO ikz=ikzs,ikze
! ip - 1 moments
moments_e(ips-1,ij,ikr,ikz,updatetlevel) = buff_rcv_L(-1,ij,ikr,ikz)
! ip - 2 moments
moments_e(ips-2,ij,ikr,ikz,updatetlevel) = buff_rcv_L(-2,ij,ikr,ikz)
END DO
END DO
END DO
!! Set up data to send to right neighbors
DO ij=ijs_e,ije_e
DO ikr=ikrs,ikre
DO ikz=ikzs,ikze
! send to right ips - 1 moments
buff_snd_R(-1,ij,ikr,ikz) = moments_e(ipe,ij,ikr,ikz,updatetlevel)
! send to right ips - 2 moments
buff_snd_R(-2,ij,ikr,ikz) = moments_e(ipe-1,ij,ikr,ikz,updatetlevel)
END DO
END DO
END DO
! Exchange data with right neighbor
CALL mpi_sendrecv(buff_snd_R, count, MPI_DOUBLE_PRECISION, dest, 0, &
buff_rcv_R, count, MPI_DOUBLE_PRECISION, source, 0, &
commx, status, ierr)
! Write received data
DO ij=ijs_e,ije_e
DO ikr=ikrs,ikre
DO ikz=ikzs,ikze
! ipe + 1 moments
moments_e(ipe+1,ij,ikr,ikz,updatetlevel) = buff_rcv_R(1,ij,ikr,ikz)
! ip + 2 moments
moments_e(ipe+2,ij,ikr,ikz,updatetlevel) = buff_rcv_R(2,ij,ikr,ikz)
END DO
END DO
END DO
END SUBROUTINE update_ghosts_p
\ No newline at end of file
module ghosts
USE basic
USE fields, ONLY : moments_e, moments_i
USE grid
USE time_integration
IMPLICIT NONE
INTEGER :: status(MPI_STATUS_SIZE), source, dest, count, ipg
PUBLIC :: update_ghosts
CONTAINS
SUBROUTINE update_ghosts
CALL cpu_time(t0_comm)
IF (num_procs_p .GT. 1) THEN ! Do it only if we share the p
CALL MPI_BARRIER(MPI_COMM_WORLD,ierr)
CALL update_ghosts_p_e
CALL MPI_BARRIER(MPI_COMM_WORLD,ierr)
CALL update_ghosts_p_i
ENDIF
CALL cpu_time(t1_comm)
tc_comm = tc_comm + (t1_comm - t0_comm)
END SUBROUTINE update_ghosts
!Communicate p+1, p+2 moments to left neighboor and p-1, p-2 moments to right one
SUBROUTINE update_ghosts_p_e
IMPLICIT NONE
count = (ijeg_e-ijsg_e+1)*(ikre-ikrs+1)*(ikze-ikzs+1)
!!!!!!!!!!! Send ghost to right neighbour !!!!!!!!!!!!!!!!!!!!!!
CALL mpi_sendrecv(moments_e(ipe_e ,ijsg_e:ijeg_e,ikrs:ikre,ikzs:ikze,updatetlevel), count, MPI_DOUBLE_COMPLEX, nbr_R, 10, &
moments_e(ips_e-1,ijsg_e:ijeg_e,ikrs:ikre,ikzs:ikze,updatetlevel), count, MPI_DOUBLE_COMPLEX, nbr_L, 10, &
comm0, status, ierr)
CALL mpi_sendrecv(moments_e(ipe_e-1,ijsg_e:ijeg_e,ikrs:ikre,ikzs:ikze,updatetlevel), count, MPI_DOUBLE_COMPLEX, nbr_R, 11, &
moments_e(ips_e-2,ijsg_e:ijeg_e,ikrs:ikre,ikzs:ikze,updatetlevel), count, MPI_DOUBLE_COMPLEX, nbr_L, 11, &
comm0, status, ierr)
!!!!!!!!!!! Send ghost to left neighbour !!!!!!!!!!!!!!!!!!!!!!
CALL mpi_sendrecv(moments_e(ips_e ,ijsg_e:ijeg_e,ikrs:ikre,ikzs:ikze,updatetlevel), count, MPI_DOUBLE_COMPLEX, nbr_L, 12, &
moments_e(ipe_e+1,ijsg_e:ijeg_e,ikrs:ikre,ikzs:ikze,updatetlevel), count, MPI_DOUBLE_COMPLEX, nbr_R, 12, &
comm0, status, ierr)
CALL mpi_sendrecv(moments_e(ips_e+1,ijsg_e:ijeg_e,ikrs:ikre,ikzs:ikze,updatetlevel), count, MPI_DOUBLE_COMPLEX, nbr_L, 13, &
moments_e(ipe_e+2,ijsg_e:ijeg_e,ikrs:ikre,ikzs:ikze,updatetlevel), count, MPI_DOUBLE_COMPLEX, nbr_R, 13, &
comm0, status, ierr)
END SUBROUTINE update_ghosts_p_e
!Communicate p+1, p+2 moments to left neighboor and p-1, p-2 moments to right one
SUBROUTINE update_ghosts_p_i
IMPLICIT NONE
count = (ijeg_i-ijsg_i+1)*(ikre-ikrs+1)*(ikze-ikzs+1) ! Number of elements sent
!!!!!!!!!!! Send ghost to right neighbour !!!!!!!!!!!!!!!!!!!!!!
CALL mpi_sendrecv(moments_i(ipe_i ,ijsg_i:ijeg_i,ikrs:ikre,ikzs:ikze,updatetlevel), count, MPI_DOUBLE_COMPLEX, nbr_R, 14, &
moments_i(ips_i-1,ijsg_i:ijeg_i,ikrs:ikre,ikzs:ikze,updatetlevel), count, MPI_DOUBLE_COMPLEX, nbr_L, 14, &
comm0, status, ierr)
CALL mpi_sendrecv(moments_i(ipe_i-1,ijsg_i:ijeg_i,ikrs:ikre,ikzs:ikze,updatetlevel), count, MPI_DOUBLE_COMPLEX, nbr_R, 15, &
moments_i(ips_i-2,ijsg_i:ijeg_i,ikrs:ikre,ikzs:ikze,updatetlevel), count, MPI_DOUBLE_COMPLEX, nbr_L, 15, &
comm0, status, ierr)
!!!!!!!!!!! Send ghost to left neighbour !!!!!!!!!!!!!!!!!!!!!!
CALL mpi_cart_shift(comm0, 0, -1, source , dest , ierr)
CALL mpi_sendrecv(moments_i(ips_i ,ijsg_i:ijeg_i,ikrs:ikre,ikzs:ikze,updatetlevel), count, MPI_DOUBLE_COMPLEX, nbr_L, 16, &
moments_i(ipe_i+1,ijsg_i:ijeg_i,ikrs:ikre,ikzs:ikze,updatetlevel), count, MPI_DOUBLE_COMPLEX, nbr_R, 16, &
comm0, status, ierr)
CALL mpi_sendrecv(moments_i(ips_i+1,ijsg_i:ijeg_i,ikrs:ikre,ikzs:ikze,updatetlevel), count, MPI_DOUBLE_COMPLEX, nbr_L, 17, &
moments_i(ipe_i+2,ijsg_i:ijeg_i,ikrs:ikre,ikzs:ikze,updatetlevel), count, MPI_DOUBLE_COMPLEX, nbr_R, 17, &
comm0, status, ierr)
END SUBROUTINE update_ghosts_p_i
END MODULE ghosts
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment