MPISEND_RECV call to exchange interface infromation between connected blocks.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=wp), | dimension(-2:dims%imx+2,-2:dims%jmx+2,-2:dims%kmx+2,1:dims%n_var) | :: | qp | Store primitive variable at cell center |
||
type(controltype), | intent(in) | :: | control | Control parameters |
||
type(boundarytype), | intent(in) | :: | bc | boundary conditions and fixed values |
||
type(extent), | intent(in) | :: | dims | Extent of the domain:imx,jmx,kmx |
subroutine apply_interface(qp, control, bc, dims)
!< MPISEND_RECV call to exchange interface infromation between
!< connected blocks.
implicit none
type(controltype), intent(in) :: control
!< Control parameters
type(extent), intent(in) :: dims
!< Extent of the domain:imx,jmx,kmx
real(wp), dimension(-2:dims%imx+2,-2:dims%jmx+2,-2:dims%kmx+2,1:dims%n_var) :: qp
!< Store primitive variable at cell center
type(boundarytype), intent(in) :: bc
!< boundary conditions and fixed values
integer:: i,j,k,n,l
integer:: status(MPI_STATUS_SIZE)
integer:: ierr
integer:: tag=1
integer:: count=0
integer :: imx, jmx, kmx, n_var
imx = dims%imx
jmx = dims%jmx
kmx = dims%kmx
n_var = control%n_var
!----------------------------------------------------------
! call pattern is change for first block = 0
! to avoid O-Grid infinite loop for mpi communication call
!-----------------------------------------------------------
if(mod(control%process_id,2)==0)then
!--- IMAX ---!
if(bc%imax_id>=0)then
!collect data
count=0
do n=1,n_var
do l=1,layers
do k=1,kmx-1
do j=1,jmx-1
count=count+1
imax_send_buf(count) = qp(imx-l,j,k,n)
end do
end do
end do
end do
call MPI_SENDRECV(imax_send_buf,ibuf_size, MPI_DOUBLE_PRECISION, bc%imax_id,tag,&
imax_recv_buf,ibuf_size, MPI_DOUBLE_PRECISION, bc%imax_id,tag,&
MPI_COMM_WORLD,status,ierr)
! redistribute data
if(bc%dir_switch(2)==0)then
count=0
do n=1,n_var
do l=1,layers
do k=Pklo(2),Pkhi(2),PkDir(2)
do j=Pjlo(2),Pjhi(2),PjDir(2)
count=count+1
qp(imx+l-1,j,k,n) = imax_recv_buf(count)
end do
end do
end do
end do
else
count=0
do n=1,n_var
do l=1,layers
do j=Pjlo(2),Pjhi(2),Pjdir(2)
do k=Pklo(2),Pkhi(2),PkDir(2)
count=count+1
qp(imx+l-1,j,k,n) = imax_recv_buf(count)
end do
end do
end do
end do
end if
end if
!--- IMIN ---!
DebugCall('apply_interface')
if(bc%imin_id>=0)then
!collect data
count=0
do n=1,n_var
do l=1,layers
do k=1,kmx-1
do j=1,jmx-1
count=count+1
imin_send_buf(count) = qp(l,j,k,n)
end do
end do
end do
end do
call MPI_SENDRECV(imin_send_buf,ibuf_size, MPI_DOUBLE_PRECISION, bc%imin_id,tag,&
imin_recv_buf,ibuf_size, MPI_DOUBLE_PRECISION, bc%imin_id,tag,&
MPI_COMM_WORLD,status,ierr)
! redistribute data
if(bc%dir_switch(1)==0)then
count=0
do n=1,n_var
do l=1,layers
do k=Pklo(1),Pkhi(1),PkDir(1)
do j=Pjlo(1),Pjhi(1),PjDir(1)
count=count+1
qp(1-l,j,k,n) = imin_recv_buf(count)
end do
end do
end do
end do
else
count=0
do n=1,n_var
do l=1,layers
do j=Pjlo(1),Pjhi(1),PjDir(1)
do k=Pklo(1),Pkhi(1),PkDir(1)
count=count+1
qp(1-l,j,k,n) = imin_recv_buf(count)
end do
end do
end do
end do
end if
end if
else
!--- IMIN ---!
DebugCall('apply_interface')
if(bc%imin_id>=0)then
!collect data
count=0
do n=1,n_var
do l=1,layers
do k=1,kmx-1
do j=1,jmx-1
count=count+1
imin_send_buf(count) = qp(l,j,k,n)
end do
end do
end do
end do
call MPI_SENDRECV(imin_send_buf,ibuf_size, MPI_DOUBLE_PRECISION, bc%imin_id,tag,&
imin_recv_buf,ibuf_size, MPI_DOUBLE_PRECISION, bc%imin_id,tag,&
MPI_COMM_WORLD,status,ierr)
! redistribute data
if(bc%dir_switch(1)==0)then
count=0
do n=1,n_var
do l=1,layers
do k=Pklo(1),Pkhi(1),PkDir(1)
do j=Pjlo(1),Pjhi(1),PjDir(1)
count=count+1
qp(1-l,j,k,n) = imin_recv_buf(count)
end do
end do
end do
end do
else
count=0
do n=1,n_var
do l=1,layers
do j=Pjlo(1),Pjhi(1),PjDir(1)
do k=Pklo(1),Pkhi(1),PkDir(1)
count=count+1
qp(1-l,j,k,n) = imin_recv_buf(count)
end do
end do
end do
end do
end if
end if
!--- IMAX ---!
if(bc%imax_id>=0)then
!collect data
count=0
do n=1,n_var
do l=1,layers
do k=1,kmx-1
do j=1,jmx-1
count=count+1
imax_send_buf(count) = qp(imx-l,j,k,n)
end do
end do
end do
end do
call MPI_SENDRECV(imax_send_buf,ibuf_size, MPI_DOUBLE_PRECISION, bc%imax_id,tag,&
imax_recv_buf,ibuf_size, MPI_DOUBLE_PRECISION, bc%imax_id,tag,&
MPI_COMM_WORLD,status,ierr)
! redistribute data
if(bc%dir_switch(2)==0)then
count=0
do n=1,n_var
do l=1,layers
do k=Pklo(2),Pkhi(2),PkDir(2)
do j=Pjlo(2),Pjhi(2),PjDir(2)
count=count+1
qp(imx+l-1,j,k,n) = imax_recv_buf(count)
end do
end do
end do
end do
else
count=0
do n=1,n_var
do l=1,layers
do j=Pjlo(2),Pjhi(2),Pjdir(2)
do k=Pklo(2),Pkhi(2),PkDir(2)
count=count+1
qp(imx+l-1,j,k,n) = imax_recv_buf(count)
end do
end do
end do
end do
end if
end if
end if
!--- JMIN ---!
if(bc%jmin_id>=0)then
!collect data
count=0
do n=1,n_var
do l=1,layers
do k=1,kmx-1
do i=1,imx-1
count=count+1
jmin_send_buf(count) = qp(i,l,k,n)
end do
end do
end do
end do
call MPI_SENDRECV(jmin_send_buf,jbuf_size, MPI_DOUBLE_PRECISION, bc%jmin_id,tag,&
jmin_recv_buf,jbuf_size, MPI_DOUBLE_PRECISION, bc%jmin_id,tag,&
MPI_COMM_WORLD,status,ierr)
! redistribute data
if(bc%dir_switch(3)==0)then
count=0
do n=1,n_var
do l=1,layers
do k=Pklo(3),Pkhi(3),PkDir(3)
do i=Pilo(3),Pihi(3),PiDir(3)
count=count+1
qp(i,1-l,k,n) = jmin_recv_buf(count)
end do
end do
end do
end do
else
count=0
do n=1,n_var
do l=1,layers
do i=Pilo(3),Pihi(3),PiDir(3)
do k=Pklo(3),Pkhi(3),PkDir(3)
count=count+1
qp(i,1-l,k,n) = jmin_recv_buf(count)
end do
end do
end do
end do
end if
end if
!--- JMAX ---!
if(bc%jmax_id>=0)then
!collect data
count=0
do n=1,n_var
do l=1,layers
do k=1,kmx-1
do i=1,imx-1
count=count+1
jmax_send_buf(count) = qp(i,jmx-l,k,n)
end do
end do
end do
end do
call MPI_SENDRECV(jmax_send_buf,jbuf_size, MPI_DOUBLE_PRECISION, bc%jmax_id,tag,&
jmax_recv_buf,jbuf_size, MPI_DOUBLE_PRECISION, bc%jmax_id,tag,&
MPI_COMM_WORLD,status,ierr)
! redistribute data
if(bc%dir_switch(4)==0)then
count=0
do n=1,n_var
do l=1,layers
do k=Pklo(4),Pkhi(4),PkDir(4)
do i=Pilo(4),Pihi(4),PiDir(4)
count=count+1
qp(i,jmx+l-1,k,n) = jmax_recv_buf(count)
end do
end do
end do
end do
else
count=0
do n=1,n_var
do l=1,layers
do i=Pilo(4),Pihi(4),PiDir(4)
do k=Pklo(4),Pkhi(4),PkDir(4)
count=count+1
qp(i,jmx+l-1,k,n) = jmax_recv_buf(count)
end do
end do
end do
end do
end if
end if
!--- KMIN ---!
if(bc%kmin_id>=0)then
!collect data
count=0
do n=1,n_var
do l=1,layers
do j=1,jmx-1
do i=1,imx-1
count=count+1
kmin_send_buf(count) = qp(i,j,l,n)
end do
end do
end do
end do
call MPI_SENDRECV(kmin_send_buf,kbuf_size, MPI_DOUBLE_PRECISION, bc%kmin_id,tag,&
kmin_recv_buf,kbuf_size, MPI_DOUBLE_PRECISION, bc%kmin_id,tag,&
MPI_COMM_WORLD,status,ierr)
! redistribute data
if(bc%dir_switch(5)==0)then
count=0
do n=1,n_var
do l=1,layers
do j=Pjlo(5),Pjhi(5),PjDir(5)
do i=Pilo(5),Pihi(5),PiDir(5)
count=count+1
qp(i,j,1-l,n) = kmin_recv_buf(count)
end do
end do
end do
end do
else
count=0
do n=1,n_var
do l=1,layers
do i=Pilo(5),Pihi(5),PiDir(5)
do j=Pjlo(5),Pjhi(5),PjDir(5)
count=count+1
qp(i,j,1-l,n) = kmin_recv_buf(count)
end do
end do
end do
end do
end if
end if
!--- KMAX ---!
if(bc%kmax_id>=0)then
!collect data
count=0
do n=1,n_var
do l=1,layers
do j=1,jmx-1
do i=1,imx-1
count=count+1
kmax_send_buf(count) = qp(i,j,kmx-l,n)
end do
end do
end do
end do
call MPI_SENDRECV(kmax_send_buf,kbuf_size, MPI_DOUBLE_PRECISION, bc%kmax_id,tag,&
kmax_recv_buf,kbuf_size, MPI_DOUBLE_PRECISION, bc%kmax_id,tag,&
MPI_COMM_WORLD,status,ierr)
! redistribute data
if(bc%dir_switch(6)==0)then
count=0
do n=1,n_var
do l=1,layers
do j=Pjlo(6),Pjhi(6),PjDir(6)
do i=Pilo(6),Pihi(6),PiDir(6)
count=count+1
qp(i,j,kmx+l-1,n) = kmax_recv_buf(count)
end do
end do
end do
end do
else
count=0
do n=1,n_var
do l=1,layers
do i=Pilo(6),Pihi(6),PiDir(6)
do j=Pjlo(6),Pjhi(6),PjDir(6)
count=count+1
qp(i,j,kmx+l-1,n) = kmax_recv_buf(count)
end do
end do
end do
end do
end if
end if
call apply_periodic_bc(qp, control, bc, dims)
end subroutine apply_interface