![]() |
||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||
![]() |
||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||
![]() |
RefNUMA – библиотека для организации виртуальной общей памяти в программах, использующих MPI.А. О. Лацис
Пример 7. Модельное приложение – сеточный аналог двумерной задачи для уравнения Пуассона. Описание некоторых дополнительных функций и макросов.Накопленные при рассмотрении предыдущих примеров знания о базовых конструкциях RefNUMA позволяют нам перейти от тривиальных примеров к рассмотрению хотя и очень простого, но законченного модельного приложения. Будем рассматривать параллельную реализацию решения сеточного аналога двумерной задачи Дирихле для уравнения Пуассона на равномерной прямоугольной сетке. Решение сеточной задачи будем искать методом Якоби при заданном числе итераций. Текст последовательной программы, параллельную реализацию которой нам надлежит выполнить, приводится ниже. Обратим внимание на то, что это – именно приводимая с целью объяснения алгоритма последовательная программа. Она оформлена как параллельная, чтобы ее можно было запускать на суперкомпьютере командой mpirun, но по существу является последовательной, то есть должна запускаться только с числом процессов, равным 1. #include <stdio.h> #include <stdlib.h> #include <mpi.h> static void usage( void ) { printf( "Usage: progrev1 <nrows> <ncols> <niter>\n" ); exit( -1 ); } int main( int argc, char *argv[] ) { int mx, my, mx1, my1, niter; void *ra; int i, j, k; FILE *fp; double rdx2, rdy2, beta, t; /***/ MPI_Init( &argc, &argv ); if ( argc != 4 ) usage(); rdx2 = 1.0; rdy2 = 2.0; beta = 0.25; mx = (int)atol( argv[1] ) + 2; my = (int)atol( argv[2] ) + 2; niter = (int)atol( argv[3] ); { double (*f)[my] = (typeof(f))malloc( mx*my*sizeof(double) ); double (*newf)[my] = (typeof(newf))malloc( mx*my*sizeof(double) ); double (*r)[my] = (typeof(r))malloc( mx*my*sizeof(double) ); if ( (!f) || (!newf) || (!r ) ) { printf( "No memory\n" ); exit( -1 ); } /***/ for ( i = 0; i < mx; i++ ) { for ( j = 0; j < my; j++ ) { f[i][j] = newf[i][j] = 0.0; r[i][j] = 0.5; if ( (i == 0) || (i == (mx-1)) || (j == 0) || (j == (my-1)) ) f[i][j] = newf[i][j] = 1.0; } } mx1 = mx - 1; my1 = my - 1; t = MPI_Wtime(); for ( k = 0; k < niter; k++ ) { 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; } } 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 ); fp = fopen( "progrev1.dat", "w" ); for ( i = 1; i < mx1; i++ ) { fwrite( newf[i]+1, sizeof(double), mx-2, fp ); } fclose( fp ); } MPI_Finalize(); return 0; } Для параллельной реализации данной задачи при помощи RefNUMA, сделаем массивы f, newf и r блочно распределенными разделяемыми, то есть распределим их по процессам параллельной программы блоками строк примерно равного размера. Каждый процесс будет обрабатывать свои блоки строк массивов f и newf, то есть присваивать элементам соответствующих строк новые значения. Для этого процессу потребуется доступ на чтение к краям соседних блоков строк, принадлежащих другим процессам, значительно менее интенсивный, чем доступ на чтение и запись к своему блоку. Таким образом, дисциплина доступа к данным в параллельном варианте рассматриваемой модельной задаче полностью соответствует рассмотренным выше рекомендациям по применению RefNUMA, и позволяет вполне естественным образом использовать конструкции параллельного и поочередного циклов. Текст параллельной программы, использующей RefNUMA, приводится ниже. #include <stdio.h> #include <stdlib.h> #include <mpi.h> #include <coarray.h> #include <shared.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] ); { COARRAY_Create_vectorized( double, f, mx, my, 1.0, fo ); COARRAY_Create_vectorized( double, newf, mx, my, 1.0, newfo ); COARRAY_Create_vectorized( double, r, mx, my, 1.0, ro ); if ( (!f) || (!newf) || (!r) ) { printf( "No memory\n" ); exit( -1 ); } /***/ PARALLEL_DO( i, 0, mx, 1, mx ) { 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, (void*)f, (void*)newf, (void*)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" ); PARALLEL_DO( i, 1, mx-1, 1, mx ) { 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**)fo, &im, &j ); printf( "In node %d for f: requested %ld used %ld\n", my_node, im, j ); coarray_report( (void**)newfo, &im, &j ); printf( "In node %d for newf: requested %ld used %ld\n", my_node, im, j ); coarray_report( (void**)ro, &im, &j ); printf( "In node %d for r: requested %ld used %ld\n", my_node, im, j ); fflush( stdout ); NODE_BY_NODE_END } 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; /***/ double **f = (typeof(f))pf; double **newf = (typeof(newf))pnewf; double **r = (typeof(r))pr; mx1 = mx - 1; my1 = my - 1; PARALLEL_DO( i, 1, mx1, 1, mx ) { 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; } } } Обратим внимание на то, что в данной программе никак не контролируется и не обеспечивается делимость нацело количества строк в массивах общей памяти на количество процессов параллельной программы. Это не ошибка – макросы COARRAY_Create_vectorized() и PARALLEL_DO() все сделают правильно в любом случае. Однако, возникает совершенно естественный вопрос о том, что делать программисту, если сложность его программы не позволяет ему ограничиться использованием макроса PARALLEL_DO(). В этом случае программисту придется использовать ряд системных функций и макросов для «ручного» управления распределением работы и данных. Описания этих функций и макросов приводятся ниже. В приведенном только что примере программы эти функции и макросы используются неявно, изнутри COARRAY_Create_vectorized() и PARALLEL_DO(), но ничто не мешает программисту, при необходимости, использовать их явно. Вычисление размера блока разделяемого массива. Функция extern long shared_blocksize( long larr, int n_nodes )возвращает число строк в блоке разделяемого массива общей длиной larr строк для распределения его по n_nodes процессам. Это число строк предполагается использовать впоследствии для указания (после перевода в байты) в качестве значения первого аргумента при вызове coarray_create(). Если число строк массива не делится нацело на число процессов, число строк в блоке берется с минимально необходимым запасом. Определение принадлежности строки разделяемого массива к своей секции (блоку). Функция extern int shared_is_local( long larr, int n_nodes, int my_node, long i )возвращает единицу, если строка i разделяемого массива длиной larr строк, распределенного по n_nodes процессам, попадает в блок номер my_node, или, что есть то же самое, попадает в свою секцию (блок) в процессе номер my_node. Иначе функция возвращает ноль. Предполагается распределение массива по процессам блоками строк такого размера, какой вернуло бы обращение к shared_blocksize( larr, n_nodes ). Вычисление границ изменения индекса при обработке разделяемого массива в параллельном цикле. Для определения нижней и верхней границ цикла PARALLEL_DO() RefNUMA использует функции shared_from() и shared_upto(), соответственно. При необходимости программист может использовать эти функции самостоятельно. Каждая функция имеет четыре аргумента: shared_from( larr, n_nodes, my_node, kbeg )и shared_upto( larr, n_nodes, my_node, kend)причем
Функция shared_from() возвращает нижнюю, а функция shared_upto() – верхнюю границу диапазона изменения переменной параллельного цикла для данного процесса и данной длины массива. Или, что есть то же самое, нижнюю и верхнюю границу пересечения двух диапазонов: диапазона номеров строк массива заданной длины, своих для указанного процесса, и диапазона изменения переменной параллельного цикла. Если диапазон индексов своих для данного процесса строк не пересекается с диапазоном изменения переменной цикла, то пересечение пусто: данному процессу не досталось ни одного витка параллельного цикла. В этом случае shared_from() возвращает значение, равное числу строк в массиве, то есть превосходящее номер любой строки, а shared_upto() возвращает 0, то есть значение, не превосходящее никакого номера строки. Цикл for по диапазону значений переменной цикла с такими границами не выполнится ни разу. Сказанное выше проще всего показать на примере. Пусть в параллельной программе из 4-х процессов имеется массив общей памяти из 20 строк. Очевидно, строки эти распределены по процессам так:
Пусть в программе имеется параллельный цикл вида: PARALLEL_DO( i, 7, 13, 1, 20 ) Это означает, что для значений i от 7 (включительно) до 13 (не включительно), с шагом 1, каждый i-й виток цикла должен выполниться в том процессе, которому принадлежит i-я строка массива длиной 20. Легко видеть, что строки с 7-й по 9-ю принадлежат первому процессу, а строки с 10-й по 12-ю включительно – второму. Таким образом: — нулевой процесс в этом параллельном цикле не делает ничего, — первый процесс выполняет витки цикла с номерами 7, 8 и 9, — второй процесс выполняет витки цикла с номерами 10, 11 и 12, — третий процесс не делает ничего. Именно эти значения границ диапазонов и можно выяснить, обратившись к функциям shared_from(20, 4, N, 7 )и shared_upto(20, 4, N, 13 )с разными значениями N (номера процесса). Возвращаемые этими двумя функциями значения составят:
Таким образом, PARALLEL_DO( i, ib, ie, is, leng )– это в точности то же самое, что и: for( i = shared_from( leng, coarray_n_nodes(), coarray_my_node(), ib ); i < shared_upto( leng, coarray_n_nodes(), coarray_my_node(), ie ); i += is ) Несколько менее эффективный, но по смыслу эквивалентный вариант такого же распределения витков цикла по процессам: for ( i = ib; i < ie; i += is ) { if ( shared_is_local( leng, coarray_n_nodes(), coarray_my_node(), i ) ) { ……… Аргументы shared_from() и shared_upto() позволяют выяснить границы диапазона изменения переменной параллельного цикла для произвольного номера процесса и даже произвольного числа процессов в параллельной программе (напомним, что этот диапазон является также пересечением диапазона номеров строк массива заданной длины, своих для данного процесса, и диапазона изменения переменной параллельного цикла). В то же время, чаще всего (например, в реализации макроса PARALLEL_DO() ) требуется выяснить границы указанного диапазона для реально имеющегося числа процессов в данной параллельной программе, и для данного процесса. Если программисту требуется сделать это в явном виде, он может воспользоваться макросами SHARED_START() и SHARED_FINISH(). Макросы SHARED_START(ib, mx) и SHARED_FINISH(ie,mx) являются сокращенной формой записи shared_from() и shared_upto(), соответственно, при числе процессов, равном реально имеющемуся числу процессов данной параллельной программы, и номере процесса, равном собственному номеру данного процесса. ib и ie суть границы диапазона номеров строк, из которого требуется вырезать диапазон строк, своих для данного процесса. mx – длина массива в строках. Оптимизация барьеров в случае простых шаблонов обращения к чужим данным. В рассматриваемой программе доступ процессов к чужим данным происходит только на чтение: ни один процесс ничего не присваивает ни одному элементу распределенного массива за пределами своего блока. Зная это, программист может (правда, исключительно на свой страх и риск) включить специальный режим несколько более быстрого, чем обычно, выполнения барьеров. Корректность работы в случае, когда программист ошибся, и запись в чужие блоки все же происходит, не гарантируется. Режим ускоренного выполнения барьеров включается присваиванием ненулевого значения специальной системной переменной remoteROM: remoteROM = 1; Описание дополнительных системных функций и макросов, не вошедших в явном виде в приведенные примеры программ, на этом закончено. Также имеет смысл обратить внимание на то, как организована запись результатов счета в файл. Используется, на первый взгляд, довольно странный прием – параллельный цикл внутри поочередного цикла. На самом деле, ничего странного в этой конструкции нет. Как мы уже отмечали выше, при первом знакомстве с конструкцией параллельного цикла, этот вид цикла не подразумевает неявного барьера в конце работы. Параллельный цикл RefNUMA, таким образом, не заставляет процессы работать параллельно, а только предоставляет им такую возможность, вычисляя для каждого процесса поручаемую этому процессу порцию работы. В данном случае нам потребовалось вычислить для каждого процесса ту часть массива, которую этот процесс должен записать в файл. Мы использовали для этого параллельный цикл, поскольку он делает только это и ничего больше. Далее, нам потребовалось сделать так, чтобы процессы выполнили порученную каждому из них работу поочередно, и мы с полным на то основанием использовали специально предусмотренную для этого конструкцию – поочередный цикл. В заключение хотелось бы сделать очень важное замечание. Может возникнуть вполне закономерный вопрос: а зачем вообще при вводе начальных данных (в нашем примере его, правда, нет), а также при выводе результатов (в нашем примере он есть) как-то делить работу между процессами? Если каждый процесс имеет прямой доступ ко всему массиву, почему бы не поручить начальный ввод и заключительный вывод какому-то одному процессу – нулевому, например? Будь наша программа написана на OpenMP, очевидно, ввод и вывод в ней были бы организованы именно так. В нашем же случае мы, к сожалению, не можем себе позволить такого упрощения программы. Как уже говорилось выше, и как подробнее объясняется в соответствующем Приложении, массовое обращение процесса к чужим данным приводит не только к замедлению работы программы (что в случае начального ввода и заключительного вывода не имеет большого значения), но и к быстрому росту затрат физической памяти. В том числе – дополнительной системной памяти, резервируемой обращением к COARRAY_MEM_ALLOC(). Если программист поленится распределить работу по начальному вводу и заключительному выводу между процессами более или менее разумно, то программа с очень большой вероятностью просто не сможет выполняться, и завершится аварийно. ◄ Пример 6 Пример 8 ► |
![]() |
||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||
Тел. +7(499)220-79-72; E-mail: inform@kiam.ru |