Запуск MPI-задач

Запуск MPI-задач #

Для запуска MPI-задач можно использовать команду mpirun. По умолчанию задача будет запущена с заказанным количеством процессов, информация об этом будет взята из системы очередей Slurm (параметр --ntasks). Если нужно запустить меньшее количество процессов, то можно использовать параметр -n np.

Например:

mpirun -n 8 ./mpiprog8.exe

Intel MPI #

Подключите модули Intel и используйте компиляторы mpiicc, mpicpc и mpiifort:

module load intel impi
mpiicc -o main main.c

Open MPI #

Вы можете использовать Open MPI с компиляторами gcc, для этого включите необходимые модули и используйте компиляторы mpicc, mpicxx и mpifort:

# компилятор gcc-12
module load gnu12 openmpi4
# компилятор gcc-9
# module load gnu9 openmpi4
# компилятор gcc-7
# module load gnu openmpi
mpicc -o main main.c

Запуск гибридных задач MPI+OpenMP #

Рассмотрим наглядный пример гибридной программы xthi. Код программы можно найти, например, на сайте NERSC.

#define _GNU_SOURCE

#include <stdio.h>
#include <unistd.h>
#include <string.h>
#include <sched.h>
#include <mpi.h>
#include <omp.h>
#include <sys/syscall.h>

/* Borrowed from util-linux-2.13-pre7/schedutils/taskset.c */
static char *cpuset_to_cstr(cpu_set_t *mask, char *str)
{
  char *ptr = str;
  int i, j, entry_made = 0;
  for (i = 0; i < CPU_SETSIZE; i++) {
    if (CPU_ISSET(i, mask)) {
      int run = 0;
      entry_made = 1;
      for (j = i + 1; j < CPU_SETSIZE; j++) {
        if (CPU_ISSET(j, mask)) run++;
        else break;
      }
      if (!run)
        sprintf(ptr, "%d,", i);
      else if (run == 1) {
        sprintf(ptr, "%d,%d,", i, i + 1);
        i++;
      } else {
        sprintf(ptr, "%d-%d,", i, i + run);
        i += run;
      }
      while (*ptr != 0) ptr++;
    }
  }
  ptr -= entry_made;
  *ptr = 0;
  return(str);
}

int main(int argc, char *argv[])
{
  int rank, thread;
  cpu_set_t coremask;
  char clbuf[7 * CPU_SETSIZE], hnbuf[64];

  MPI_Init(&argc, &argv);
  MPI_Comm_rank(MPI_COMM_WORLD, &rank);
  memset(clbuf, 0, sizeof(clbuf));
  memset(hnbuf, 0, sizeof(hnbuf));
  (void)gethostname(hnbuf, sizeof(hnbuf));
  #pragma omp parallel private(thread, coremask, clbuf)
  {
    thread = omp_get_thread_num();
    (void)sched_getaffinity(0, sizeof(coremask), &coremask);
    cpuset_to_cstr(&coremask, clbuf);
    #pragma omp barrier
    printf("Hello from rank %d, thread %d, on %s. (core affinity = %s)\n",
            rank, thread, hnbuf, clbuf);
  }
/*  sleep(60); */
  MPI_Finalize();
  return(0);
}

Для компиляции с помощью Intel MPI выполните следующие команды:

# загружаем модули Intel MPI
module load intel impi
# используем компилятор mpiicc и добавляем параметр -qopenmp
mpiicc -qopenmp xthi.c -o xthi

Для компиляции с помощью gcc и Open MPI выполните следующие команды:

# загружаем модули gnu и openmpi
module load gnu12 openmpi4
# используем компилятор mpicc и добавляем параметр -fopenmp
mpicc -fopenmp xthi.c -o xthi

Для запуска задачи воспользуемся следующим скриптом test-xthi.sh:

#!/bin/bash

# A hybrid MPI+OpenMP example xthi
#SBATCH --partition=short
#SBATCH --job-name=xthi
#SBATCH --output=xthi.out
#SBATCH --time=5

####### 4 MPI ranks
#SBATCH --ntasks=4

####### 8 OMP threads per MPI rank
#SBATCH --cpus-per-task=8

## load MPI libs if needed
## for Intel MPI:
# module load intel impi

## for Open MPI:
# module load gnu12 openmpi4


## explicitly set number of OMP threads to the value of --cpus-per-task
## (required for mpich and mvapich2)
# export OMP_NUM_THREADS=$SLURM_CPUS_PER_TASK

## run with mpirun
mpirun ./xthi

## Open MPI requires additional parameter
# mpirun --map-by slot:PE=$SLURM_CPUS_PER_TASK ./xthi

Система очередей сохраняет текущие переменные окружения, поэтому в самом скрипте активировать модули для MPI не обязательно, если они уже активированы в текущей сессии.

Переменную OMP_NUM_THREADS можно явно приравнять значению $SLURM_CPUS_PER_TASK. В некоторых случаях в этом нет необходимости, так как каждый MPI процесс будет запущен на заданном количестве ядер, и OpenMP автоматически выберет нужное количество потоков.

При использовании Open MPI нужно добавить параметр --map-by slot:PE=$SLURM_CPUS_PER_TASK, иначе независимо от значения OMP_NUM_THREADS все нити OpenMP будут выполняться на одном ядре.

Ставим задачу в очередь:

sbatch test-xthi.sh

Результат работы будет в файле xthi.out, пример:

Hello from rank 3, thread 0, on n38. (core affinity = 0-7)
Hello from rank 3, thread 1, on n38. (core affinity = 0-7)
Hello from rank 3, thread 7, on n38. (core affinity = 0-7)
Hello from rank 3, thread 3, on n38. (core affinity = 0-7)
Hello from rank 3, thread 4, on n38. (core affinity = 0-7)
Hello from rank 3, thread 5, on n38. (core affinity = 0-7)
Hello from rank 3, thread 2, on n38. (core affinity = 0-7)
Hello from rank 3, thread 6, on n38. (core affinity = 0-7)
Hello from rank 2, thread 0, on n37. (core affinity = 4-11)
Hello from rank 2, thread 2, on n37. (core affinity = 4-11)
Hello from rank 2, thread 3, on n37. (core affinity = 4-11)
Hello from rank 2, thread 4, on n37. (core affinity = 4-11)
Hello from rank 2, thread 5, on n37. (core affinity = 4-11)
Hello from rank 2, thread 1, on n37. (core affinity = 4-11)
Hello from rank 2, thread 6, on n37. (core affinity = 4-11)
Hello from rank 2, thread 7, on n37. (core affinity = 4-11)
Hello from rank 0, thread 0, on n37. (core affinity = 12-19)
Hello from rank 0, thread 2, on n37. (core affinity = 12-19)
Hello from rank 0, thread 1, on n37. (core affinity = 12-19)
Hello from rank 0, thread 5, on n37. (core affinity = 12-19)
Hello from rank 0, thread 3, on n37. (core affinity = 12-19)
Hello from rank 0, thread 7, on n37. (core affinity = 12-19)
Hello from rank 0, thread 6, on n37. (core affinity = 12-19)
Hello from rank 0, thread 4, on n37. (core affinity = 12-19)
Hello from rank 1, thread 3, on n37. (core affinity = 0-3,20-23)
Hello from rank 1, thread 2, on n37. (core affinity = 0-3,20-23)
Hello from rank 1, thread 5, on n37. (core affinity = 0-3,20-23)
Hello from rank 1, thread 0, on n37. (core affinity = 0-3,20-23)
Hello from rank 1, thread 7, on n37. (core affinity = 0-3,20-23)
Hello from rank 1, thread 6, on n37. (core affinity = 0-3,20-23)
Hello from rank 1, thread 1, on n37. (core affinity = 0-3,20-23)
Hello from rank 1, thread 4, on n37. (core affinity = 0-3,20-23)

Задача была запущена на двух узлах (n37 и n38), на первом запущены 3 MPI процесса по 8 ядер, на втором – один.

Обратите внимание на привязку к ядрам core affinity. Для одного из MPI-процессов (rank 1) эта привязка не оптимальна (core affinity = 0-3,20-23), так как задействована часть ядер одного процессора и часть ядер второго процессора (в очереди short на каждом узле два процессора с 12 ядрами (0-11 и 12-23). Сравните с привязкой к ядрам при использовании дополнительного ограничения на 2 MPI-процесса на узел:

sbatch --ntasks-per-node=2 test-xthi.sh

Дополнительные опции команд sbatch, srun, mpirun могут влиять на конфигурацию выделенных ядер и на привязку процессов к ядрам. Обратитесь к документации Slurm, Open MPI, Intel MPI за дополнительной информацией.