Programming Taskbook


E-mail:

Пароль:

Регистрация пользователя   Восстановление пароля

 

ЮФУ SMBU

Электронный задачник по программированию

©  М. Э. Абрамян (Южный федеральный университет, Университет МГУ-ППИ в Шэньчжэне), 1998–2024

 

PT for MPI-2 | Решения | Отладка параллельных программ

PrevNext


Отладка параллельных программ

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

Отмеченные особенности задачника могут оказаться полезными не только при решении учебных задач по параллельному MPI-программированию, но и при разработке параллельных программ общего назначения. Поэтому, наряду с «обычными» группами заданий, позволяющими ознакомиться с большинством возможностей библиотеки MPI, в задачник Programming Taskbook for MPI-2 включена особая группа MPIDebug, предназначенная не для выполнения определенных заданий, а для отладки любых параллельных программ. Данная группа состоит из 36 «заданий», причем задание с номером N обеспечивает запуск параллельной программы с N процессами. Таким образом, если создать программу-заготовку для одного из заданий группы MPIDebug, то при запуске этой программы из интегрированной среды она автоматически будет выполняться в параллельном режиме с требуемым числом процессов.

Воспользуемся группой MPIDebug для разработки программы, реализующей один из параллельных алгоритмов умножения матрицы A на вектор b — так называемый самопланирующий алгоритм, в котором главный процесс координирует работу подчиненных процессов, пересылая им требуемые наборы данных и получая от них результаты. Исходные данные вводятся в главном процессе. Затем главный процесс пересылает во все подчиненные процессы вектор b и по одной из строк матрицы A. После этого запускается цикл, в котором главный процесс принимает от подчиненных процессов результаты умножения строки матрицы на вектор и пересылает им оставшиеся строки матрицы A. Цикл завершится, когда главный процесс перешлет подчиненным процессам все строки матрицы A и получит от них все элементы результирующего вектора Ab.

Следует заметить, что при реализации данного алгоритма нельзя заранее определить, какие именно строки матрицы будет обрабатывать тот или иной подчиненный процесс. Если, к примеру, программа использует 5 процессов, то строка матрицы номер 5 будет переслана тому процессу (из числа подчиненных процессов с рангами от 1 до 4), который первым завершит умножение ранее посланной ему строки и вернет результат главному процессу.

Для тестирования алгоритма будем обрабатывать матрицу A порядка N = 20, каждая строка которой содержит одинаковые элементы, равные порядковому номеру строки: aij = i для i, j = 1, …, N. Вектор b того же размера N будет содержать одинаковые элементы, равные 0.01. Таким образом, элементы вектора c = Ab будут иметь вид ci = 0.2i, i = 1, …, N:

с = (0.2, 0.4, 0.6, …, 3.8, 4.0).

Для правильной работы алгоритма необходимо, чтобы количество подчиненных процессов не превосходило порядок матрицы N, так как на начальном этапе алгоритма каждому подчиненному процессу посылается одна из строк данной матрицы. Выберем для определенности число процессов нашей программы равным 10; для того чтобы обеспечить автоматический запуск параллельной программы с этим числом процессов, будем использовать задание MPIDebug10.

При запуске программы-заготовки, созданной для этого задания, окно задачника будет иметь следующиий вид:

Вместо исходных и результирующих данных в этом окне выводится описание тех средств задачника, которые предназначены для вывода различной отладочной информации.

Дополним созданную заготовку фрагментом, в котором описываются необходимые переменные, выполняется формирование в главном процессе исходной матрицы A и вектора b, а также осуществляется пересылка вектора b во все подчиненные процессы:

void Solve()
{
  Task("MPIDebug10");
  int flag;
  MPI_Initialized(&flag);
  if (flag == 0)
    return;
  int rank, size;
  MPI_Comm_size(MPI_COMM_WORLD, &size);
  MPI_Comm_rank(MPI_COMM_WORLD, &rank);
  const int n = 20;
  double a[n][n], b[n], c[n], d;
  int cur_row;
  MPI_Status s;
  // генерация исходных данных в главном процессе
  if (rank == 0)
    for (int i = 0; i < n; ++i)
    {
      b[i] = 0.01;
      for (int j = 0; j < n; ++j)
        a[i][j] = i + 1;
    }
  // пересылка всем процессам вектора b
  MPI_Bcast(b, n, MPI_DOUBLE, 0, MPI_COMM_WORLD);
  ShowLine("Vector b: ..., ", b[n-1]);

При запуске данного варианта программы в окне задачника появится раздел отладки:

Поскольку в новом варианте программы используются вызовы функции Show (или ShowLine), в окне задачника автоматически отображается раздел отладки с информацией, которая была выведена с помощью этих функций. В нашем случае каждый из процессов (включая главный) вывел информацию о том, что он успешно получил вектор b (для краткости выведено только значение последнего элемента массива b с индексом n – 1; чтобы подчеркнуть, что выведен последний элемент, перед ним указано многоточие).

Так как в основных разделах окна задачника в нашем случае содержится лишь справочная информация, удобно скрыть все разделы окна, кроме раздела отладки. Если мы уже находимся в окне задачника, то для этого достаточно нажать комбинацию Ctrl+пробел. Однако еще удобнее автоматически скрывать ненужные разделы при выводе окна задачника на экран. Для этого достаточно вызвать в любом месте программы вспомогательную функцию HideTask. Если добавить вызов этой функции в конец нашей программы и повторно запустить ее на выполнение, то окно задачника будет содержать только раздел отладки:

Вернемся к реализации алгоритма умножения матрицы на вектор и приведем его завершающую часть:

if (rank == 0)
{
  // пересылка подчиненным процессам начальных строк матрицы A
  for (int i = 1; i < size; ++i)
  {
    MPI_Send(a[i-1], n, MPI_DOUBLE, i, i - 1, MPI_COMM_WORLD);
    Show("Sending initial data: dest=", i);
    Show("tag=", i - 1);
    ShowLine("-> ..., ", a[i - 1][0]);
  }
  cur_row = size - 2;
  // получение элементов произведения и пересылка остальных строк
  for (int i = 0; i < n; ++i)
  {
    MPI_Recv(&d, 1, MPI_DOUBLE, MPI_ANY_SOURCE, MPI_ANY_TAG,
      MPI_COMM_WORLD, &s);
    Show("Receiving result: source=", s.MPI_SOURCE);
    Show("tag=", s.MPI_TAG, 2);
    Show("<- ", d);
    c[s.MPI_TAG] = d;
    ++cur_row;
    if (cur_row < n)
      MPI_Send(a[cur_row], n, MPI_DOUBLE, s.MPI_SOURCE,
        cur_row, MPI_COMM_WORLD);
    else
      MPI_Send(&d, 1, MPI_DOUBLE, s.MPI_SOURCE, cur_row,
        MPI_COMM_WORLD);
    Show("| Sending data: dest=", s.MPI_SOURCE);
    ShowLine("tag=", cur_row, 2);
  }
  // вывод полученного вектора Ab
  SetPrecision(1);
  ShowLine("Resulting vector Ab:");
  for (int i = 0; i < n; ++i)
    Show(c[i], 3);
}
else
  // обработка данных в подчиненных процессах
  while (true)
  {
    MPI_Recv(c, n, MPI_DOUBLE, 0, MPI_ANY_TAG, MPI_COMM_WORLD,
      &s);
    Show("Receiving data: tag=", s.MPI_TAG, 2);
    if (s.MPI_TAG >= n)
      break;
    d = 0;
    for (int i = 0; i < n; ++i)
      d += c[i] * b[i];
    MPI_Send(&d, 1, MPI_DOUBLE, 0, s.MPI_TAG, MPI_COMM_WORLD);
    Show("| Sending result: tag=", s.MPI_TAG, 2);
    ShowLine("-> ", d);
  }

Данная часть состоит из двух больших фрагментов, первый из которых должен выполняться в главном процессе, а второй — в каждом из подчиненных процессов.

Главный процесс вначале пересылает каждому подчиненному процессу по одной из строк исходной матрицы и инициализирует переменную cur_row, в которой содержится индекс последней обработанной, т. е. пересланной в подчиненный процесс, строки матрицы (строки индексируются от нуля). После этого запускается основной цикл, в котором главный процесс получает от подчиненных процессов результирующие элементы произведения и пересылает им очередные строки матрицы. Поскольку порядок получения данных от подчиненных процессов не определен, в функции MPI_Recv используется параметр MPI_ANY_SOURCE, обеспечивающий получение данных от любого пославшего их процесса. В метке сообщения (tag) подчиненный процесс передает главному процессу дополнительную информацию — индекс вычисленного элемента произведения. Для того чтобы главный процесс мог принимать сообщения с любыми метками, параметр msgtag функции MPI_Recv полагается равным MPI_ANY_TAG. Информацию о ранге процесса, пославшего сообщение, и о метке сообщения главный процесс получает с помощью записи s типа MPI_Status, обращаясь к ее полям MPI_SOURCE и MPI_TAG соответственно. После получения сообщения от подчиненного процесса главный процесс отправляет ему новое сообщение, которое либо содержит очередную необработанную строку матрицы (если такие строки еще остались), либо специальный признак завершения. Если посылается очередная строка, то метка посылаемого сообщения полагается равной индексу строки cur_row; если же посылается признак завершения, то метка содержит число, которое больше максимального индекса строки матрицы. После завершения цикла for в главном процессе выводятся все элементы найденного произведения.

В подчиненном процессе запускается цикл, в котором вначале выполняется прием очередной строки матрицы из главного процесса (ее индекс передается в метке сообщения), затем полученная строка умножается на ранее полученный вектор, и наконец, результат умножения пересылается главному процессу. Если же при приеме очередного сообщения от главного процесса получена метка, превышающая номер максимального индекса строки матрицы, то происходит немедленный выход из цикла, и подчиненный процесс завершает свою работу.

Все этапы описанного выше алгоритма сопровождаются отладочной печатью. Для этого используется функция Show и ее вариант ShowLine (обеспечивающий дополнительный переход на новую строку в разделе отладки после вывода данных). Перед выводом очередного числового элемента, как правило, выводится предваряющий его строковый комментарий, а в некоторых случаях после числового элемента дополнительно указывается количество экранных позиций, которые следует отвести для его вывода (это позволяет обеспечить выравнивание выводимых данных по вертикали). Перед выводом в главном процессе элементов результирующего произведения вызывается функция SetPrecision(1), обеспечивающая вывод всех последующих данных вещественного типа с одним знаком после десятичного разделителя (заметим, что по умолчанию число дробных знаков полагается равным 2, однако при использовании этого значения все результирующие элементы не удалось бы вывести в одной строке, так как ширина области выводимых данных в разделе отладки фиксирована и равна 78 позициям).

При оформлении вывода было учтено, что данные, выводимые функцией Show или ShowLine, отделяются одним пробелом от ранее выведенных данных, находящихся на той же экранной строке.

Приведем отладочные данные, которые были получены при запуске программы (напомним, что перед первым символом «|» указывается ранг процесса, который вывел соответствующую строку отладочной информации, а перед символом «>» — порядковый номер выведенной строки; при этом строки, связанные с различными процессами, нумеруются независимо).

 0|  1>  Vector b: ..., 0.01
 0|  2>  Sending initial data: dest=1 tag=0 -> ..., 1.00
 0|  3>  Sending initial data: dest=2 tag=1 -> ..., 2.00
 0|  4>  Sending initial data: dest=3 tag=2 -> ..., 3.00
 0|  5>  Sending initial data: dest=4 tag=3 -> ..., 4.00
 0|  6>  Sending initial data: dest=5 tag=4 -> ..., 5.00
 0|  7>  Sending initial data: dest=6 tag=5 -> ..., 6.00
 0|  8>  Sending initial data: dest=7 tag=6 -> ..., 7.00
 0|  9>  Sending initial data: dest=8 tag=7 -> ..., 8.00
 0| 10>  Sending initial data: dest=9 tag=8 -> ..., 9.00
 0| 11>  Receiving result: source=1 tag= 0 <- 0.20 | Sending data: dest=1 tag= 9
 0| 12>  Receiving result: source=2 tag= 1 <- 0.40 | Sending data: dest=2 tag=10
 0| 13>  Receiving result: source=3 tag= 2 <- 0.60 | Sending data: dest=3 tag=11
 0| 14>  Receiving result: source=1 tag= 9 <- 2.00 | Sending data: dest=1 tag=12
 0| 15>  Receiving result: source=2 tag=10 <- 2.20 | Sending data: dest=2 tag=13
 0| 16>  Receiving result: source=3 tag=11 <- 2.40 | Sending data: dest=3 tag=14
 0| 17>  Receiving result: source=2 tag=13 <- 2.80 | Sending data: dest=2 tag=15
 0| 18>  Receiving result: source=1 tag=12 <- 2.60 | Sending data: dest=1 tag=16
 0| 19>  Receiving result: source=3 tag=14 <- 3.00 | Sending data: dest=3 tag=17
 0| 20>  Receiving result: source=1 tag=16 <- 3.40 | Sending data: dest=1 tag=18
 0| 21>  Receiving result: source=2 tag=15 <- 3.20 | Sending data: dest=2 tag=19
 0| 22>  Receiving result: source=3 tag=17 <- 3.60 | Sending data: dest=3 tag=20
 0| 23>  Receiving result: source=1 tag=18 <- 3.80 | Sending data: dest=1 tag=21
 0| 24>  Receiving result: source=2 tag=19 <- 4.00 | Sending data: dest=2 tag=22
 0| 25>  Receiving result: source=6 tag= 5 <- 1.20 | Sending data: dest=6 tag=23
 0| 26>  Receiving result: source=7 tag= 6 <- 1.40 | Sending data: dest=7 tag=24
 0| 27>  Receiving result: source=8 tag= 7 <- 1.60 | Sending data: dest=8 tag=25
 0| 28>  Receiving result: source=9 tag= 8 <- 1.80 | Sending data: dest=9 tag=26
 0| 29>  Receiving result: source=5 tag= 4 <- 1.00 | Sending data: dest=5 tag=27
 0| 30>  Receiving result: source=4 tag= 3 <- 0.80 | Sending data: dest=4 tag=28
 0| 31>  Resulting vector Ab:
 0| 32>  0.2 0.4 0.6 0.8 1.0 1.2 1.4 1.6 1.8 2.0 2.2 2.4 2.6 2.8 3.0 3.2 3.4 3.6 3.8 4.0
 1|  1>  Vector b: ..., 0.01
 1|  2>  Receiving data: tag= 0 | Sending result: tag= 0 -> 0.20
 1|  3>  Receiving data: tag= 9 | Sending result: tag= 9 -> 2.00
 1|  4>  Receiving data: tag=12 | Sending result: tag=12 -> 2.60
 1|  5>  Receiving data: tag=16 | Sending result: tag=16 -> 3.40
 1|  6>  Receiving data: tag=18 | Sending result: tag=18 -> 3.80
 1|  7>  Receiving data: tag=21
 2|  1>  Vector b: ..., 0.01
 2|  2>  Receiving data: tag= 1 | Sending result: tag= 1 -> 0.40
 2|  3>  Receiving data: tag=10 | Sending result: tag=10 -> 2.20
 2|  4>  Receiving data: tag=13 | Sending result: tag=13 -> 2.80
 2|  5>  Receiving data: tag=15 | Sending result: tag=15 -> 3.20
 2|  6>  Receiving data: tag=19 | Sending result: tag=19 -> 4.00
 2|  7>  Receiving data: tag=22
 3|  1>  Vector b: ..., 0.01
 3|  2>  Receiving data: tag= 2 | Sending result: tag= 2 -> 0.60
 3|  3>  Receiving data: tag=11 | Sending result: tag=11 -> 2.40
 3|  4>  Receiving data: tag=14 | Sending result: tag=14 -> 3.00
 3|  5>  Receiving data: tag=17 | Sending result: tag=17 -> 3.60
 3|  6>  Receiving data: tag=20
 4|  1>  Vector b: ..., 0.01
 4|  2>  Receiving data: tag= 3 | Sending result: tag= 3 -> 0.80
 4|  3>  Receiving data: tag=28
 5|  1>  Vector b: ..., 0.01
 5|  2>  Receiving data: tag= 4 | Sending result: tag= 4 -> 1.00
 5|  3>  Receiving data: tag=27
 6|  1>  Vector b: ..., 0.01
 6|  2>  Receiving data: tag= 5 | Sending result: tag= 5 -> 1.20
 6|  3>  Receiving data: tag=23
 7|  1>  Vector b: ..., 0.01
 7|  2>  Receiving data: tag= 6 | Sending result: tag= 6 -> 1.40
 7|  3>  Receiving data: tag=24
 8|  1>  Vector b: ..., 0.01
 8|  2>  Receiving data: tag= 7 | Sending result: tag= 7 -> 1.60
 8|  3>  Receiving data: tag=25
 9|  1>  Vector b: ..., 0.01
 9|  2>  Receiving data: tag= 8 | Sending result: tag= 8 -> 1.80
 9|  3>  Receiving data: tag=26

Из приведенного текста видно, что в ходе выполнения алгоритма процессы рангов 1 и 2 обработали по 5 строк исходной матрицы, процесс ранга 3 обработал 4 строки, а остальные 6 подчиненных процессов — по 1 строке. В последней, 32-й строке, связанной с главным процессом, содержатся элементы найденного произведения.

Заметим, что текст, который отображается в разделе отладки окна задачника, можно скопировать в буфер обмена Windows; для этого достаточно либо нажать стандартную комбинацию клавиш Ctrl+C, либо вызвать контекстное меню раздела отладки (нажав правую кнопку мыши) и выполнить его команду «Скопировать данные из раздела отладки в буфер».

Изменив имя задания в функции Task на «MPIDebug5», мы сможем протестировать разработанный алгоритм в параллельном режиме, использующем 5 процессов. Приведем отладочную информацию, которая была выведена при таком варианте запуска программы:

 0|  1>  Vector b: ..., 0.01
 0|  2>  Sending initial data: dest=1 tag=0 -> ..., 1.00
 0|  3>  Sending initial data: dest=2 tag=1 -> ..., 2.00
 0|  4>  Sending initial data: dest=3 tag=2 -> ..., 3.00
 0|  5>  Sending initial data: dest=4 tag=3 -> ..., 4.00
 0|  6>  Receiving result: source=4 tag= 3 <- 0.80 | Sending data: dest=4 tag= 4
 0|  7>  Receiving result: source=2 tag= 1 <- 0.40 | Sending data: dest=2 tag= 5
 0|  8>  Receiving result: source=3 tag= 2 <- 0.60 | Sending data: dest=3 tag= 6
 0|  9>  Receiving result: source=4 tag= 4 <- 1.00 | Sending data: dest=4 tag= 7
 0| 10>  Receiving result: source=1 tag= 0 <- 0.20 | Sending data: dest=1 tag= 8
 0| 11>  Receiving result: source=2 tag= 5 <- 1.20 | Sending data: dest=2 tag= 9
 0| 12>  Receiving result: source=3 tag= 6 <- 1.40 | Sending data: dest=3 tag=10
 0| 13>  Receiving result: source=4 tag= 7 <- 1.60 | Sending data: dest=4 tag=11
 0| 14>  Receiving result: source=1 tag= 8 <- 1.80 | Sending data: dest=1 tag=12
 0| 15>  Receiving result: source=2 tag= 9 <- 2.00 | Sending data: dest=2 tag=13
 0| 16>  Receiving result: source=3 tag=10 <- 2.20 | Sending data: dest=3 tag=14
 0| 17>  Receiving result: source=4 tag=11 <- 2.40 | Sending data: dest=4 tag=15
 0| 18>  Receiving result: source=2 tag=13 <- 2.80 | Sending data: dest=2 tag=16
 0| 19>  Receiving result: source=4 tag=15 <- 3.20 | Sending data: dest=4 tag=17
 0| 20>  Receiving result: source=2 tag=16 <- 3.40 | Sending data: dest=2 tag=18
 0| 21>  Receiving result: source=3 tag=14 <- 3.00 | Sending data: dest=3 tag=19
 0| 22>  Receiving result: source=1 tag=12 <- 2.60 | Sending data: dest=1 tag=20
 0| 23>  Receiving result: source=4 tag=17 <- 3.60 | Sending data: dest=4 tag=21
 0| 24>  Receiving result: source=2 tag=18 <- 3.80 | Sending data: dest=2 tag=22
 0| 25>  Receiving result: source=3 tag=19 <- 4.00 | Sending data: dest=3 tag=23
 0| 26>  Resulting vector Ab:
 0| 27>  0.2 0.4 0.6 0.8 1.0 1.2 1.4 1.6 1.8 2.0 2.2 2.4 2.6 2.8 3.0 3.2 3.4 3.6 3.8 4.0
 1|  1>  Vector b: ..., 0.01
 1|  2>  Receiving data: tag= 0 | Sending result: tag= 0 -> 0.20
 1|  3>  Receiving data: tag= 8 | Sending result: tag= 8 -> 1.80
 1|  4>  Receiving data: tag=12 | Sending result: tag=12 -> 2.60
 1|  5>  Receiving data: tag=20
 2|  1>  Vector b: ..., 0.01
 2|  2>  Receiving data: tag= 1 | Sending result: tag= 1 -> 0.40
 2|  3>  Receiving data: tag= 5 | Sending result: tag= 5 -> 1.20
 2|  4>  Receiving data: tag= 9 | Sending result: tag= 9 -> 2.00
 2|  5>  Receiving data: tag=13 | Sending result: tag=13 -> 2.80
 2|  6>  Receiving data: tag=16 | Sending result: tag=16 -> 3.40
 2|  7>  Receiving data: tag=18 | Sending result: tag=18 -> 3.80
 2|  8>  Receiving data: tag=22
 3|  1>  Vector b: ..., 0.01
 3|  2>  Receiving data: tag= 2 | Sending result: tag= 2 -> 0.60
 3|  3>  Receiving data: tag= 6 | Sending result: tag= 6 -> 1.40
 3|  4>  Receiving data: tag=10 | Sending result: tag=10 -> 2.20
 3|  5>  Receiving data: tag=14 | Sending result: tag=14 -> 3.00
 3|  6>  Receiving data: tag=19 | Sending result: tag=19 -> 4.00
 3|  7>  Receiving data: tag=23
 4|  1>  Vector b: ..., 0.01
 4|  2>  Receiving data: tag= 3 | Sending result: tag= 3 -> 0.80
 4|  3>  Receiving data: tag= 4 | Sending result: tag= 4 -> 1.00
 4|  4>  Receiving data: tag= 7 | Sending result: tag= 7 -> 1.60
 4|  5>  Receiving data: tag=11 | Sending result: tag=11 -> 2.40
 4|  6>  Receiving data: tag=15 | Sending result: tag=15 -> 3.20
 4|  7>  Receiving data: tag=17 | Sending result: tag=17 -> 3.60
 4|  8>  Receiving data: tag=21

PrevNext

 

Рейтинг@Mail.ru

Разработка сайта:
М. Э. Абрамян, В. Н. Брагилевский

Последнее обновление:
01.01.2024