excjs Subroutine

public subroutine excjs(a, sx, ex, sy, ey, sz, ez, ih, jh)

Arguments

Type IntentOptional Attributes Name
real :: a(sx-ih:ex+ih,sy-jh:ey+jh,sz:ez)
integer :: sx
integer :: ex
integer :: sy
integer :: ey
integer :: sz
integer :: ez
integer :: ih
integer :: jh

Calls

proc~~excjs~~CallsGraph proc~excjs modmpi::excjs mpi_sendrecv mpi_sendrecv proc~excjs->mpi_sendrecv

Called by

proc~~excjs~~CalledByGraph proc~excjs modmpi::excjs proc~bcpup modboundary::bcpup proc~bcpup->proc~excjs proc~closurebc modboundary::closurebc proc~closurebc->proc~excjs proc~cyclichj modboundary::cyclichj proc~cyclichj->proc~excjs proc~cyclicmj modboundary::cyclicmj proc~cyclicmj->proc~excjs proc~cyclicqj modboundary::cyclicqj proc~cyclicqj->proc~excjs proc~cyclicsj modboundary::cyclicsj proc~cyclicsj->proc~excjs proc~readinletfile modinlet::readinletfile proc~readinletfile->proc~excjs proc~tkestatsdump modstatsdump::tkestatsdump proc~tkestatsdump->proc~excjs proc~boundary modboundary::boundary proc~boundary->proc~cyclichj proc~boundary->proc~cyclicmj proc~boundary->proc~cyclicqj proc~boundary->proc~cyclicsj proc~inletgen modinlet::inletgen proc~boundary->proc~inletgen proc~inletgennotemp modinlet::inletgennotemp proc~boundary->proc~inletgennotemp proc~closure modsubgrid::closure proc~closure->proc~closurebc proc~fillps modpois::fillps proc~fillps->proc~bcpup proc~inletgen->proc~readinletfile proc~inletgennotemp->proc~readinletfile proc~readinitfiles modstartup::readinitfiles proc~readinitfiles->proc~readinletfile proc~readinitfiles->proc~boundary proc~statsdump modstatsdump::statsdump proc~statsdump->proc~tkestatsdump proc~poisson modpois::poisson proc~poisson->proc~fillps proc~startup modstartup::startup proc~startup->proc~readinitfiles proc~subgrid modsubgrid::subgrid proc~subgrid->proc~closure program~dalesurban DALESURBAN program~dalesurban->proc~boundary program~dalesurban->proc~statsdump program~dalesurban->proc~poisson program~dalesurban->proc~startup program~dalesurban->proc~subgrid

Contents

Source Code


Source Code

  subroutine excjs( a, sx, ex, sy, ey, sz,ez,ih,jh)
    implicit none
  integer sx, ex, sy, ey, sz, ez, ih, jh
  real a(sx-ih:ex+ih, sy-jh:ey+jh, sz:ez)
  integer status(MPI_STATUS_SIZE), iiget
  integer ii, i, j, k
  real,allocatable, dimension(:) :: buffj1,buffj2,buffj3,buffj4
  iiget = jh*(ex - sx + 1 + 2*ih)*(ez - sz + 1)

  allocate( buffj1(iiget),&
            buffj2(iiget),&
            buffj3(iiget),&
            buffj4(iiget))

  if(nbrtop/=MPI_PROC_NULL)then
    ii = 0
    do j=1,jh
    do k=sz,ez
    do i=sx-ih,ex+ih
      ii = ii + 1
      buffj1(ii) = a(i,ey-j+1,k) ! tg3315 buffj1 is je-jhc ghost cells on myid
    enddo
    enddo
    enddo
  endif

  call MPI_SENDRECV(  buffj1,  ii    , MY_REAL, nbrtop, 4, &
                           buffj2,  iiget , MY_REAL, nbrbottom,  4, &
                           comm3d, status, mpierr )

  ! tg3315 sends this to nbrtop and pulls buffj2 from nbrbottom! send and receive process that is good for executing chain shifts.

  if(nbrbottom/=MPI_PROC_NULL)then
    ii = 0
    do j=1,jh
    do k=sz,ez
    do i=sx-ih,ex+ih
      ii = ii + 1
      a(i,sy-j,k) = buffj2(ii) ! set the previous ghost cells to buffj2 (last cells of nbrbottom I think)
    enddo
    enddo
    enddo
  endif

!   call barrou()

  ! repeats this process for other way round

  if(nbrbottom/=MPI_PROC_NULL)then
    ii = 0
    do j=1,jh
    do k=sz,ez
    do i=sx-ih,ex+ih
      ii = ii + 1
      buffj3(ii) = a(i,sy+j-1,k)
    enddo
    enddo
    enddo
  endif
  call MPI_SENDRECV(  buffj3,  ii    , MY_REAL, nbrbottom,  5, &
                          buffj4,  iiget , MY_REAL, nbrtop, 5, &
                          comm3d, status, mpierr )
  if(nbrtop/=MPI_PROC_NULL)then
    ii = 0
    do j=1,jh
    do k=sz,ez
    do i=sx-ih,ex+ih
      ii = ii + 1
      a(i,ey+j,k) = buffj4(ii)
    enddo
    enddo
    enddo
  endif

!   call barrou()
  deallocate (buffj1,buffj2,buffj3,buffj4)

  return
  end subroutine excjs