Средства измерения производительности программ на Blue Gene/P

Доступ к внутреннему счетчику производительности (Universal Performance Counter (UPC) Unit) вычислительного узла суперкомпьютера Blue Gene/P возможен через системные программные интерфейсы (system programming interfaces, SPIs).

Перед применением этих функций внимательно ознакомьтесь с доступной документацией и заголовочными файлами:

Ниже приведен пример программы, написанной на основе указанных материалов. Обратите внимание, что компилировать код, включающий заголовочные файлы <spi/UPC.h> и <spi/UPC_Events.h>, необходимо командой gcc. Указанное приложение допустимо запускать только в SMP-режиме, либо синхронизировать доступ к UPC-процедурам с помощью локального коммуникатора (этот вопрос мы здесь не рассматриваем). Программа выдает результаты лишь для нулевого процессорного ядра. Актуальная версия кода (она может отличаться от приведенной ниже) доступна на системе в папке

/home/pozdneev/test/upc_wrappers

При запуске программы с командным файлом

# @ job_type = bluegene
# @ class = large
# @ output = $(jobid).out
# @ error = $(jobid).err
# @ wall_clock_limit = 00:05:00
# @ bg_size = 128
# @ queue
/bgsys/drivers/ppcfloor/bin/mpirun \
        -exe /gpfs/data/pozdneev/test/upc_wrappers/mpi-omp \
        -mode SMP \
        -env "OMP_NUM_THREADS=2"
в стандартном выводе в будет напечатано следующее:
Threads:        2
Flop:   1e+07
Time:   0.0502582
Flops:  1.98973e+08
MFlops: 198.973

Файл mpi-omp.c:

/*
 * mpi-omp.c
 *
 * Example of BG/P UPC wrappers usage in MPI+OpenMP application
 * /bgsys/drivers/ppcfloor/comm/bin/mpixlc_r -qsmp=omp \
 *          -O5 -qarch=450d -qtune=450 \
 *          -o mpi-omp mpi-omp.c upc_wrappers.o
 *
 * Alexander Pozdneev
 * http://hpc.cmc.msu.ru
 *
 * Based on the
 * - http://www.nccs.gov/wp-content/training/2008_bluegene/BobWalkup_BGP_UPC.pdf
 * - https://wiki.alcf.anl.gov/index.php/Performance
 *
 * See also:
 * - https://wiki.alcf.anl.gov/images/9/95/Blue_GeneP_UPC_User_Guide_3.0.pdf
 * - /bgsys/drivers/ppcfloor/arch/include/spi/UPC.h
 * - /bgsys/drivers/ppcfloor/arch/include/spi/UPC_Events.h
 * - IBM RedPaper REDP-4256-01 http://www.redbooks.ibm.com
 */


#include <stdio.h>
#include <stdlib.h>
#include <mpi.h>
#include <omp.h>

#include "upc_wrappers.h"

#define N (10 * 1000  * 1000)

int main(int argc, char **argv)
{
    MPI_Comm comm = MPI_COMM_WORLD;
    int rank;
    const int root = 0;

    int i;
    double *x, *y;
    double a;
    double t0, t1;

    MPI_Init(&argc, &argv);
    MPI_Comm_rank(comm, &rank);

    x = malloc(N * sizeof(*x));
    y = malloc(N * sizeof(*y));

#pragma omp parallel for schedule(static)
    for (i = 0; i < N; ++i)
    {
        x[i] = i;
        y[i] = N - i;
    }

    a = random() / (double)RAND_MAX;

    UPC_Wrapper_Init();
    UPC_Wrapper_Start_MODE_0();
    t0 = MPI_Wtime();

#pragma omp parallel for schedule(static)
    for (i = 0; i < N; ++i)
        y[i] += a * x[i];

    t1 = MPI_Wtime();
    UPC_Wrapper_Stop();

    if (rank == root)
    {
        double flop = UPC_Wrapper_Get_PU0_Flop();
        double time = t1 - t0;

        printf("Threads:\t%d\nFlop:\t%g\nTime:\t%g\nFlops:\t%g\nMFlops:\t%g\n",
                omp_get_max_threads(),
                flop, time, flop / time, flop / (time * 1000.0 * 1000.0));
    }

    free(x);
    free(y);

    MPI_Finalize();

    return 0;
}

Файл upc_wrappers.h:

/*
 * upc_wrappers.h
 *
 * Wrappers for BG/P Universal Performance Counter (UPC) Unit routines
 */


#ifndef UPC_WRAPPERS_H
#define UPC_WRAPPERS_H

#ifdef __cplusplus
extern "C" {
#endif

void UPC_Wrapper_Init(void);
void UPC_Wrapper_Start_MODE_0(void);
void UPC_Wrapper_Stop(void);
double UPC_Wrapper_Get_PU0_Flop(void);

#ifdef __cplusplus
} /* extern "C" */
#endif

#endif /* UPC_WRAPPERS_H */

Файл upc_wrappers.c:

/*
 * upc_wrappers.c
 *
 * Wrappers for BG/P Universal Performance Counter (UPC) Unit routines
 * gcc -c -I/bgsys/drivers/ppcfloor/arch/include -Wall -Wextra upc_wrappers.c
 */


#include <spi/UPC.h>
#include <spi/UPC_Events.h>

#include "upc_wrappers.h"

typedef struct Counter_Struct_s {
    int32_t rank;                       /* rank */
    int32_t core;                       /* core */
    int32_t upc_number;                 /* UPC number */
    int32_t number_processes_per_upc;   /* number of processes per UPC unit */
    BGP_UPC_Mode_t mode;                /* user mode */
    int32_t number_of_counters;         /* number of counter values returned */
    char location[24];                  /* location */
    int64_t elapsed_time;               /* elapsed time */
    uint32_t reserved_1;                /* reserved for alignment */
    uint32_t reserved_2;                /* reserved for alignment */
    int64_t values[256];                /* counter values */
} Counter_Struct;

static Counter_Struct DATA;

void UPC_Wrapper_Init(void)
{
    BGP_UPC_Initialize();
}

void UPC_Wrapper_Start_MODE_0(void)
{
    BGP_UPC_Mode_t user_mode = BGP_UPC_MODE_0;
    BGP_UPC_Event_Edge_t event_edge = BGP_UPC_CFG_EDGE_DEFAULT;

    BGP_UPC_Initialize_Counter_Config(user_mode, event_edge);
    BGP_UPC_Zero_Counter_Values();
    BGP_UPC_Start(0);
}

void UPC_Wrapper_Stop(void)
{
    BGP_UPC_Stop();
    BGP_UPC_Read_Counter_Values(&DATA, sizeof(DATA), BGP_UPC_READ_EXCLUSIVE);
}

/* See IBM RedPaper REDP-4256-01 http://www.redbooks.ibm.com */
double UPC_Wrapper_Get_PU0_Flop(void)
{
    return
        +  1.0 * DATA.values[22] /* 0th core fadd, fsub */
        +  1.0 * DATA.values[23] /* 0th core fmul */
        +  2.0 * DATA.values[24] /* 0th core fmadd family */
        + 13.0 * DATA.values[25] /* 0th core fdiv */
        +  1.0 * DATA.values[26] /* 0th core other fp non-storage ops */
        +  2.0 * DATA.values[27] /* 0th core SIMD fpadd, fpsub */
        +  2.0 * DATA.values[28] /* 0th core SIMD fpmul */
        +  4.0 * DATA.values[29] /* 0th core SIMD fpmadd family */
        +  2.0 * DATA.values[30] /* 0th core SIMD fp non-storage ops */
        ;
}

Для компиляции воспользуйтесь скриптом make.sh

#!/bin/sh

gcc -c -I/bgsys/drivers/ppcfloor/arch/include -Wall -Wextra upc_wrappers.c

/bgsys/drivers/ppcfloor/comm/bin/mpixlc_r -qsmp=omp \
        -O5 -qarch=450d -qtune=450 \
        -o mpi-omp mpi-omp.c upc_wrappers.o