flow_gather_owned_scalar_root Subroutine

public subroutine flow_gather_owned_scalar_root(flow, field, root_field)

Uses

  • proc~~flow_gather_owned_scalar_root~~UsesGraph proc~flow_gather_owned_scalar_root mod_mpi_flow::flow_gather_owned_scalar_root module~mod_profiler mod_profiler proc~flow_gather_owned_scalar_root->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

Gathers owned scalar cell values to rank 0 only.

Arguments

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

Calls

proc~~flow_gather_owned_scalar_root~~CallsGraph proc~flow_gather_owned_scalar_root mod_mpi_flow::flow_gather_owned_scalar_root mpi_gatherv mpi_gatherv proc~flow_gather_owned_scalar_root->mpi_gatherv proc~check_mpi~2 mod_mpi_flow::check_mpi proc~flow_gather_owned_scalar_root->proc~check_mpi~2 proc~profiler_start mod_profiler::profiler_start proc~flow_gather_owned_scalar_root->proc~profiler_start proc~profiler_stop mod_profiler::profiler_stop proc~flow_gather_owned_scalar_root->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

Source Code

   subroutine flow_gather_owned_scalar_root(flow, field, root_field)
      use mod_profiler, only : profiler_start, profiler_stop
      type(flow_mpi_t), intent(inout) :: flow
      real(rk), intent(in) :: field(:)
      real(rk), intent(inout) :: root_field(:)
      integer :: ierr, r, first

      flow%gather_sendbuf = field(flow%first_cell:flow%last_cell)

      call profiler_start('MPI_Communication')
      call MPI_Gatherv(flow%gather_sendbuf, flow%nlocal, MPI_DOUBLE_PRECISION, &
                       flow%gather_recvbuf, flow%gather_counts, flow%gather_displs, &
                       MPI_DOUBLE_PRECISION, 0, flow%comm, ierr)
      call check_mpi(ierr, 'MPI_Gatherv owned scalar root')
      call profiler_stop('MPI_Communication')

      if (flow%rank == 0) then
         root_field = zero
         do r = 1, flow%nprocs
            first = flow%gather_firsts(r)
            root_field(first:first + flow%gather_counts(r) - 1) = &
               flow%gather_recvbuf(flow%gather_displs(r) + 1: &
                                   flow%gather_displs(r) + flow%gather_counts(r))
         end do
      end if
   end subroutine flow_gather_owned_scalar_root