![]() |
||||||||||||||||||||||||||||||||||
![]() |
||||||||||||||||||||||||||||||||||
![]() |
Инженерная методика адаптации приложения к гибридному кластеру с ускорителями на ПЛИСЧасть 5. Реализация итерационного шага вычислений с использованием буферных массивов.С. С. Андреев, С. А. Дбар, А. О. Лацис, Е. А. Плоткина
Прежде, чем переходить к разбиению обрабатываемых данных на мелкие блоки, попробуем записать саму обработку более удобным и эффективным способом. Как мы видели выше, использование двух массивов значений температуры, меняющихся ролями на каждой итерации, довольно неудобно, не говоря уже о двойном расходе дефицитной внутренней памяти сопроцессора. Кроме того, из самых общих соображений о том, как устроены эффективные сопроцессоры, реализуемые в ПЛИС, довольно очевидно, что функция itstep() у нас написана довольно неэффективным способом. В самом деле, все старые значения температуры выбираются из единого массива f. Как мы узнаем позже, главный источник ускорения обработки в ПЛИС - это усовершенствование именно хранения и перемещения данных. Если мы учтем также, что каждый массив, в терминах сопроцессора, это, как правило, отдельное электронное устройство, работающее независимо от других и параллельно с ними, то нам захочется переписать эту функцию в несколько другой форме. Перебирая строки массива старых значений температуры, мы используем каждую строку трижды. Сначала - как "верхнюю", или следующую, потом как "среднюю", или текущую, и, наконец, как "нижнюю", или предыдущую. Выделив для "средней" и "нижней" строк отдельные одномерные массивы, мы не только создадим возможность для последующей более эффективной реализации алгоритма в ПЛИС, но и избавимся от необходимости иметь два меняющихся местами массива температур. В самом деле, два массива температур нужны были нам для того, чтобы, согласно итерационному процессу Якоби, использовать для вычисления нового значения температуры в каждой ячейке сетки старые, еще не обновленные, значения всех соседей. Легко видеть, что в приведенном ниже варианте функции itstep() это условие соблюдено, причем с использованием единственного массива температур: #include <stdio.h> #include "dimension.h" #define MAXMY 10000 int itstep_test( int my_number, int n_of_nodes, int mx, int my, int stw ) { if ( my > MAXMY ) { fprintf( stderr, " my (%d) greater than maxmy (%d)\n", my, MAXMY ); return 1; } else return 0; } void itstep( int mx, int my, void *pf, void *pr, double rdx2, double rdy2, double beta, int niter ) { DIM2( double, f, my ) = (typeof(f))pf; DIM2( double, r, my ) = (typeof(r))pr; static double lower[MAXMY]; static double middle[MAXMY]; int i, j, k, mx1, my1; /***/ mx1 = mx - 1; my1 = my - 1; for ( k = 0; k < niter; k++ ) { for ( i = 0; i < mx1; i++ ) { for ( j = 0; j < my; j++ ) { lower[j] = middle[j]; middle[j] = f[i][j]; } if ( i ) { for ( j = 1; j < my1; j++ ) { f[i][j] = ((lower[j]+f[i+1][j])*rdx2+(middle[j-1]+middle[j+1])*rdy2-r[i][j]) *beta; } } } } } Обратим внимание на то, что к функции itstep() добавилась новая функция itstep_test(). Эта функция будет вызвана вызывающей программой однократно, в начале работы, и должна выполнить все однократные начальные проверки и/или настройки, необходимые для дальнейшего функционирования itstep() и, в перспективе, реализующего ее сопроцессора. В данном случае эта функция выполняет единственную проверку - не больше ли размер строки обрабатываемого массива значений температуры, чем размер буферных массивов для "нижней" и "средней" строки, используемых в itstep(). Можно было бы потребовать, чтобы размер строки массива обрабатываемых данных не сравнивался со статически заданной константой, а использовался для динамического выделения памяти именно такого размера под буферные массивы. Однако, не надо забывать и о том, что мы постепенно готовимся к аппаратной реализации itstep() в качестве сопроцессора, а внутри сопроцессора, как легко догадаться, никакого динамического выделения памяти не бывает. Отметим также, что к оформлению проверок корректности, даже тривиальных, в виде отдельных функций, не следует относиться пренебрежительно. По мере подготовки к аппаратной реализации сопроцессора мы будем сталкиваться с постоянно растущим набором условий его корректного использования, которые очень легко сначала - забыть проверить, затем - невольно нарушить, и, наконец, долго и безуспешно искать ошибку там, где ее не было и нет. С другой стороны, засорение текста программы более высокого уровня (в нашем случае - main()) хаотичным набором проверок корректности использования гораздо более низкоуровневых функций - тоже не очень хорошая идея. Ведь, по мере развития и совершенствования этих самых низкоуровневых функций, конкретный набор проверок корректности их использования также может меняться! Поэтому, элементарные принципы правильной программной инженерии подсказывают нам выделять эти проверки в специальные функции, а текст этих функций располагать в файлах именно с тем исходным текстом, корректность использования которого проверяется. Далее приводится модифицированный, согласно описанным только что изменениям в вызываемых функциях, текст функции main(): #include <mpi.h> #include <stdio.h> #include <stdlib.h> #include <string.h> #include "vsda.h" #include "dimension.h" // Change here the number of steps, the cell geometry, etc #define STEPITER 1000 #define delx 0.5 #define dely 0.25 #define NITER 10002 #define STW 3 // end change here. /***/ extern void itstep( int mx, int my, void *f, void *r, double rdx2, double rdy2, double beta, int niters ); extern int itstep_test( int my_number, int n_of_nodes, int mx, int my, int niter ); int main( int argc, char **argv ) { int i, j, k, n, mx1, mx2, my1, my_number, n_of_nodes, totalmx, partmx, leftmx, mx, my, niter, stw; FILE *fp; double howlong; double rdx2, rdy2, beta; VSDA vsda; long *distr_out; MPI_Request rq[5]; /***/ MPI_Init( &argc, &argv ); MPI_Comm_rank( MPI_COMM_WORLD, &my_number ); MPI_Comm_size( MPI_COMM_WORLD, &n_of_nodes ); if ( argc != 3 ) { if ( !my_number ) { fprintf( stderr, "Usage: %s <nrows> <ncolumns>\n", argv[0] ); fflush( stderr ); } return (-1); } totalmx = mx = (int)atol( argv[1] ); my = (int)atol( argv[2] ); if ( my < 1 ) { if ( !my_number ) { fprintf( stderr, "Number of columns (%d) should be positive\n", my ); fflush( stderr ); } return (-1); } distr_out = (long*)malloc( n_of_nodes*sizeof( *distr_out ) ); if ( !distr_out ) { fprintf( stderr, "No memory\n" ); fflush( stderr ); return( -1 ); } totalmx += 2; my += 2; stw = STW; niter = NITER; /* Make NITER divisible by STW: */ niter -= (niter%stw); if ( itstep_test( my_number, n_of_nodes, mx, my, stw ) ) { fprintf( stderr, "itstep_test not passed, exiting\n" ); MPI_Finalize(); return -1; } /* Here we know the array width, so make the arrays themselves: */ { DIM2( double, f, my ); DIM2( double, r, my ); /* Make the distributed array, discover the local portion length: */ VSDA_Create( MPI_COMM_WORLD, VSDA_BLOCKED, 0, (long)totalmx, (long)sizeof(*f), (long)stw, NULL, distr_out, &vsda ); mx = distr_out[my_number]; if ( n_of_nodes > 1 ) { if ( my_number > 0 ) mx += stw; if ( my_number < (n_of_nodes-1) ) mx += stw; } f = (typeof(f))malloc( mx*sizeof(*f) ); r = (typeof(r))malloc( mx*sizeof(*r) ); if ( (!f) || (!r) ) { fprintf( stderr, "Cannot allocate, exiting\n" ); exit( -1 ); } rdx2 = 1./delx/delx; rdy2 = 1./dely/dely; beta = 1.0/(2.0*(rdx2+rdy2)); if ( !my_number ) { printf( "Solving heat conduction task on %d by %d grid by %d processors\n", totalmx-2, my-2, n_of_nodes ); fflush( stdout ); } for ( i = 0; i < mx; i++ ) { for ( j = 0; j < my; j++ ) { if ( ((i==0) && (my_number==0)) || (j==0) || ((i==(mx-1)) && (my_number==(n_of_nodes-1))) || (j==(my-1)) ) f[i][j] = 1.0; else f[i][j] = 0.0; r[i][j] = 0.0; } } /* Iteration loop: */ howlong = MPI_Wtime(); for ( n = 0; n < niter/stw; n++ ) { if ( !my_number ) { if ( !(n%STEPITER) ) printf( "Iteration %d\n", n ); } /* Step of calculation starts here: */ itstep( mx, my, f, r, rdx2, rdy2, beta, stw ); /* Do all the transfers: */ VSDA_Update_begin( &vsda, rq, f ); VSDA_Update_end( rq ); } if ( !my_number ) { printf( "Elapsed time: %f sec\n", (float)(MPI_Wtime()-howlong) ); fp = fopen( "output_par_3.dat", "w" ); fclose( fp ); } for ( j = 0; j < n_of_nodes; j++ ) { MPI_Barrier( MPI_COMM_WORLD ); if ( j == my_number ) { int istart, isize; istart = 1; isize = mx-1; if ( n_of_nodes > 1 ) { if ( my_number > 0 ) istart = stw; if ( my_number < (n_of_nodes-1) ) isize = mx-stw; } fp = fopen( "output_par_3.dat", "a" ); for ( i = istart; i < isize; i++ ) { fwrite( &(f[i][1]), my-2, sizeof(f[0][0]), fp ); } fclose( fp ); } } } MPI_Finalize(); return 0; } По сравнению с предыдущим вариантом, исчез массив newf, а также требование, чтобы значение STW было нечетным. Появилось обращение к itstep_test(). Сохранившееся требование делимости числа итераций на STW, в принципе, не так уж сложно снять, но заниматься этим в этом документе мы не будем. ◄ Часть 4 Часть 6 ► |
![]() |
||||||||||||||||||||||||||||||||
Тел. +7(499)220-79-72; E-mail: inform@kiam.ru |