radiation_mpi_initialize Subroutine

public subroutine radiation_mpi_initialize(rad, comm_parent, n_tasks)

Initializes the radiation MPI context by duplicating a parent communicator.

Arguments

Type IntentOptional Attributes Name
type(radiation_mpi_t), intent(inout) :: rad

The radiation context to initialize.

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

The parent communicator (usually MPI_COMM_WORLD).

integer, intent(in), optional :: n_tasks

Optional total number of tasks to decompose immediately.


Calls

proc~~radiation_mpi_initialize~~CallsGraph proc~radiation_mpi_initialize mod_mpi_radiation::radiation_mpi_initialize mpi_comm_dup mpi_comm_dup proc~radiation_mpi_initialize->mpi_comm_dup mpi_comm_rank mpi_comm_rank proc~radiation_mpi_initialize->mpi_comm_rank mpi_comm_size mpi_comm_size proc~radiation_mpi_initialize->mpi_comm_size proc~check_mpi~4 mod_mpi_radiation::check_mpi proc~radiation_mpi_initialize->proc~check_mpi~4 proc~radiation_mpi_finalize mod_mpi_radiation::radiation_mpi_finalize proc~radiation_mpi_initialize->proc~radiation_mpi_finalize proc~radiation_task_bounds mod_mpi_radiation::radiation_task_bounds proc~radiation_mpi_initialize->proc~radiation_task_bounds proc~fatal_error mod_kinds::fatal_error proc~check_mpi~4->proc~fatal_error proc~radiation_mpi_finalize->proc~check_mpi~4 mpi_comm_free mpi_comm_free proc~radiation_mpi_finalize->mpi_comm_free proc~radiation_task_bounds->proc~fatal_error

Called by

proc~~radiation_mpi_initialize~~CalledByGraph proc~radiation_mpi_initialize mod_mpi_radiation::radiation_mpi_initialize program~lowmach_react_hex lowmach_react_hex program~lowmach_react_hex->proc~radiation_mpi_initialize

Source Code

   subroutine radiation_mpi_initialize(rad, comm_parent, n_tasks)
      type(radiation_mpi_t), intent(inout) :: rad
      type(MPI_Comm), intent(in) :: comm_parent
      integer, intent(in), optional :: n_tasks

      integer :: ierr

      call radiation_mpi_finalize(rad)

      call MPI_Comm_dup(comm_parent, rad%comm, ierr)
      call check_mpi(ierr, 'MPI_Comm_dup radiation')

      call MPI_Comm_rank(rad%comm, rad%rank, ierr)
      call check_mpi(ierr, 'MPI_Comm_rank radiation')

      call MPI_Comm_size(rad%comm, rad%nprocs, ierr)
      call check_mpi(ierr, 'MPI_Comm_size radiation')

      if (present(n_tasks)) then
         call radiation_task_bounds(rad, n_tasks)
      end if
   end subroutine radiation_mpi_initialize