La bibliothèque PETSc
Description
PETSc est une bibliothèque de fonctions écrites en C (mais offrant également une interface en Fortran) permettant de manipuler des vecteurs et matrices denses ou creuses et de résoudre les systèmes linéaires correspondants (souvent liées aux équations aux dérivées partielles) avec des solveurs directs ou itératifs. PETSc propose aussi des solveurs non-linéaires, des méthodes de résolution d'équations différentielles et beaucoup d'autres outils. Enfin, PETSc permet aussi d'utiliser des bibliothèques externes.
PETSc est disponible pour CPU et pour GPU.
Versions installées et environnements
Diverses versions sont installées sur Jean Zay et sont accessibles par la commande module. Elles peuvent être listées avec la commande suivante :
module avail petsc
Des informations sur une version particulière peuvent être obtenues via la commande :
module show petsc-xxx
Par exemple :
module show petsc/3.21.4-mpi-complex-int64--------------------------------------------------------------------------/lustre/fshomisc/sup/hpe/pub/modules-idris-env4/modulefiles/linux-rhel9-skylake_avx512/petsc/3.21.4-mpi-complex-int64:
module-whatis {PETSc is a suite of data structures and routines for the scalable (parallel) solution of scientific applications modeled by partial differential equations. }prereq intel-compilers/19.1.3 intel-compilers/2021.9.0conflict petsc
Available software environment(s):- intel-compilers/19.1.3 intel-mpi/2019.9- intel-compilers/2021.9.0 intel-mpi/2021.9
If you want to use this module with another software environment, please contact the support team.--------------------------------------------------------------------------
Cette installation est disponible dans deux environnement différents, avec le premier (intel-compilers/19.1.3 intel-mpi/2019.9 ici) venant par défaut. Si vous désirez un autre environnement, il faut le charger au préalable. Par exemple :
module load intel-compilers/2021.9.0 intel-mpi/2019.9 petsc
Les options utilisées pour la compilation de PETSc sont accessibles dans le fichier $PETSC_DIR/$PETSC_ARCH/include/petscconfiginfo.h et les options installées dans $PETSC_DIR/$PETSC_ARCH/include/petscconf.h. Les variables $PETSC_DIR et $PETSC_ARCH sont automatiquement valorisées lors du chargement de PETSc via la commande module.
Il est possible d'ajouter d'autres versons de PETSc. Pour cela, veuillez contacter l'équipe Support aux Utilisateurs.
Documentation
📝 La documentation sur le site de référence.
Exemples d'utilisations
Dans cet exemple, nous utilisons un exercice basé sur la sous-classe DMDA de DM, qui gère la communication entre les structures algébriques de PETSc et les structures de données maillées. Cet exercice s'inscrit dans le cadre de notre formation PETSc. Vous pouvez consulter davantage de détails dans le support de cours accessible à l'adresse suivante : Formation PETSc.
Code source en C et fichier de compilation:dmda.c:
#include <petscdmda.h>
#include <petscksp.h>
#include <petsctime.h>
int main(int argc, char **argv)
{
const PetscInt nx = 10, ny = 10, stencil_size = 5;
PetscInt i, j, its;
Mat A;
Vec x, b;
KSP solver;
const PetscReal rtol = 1.e-8;
KSPConvergedReason reason;
PetscReal errorNorm, rnorm;
PetscLogDouble t1, t2;
DM dm;
DMDALocalInfo info;
const PetscInt stencilWidth = 1;
MatStencil row, col5[stencil_size];
PetscScalar hx2, hy2, coef, coef5[stencil_size];
PetscScalar **bgrid;
PetscFunctionBeginUser;
PetscCall(PetscInitialize(&argc, &argv, NULL, NULL));
// Create the 2-D DMDA object
PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_STAR,
nx, ny, PETSC_DECIDE, PETSC_DECIDE, 1, stencilWidth,
NULL, NULL, &dm));
PetscCall(DMSetFromOptions(dm));
PetscCall(DMSetUp(dm));
// View the DMDA object
PetscCall(DMView(dm, PETSC_VIEWER_STDOUT_WORLD));
// Create the A matrix from the DMDA object
PetscCall(DMCreateMatrix(dm, &A));
// Retrieve local information from the DMDA object
PetscCall(DMDAGetLocalInfo(dm, &info));
hx2 = 1.0 / ((info.mx - 1) * (info.mx - 1));
hy2 = 1.0 / ((info.my - 1) * (info.my - 1));
coef = 1.0;
coef5[0] = 2.0 / hx2 + 2.0 / hy2;
coef5[1] = -1.0 / hx2; coef5[2] = -1.0 / hx2;
coef5[3] = -1.0 / hy2; coef5[4] = -1.0 / hy2;
// Loop on the grid points
for (j = info.ys; j < info.ys + info.ym; j++) {
for (i = info.xs; i < info.xs + info.xm; i++) {
row.i = i; row.j = j; row.c = 0;
if (i == 0 || i == (info.mx - 1) || j == 0 || j == (info.my - 1)) {
// Set matrix values to enforce boundary conditions (homogeneous Dirichlet conditions)
PetscCall(MatSetValuesStencil(A, 1, &row, 1, &row, &coef, INSERT_VALUES));
} else {
// Set matrix values fo interior points
col5[0].i = i; col5[0].j = j; col5[0].c = 0;
col5[1].i = i - 1; col5[1].j = j; col5[1].c = 0;
col5[2].i = i + 1; col5[2].j = j; col5[2].c = 0;
col5[3].i = i; col5[3].j = j - 1; col5[3].c = 0;
col5[4].i = i; col5[4].j = j + 1; col5[4].c = 0;
PetscCall(MatSetValuesStencil(A, 1, &row, stencil_size, col5, coef5, INSERT_VALUES));
}
}
}
PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
// PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD));
// Create global vectors b and x from the DMDA object
PetscCall(DMCreateGlobalVector(dm, &b));
PetscCall(DMCreateGlobalVector(dm, &x));
PetscCall(VecSet(b, 0.0));
PetscCall(VecSetRandom(x, NULL));
PetscCall(KSPCreate(PETSC_COMM_WORLD, &solver));
PetscCall(KSPSetOperators(solver, A, A));
PetscCall(KSPSetTolerances(solver, rtol, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
PetscCall(KSPSetInitialGuessNonzero(solver, PETSC_TRUE));
PetscCall(KSPSetFromOptions(solver));
PetscCall(KSPSetUp(solver));
PetscCall(KSPView(solver, PETSC_VIEWER_STDOUT_WORLD));
PetscCall(PetscTime(&t1));
PetscCall(KSPSolve(solver, b, x));
PetscCall(PetscTime(&t2));
PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Elapse time: %lf s\n", t2 - t1));
PetscCall(KSPGetConvergedReason(solver, &reason));
PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Converged reason: %s\n", KSPConvergedReasons[reason]));
PetscCall(KSPGetIterationNumber(solver, &its));
PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Iterations number: %d\n", its));
PetscCall(KSPGetResidualNorm(solver, &rnorm));
PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Residual norm: %g\n", rnorm));
PetscCall(VecNorm(x, NORM_2, &errorNorm));
PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error norm (NORM_2): %g\n", errorNorm));
PetscCall(VecNorm(x, NORM_INFINITY, &errorNorm));
PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error norm (NORM_INFINITY): %g\n", errorNorm));
PetscCall(KSPDestroy(&solver));
PetscCall(MatDestroy(&A));
PetscCall(VecDestroy(&x));
PetscCall(VecDestroy(&b));
PetscCall(DMDestroy(&dm));
PetscCall(PetscFinalize());
return 0;
}
}
Makefile:
ALL: dmda
SOURCES = dmda.c
OBJ = $(SOURCES:.c=.o)
EXE = dmda.exe
CLEANFILES = ${OBJ} ${EXE}
include ${PETSC_DIR}/lib/petsc/conf/variables
include ${PETSC_DIR}/lib/petsc/conf/rules
dmda: ${OBJ}
${CLINKER} -o ${EXE} ${OBJ} ${PETSC_LIB}
Exemple CPU
Pour le cas CPU, la compilation s’effectue après le chargement du module petsc/3.20.1-mpi, par exemple. L’exécution du job se fait ensuite à l’aide du script de soumission suivant :
#!/bin/bash
#SBATCH --job-name=dmda_cpu # Nom du job
#SBATCH --nodes=1
#SBATCH --ntasks=4 # Nombre total de processus MPI
#SBATCH --hint=nomultithread # 1 processus MPI par coeur physique (pas d'hyperthreading)
#SBATCH --time=00:10:00 # Temps d’exécution maximum demande (HH:MM:SS)
#SBATCH --output=dmda_cpu_%j.out # Nom du travail
#SBATCH --error=dmda_cpu_%j.out # Nom du fichier d'erreur (ici commun avec la sortie)
#SBATCH --qos=qos_cpu-dev
cd ${SLURM_SUBMIT_DIR}
module purge
# Utilise la meme version que pour la compilation du code
module load petsc/3.20.1-mpi
module list
view_args="-log_view" # to avoid slowing down due to -log_view_gpu_time
cpu_args="-dm_mat_type mpiaij -dm_vec_type mpi"
set -x
srun ./dmda.exe -da_grid_x 1000 -da_grid_y 1000 -ksp_type cg $cpu_args $view_args
Exemple GPU
La version GPU de PETSc, utilisant CUDA, est installée sur toutes nos partitions GPU : V100, A100, et H100. Les modules correspondants sont disponibles sous le nom petsc/3.22.1-mpi-cuda (pour les partitions A100 et H100, faire au préalable un module load arch/a100 ou module load arch/h100, respectivement).
Afin de générer un exécutable compatible avec la partition H100, pensez à compiler votre code PETSc sur la partition compil_H100 en effectuant un module load arch/h100.
Sur Jean-Zay, des bibliothèques OpenMPI CUDA-aware prennent en charge le GPUDirect, qui permet d’optimiser les communications entre GPU. Pour plus d'informations, veuillez consulter la documentation MPI CUDA-aware et GPUDirect. Dans cet exemple, nous exploitons GPUDirect pour une performance optimale dans un cas multi-GPU. Cela nécessite d’adapter le code dmda.c présenté dans la section précédente en effectuant les modifications nécessaires suivantes :
Ajout de l’étape de binding, inclure init_cuda.c dans les SOURCES du fichier Makefile :
init_cuda.c:
#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);
}
}
Initialisation de CUDA avant MPI :
dmda_gpu.c:
…
/* On initialise CUDA... */
initialisation_cuda();
PetscFunctionBeginUser;
PetscCall(PetscInitialize(&argc, &argv, NULL, NULL));
…
Voici un script de soumission adapté pour exécuter cet exercice sur les GPU A100 :
#!/bin/bash
#SBATCH --job-name=dmda_gpu # nom du job
#SBATCH --nodes=1
#SBATCH --ntasks-per-node=2 # nombre de tâche MPI par noeud
#SBATCH --gres=gpu:2 # nombre de GPU à réserver par nœud
#SBATCH --cpus-per-task=8 # nombre de coeurs à réserver par tâche
#SBATCH --hint=nomultithread # 1 processus MPI par coeur physique (pas d'hyperthreading)
#SBATCH --time=00:10:00 # Temps d’exécution maximum demande (HH:MM:SS)
#SBATCH --output=dmda_gpu_a100%j.out # Nom du travail
#SBATCH --error=dmda_gpu_a100%j.out # Nom du fichier d'erreur (ici commun avec la sortie)
#SBATCH --qos=qos_gpu_a100-dev
#SBATCH -C a100
cd ${SLURM_SUBMIT_DIR}
module purge
module load arch/a100
# Utilise la meme version que pour la compilation du code
module load petsc/3.22.1-mpi-cuda
module list
view_args="-log_view" # to avoid slowing down due to -log_view_gpu_time
if [ $SLURM_NTASKS == 1]
then
### mono-GPU :
gpu_args="-dm_mat_type aijcusparse -dm_vec_type cuda"
else
### multi-GPU :
gpu_args="-dm_mat_type mpiaijcusparse -dm_vec_type mpicuda"
fi
set -x
srun ./dmda.exe -da_grid_x 1000 -da_grid_y 1000 -ksp_type cg -pc_type hypre $gpu_args $view_args