flow_exchange_face_scalar Subroutine

public subroutine flow_exchange_face_scalar(flow, face_field)

Uses

  • proc~~flow_exchange_face_scalar~~UsesGraph proc~flow_exchange_face_scalar mod_mpi_flow::flow_exchange_face_scalar module~mod_profiler mod_profiler proc~flow_exchange_face_scalar->module~mod_profiler iso_fortran_env iso_fortran_env module~mod_profiler->iso_fortran_env module~mod_kinds mod_kinds module~mod_profiler->module~mod_kinds mpi_f08 mpi_f08 module~mod_profiler->mpi_f08 module~mod_kinds->iso_fortran_env

Exchanges owner-computed face scalar values to ranks owning the neighbor cell.

Arguments

Type IntentOptional Attributes Name
type(flow_mpi_t), intent(inout) :: flow
real(kind=rk), intent(inout) :: face_field(:)

Calls

proc~~flow_exchange_face_scalar~~CallsGraph proc~flow_exchange_face_scalar mod_mpi_flow::flow_exchange_face_scalar mpi_irecv mpi_irecv proc~flow_exchange_face_scalar->mpi_irecv mpi_isend mpi_isend proc~flow_exchange_face_scalar->mpi_isend mpi_waitall mpi_waitall proc~flow_exchange_face_scalar->mpi_waitall proc~check_mpi~2 mod_mpi_flow::check_mpi proc~flow_exchange_face_scalar->proc~check_mpi~2 proc~profiler_start mod_profiler::profiler_start proc~flow_exchange_face_scalar->proc~profiler_start proc~profiler_stop mod_profiler::profiler_stop proc~flow_exchange_face_scalar->proc~profiler_stop proc~fatal_error mod_kinds::fatal_error proc~check_mpi~2->proc~fatal_error mpi_wtime mpi_wtime proc~profiler_start->mpi_wtime proc~find_or_create_timer mod_profiler::find_or_create_timer proc~profiler_start->proc~find_or_create_timer proc~profiler_stop->mpi_wtime proc~profiler_stop->proc~find_or_create_timer proc~record_edge mod_profiler::record_edge proc~profiler_stop->proc~record_edge

Called by

proc~~flow_exchange_face_scalar~~CalledByGraph proc~flow_exchange_face_scalar mod_mpi_flow::flow_exchange_face_scalar proc~advance_projection_step mod_flow_projection::advance_projection_step proc~advance_projection_step->proc~flow_exchange_face_scalar program~lowmach_react_hex lowmach_react_hex program~lowmach_react_hex->proc~advance_projection_step

Source Code

   subroutine flow_exchange_face_scalar(flow, face_field)
      use mod_profiler, only : profiler_start, profiler_stop
      type(flow_mpi_t), intent(inout) :: flow
      real(rk), intent(inout) :: face_field(:)
      integer, parameter :: face_halo_tag = 9283
      integer :: i, j, nreq, ierr, offset, count

      if (.not. allocated(flow%face_sendbuf)) return
      if (flow%nface_recv_ranks + flow%nface_send_ranks == 0) return

      do i = 1, flow%nface_send_ranks
         offset = flow%face_send_displs(i)
         count = flow%face_send_counts(i)
         do j = 1, count
            flow%face_sendbuf(offset + j) = face_field(flow%face_send_faces(offset + j))
         end do
      end do

      call profiler_start('MPI_Communication')
      nreq = 0
      do i = 1, flow%nface_recv_ranks
         offset = flow%face_recv_displs(i)
         count = flow%face_recv_counts(i)
         nreq = nreq + 1
         call MPI_Irecv(flow%face_recvbuf(offset + 1), count, MPI_DOUBLE_PRECISION, &
                        flow%face_recv_ranks(i), face_halo_tag, flow%comm, flow%face_requests(nreq), ierr)
         call check_mpi(ierr, 'face halo irecv')
      end do
      do i = 1, flow%nface_send_ranks
         offset = flow%face_send_displs(i)
         count = flow%face_send_counts(i)
         nreq = nreq + 1
         call MPI_Isend(flow%face_sendbuf(offset + 1), count, MPI_DOUBLE_PRECISION, &
                        flow%face_send_ranks(i), face_halo_tag, flow%comm, flow%face_requests(nreq), ierr)
         call check_mpi(ierr, 'face halo isend')
      end do
      call MPI_Waitall(nreq, flow%face_requests(1:nreq), MPI_STATUSES_IGNORE, ierr)
      call check_mpi(ierr, 'face halo waitall')
      call profiler_stop('MPI_Communication')

      do i = 1, size(flow%face_recv_faces)
         face_field(flow%face_recv_faces(i)) = flow%face_recvbuf(i)
      end do
   end subroutine flow_exchange_face_scalar