Compilation avec MPI CUDA-aware et GPUDirect
Pour une performance optimale, des bibliothĂšques OpenMPI CUDA-aware supportant le GPUDirect sont disponibles sur Jean Zay.
Ces bibliothÚques MPI permettent d'effectuer des communications utilisant des buffers d'envoi et de réception alloués sur la mémoire du GPU. Grùce au support de GPUDirect, ces transferts se font directement de GPU à GPU sans recopie intermédiaire sur la mémoire du CPU, lorsque c'est possible.
Le GPUDirect n'est pas utilisable sur Jean Zay pour les codes utilisant plusieurs GPU par processus MPI.
Compilation du codeâ
Il est nécessaire de compiler le code en utilisant l'une des bibliothÚques OpenMPI CUDA-aware disponibles sur Jean Zay.
AprÚs avoir chargé le compilateur que vous souhaitez utiliser, vous devez charger l'un des modules suivants :
$ module avail openmpi/*-cuda
----- /lustre/fshomisc/sup/hpe/pub/modules-idris-env4/modulefiles/linux-rhel9-skylake_avx512 ------
openmpi/3.1.4-cuda openmpi/4.0.2-cuda openmpi/4.1.0-cuda openmpi/4.1.5-cuda
openmpi/3.1.6-cuda openmpi/4.0.4-cuda openmpi/4.1.1-cuda openmpi/4.1.6-cuda
openmpi/4.0.1-cuda openmpi/4.0.5-cuda openmpi/4.1.4-cuda openmpi/4.1.8-cuda
$ module load openmpi/4.0.4-cuda
Si OpenMPI n'est pas disponible pour le compilateur désiré, un message d'erreur sera affiché. N'hésitez pas à contacter l'assistance à assist@idris.fr pour demander une nouvelle installation si nécessaire.
Pour connaĂźtre la liste des compilateurs pour lesquels une version
donnée d'OpenMPI CUDA-aware est disponible, vous pouvez utiliser la commande
module show openmpi/<version>-cuda. Par exemple :
$ module show openmpi/4.1.8-cuda
------------------------------------------------------------------------------------
/lustre/fshomisc/sup/hpe/pub/modules-idris-env4/modulefiles/linux-rhel9-skylake_avx512/openmpi/4.1.8-cuda:
module-whatis {An open source Message Passing Interface implementation.}
prereq intel-compilers/2021.9.0 nvidia-compilers/25.1 gcc/14.2.0 intel-oneapi-compilers/2023.1.0
conflict openmpi
conflict intel-mpi
Available software environment(s):
- intel-compilers/2021.9.0 cuda/12.8.0
- nvidia-compilers/25.1 cuda/12.8.0
- gcc/14.2.0 cuda/12.8.0
- intel-oneapi-compilers/2023.1.0 cuda/12.8.0
If you want to use this module with another software environment,
please contact the support team.
------------------------------------------------------------------------------------
La rubrique "Available software environment(s):" vous indique les compilateurs compatibles avec la bibliothÚque indiquée.
La compilation se fait en utilisant les wrappers d'OpenMPI. Par exemple :
mpifort source.f90
mpicc source.c
mpic++ source.C
Aucune option particuliÚre n'est nécessaire pour la compilation, vous pouvez vous référer aux pages de documentation sur les compilateurs NVIDIA/PGI et sur la compilation de code OpenACC pour plus d'information sur la compilation de codes utilisant les GPU.
Adaptation du codeâ
Sur les partitions V100 ou A100 de Jean Zay, l'utilisation de la fonctionnalité MPI CUDA-aware GPUDirect impose de respecter un ordre d'initialisation bien précis pour CUDA ou OpenACC et MPI dans le code :
- initialisation de CUDA ou OpenACC ;
- choix du GPU que chaque processus MPI doit utiliser (étape de binding) ;
- initialisation de MPI.
Si cet ordre d'initialisation n'est pas respecté,
l'exécution de votre code risque de planter avec l'erreur suivante :
CUDA failure: cuCtxGetDevice().
Cette adaptation est indispensable si vous travaillez sur une des partitions dépendant du réseau d'interconnexion OmniPAth, soient les partitions V100 et A100 de Jean Zay.
Vous pouvez aussi l'appliquer si vous travaillez sur la partition H100 (sur réseau InfiniBand) mais elle n'est pas nécessaire.
Vous trouverez ci-dessous un exemple de sous-routine en Fortran et en C permettant d'initialiser OpenACC avant d'initialiser MPI ainsi qu'un exemple CUDA.
Exemple OpenACCâ
Cet exemple ne fonctionne que dans le cas oĂč on alloue exactement un processus MPI par GPU.
Version Fortranâ
#ifdef _OPENACC
subroutine initialisation_openacc
USE openacc
character(len=6) :: local_rank_env
integer :: local_rank_env_status, local_rank
! Initialisation d'OpenACC
!$acc init
! Récupération du rang local du processus via la variable d'environnement
! positionnée par Slurm, l'utilisation de MPI_Comm_rank n'étant pas encore
! possible puisque cette routine est utilisée AVANT l'initialisation de MPI
call get_environment_variable(name="SLURM_LOCALID", value=local_rank_env, status=local_rank_env_status)
if (local_rank_env_status == 0) then
read(local_rank_env, *) local_rank
! Définition du GPU à utiliser via OpenACC
call acc_set_device_num(local_rank, acc_get_device_type())
else
print *, "Erreur : impossible de déterminer le rang local du processus"
stop 1
end if
end subroutine initialisation_openacc
#endif
Exemple d'utilisation :
! On initialise OpenACC...
#ifdef _OPENACC
call initialisation_openacc
#endif
! ... avant d'initialiser MPI
call mpi_init(code)
Version Câ
#ifdef _OPENACC
void initialisation_openacc()
{
char* local_rank_env;
int local_rank;
/* Initialisation d'OpenACC */
#pragma acc init
/* Récupération du rang local du processus via la variable d'environnement
positionnée par Slurm, l'utilisation de MPI_Comm_rank n'étant pas encore
possible puisque cette routine est utilisée AVANT l'initialisation de MPI */
local_rank_env = getenv("SLURM_LOCALID");
if (local_rank_env) {
local_rank = atoi(local_rank_env);
/* Définition du GPU à utiliser via OpenACC */
acc_set_device_num(local_rank, acc_get_device_type());
} else {
printf("Erreur : impossible de déterminer le rang local du processus\n");
exit(1);
}
}
#endif
Exemple d'utilisation :
#ifdef _OPENACC
/* On initialise OpenACC... */
initialisation_openacc();
#endif
/* ... avant d'initialiser MPI */
MPI_Init(&argc, &argv);
Exemple CUDAâ
Cet exemple ne fonctionne que dans le cas oĂč on alloue exactement un processus MPI par GPU.
#include <cuda.h>
#include <cuda_runtime.h>
void initialisation_cuda()
{
char* local_rank_env;
int local_rank;
cudaError_t cudaRet;
/* Récupération du rang local du processus via la variable d'environnement
positionnée par Slurm, l'utilisation de MPI_Comm_rank n'étant pas encore
possible puisque cette routine est utilisée AVANT l'initialisation de MPI */
local_rank_env = getenv("SLURM_LOCALID");
if (local_rank_env) {
local_rank = atoi(local_rank_env);
/* Définition du GPU à utiliser pour chaque processus MPI */
cudaRet = cudaSetDevice(local_rank);
if(cudaRet != CUDA_SUCCESS) {
printf("Erreur: cudaSetDevice a échoué\n");
exit(1);
}
} else {
printf("Erreur : impossible de déterminer le rang local du processus\n");
exit(1);
}
}
Exemple d'utilisation :
/* On initialise CUDA... */
initialisation_cuda();
/* ... avant d'initialiser MPI */
MPI_Init(&argc, &argv);
ExĂ©cution du codeâ
Lors de l'exécution, vous devez vous assurez de charger avec la
commande module la mĂȘme bibliothĂšque MPI que celle qui a Ă©tĂ© utilisĂ©e
pour la compilation du code, puis de bien utiliser la commande srun
pour démarrer celui-ci.
Le support de CUDA-aware et GPUDirect est alors activé par défaut sans opération supplémentaire.
Pour plus d'informations sur la soumission des travaux multi-GPU utilisant MPI, consultez la page sur l'exécution d'un travail multi-GPU MPI CUDA-aware et GPUDirect en batch.