Запуск MPI-задач #
Короткий совет: заменитеmpiexec
наmpirun
.
На странице со статусом работы InfiniBand есть дополнительная информация об известных проблемах.
Для запуска MPI-задач можно использовать команду mpirun
.
По умолчанию задача будет запущена на заказанном количестве процессов, информация об этом будет взята из системы очередей Slurm.
Если нужно запустить меньшее количество процессов, то можно использовать параметр -n np
.
например:
mpirun -n 8 ./mpiprog8.exe
Intel MPI #
При использовании Intel MPI команда mpiexec
работать не будет, используйте mpirun
.
OpenMPI #
Вы можете использовать OpenMPI с компиляторами gcc, для этого включите модуль openmpi:
module load gnu openmpi
Запуск гибридных задач MPI+OpenMPI #
Рассмотрим следующий пример простой гибридной задачи hybrid.c
:
#include <stdio.h>
#include <omp.h>
#include "mpi.h"
int main(int argc, char *argv[]) {
int numprocs, rank, namelen;
char processor_name[MPI_MAX_PROCESSOR_NAME];
int iam = 0, np = 1;
MPI_Init(&argc, &argv);
MPI_Comm_size(MPI_COMM_WORLD, &numprocs);
MPI_Comm_rank(MPI_COMM_WORLD, &rank);
MPI_Get_processor_name(processor_name, &namelen);
#pragma omp parallel default(shared) private(iam, np)
{
np = omp_get_num_threads();
iam = omp_get_thread_num();
printf("Node: %12s rank: %2d/%2d thread_id: %2d/%2d\n",
processor_name, rank, numprocs, iam, np);
}
MPI_Finalize();
}
Для компиляции с помощью Intel MPI выполните следующие команды:
. /opt/intel/composer_xe_2015/bin/compilervars.sh intel64
. /opt/intel/impi_5.0.3/bin64/mpivars.sh
mpiicc -qopenmp hybrid.c -o hybrid
Для компиляции с помощью gcc + OpenMPI выполните следующие команды:
module load gnu openmpi
mpicc -fopenmp hybrid.c -o hybrid
Для запуска задачи воспользуемся следующим скриптом test-hybrid.sh
:
#!/bin/bash
# A hybrid MPI+OpenMPI example
#SBATCH --job-name=hybrid
#SBATCH --output=hybrid.out
#SBATCH --partition=x6core
#SBATCH --time=5
####### 4 MPI ranks
#SBATCH --ntasks=4
####### 6 OMP threads per MPI rank
#SBATCH --cpus-per-task=6
## load MPI libs if needed
## for Intel:
# . /opt/intel/composer_xe_2015/bin/compilervars.sh intel64
# . /opt/intel/impi_5.0.3/bin64/mpivars.sh
## for OpenMPI:
# module load gnu openmpi
## set number of OMP thread to the value of --cpus-per-task
export OMP_NUM_THREADS=$SLURM_CPUS_PER_TASK
mpirun ./hybrid
Система очередей передает текущие переменные окружения, поэтому в самом скрипте активировать профиль для MPI не обязательно, если он уже активирован у текущего пользователя.
Переменную OMP_NUM_THREADS
лучше явно приравнять значению $SLURM_CPUS_PER_TASK
. В некоторых случаях в этом нет необходимости, так как каждый MPI процесс будет запущен на ограниченном количестве ядер, и OpenMP автоматически выберет нужное количество потоков.
Запускаем задачу в очередь:
sbatch test-hybrid.sh
Результат работы будет в файле hybrid.out
, пример:
Node: cl1n031 rank: 1/ 4 thread_id: 0/ 6
Node: cl1n031 rank: 1/ 4 thread_id: 1/ 6
Node: cl1n031 rank: 1/ 4 thread_id: 5/ 6
Node: cl1n031 rank: 1/ 4 thread_id: 3/ 6
Node: cl1n031 rank: 0/ 4 thread_id: 5/ 6
Node: cl1n031 rank: 0/ 4 thread_id: 0/ 6
Node: cl1n031 rank: 0/ 4 thread_id: 4/ 6
Node: cl1n031 rank: 1/ 4 thread_id: 2/ 6
Node: cl1n031 rank: 1/ 4 thread_id: 4/ 6
Node: cl1n031 rank: 0/ 4 thread_id: 1/ 6
Node: cl1n031 rank: 0/ 4 thread_id: 2/ 6
Node: cl1n031 rank: 0/ 4 thread_id: 3/ 6
Node: cl1n032 rank: 2/ 4 thread_id: 0/ 6
Node: cl1n032 rank: 2/ 4 thread_id: 2/ 6
Node: cl1n032 rank: 2/ 4 thread_id: 1/ 6
Node: cl1n032 rank: 2/ 4 thread_id: 4/ 6
Node: cl1n032 rank: 2/ 4 thread_id: 5/ 6
Node: cl1n032 rank: 2/ 4 thread_id: 3/ 6
Node: cl1n032 rank: 3/ 4 thread_id: 1/ 6
Node: cl1n032 rank: 3/ 4 thread_id: 0/ 6
Node: cl1n032 rank: 3/ 4 thread_id: 5/ 6
Node: cl1n032 rank: 3/ 4 thread_id: 4/ 6
Node: cl1n032 rank: 3/ 4 thread_id: 3/ 6
Node: cl1n032 rank: 3/ 4 thread_id: 2/ 6
Задача была запущена на двух узлах (cl1n031
и cl1n032
), на каждом узле были запущены две копии MPI-программы, и в каждой используется 6 OpenMP потоков.