Новости
Вычислительные ресурсы
Доступ к ресурсам
Регистрация пользователя
Документация
Исследования
Вопросы - ответы
Контакты

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

А. О. Лацис

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

Гибридный вариант главной программы рассмотренного в предыдущем разделе приложения приводится ниже.

#include <mpi.h>
#include <stdio.h>
#include <stdlib.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;
#ifdef FORACCEL
    long accelstart, accelfinish, accelsize;
#endif
/***/
    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;
#ifdef FORACCEL
      void *af, *anewf, *ar;
#endif
      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 );
       }
#ifdef FORACCEL
      accelstart = SHARED_START( 1, mx );
      accelfinish = SHARED_FINISH( mx-1, mx );
      accelsize = accelfinish-accelstart;
      if ( accelsize <= 0 )
       {
        accelsize = 0;
       }
      else
       {
        accelstart--;
        accelfinish++;
        accelsize += 2;
       }
//      printf( "Process %d: accelstart=%ld accelfinish=%ld accelsize=%ld\n", my_node,
                                                  accelstart, accelfinish, accelsize );
      if ( accelsize )
       {
        af = calloc_accel( accelsize, my*sizeof(double) );
        anewf = calloc_accel( accelsize, my*sizeof(double) );
        ar = calloc_accel( accelsize, my*sizeof(double) );
        if ( (!af) || (!anewf) || (!ar) )
         {
          printf( "No device memory\n" );
          exit( -1 );
         }
       }
#endif
/***/
#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;
         }
       }
#ifdef FORACCEL
      shared_to_accel( af, (void **)f, mx, accelstart, accelsize,
                                                          my*sizeof(double), REALLYCOPY );
      shared_to_accel( anewf, (void **)newf, mx, accelstart, accelsize,
                                                          my*sizeof(double), REALLYCOPY );
      shared_to_accel( ar, (void **)r, mx, accelstart, accelsize,
                                                          my*sizeof(double), REALLYCOPY );
#endif
      coarray_barrier();
      t = MPI_Wtime();
      for ( k = 0; k < niter; k++ )
       {
#ifdef FORACCEL
	if ( accelsize )
	 {
	  shared_to_accel( af, (void **)f, mx, accelstart, accelsize,
                                           my*sizeof(double),  JUSTUPDATE );
	  itstep( accelsize, my, af, anewf, ar, rdx2, rdy2, beta );
	  shared_from_accel( (void **)newf, mx, accelstart, anewf, accelsize,
                                                my*sizeof(double),  JUSTUPDATE );
         }
#else
	itstep( mx, my, f, newf, r, rdx2, rdy2, beta );
#endif
        coarray_barrier();
        ra = (void*)f;
        f = newf;
        newf = (typeof(newf))ra;
#ifdef FORACCEL
        ra = (void*)af;
        af = anewf;
        anewf = (typeof(anewf))ra;
#endif
       }
      t = MPI_Wtime() - t;
      printf( "Time: %f speed %f flops\n", t,
                  (((double)(mx-2))*((double)(my-2))*7.0*niter)/t );
#ifdef FORACCEL
      shared_from_accel( (void **)f, mx, accelstart+1, ((double*)af)+my, accelsize-2,
                                        my*sizeof(double), REALLYCOPY );
#endif
      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;
}
#ifndef FORACCEL
    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
}
#endif

По сравнению с рассмотренными ранее примерами, в этом тексте интерес представляют впервые появившиеся функции shared_to_accel() и shared_from_accel(). Это стандартные функции копирования данных между процессорной и сопроцессорной памятью, очень похожие на уже рассмотренные выше функции vectorized_to_accel() и vectorized_from_accel(). Отличие состоит в отсутствии аргумента, задающего длину строки массива общей памяти. При использовании этих функций следует обратить внимание на правильность задания значения аргумента «длина элемента массива в байтах». Для функций vectorized…() под элементом понимался элемент двумерного массива, то есть, в нашем примере, элемент имееет тип double, размер элемента равен sizeof(double). Для функций shared…() под элементом понимается элемент одномерного массива, то есть, в нашем примере, элемент имеет тип element, размер элемента равен sizeof(element), или, что есть то же самое, my*sizeof(double).

◄ Пример 9 Приложение 1 ►
 
 
 
 
 
 
 
 
  Тел. +7(499)220-79-72; E-mail: inform@kiam.ru