RefNUMA – библиотека для организации виртуальной общей памяти в программах, использующих MPI.

А. О. Лацис

Пример 9. То же модельное приложение на простых, не векторизованных разделяемых массивах.

Вводя понятие векторизованного, блочно распределенного разделяемого массива, мы отмечали тот очевидный факт, что такой массив не может иметь размерность, меньшую, чем 2, а его основная «деталь» – вектор указателей на строки – занимает место в памяти, иногда – довольно заметное. В то же время, многие программы на С традиционно пишутся с использованием одномерных массивов. В этом разделе мы рассмотрим специальные возможности RefNUMA по организации одномерных массивов. Поскольку в С, строго говоря, все массивы одномерные («многомерным массивом» для краткости называют одномерный массив массивов меньшей размерности), описываемые ниже возможности правильнее называть возможностями организации произвольных массивов.

Недостаток этого аппарата – в некотором снижении быстродействия. Достоинство, помимо универсальности, в том, что не требуется использовать дополнительную память, как в случае векторизованного массива.

Пусть указатель на массив указателей на секции секционированного массива объявлен как double **f, и в каждой секции предполагается хранить nb элементов типа double. Очевидно, если мы захотим использовать этот секционированный массив как одномерный, блочно распределенный разделяемый массив, обращение к i-му элементу этого массива должно иметь вид: f[i/nb][i%nb]. Договоримся записывать в нашей программе обращение к i-му элементу этого массива в виде F[i]. Если бы теперь какая-нибудь специальная программа обработала перед трансляцией текст нашей программы, заменив все вхождения F[i] на f[i/nb][i%nb], в результате такой обработки получилась бы правильная программа работы с секционированным массивом как со сплошным одномерным. Конечно, на практике все не так примитивно, но общая идея понятна. При выполнении такой программы, необходимость делить значение индекса на nb, а также брать остаток от этого деления, несколько замедлит работу программы, по сравнению с использованием «настоящего» одномерного массива F. Однако, в случае, если nb – переменная, хранимая в регистре, замедление это будет вполне терпимым.

Такие программы предварительной обработки исходных текстов программ, как правило, предназначенные для систематической замены одних несложных фрагментов текста на другие, применяются довольно широко. Называются они препроцессорами исходного текста, или просто препроцессорами. Наиболее широко известен повсеместно используемый стандартный препроцессор С, предназначенный для раскрытия в тексте программы конструкций #include, #define и вызовов соответствующих макросов.

Препроцессоры бывают общего назначения, предназначенные для задания более или менее произвольных макросов, и раскрытия их вызовов в более или менее произвольных текстах, и специальные, нацеленные на одно, совершенно конкретное, преобразование текстов программ на совершенно конкретном языке. Для решения поставленной выше задачи, грубо говоря, систематической замены F[i] на f[i/nb][i%nb], можно попробовать воспользоваться каким-нибудь препроцессором общего назначения, а можно написать специальный препроцессор, предназначенный для решения этой и только этой задачи. Стандартный препроцессор С решает эту задачу «почти наполовину» - он мог бы легко обеспечить необходимую замену для F(i), но не для F[i]. Решено было разработать специальный препроцессор, и сделать его входной язык похожим на язык директив OpenMP.

Исходный текст процессорного (не гибридного) варианта модельного приложения, написанного в расчете на предварительную обработку препроцессором, приводится ниже. Директива управления препроцессированием должна располагаться в отдельной строке, и всегда начинается с #pragma refnuma

#include <stdio.h>
#include <stdlib.h>
#include <mpi.h>
    COARRAY_MEM_ALLOC( 2000000000l );
    static void usage( void )
{
    printf( "Usage: progrev_coarr <nrows> <ncols> <niter>\n" );
    exit( -1 );
}
    void itstep( int mx, int my, void *pf, void *pnewf, void *pr,
                       double rdx2, double rdy2, double beta );
    int main( int argc, char *argv[] )
{
    int mx, my, mx1, my1, niter, n_nodes, my_node;
    void *ra;
    long im, j, k;
    register int i;
    FILE *fp;
    double rdx2, rdy2, beta, t;
/***/
    COARRAY_Init( &argc, &argv );
    if ( argc != 4 ) usage();
    rdx2 = 1.0;
    rdy2 = 1.0;
    beta = 0.25;
    my_node = coarray_my_node();
    n_nodes = coarray_n_nodes();
    mx = (int)atol( argv[1] ) + 2;
    my = (int)atol( argv[2] ) + 2;
    niter = (int)atol( argv[3] );
     {
      typedef double element[my];
      register long nblock;
      nblock = shared_blocksize( (long)mx, n_nodes );
      COARRAY_Create( element, fo, mx, 1.0 );
      COARRAY_Create( element, newfo, mx, 1.0 );
      COARRAY_Create( element, ro, mx, 1.0 );
#pragma refnuma shared f    coarray fo blocksize nblock
#pragma refnuma shared newf coarray newfo blocksize nblock
#pragma refnuma shared r    coarray ro blocksize nblock
      if ( (!fo) || (!newfo) || (!ro) )
       {
        printf( "No memory\n" );
        exit( -1 );
       }
/***/
#pragma refnuma parallel for arraysize mx
      for ( i = 0; i < mx; i++ )
       {
        for ( j = 0; j < my; j++ )
         {
          f[i][j] = newf[i][j] = 0.0;
          r[i][j] = 0.0;
          if ( (i == 0) || (i == (mx-1)) || (j == 0) || (j == (my-1)) )
                                                         f[i][j] = newf[i][j] = 1.0;
         }
       }
      coarray_barrier();
      t = MPI_Wtime();
      for ( k = 0; k < niter; k++ )
       {
        itstep( mx, my, f, newf, r, rdx2, rdy2, beta );
        coarray_barrier();
        ra = (void*)f;
        f = newf;
        newf = (typeof(newf))ra;
       }
      t = MPI_Wtime() - t;
      printf( "Time: %f speed %f flops\n", t,
                   (((double)(mx-2))*((double)(my-2))*7.0*niter)/t );
      if ( !my_node )
       {
        fp = fopen( "progrev_coarr.dat", "w" );
        fclose( fp );
       }
      NODE_BY_NODE_BEGIN( k, 0, n_nodes )
      fp = fopen( "progrev_coarr.dat", "a" );
#pragma refnuma parallel for arraysize mx
      for ( i = 1; i < (mx-1); i++ )
       {
        fwrite( f[i]+1, sizeof(double), my-2, fp );
//	for ( j = 1; j < my1; j++ ) fprintf( fp, "%f  ", f[i][j] ); fprintf( fp, "\n" );
       }
      fclose( fp );
      coarray_report( (void**)f, &im, &j );
      printf( "In node %d for f: requested %ld used %ld\n", my_node, im, j );
      coarray_report( (void**)newf, &im, &j );
      printf( "In node %d for newf: requested %ld used %ld\n", my_node, im, j );
      coarray_report( (void**)r, &im, &j );
      printf( "In node %d for r: requested %ld used %ld\n", my_node, im, j );
      fflush( stdout );
      NODE_BY_NODE_END
#pragma refnuma shared f off
#pragma refnuma shared newf off
#pragma refnuma shared r off
     }
    COARRAY_Finalize();
    return 0;
}
    void itstep( int mx, int my, void *pf, void *pnewf, void *pr,
                        double rdx2, double rdy2, double beta )
{
    int mx1, my1;
    register int i, j;
/***/
    typedef double element[my];
    register long nblock;
    nblock = shared_blocksize( (long)mx, coarray_n_nodes() );
    COARRAY_Create_from_arg( element, fo, pf );
    COARRAY_Create_from_arg( element, newfo, pnewf );
    COARRAY_Create_from_arg( element, ro, pr );
#pragma refnuma shared f    coarray fo blocksize nblock
#pragma refnuma shared newf coarray newfo blocksize nblock
#pragma refnuma shared r    coarray ro blocksize nblock
    mx1 = mx - 1;
    my1 = my - 1;
#pragma refnuma parallel for arraysize mx
    for ( i = 1; i < mx1; i++ )
     {
      for ( j = 1; j < my1; j++ )
       {
        newf[i][j] =
                 ((f[i-1][j]+f[i+1][j])*rdx2+(f[i][j-1]+f[i][j+1])*rdy2-r[i][j])*beta;
       }
     }
#pragma refnuma shared f off
#pragma refnuma shared newf off
#pragma refnuma shared r off
}

Замечание.

В приведенных ранее примерах, в особенности – при рассмотрении дополнительных функций и макросов, описывались возможности выяснения тех или иных диапазонов номеров строк векторизованного массива общей памяти, а также принадлежности заданной строки векторизованного массива общей памяти тому или иному процессу параллельной программы. На базе именно этих возможностей строится, в частности, понятие параллельного цикла. Все эти возможности вовсе не являются специфичными именно для векторизованных массивов общей памяти. Каждая из них в равной мере относится и к простым разделяемым массивам общей памяти, рассматриваемым в настоящем примере. Конечно, при замене векторизованных массивов общей памяти на простые разделяемые массивы следует во всех соответствующих определениях и пояснениях заменить упоминание строки (векторизованного массива общей памяти) на упоминание элемента (простого, одномерного разделяемого массива общей памяти). В последующем изложении мы не будем заново разъяснять понятия параллельного цикла, смысл макросов SHARED_START() и SHARED_FINISH(), и прочих подобных конструкций RefNUMA.

Приступим к разбору текста программы.

Вплоть до строки

typedef double element[my];

приведенный выше текст ничем не отличается от вариантов, рассмотренных ранее. Необходимость объявления element как типа «массив my элементов типа double» объясняется тем, что мы собираемся рассматривать массивы f, newf и r как одномерные. По смыслу задачи эти массивы являются двумерными массивами элементов типа double. Если рассматривать их как одномерные массивы, то это будут одномерные массивы строк, по my значений типа double в строке. Тип именно таких строк мы и объявили как element.

Переменная register long nblock; используется далее в качестве размера секции секционированного массива. Говоря далее о размере, или длине, секционированного массива, будем понимать иметь в виду размер, или длину, исключительно и только в element'ах (в строках). Переменную nblock важно объявить явным образом как register. Если транслятор не поместит эту переменную в регистр, то конструкция вида f[i/nblock][i%nblock] будет работать недопустимо медленно. Если транслятор поместит эту переменную в регистр сам, без явного указания со стороны программиста, а потом, при внесении в программу тех или иных мелких изменений, «передумает», программа внезапно и необъяснимо замедлится, возможно, в разы. Поэтому слово «register» в объявлении этой переменной надо написать обязательно, а если транслятор будут выдавать связанные с этим словом предупредительные диагностические сообщения, то от них обязательно надо избавиться.

Макрос COARRAY_Create() похож на рассмотренный нами выше макрос COARRAY_Create_vectorized(), но создает не векторизованный распределенный массив, а обычный секционированный массив элементов указанного типа, заданной длины, с заданным коэффициентом запаса по буферной памяти.

Чтобы связать с секционированным массивом простой, блочно распределенный, разделяемый массив, необходимо использовать директиву #pragma refnuma shared.

Так, директива

#pragma refnuma shared f coarray fo blocksize nblock
означает, что далее по тексту все упоминания f[i] будут заменяться на fo[i/nblock][i%nblock], а все упоминания f без индексных скобок будут заменяться на fo (это необходимо для передачи разделяемого массива в функцию в качестве аргумента). То есть программу можно писать так, как если бы в ее начале был объявлен одномерный разделяемый, блочно распределенный массив общей памяти f.

Область действия рассмотренной директивы завершается директивой

#pragma refnuma shared f off

Эту директиву рекомендуется размещать перед закрывающей фигурной скобкой блока, в котором объявлен соответствующий массив.

Директива

#pragma refnuma parallel for arraysize mx
превращает ближайший следующий за ней цикл for, в данном случае
for ( i = 0; i < mx; i++ )
в параллельный цикл, работающий в точности, как:
PARALLEL_DO( i, 0, mx, 1, mx );

Теперь рассмотрим на примере функции itstep() порядок передачи простого разделяемого массива в функцию в качестве аргумента.

При вызове функции itstep() ей в качестве аргументов передаются массивы f, newf и r. Соответствующие формальные параметры в самой функции объявлены как void*.

Внутри функции снова объявляется тип element. Поскольку тип этот является динамическим (число double'ов в element'е становится известно только во время выполнения программы), тип приходится объявлять в блоке, и за пределами блока это объявление не действует. Зная значение my, которое также передается в itstep() в качестве аргумента, тип element легко объявить таким, чтобы он совпадал с типом element из вызывающей программы. Чтобы превратить переданный в функцию указатель void* в указатель на секционированный массив, используется специальный макрос COARRAY_Create_from_arg(). Все остальное делается так же, как в главной программе.

◄ Пример 8 Пример 10 ►
 
 
 
 
 
 
 
 
  Тел. +7(499)220-79-72; E-mail: inform@kiam.ru