|
Система параллельного программирования RefNUMA.Руководство программиста.А. О. Лацис Введение.Отсутствие полноценной, адресуемой простым машинным указателем, общей памяти, разделяемой ветвями параллельной программы, было и остается важнейшим источником сложности разработки параллельных программ для вычислительных систем с массовым параллелизмом. Описываемая ниже система параллельного программирования RefNUMA призвана восполнить этот пробел применительно к достаточно широкому кругу параллельных приложений. Аппаратная реализация общей памяти в системах массового параллелизма ограничивается обычно пределами одного вычислительного узла. Предоставление программисту возможности доступа к общей памяти, разделяемой процессами на разных узлах, возможно путем реализации системы зеркал – копий фрагментов «чужой» памяти, поддерживаемых ОС в каждом узле. Приведение всех копий в когерентное состояние происходит обычно во время барьерной синхронизации. Наиболее широкую известность получили 2 системы, основанные на технологии такой зеркалированнной, или рефлективной, общей памяти. Система vSMP фирмы ScaleMP [1] предоставляет рефлективную общую память уже на уровне ОС. С использованием такой системы на кластере, связанном сетью Ethernet, можно запустить единую многопроцессорную ОС, покрывающую все узлы, так, как если бы все процессоры кластера входили в единый многопроцессорный сервер. Параллельные приложения под управлением такой ОС разрабатываются на OpenMP, и так же, как и сама ОС, работают на всех (или некоторых) вычислительных узлах кластера как на едином сервере. Система Cluster OpenMP фирмы Intel [2] построена по тому же принципу, но предоставляет рефлективную общую память на уровне отдельной параллельной программы. Программа разрабатывается с использованием немного модифицированного варианта OpenMP, а запускается на кластере как параллельная программа для кластера, но самой программе во время выполнения «кажется», что она представляет собой единый процесс, работающий на едином многопроцессорном сервере с общей памятью. Не секрет, что обе системы на практике демонстрируют довольно низкую эффективность. Параллельные приложения, написанные для этих систем, часто плохо работают и еще хуже масштабируются. Чтобы понять, почему это так, рассмотрим тривиальную, но формально корректную реализацию рефлективной памяти. Рассуждения для краткости проведем на примере единственного разделяемого всеми процессами массива. Пусть каждый процесс хранит в своей локальной памяти два экземпляра всего разделяемого (общего) массива: оригинал и копию. По мере работы процесса содержимое оригинала, естественно, меняется. Во время барьерной синхронизации, которая неизбежно задается программистом явно во всех программах, работающих с общей памятью, каждый процесс сравнивает оригинал с копией, а информацию обо всех обнаруженных различиях рассылает оригиналам всех остальных процессов, попутно обновляя копию. На следующем барьере взаимная рассылка накопившихся со времени прошлого барьера изменений повторяется, и т. д. Вполне очевидно, что такая тривиальная реализация вопиюще неэффективна по трем причинам:
Все предлагаемые для практического применения системы рефлективной памяти используют те или иные алгоритмы (довольно сложные) для сокращения накладных расходов по каждому из трех пунктов, то есть хранят не все, сравнивают не все и рассылают не всем. Обнаруживаемая на практике низкая эффективность работы означает, что или достигнутая «экономия» недостаточна, или использованная для реализации коммуникационная сеть недостаточно хороша в получающихся при таком подходе режимах работы. Описываемая ниже система RefNUMA призвана восполнить оба этих недостатка. Во-первых, ее реализация базируется на библиотеке shmem [3], обычно весьма эффективной в режимах нерегулярной рассылки коротких сообщений. Во-вторых, главным источником повышения эффективности рефлективной общей памяти остается сокращение накладных расходов, то есть совершенствование самого алгоритма обеспечения когерентности «зеркал». Если давно разработанные и хорошо изученные алгоритмы не обеспечивают такого сокращения в должной мере, то для решения проблемы остаются всего 2 пути. Первый путь заключается в том, чтобы потребовать от программиста указания в программе, в том или ином виде, дополнительной информации о характере доступа процессов в общую память, чтобы система могла выполнить более глубокую оптимизацию хранения и рассылок. Например, некоторый процесс может «обещать» вообще не трогать некоторую часть общего массива, или не изменять ее. По вполне очевидным причинам, идти на такое требование к программисту не желательно. Второй путь заключается в том, чтобы выбрать достаточно широкий класс задач, отличающихся некоторым «разумным», с точки зрения характера доступа в общую память, поведением, и оптимизировать систему под именно такое «разумное» (в среднем) поведение прикладных программ. Формальная работоспособность при этом обеспечивается для любых программ, но степень эффективности (точнее – степень отклонения эффективности от качественной реализации с помощью MPI) будет зависеть от степени «разумности» поведения конкретной программы. Крайне важно при этом добиться того, чтобы выбранный критерий «разумности» не только допускал действительно глубокую оптимизацию работы системы, но и был понятен прикладному программисту, а степень соответствия конкретной программы этому критерию легко оценивалась бы программистом из самых общих соображений. Описываемая ниже система RefNUMA идет в реализации рефлективной общей памяти именно по этому пути. 1. Процессы и адресные пространства.Система RefNUMA, как и Cluster OpenMP, предоставляет возможность использования разными процессами единой общей памяти в рамках конкретной программы. Сама программа при этом напоминает скорее программу, написанную с использованием MPI, нежели программу на OpenMP. Программа запускается командой mpirun и состоит из заданного пользователем, фиксированного количества параллельных ветвей. Каждая ветвь – отдельный, полноценный процесс, со своим, отдельным от других, комплектом переменных и массивов. Массивы общей памяти, доступные всем процессам, требуется создавать в программе специально, обращаясь к соответствующим системным функциям. Для синхронизации доступа к таким массивам процессы используют операцию барьерной синхронизации. Изменения, вносимые тем или иным процессом в массивы общей памяти, не сразу становятся известны другим процессам. Доведение этих изменений «до сведения» других процессов обычно происходит во время выполнения барьера. Сразу после барьера все процессы видят в каждом из массивов общей памяти одинаковое содержимое. Всеми этими качествами система RefNUMA напоминает такие языки параллельного программирования, как UPC [4] и Co-array Fortran [5]. Однако, оформлена она как библиотека, а не как специализированный язык. Форма записи, как мы увидим ниже, специально приспособлена для наиболее безболезненного использования библиотеки совместно с OpenMP. В настоящее время RefNUMA может быть использована в программах на С (С++). 2. Секционированные и разделяемые массивы.В основе предлагаемой модели общей памяти лежит понятие NUMA (Non-Uniform Memory Access). В своем наиболее общем виде это понятие сводится к следующему. Общая, доступная всем процессам, память делится на блоки, каждый из которых реально расположен в некотором процессе (является для него своим), а прочим процессам доступен как чужой. Обращения к данным в чужом блоке выполняются медленнее (иногда – в десятки и сотни раз), чем к своим данным. Программа будет работать эффективно при условии, что обращений к чужим данным в ней гораздо меньше, чем к своим. Система RefNUMA реализует довольно своеобразный вариант этой концепции. С одной стороны, сами по себе обращения к чужим данным в этой системе выполняются так же быстро, как и к своим (за исключением первого за время работы программы обращения по конкретному чужому адресу, которое может быть очень медленным). С другой стороны, чем больше становится объем чужих данных, к которым программа успела обратиться за время своей работы, тем медленнее будет впоследствии выполняться барьерная синхронизация, и тем больше система вынуждена будет использовать внутренней буферной памяти для организации «зеркал». Интенсивно обращающаяся по все новым и новым адресам чужой памяти программа постепенно обнаружит, что барьерная синхронизация стала очень медленной, и, быть может, даже не сможет в какой-то момент продолжать работу из-за нехватки буферной памяти. Если же программа работает в основном со своей частью общей памяти, обращаясь лишь к ограниченному объему одних и тех же чужих данных, то все будет хорошо: нарастание времени барьерной синхронизации и объема используемой буферной памяти в некоторый момент времени прекратится на приемлемом уровне. Конкретно, модель общей памяти RefNUMA основана на понятии секционированного массива (co-array). Секционированный массив состоит из секций равного размера, по числу процессов в параллельной программе. Секция номер N располагается в процессе номер N (является для него своей), остальные секции для этого процесса являются чужими. Функция создания секционированного массива возвращает массив указателей на секции, в порядке их номеров. Например, если этот массив называется F, а каждая секция есть двумерный массив, то запись F[N][I][J] в программе означает обращение к элементу с индексами (I,J) из секции номер N. Весь массив, таким образом, оказывается трехмерным, но первое измерение – особое. Индекс по этому измерению соответствует номеру процесса, для которого данная секция – своя. Здесь просматривается прямая аналогия с Co-array Fortran – только там «процессная размерность» последняя, а не первая, и указывается, в отличие от «обычных» размерностей, в не свойственных Фортрану квадратных скобках. Обращение к элементам своей секции, как мы уже говорили выше, «бесплатно» для процесса во всех отношениях. Факт такого обращения не увеличивает время выполнения последующих барьеров, и не требует каких-либо дополнительных системных ресурсов. Иначе обстоит дело, когда процесс обращается к элементу чужой для него секции секционированного массива. Каждое такое обращение к новым элементам массива, к которым данный процесс еще не обращался, увеличивает как время последующих барьеров, так и потребности в системной памяти для их обслуживания. Если программист готов писать программу непосредственно в терминах секционированных массивов, то сказанное выше – это практически все, что ему (программисту) необходимо знать о системе RefNUMA. Для многих программ, тем не менее, непосредственное использование секционированных массивов неудобно. Например, для программиста может показаться естественным организовать одномерный массив общей памяти, блочно распределенный по процессам, наподобие секционированного массива, но с «плавным» переходом индексации с блока на блок. В таком массиве секция номер 0 состоит, например, из элементов с 0-го по 99-й, секция номер 1 – из элементов с 100-го по 199-й, и так далее. Специальной «процессной» размерности в этом случае, как мы видим, нет. Зато первая из «содержательных» размерностей (в нашем случае - единственная) блочно распределена по процессам. Здесь мы видим прямую аналогию с разделяемыми массивами UPC. RefNUMA реализует такие блочно распределенные разделяемые (shared) массивы как «надстройку» над секционированными массивами. Создав секционированный массив, программа может впоследствии рассматривать его либо как секционированный, либо как разделяемый (либо то так, то эдак, что, по понятным причинам, не рекомендуется). Для обращения к массиву общей памяти как к разделяемому используются специальные системные макросы. Из сказанного в этом разделе становится ясно, для каких программ годится предлагаемая система. Пусть параллельную программу удается написать в терминах блочного распределения обработки, так, чтобы каждый процесс обращался преимущественно к своим блокам некоторых распределенных массивов общей памяти, а данные из чужих блоков использовал сравнительно редко, и почти всегда примерно одни и те же. Такая программа, скорее всего, будет работать на RefNUMA почти так же хорошо, как если бы ее реализовали с использованием MPI. При этом совершенно не имеет значения ни регулярность расположения чужих данных в памяти, ни «точность попадания» локальной обработки в блок. Указанное выше «разумное поведение» должно соблюдаться в среднем за время работы программы, отдельные отклонения допускаются. Категорически не годятся для реализации с использованием нынешней версии RefNUMA программы, в которых обращения каждого из процессов к массивам общей памяти не локализовано, либо постоянно меняется со временем. Например, программы расчетов на адаптивных сетках. Впоследствии планируется реализация версии, допускающей сравнительно редкие изменения шаблона обращения к общей памяти в моменты, явно указываемые в программе программистом. Например, программы расчетов на адаптивных сетках, перестраиваемых в редкие, четко определенные в программе моменты времени. 3. Пример использования.Рассмотрим пример параллельной реализации модельной задачи – решение сеточного аналога задачи Дирихле итерационным методом Якоби с фиксированным числом итераций. Текст последовательного варианта программы приводится ниже. #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)) ) Программа использует инициализацию и терминацию MPI, чтобы ее можно было запускать на кластере, а также функцию MPI_Wtime() для измерения времени, но является по смыслу последовательной. В ней нет обменов данными между процессами, и предполагается, что процесс запущен всего один. Рассмотрим параллельную реализацию этой программы с использованием RefNUMA. #include <stdio.h> #include <stdlib.h> #include <mpi.h> #include "coarray.h" #include "shared.h" COARRAY_MEM_ALLOC( 100000000l ); static void usage( void ) { printf( "Usage: progrev_coarr <nrows> <ncols> <niter>\n" ); exit( -1 ); } int main( int argc, char *argv[] ) { int mx, my, mx1, my1, niter, n_nodes, my_node; void *ra; long im, j, k, nblock; 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]; element **f, **newf, **r; nblock = shared_blocksize( (long)mx, n_nodes ); f = (element**)coarray_create( nblock*sizeof(element), nblock*sizeof(element) ); newf = (element**)coarray_create( nblock*sizeof(element), nblock*sizeof(element) ); r = (element**)coarray_create( nblock*sizeof(element), nblock*sizeof(element) ); if ( (!f) || (!newf) || (!r) ) { printf( "No memory\n" ); exit( -1 ); } #define F(i) shared(f, nblock, i) #define NEWF(i) shared(newf, nblock, i) #define R(i) shared(r, nblock, i) #define IFROM( i ) shared_from( mx, n_nodes, my_node, (i) ) #define IUPTO( i ) shared_upto( mx, n_nodes, my_node, (i) ) /***/ for ( i = IFROM( 0 ); i < IUPTO( 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)) ) Файл заголовков "coarray.h" - это главный заголовочный файл библиотеки RefNUMA. Файл заголовков "shared.h" необходим, если предполагается использовать не только секционированные, но и разделяемые массивы, и некоторый дополнительный сервис, связанный с ними. Строка COARRAY_MEM_ALLOC( 100000000l ); предписывает системе выделить указанный в скобках объем памяти (в байтах) в каждом процессе для реализации массивов общей памяти. В этот объем должны вписаться (суммарно) все локальные (свои) секции всех массивов общей памяти, причем те их части, к которым будут обращаться другие процессы, будут представлены в двух экземплярах. Также сюда входят копии тех частей чужих секций, к которым будет обращаться данный процесс, и некоторый не очень большой «довесок». В процессе работы программы исчерпание этого буфера памяти проверяется. Как и в случае последовательного варианта, программа лишь формально использует MPI, только для инициализации, терминации и измерения времени. Для корректной инициализации и терминации всех используемых внутри RefNUMA коммуникационных библиотек, включая MPI, независимо от конкретного варианта параллельной системы программирования, вместо функций MPI_Init() и MPI_Finalize() следует использовать макросы COARRAY_Init() и COARRAY_Finalize() Строки my_node = coarray_my_node(); n_nodes = coarray_n_nodes(); означают извлечение из системы «стандартных координат» процесса – собственного номера и общего числа процессов. Они совпадают с соответствующими «координатами» процесса в MPI-коммуникаторе MPI_COMM_WORLD. В приведенном ниже фрагменте программы создаются три секционированных массива общей памяти – f, newf и r. Элементами массива будут строки сетки и правых частей, размер каждой строки – my элементов типа double typedef double element[my]; element **f, **newf, **r; nblock = shared_blocksize( (long)mx, n_nodes ); f = (element**)coarray_create( nblock*sizeof(element), nblock*sizeof(element) ); newf = (element**)coarray_create( nblock*sizeof(element), nblock*sizeof(element) ); r = (element**)coarray_create( nblock*sizeof(element), nblock*sizeof(element) ); if ( (!f) || (!newf) || (!r) ) { printf( "No memory\n" ); exit( -1 ); } Системная функция shared_blocksize() вычисляет подходящий размер блока (секции) для разделяемого массива длиной mx элементов, при распределении его по n_nodes процессам. Затем создаются секционированные массивы общей памяти, каждая секция которых содержит вычисленное таким образом число элементов (строк по my значений double в каждой). За счет того, что mx может не делиться нацело на n_nodes, размер секционированного массива получится с некоторым «запасом». Функция создания секционированного массива общей памяти coarray_create() реализует коллективную операцию, все процессы должны вызвать ее, создавая один и тот же массив, логически одновременно, и с одними и теми же значениями аргументов. Требуемый размер секции в байтах указывается первым аргументом при обращении к coarray_create(). Вторым аргументом указывается размер дополнительной памяти для системных нужд (для организации «зеркал» создаваемого массива). Необходимый объем этой дополнительной системной памяти трудно вычислить заранее. Как уже неоднократно говорилось выше, он определяется объемом адресуемой процессами задачи «чужой» памяти. Обычно имеет смысл при первом запуске программы заказать под системные нужды столько же памяти, каков размер секции массива, а в конце программы спросить у системы, сколько ее реально было израсходовано, чтобы при последующих запусках программы, если надо, «уменьшить аппетиты». Поскольку память под системные нужды выделяется в нескольких разных местах страницами по 4К байт, при очень маленьких, «игрушечных» размерах секции массива, необходимых иногда на ранних стадиях отладки, системной памяти может потребоваться в несколько раз больше, чем размер секции. При реальной работе ее обычно требуется во много раз меньше ее (секции) размера. Строки #define F(i) shared(f, nblock, i) #define NEWF(i) shared(newf, nblock, i) #define R(i) shared(r, nblock, i) #define IFROM( i ) shared_from( mx, n_nodes, my_node, (i) ) #define IUPTO( i ) shared_upto( mx, n_nodes, my_node, (i) ) сосредотачивают в себе практически весь сервис по использованию разделяемых массивов в данной программе. Первые 3 строки используются ниже для адресации секционированных массивов как разделяемых. Проведем подробное рассмотрение на примере массива f Системный макрос shared() из заголовочного файла «shared.h» позволяет адресовать секционированный массив, указанный в качестве его первого аргумента, как разделяемый. Размер блока указывается вторым аргументом, а индекс (считая от начала массива) – третьим. Таким образом, всюду в программе, где нам нужен j-й элемент i-й строки разделяемого массива f, распределенного блоками размера nblock, мы могли бы писать примерно так: shared( f, nblock, i )[j] = 0.33; Чтобы не писать так много лишнего, мы можем спрятать эту запись в собственный (не системный, а наш, присущий именно этой программе) макрос F, чтобы далее в программе иметь право писать просто: F(i)[j] = 0.33; Теперь рассмотрим две завершающих строки приведенного выше фрагмента программы. Использованные в них обращения к системным функциям shared_from() и shared_upto() служат для вычисления границ индексов в параллельных циклах, обрабатывающих разделяемые массивы. Пусть в последовательной программе имеется цикл вида: for ( i = istart; i < ifinish; i++ ) который используется для прохода по некоторым массивам с индексом i. Будем также предполагать, что витки этого цикла независимы. Очевидно, в параллельной программе, в которой эти массивы станут разделяемыми, витки этого цикла потребуется распределить по процессам таким образом, чтобы для каждого значения i соответствующий виток цикла исполнялся в процессе, для которого i-е элементы обрабатываемых массивов являются своими. Поскольку массивы наши распределены блочно, витки цикла тоже распределятся блочно, то есть в каждом процессе индекс i будет пробегать некоторый сплошной диапазон значений. Функции shared_from() и shared_upto() как раз и вычисляют границы этого диапазона индексов для заданного размера массива, общего числа процессов и номера данного процесса. Как и в случае с макросом shared(), мы для сокращения текста программы и повышения его наглядности «оборачиваем» обращения к функциям shared_from() и shared_upto() в свои собственные, только для данной программы предназначенные, макросы. Далее, вплоть до строки с печатью времени выполнения критического цикла и достигнутого полезного быстродействия, смысл текста программы довольно очевиден. От последовательного прототипа он отличается лишь использованием того сервиса по работе с разделяемыми массивами, который мы обсудили только что. При записи насчитанных разными процессами фрагментов сетки в общий файл с результатами расчет возникает хорошо известная проблема синхронизации доступа процессов к этому файлу. Первоначально предполагалось рекомендовать пользователям RefNUMA использовать для ввода-вывода библиотеку MPI-IO [6], которая разработана специально для решения этой проблемы. Однако, выяснилось, что в MPI-IO отсутствует форматный ввод-вывод (есть аналог fwrite(), но нет аналога fprintf()). Поскольку самому автору неоднократно приходилось, в процессе отладки RefNUMA и настоящего примера, неоднократно заменять fwrite() на fprintf(), решено было разработать собственное средство для записи в файл массивов общей памяти. Это средство можно назвать «последовательной формой параллельного цикла». Заголовком этого цикла является строка NODE_BY_NODE_BEGIN( k, 0, n_nodes ) а концом (закрывающей операторной скобкой) – строка NODE_BY_NODE_END Тело цикла, расположенное между этими двумя строками, выполняется по k от 0 до n_nodes, причем каждый его виток выполняется только в процессе номер k. Витки цикла выполняются строго последовательно: сначала – нулевой виток в нулевом процессе, затем – первый в первом и т. д. Естественно, данную конструкцию можно использовать не только для записи в файл. Выше мы уже говорили о том, что в конце выполнения программы имеет смысл поинтересоваться у системы, сколько реально дополнительной системной памяти было использовано для нужд каждого из наших массивов, чтобы впоследствии уменьшить значение второго аргумента coarray_create() до некоторых разумных величин. Для этого в конце программы используются обращения к функции coarray_report() 4. Описания функций и макросов.4.1. Работа с секционированными массивами.Для работы с секционированными массивами (и, тем самым, для работы с системой вообще) программа должна включать заголовочный файл "coarray.h". Если в программе предусмотрено хотя бы формальное использование MPI, обращение к MPI_Init() должно предшествовать любым обращениям к функциям RefNUMA. Пользователю следует постоянно помнить, что все целочисленные размеры, смещения, индексы и т. п. в системе RefNUMA имеют тип long, а не int. Тип int имеют только номера процессов и их количество. Это важно, если в функцию передается адрес возвращаемого параметра. 4.1.1. Выделение памяти для организации секционированных массивовCOARRAY_MEM_ALLOC( long N ); Выделить для размещения массивов общей памяти в данном процессе N байт памяти. Данная директива является неисполняемым макросом. Должна присутствовать в одном и только одном файле с исходным текстом программы пользователя вне тел функций, в неисполняемой части текста. В выделяемой памяти располагаются (и должны умещаться): - локальные (свои) секции всех секционированных массивов, - дополнительная системная память для реализации секционированных массивов, заказываемая указанием второго аргумента функции coarray_create() (см. ниже), - системные таблицы, размер которых для каждого массива приблизительно равен двум процентам произведения длины секции массива на число процессов. 4.1.2. Общее управлениеint coarray_my_node( void )– узнать собственный номер процесса. int coarray_n_nodes( void )– узнать общее число процессов. void coarray_barrier( void )– выполнить барьерную синхронизацию. Изменения, внесенные процессом в секционированный массив, становятся видны другим процессам в результате барьерной синхронизации. Сразу после барьера каждый секционированный массив выглядит одинаково с точки зрения любого процесса. 4.1.3. Создание секционированного массиваvoid **coarray_create( long n, long m )– создать секционированный массив с длиной секции n байтов, отведя для системных нужд дополнительно m байтов. Функция возвращает адрес массива указателей на секции созданного секционированного массива или 0, если создать его не удалось. Например: double **f; ... f = (double **)coarray_create( n, m ); f[i][j] = 0.5; // j-му элементу i-й секции массива f присвоить значение 0.5 Создание секционированного массива – коллективная операция, выполняемая во всех процессах логически одновременно и с одними и теми же значениями аргументов. Выбирая значение m (объем дополнительной системной памяти для реализации массива), следует учитывать следующее. Дополнительная системная память не требуется, если каждый процесс обращается только к элементам своей секции массива (той, номер которой равен номеру процесса). При первом обращении программы к любой странице какой-либо чужой секции (размер страницы – 4КБ) происходит выделение одной страницы дополнительной системной памяти в обратившемся процессе, и еще одной страницы – в процессе, для которого эта секция является своей. Естественно, точно подсчитывать число требуемых страниц программист не должен. Рекомендуется очень грубо оценить этот объем, выделить дополнительную системную память с запасом, а затем, в конце работы программы, узнать, сколько ее действительно потребовалось, чтобы проверить справедливость сделанных оценок. Первое обращение процесса к ранее не тронутой этим процессом чужой странице выполняется очень медленно, все остальные обращения – с нормальной скоростью обращения к памяти. Тип данных секции секционированного массива может быть любым, по усмотрению программиста, но следует учитывать скважность доступа. Например, пусть скважность доступа в реализации библиотеки равна одному машинному слову (4 байта). Пусть два разных процесса модифицируют разные байты одного и того же слова. Тогда никакой третий процесс не получит после барьера в своем «зеркале» правильного значения этого слова. Он получит либо целое слово из одного обратившегося к данному слову процесса, либо целое слово из другого, но никогда не комбинацию байтов, совместно сформированную в одном слове двумя процессами. В настоящее время скважность доступа равна двойному слову (значение double или long, 8 байтов). Работа со значениями этих типов (и с массивами из них), таким образом, вполне безопасна. 4.1.4. Опрос объема использованной дополнительной системной памятиvoid coarray_report( void **coarray, long *requested, long *used ) Для секционированного массива coarray выясняется объем запрошенной при создании (значение второго аргумента coarray_create()) и израсходованной в данный момент в данном процессе дополнительной системной памяти в байтах. Объем запрошенной памяти возвращается по указателю, заданному вторым аргументом (*requested), объем израсходованной – по указателю, заданному третьим аргументом (*used). В качестве первого аргумента следует передавать строго адрес массива указателей на секции, полученный успешным обращением к coarray_create(), а не «самодельной» копии этого массива. 4.1.5. Организация последовательной работы процессов со взаимным исключениемИсполняемые макросы NODE_BY_NODE_BEGIN( k, 0, n_nodes )и NODE_BY_NODE_END являются заголовком и, соответственно, концом специальной формы цикла по k, аналогичного циклу for ( k = 0; k < n_nodes; k++ ) Соответствующий каждому значению k виток этого цикла выполняется в процессе номер k. Витки выполняются строго последовательно, следующий – только после предыдущего. Диапазон значений k не обязан охватывать все номера процессов. В первую очередь, цикл предназначен для записи секционированных массивов в файл. 4.2. Работа с разделяемыми массивамиДля работы с разделяемыми массивами программа должна включать заголовочный файл "shared.h". Как уже говорилось выше, разделяемый массив – это не специальный вид системного объекта, а режим использования в программе секционированного массива. Чтобы использовать разделяемый массив, следует создать секционированный массив, а затем использовать для работы с ним описываемые ниже функции и макросы. 4.2.1. Вычисление размера блока разделяемого массиваФункция extern long shared_blocksize( long larr, int n_nodes ) возвращает размер блока разделяемого массива длиной larr для распределения его по n_nodes процессам. Этот размер предполагается использовать впоследствии для указания (после перевода в байты) в качестве значения первого аргумента при вызове coarray_create(). Если длина массива не делится нацело на число процессов, размер блока берется с минимально необходимым запасом. Здесь и далее под длиной разделяемого массива понимается его длина в тех элементах, тип которых выбран при объявлении в программе указателей на секции. Так, в приведенном выше примере применения указатели на секции имеют тип «массив элементов double длиной my», и именно в таких строках сетки, каждая из которых состоит из my значений типа double, и измеряется длина разделяемого массива. 4.2.2. Определение принадлежности элемента разделяемого массива к своей секции (блоку).Функция 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 ) 4.2.3. Индексация разделяемого массиваМакрос shared( f, nb, i )обеспечивает доступ к i-му элементу разделяемого массива, для реализации которого используется секционированный массив f, при размере блока разделяемого массива (секции f), равном nb. Макрос устроен так, что, если в качестве его первого и второго аргумента используются скалярные переменные, константы или очень простые выражения, то время доступа к элементу разделяемого массива при помощи макроса практически не возрастает по сравнению со временем доступа в режиме простого секционированного массива. Но это верно только в том случае, когда компилятор помещает соответствующую скалярную переменную в регистр. В противном случае время обращения к элементу разделяемого массива возрастает более чем на порядок. Причины, по которым компилятор помещает или не помещает ту или иную скалярную переменную в регистр по умолчанию, не всегда очевидны программисту. Поэтому, во избежание неожиданной катастрофической потери быстродействия, необходимо те переменные, из которых получается второй аргумент этого макроса, явно объявлять в программе как "register", и добиваться исчезновения диагностических сообщений компилятора, связанных с таким объявлением (эти сообщения обычно свидетельствуют о том, что переменную фактически поместить в регистр не удалось). 4.2.4. Вычисление границ изменения индекса при обработке разделяемого массива в параллельном циклеПусть в последовательной программе, параллельную реализацию которой требуется выполнить при помощи RefNUMA, имеется цикл: for ( k = kbeg; k < kend; k++ ) с независимыми витками. Пусть программист планирует распределить витки этого цикла по процессам таким образом, чтобы k-й виток цикла выполнялся в некотором процессе тогда и только тогда, когда k-й элемент разделяемого массива длиной larr является для этого процесса своим. Тогда реализация такого цикла в параллельной программе должна иметь вид: for ( k = shared_from( larr, n_nodes, my_node, kbeg ); k < shared_upto( larr, n_nodes, my_node, kend; k++ ) где n_nodes и my_node – общее число процессов и номер данного процесса, соответственно. Несколько менее эффективный, но по смыслу эквивалентный вариант такого же распределения витков цикла по процессам: for ( k = kbeg; k < kend; k++ ) { if ( shared_is_local( larr, n_nodes, my_node, k ) ) { ... 5. Рекомендации по использованию на К-100В настоящее время готовы реализации RefNUMA на основе Qlogic shmem и на основе shmem МВС-Экспресс. Последняя работает ощутимо быстрее. Для их использования следует выбрать вариант системы параллельного программирования "intelopenmpishmem" или "intelintelmpimvse", соответственно [7]. Текст примера и скрипт для его компиляции находятся в /common/coarray/example, а настоящий документ – в /common/coarray/doc Литература |
|
||||||||||||||||||||||||||||||||
Тел. +7(499)220-79-72; E-mail: inform@kiam.ru |