flow_exchange_cell_scalars Subroutine

public subroutine flow_exchange_cell_scalars(flow, field1, field2, field3, field4)

Exchanges two to four related cell scalars in one halo transaction.

This avoids allocating a temporary matrix at call sites such as energy and transport, while keeping the same one-message-per-neighbor behavior as flow_exchange_cell_matrix. The separate radiation communicator is not touched; this routine operates only on the hydrodynamic flow communicator.

Arguments

Type IntentOptional Attributes Name
type(flow_mpi_t), intent(inout) :: flow
real(kind=rk), intent(inout) :: field1(:)
real(kind=rk), intent(inout), optional :: field2(:)
real(kind=rk), intent(inout), optional :: field3(:)
real(kind=rk), intent(inout), optional :: field4(:)

Source Code

   subroutine flow_exchange_cell_scalars(flow, field1, field2, field3, field4)
      use mod_profiling, only : profiler_start, profiler_stop
      type(flow_mpi_t), intent(inout) :: flow
      real(rk), intent(inout) :: field1(:)
      real(rk), intent(inout), optional :: field2(:), field3(:), field4(:)
      integer, parameter :: cell_scalars_halo_tag = 9284
      integer :: ncomp, i, j, nreq, ierr, offset, count, pos, cell

      if (.not. allocated(flow%cell_sendbuf)) return
      if (flow%ncell_recv_ranks + flow%ncell_send_ranks == 0) return

      ncomp = 1
      if (present(field2)) ncomp = ncomp + 1
      if (present(field3)) ncomp = ncomp + 1
      if (present(field4)) ncomp = ncomp + 1
      if (ncomp > flow%cell_halo_max_components) then
         call fatal_error('mpi_flow', 'cell scalar halo component count exceeds cached buffer size')
      end if

      do i = 1, flow%ncell_send_ranks
         offset = flow%cell_send_displs(i)
         count = flow%cell_send_counts(i)
         do j = 1, count
            cell = flow%cell_send_cells(offset + j)
            pos = (offset + j - 1) * ncomp
            flow%cell_sendbuf(pos + 1) = field1(cell)
            if (present(field2)) flow%cell_sendbuf(pos + 2) = field2(cell)
            if (present(field3)) flow%cell_sendbuf(pos + 3) = field3(cell)
            if (present(field4)) flow%cell_sendbuf(pos + 4) = field4(cell)
         end do
      end do

      call profiler_start('MPI_Communication')
      nreq = 0
      do i = 1, flow%ncell_recv_ranks
         offset = flow%cell_recv_displs(i) * ncomp
         count = flow%cell_recv_counts(i) * ncomp
         nreq = nreq + 1
         call MPI_Irecv(flow%cell_recvbuf(offset + 1), count, MPI_DOUBLE_PRECISION, &
                        flow%cell_recv_ranks(i), cell_scalars_halo_tag, flow%comm, flow%cell_requests(nreq), ierr)
         call check_mpi(ierr, 'cell scalars halo irecv')
      end do
      do i = 1, flow%ncell_send_ranks
         offset = flow%cell_send_displs(i) * ncomp
         count = flow%cell_send_counts(i) * ncomp
         nreq = nreq + 1
         call MPI_Isend(flow%cell_sendbuf(offset + 1), count, MPI_DOUBLE_PRECISION, &
                        flow%cell_send_ranks(i), cell_scalars_halo_tag, flow%comm, flow%cell_requests(nreq), ierr)
         call check_mpi(ierr, 'cell scalars halo isend')
      end do
      call MPI_Waitall(nreq, flow%cell_requests(1:nreq), MPI_STATUSES_IGNORE, ierr)
      call check_mpi(ierr, 'cell scalars halo waitall')
      call profiler_stop('MPI_Communication')

      do i = 1, size(flow%cell_recv_cells)
         cell = flow%cell_recv_cells(i)
         pos = (i - 1) * ncomp
         field1(cell) = flow%cell_recvbuf(pos + 1)
         if (present(field2)) field2(cell) = flow%cell_recvbuf(pos + 2)
         if (present(field3)) field3(cell) = flow%cell_recvbuf(pos + 3)
         if (present(field4)) field4(cell) = flow%cell_recvbuf(pos + 4)
      end do
   end subroutine flow_exchange_cell_scalars