flow_mpi_initialize Subroutine

public subroutine flow_mpi_initialize(mesh, flow, comm_parent, max_gather_components)

Sets up domain decomposition for a given mesh.

Splits the total number of cells among available processors using a contiguous block decomposition.

Arguments

Type IntentOptional Attributes Name
type(mesh_t), intent(in) :: mesh

The mesh to decompose.

type(flow_mpi_t), intent(inout) :: flow

The MPI context to populate.

type(MPI_Comm), intent(in) :: comm_parent

The parent communicator (usually MPI_COMM_WORLD).

integer, intent(in), optional :: max_gather_components

Calls

proc~~flow_mpi_initialize~~CallsGraph proc~flow_mpi_initialize mod_mpi_flow::flow_mpi_initialize mpi_comm_dup mpi_comm_dup proc~flow_mpi_initialize->mpi_comm_dup mpi_comm_rank mpi_comm_rank proc~flow_mpi_initialize->mpi_comm_rank mpi_comm_size mpi_comm_size proc~flow_mpi_initialize->mpi_comm_size proc~check_mpi~2 mod_mpi_flow::check_mpi proc~flow_mpi_initialize->proc~check_mpi~2 proc~flow_mpi_finalize mod_mpi_flow::flow_mpi_finalize proc~flow_mpi_initialize->proc~flow_mpi_finalize proc~setup_cell_halo mod_mpi_flow::setup_cell_halo proc~flow_mpi_initialize->proc~setup_cell_halo proc~setup_cell_owners mod_mpi_flow::setup_cell_owners proc~flow_mpi_initialize->proc~setup_cell_owners proc~setup_face_halo mod_mpi_flow::setup_face_halo proc~flow_mpi_initialize->proc~setup_face_halo proc~setup_owned_faces mod_mpi_flow::setup_owned_faces proc~flow_mpi_initialize->proc~setup_owned_faces proc~setup_owned_gather mod_mpi_flow::setup_owned_gather proc~flow_mpi_initialize->proc~setup_owned_gather proc~fatal_error mod_kinds::fatal_error proc~check_mpi~2->proc~fatal_error proc~flow_mpi_finalize->proc~check_mpi~2 mpi_comm_free mpi_comm_free proc~flow_mpi_finalize->mpi_comm_free proc~mesh_neighbor_for_cell mod_mpi_flow::mesh_neighbor_for_cell proc~setup_cell_halo->proc~mesh_neighbor_for_cell proc~pack_rank_metadata mod_mpi_flow::pack_rank_metadata proc~setup_cell_halo->proc~pack_rank_metadata proc~prefix_counts mod_mpi_flow::prefix_counts proc~setup_cell_halo->proc~prefix_counts proc~setup_face_halo->proc~pack_rank_metadata proc~setup_face_halo->proc~prefix_counts proc~setup_owned_gather->proc~check_mpi~2 mpi_allgather mpi_allgather proc~setup_owned_gather->mpi_allgather proc~setup_owned_gather->proc~fatal_error

Called by

proc~~flow_mpi_initialize~~CalledByGraph proc~flow_mpi_initialize mod_mpi_flow::flow_mpi_initialize program~lowmach_react_hex lowmach_react_hex program~lowmach_react_hex->proc~flow_mpi_initialize

Source Code

   subroutine flow_mpi_initialize(mesh, flow, comm_parent, max_gather_components)
      type(mesh_t), intent(in) :: mesh
      type(flow_mpi_t), intent(inout) :: flow
      type(MPI_Comm), intent(in) :: comm_parent
      integer, intent(in), optional :: max_gather_components

      integer :: ierr, base, rem, c
      integer :: max_components

      call flow_mpi_finalize(flow)

      call MPI_Comm_dup(comm_parent, flow%comm, ierr)
      call check_mpi(ierr, 'MPI_Comm_dup flow')

      call MPI_Comm_rank(flow%comm, flow%rank, ierr)
      call check_mpi(ierr, 'MPI_Comm_rank flow')

      call MPI_Comm_size(flow%comm, flow%nprocs, ierr)
      call check_mpi(ierr, 'MPI_Comm_size flow')

      ! Contiguous cell range calculation:
      ! The total ncells are split into roughly equal blocks. If ncells is not 
      ! divisible by nprocs, the first 'rem' ranks get an extra cell to ensure 
      ! all cells are covered.
      base = mesh%ncells / flow%nprocs
      rem = mod(mesh%ncells, flow%nprocs)

      if (flow%rank < rem) then
         flow%nlocal = base + 1
         flow%first_cell = flow%rank * (base + 1) + 1
      else
         flow%nlocal = base
         flow%first_cell = rem * (base + 1) + (flow%rank - rem) * base + 1
      end if

      ! The global cell indices owned by this rank are [first_cell, last_cell].
      flow%last_cell = flow%first_cell + flow%nlocal - 1

      allocate(flow%owned(mesh%ncells))
      flow%owned = .false.

      do c = flow%first_cell, flow%last_cell
         if (c >= 1 .and. c <= mesh%ncells) flow%owned(c) = .true.
      end do

      max_components = 4
      if (present(max_gather_components)) max_components = max(max_components, max_gather_components)
      flow%cell_halo_max_components = max(1, max_components)

      call setup_owned_gather(mesh, flow, max_gather_components)
      call setup_cell_owners(mesh, flow)
      call setup_owned_faces(mesh, flow)
      call setup_cell_halo(mesh, flow)
      call setup_face_halo(mesh, flow)
   end subroutine flow_mpi_initialize