| | | |

Computação Paralela com C++

Ajude a manter o site livre, gratuito e sem propagandas. Colabore!

2.5 Sincronização

A sincronização de threads é uma forma de controlar o acesso a recursos compartilhados e de garantir que as operações sejam executadas em uma ordem específica. A sincronização é importante para evitar condições de corrida e garantir a consistência dos dados. Entretanto, a sincronização pode levar a uma diminuição do desempenho, pois threads podem ficar bloqueadas esperando que a sincronização ocorra.

As sincronizações podem ocorrer por uso de construtores específicos (por exemplo, omp critical e omp barrier) ou de forma implícita ao final de regiões paralelas. Por exemplo, ao final de uma região paralela, o master thread espera que todos os outros threads terminem suas operações antes de continuar a execução do código.

1#pragma omp parallel
2{
3 // região paralela
4
5} // sincronização

Alguns construtores, como o do laço for, também fazem sincronização implícita.

1#pragma omp parallel
2{
3 // região paralela
4
5 #pragma omp for
6 for (int i = 0; i < n; i++)
7 {
8 // laço paralelo
9
10 } // sincronização
11
12 // região paralela
13
14} // sincronização

Nestes casos, a sincronização pode ser evitada com o uso da cláusula nowait.

1#pragma omp parallel
2{
3 // região paralela
4
5 #pragma omp for nowait
6 for (int i = 0; i < n; i++)
7 {
8 // laço paralelo
9
10 } // não há sincronização
11
12 // região paralela
13
14} // sincronização
Exemplo 2.5.1.

O seguinte método OpenMP/C++ computa a norma do residuo de um sistema linear Ax=b, i.e. Axb2.

1#include <omp.h>
2#include <cmath>
3#include <eigen3/Eigen/Dense>
4
5double res_norm(const Eigen::VectorXd &x,
6 const Eigen::MatrixXd &A,
7 const Eigen::VectorXd &b)
8{
9 int n = A.rows();
10 int m = A.cols();
11 Eigen::VectorXd r = Eigen::VectorXd::Zero(n);
12 double r_norm2 = 0.0;
13
14 #pragma omp parallel
15 {
16 // r = Ax - b
17 #pragma omp for schedule(static) nowait
18 for (int i = 0; i < n; i++)
19 {
20 for (int j = 0; j < m; j++)
21 {
22 r[i] += A(i, j) * x[j];
23 }
24 r[i] -= b[i];
25 } // não há sincronização
26
27 // res = ||r||^2
28 #pragma omp for schedule(static) reduction(+:r_norm2)
29 for (int i = 0; i < n; i++)
30 {
31 r_norm2 += r[i] * r[i];
32 } // sincronização
33
34 } // sincronização
35
36 return sqrt(r_norm2);
37
38}

Verifique!

2.5.1 Barreira

Barreira é um tipo de sincronização em que as threads devem esperar até que todas tenham alcançado um ponto específico do código. Uma barreira explícita pode ser criada com a diretiva omp barrier.

Exemplo 2.5.2.(Método de Jacobi)

O seguinte código é uma implementação OpenMP/C++ do método de Jacobi para resolver um sistema linear Ax=b. Nele, a distribuição do laço é feita de forma manual e estática. A cada iteração, as threads devem esperar até que todas tenham terminado a iteração anterior. Isso é feito com a diretiva omp barrier.

Código 7: jacobi_v0.hpp
1#include <iostream>
2#include <omp.h>
3#include <eigen3/Eigen/Dense>
4
5Eigen::VectorXd jacobi(const Eigen::MatrixXd &A,
6 const Eigen::VectorXd &b)
7{
8 // params
9 int itmax = 1000;
10 double tol = 1e-8;
11
12 int n = A.rows();
13 Eigen::VectorXd x = Eigen::VectorXd::Zero(n);
14
15 double r_norm = 0.0;
16
17 #pragma omp parallel num_threads(3)
18 {
19 // distribuição
20 int nthreads = omp_get_num_threads();
21 int tid = omp_get_thread_num();
22 int chunk = n / nthreads;
23 int start = tid * chunk;
24 int end = (tid == nthreads - 1) ? n : start + chunk;
25
26 // iteração de Jacobi
27 for (int it = 0; it < itmax; ++it)
28 {
29 for (int i = start; i < end; ++i)
30 {
31 double sum = 0.0;
32 for (int j = 0; j < n; ++j)
33 {
34 if (i != j)
35 sum += A(i, j) * x(j);
36 }
37 x(i) = (b(i) - sum) / A(i, i);
38 }
39
40 // sincronização
41 #pragma omp barrier
42
43 // verificação
44 if (tid == 0)
45 {
46 r_norm = (A * x - b).norm();
47 if (r_norm < tol)
48 {
49 std::cout << "Converged in "
50 << it << " iterations."
51 << std::endl;
52 }
53 }
54
55 // sincronização
56 #pragma omp barrier
57 if (r_norm < tol)
58 break;
59 }
60 }
61
62 return x;
63
64}

Verifique!

2.5.2 Exercícios

E. 2.5.1.

Desenvolva um método OpenMP/C++ para computar a multiplicação matriz-matriz AB. O código deve fazer a distribuição das threads nas colunas da matriz B. Use a cláusula omp nowait para evitar sincronizações desnecessárias.

E. 2.5.2.

Implemente uma versão OpenMP/C++ do método de Jacobi para a solução de differenças finitas da equação de Laplace

u=f (2.4)

em um domínio retangular e com condições de contorno de Dirichlet homogêneas.

E. 2.5.3.

Implemente uma versão OpenMP/C++ do método de Euler explícito para a solução de diferenças finitas da equação do calor

ut=2u (2.5)

em um domínio retangular e com condições de contorno de Dirichlet homogêneas e uma condição inicial arbitrária.


Envie seu comentário

Aproveito para agradecer a todas/os que de forma assídua ou esporádica contribuem enviando correções, sugestões e críticas!

Opcional. Preencha seu nome para que eu possa lhe contatar.
Opcional. Preencha seu e-mail para que eu possa lhe contatar.
As informações preenchidas são enviadas por e-mail para o desenvolvedor do site e tratadas de forma privada. Consulte a política de uso de dados para mais informações.

Licença Creative Commons
Este texto é disponibilizado nos termos da Licença Creative Commons Atribuição-CompartilhaIgual 4.0 Internacional. Ícones e elementos gráficos podem estar sujeitos a condições adicionais.

Computação Paralela com C++

Ajude a manter o site livre, gratuito e sem propagandas. Colabore!

2.5 Sincronização

A sincronização de threads é uma forma de controlar o acesso a recursos compartilhados e de garantir que as operações sejam executadas em uma ordem específica. A sincronização é importante para evitar condições de corrida e garantir a consistência dos dados. Entretanto, a sincronização pode levar a uma diminuição do desempenho, pois threads podem ficar bloqueadas esperando que a sincronização ocorra.

As sincronizações podem ocorrer por uso de construtores específicos (por exemplo, omp critical e omp barrier) ou de forma implícita ao final de regiões paralelas. Por exemplo, ao final de uma região paralela, o master thread espera que todos os outros threads terminem suas operações antes de continuar a execução do código.

1#pragma omp parallel
2{
3 // região paralela
4
5} // sincronização

Alguns construtores, como o do laço for, também fazem sincronização implícita.

1#pragma omp parallel
2{
3 // região paralela
4
5 #pragma omp for
6 for (int i = 0; i < n; i++)
7 {
8 // laço paralelo
9
10 } // sincronização
11
12 // região paralela
13
14} // sincronização

Nestes casos, a sincronização pode ser evitada com o uso da cláusula nowait.

1#pragma omp parallel
2{
3 // região paralela
4
5 #pragma omp for nowait
6 for (int i = 0; i < n; i++)
7 {
8 // laço paralelo
9
10 } // não há sincronização
11
12 // região paralela
13
14} // sincronização
Exemplo 2.5.1.

O seguinte método OpenMP/C++ computa a norma do residuo de um sistema linear Ax=b, i.e. Axb2.

1#include <omp.h>
2#include <cmath>
3#include <eigen3/Eigen/Dense>
4
5double res_norm(const Eigen::VectorXd &x,
6 const Eigen::MatrixXd &A,
7 const Eigen::VectorXd &b)
8{
9 int n = A.rows();
10 int m = A.cols();
11 Eigen::VectorXd r = Eigen::VectorXd::Zero(n);
12 double r_norm2 = 0.0;
13
14 #pragma omp parallel
15 {
16 // r = Ax - b
17 #pragma omp for schedule(static) nowait
18 for (int i = 0; i < n; i++)
19 {
20 for (int j = 0; j < m; j++)
21 {
22 r[i] += A(i, j) * x[j];
23 }
24 r[i] -= b[i];
25 } // não há sincronização
26
27 // res = ||r||^2
28 #pragma omp for schedule(static) reduction(+:r_norm2)
29 for (int i = 0; i < n; i++)
30 {
31 r_norm2 += r[i] * r[i];
32 } // sincronização
33
34 } // sincronização
35
36 return sqrt(r_norm2);
37
38}

Verifique!

2.5.1 Barreira

Barreira é um tipo de sincronização em que as threads devem esperar até que todas tenham alcançado um ponto específico do código. Uma barreira explícita pode ser criada com a diretiva omp barrier.

Exemplo 2.5.2.(Método de Jacobi)

O seguinte código é uma implementação OpenMP/C++ do método de Jacobi para resolver um sistema linear Ax=b. Nele, a distribuição do laço é feita de forma manual e estática. A cada iteração, as threads devem esperar até que todas tenham terminado a iteração anterior. Isso é feito com a diretiva omp barrier.

Código 7: jacobi_v0.hpp
1#include <iostream>
2#include <omp.h>
3#include <eigen3/Eigen/Dense>
4
5Eigen::VectorXd jacobi(const Eigen::MatrixXd &A,
6 const Eigen::VectorXd &b)
7{
8 // params
9 int itmax = 1000;
10 double tol = 1e-8;
11
12 int n = A.rows();
13 Eigen::VectorXd x = Eigen::VectorXd::Zero(n);
14
15 double r_norm = 0.0;
16
17 #pragma omp parallel num_threads(3)
18 {
19 // distribuição
20 int nthreads = omp_get_num_threads();
21 int tid = omp_get_thread_num();
22 int chunk = n / nthreads;
23 int start = tid * chunk;
24 int end = (tid == nthreads - 1) ? n : start + chunk;
25
26 // iteração de Jacobi
27 for (int it = 0; it < itmax; ++it)
28 {
29 for (int i = start; i < end; ++i)
30 {
31 double sum = 0.0;
32 for (int j = 0; j < n; ++j)
33 {
34 if (i != j)
35 sum += A(i, j) * x(j);
36 }
37 x(i) = (b(i) - sum) / A(i, i);
38 }
39
40 // sincronização
41 #pragma omp barrier
42
43 // verificação
44 if (tid == 0)
45 {
46 r_norm = (A * x - b).norm();
47 if (r_norm < tol)
48 {
49 std::cout << "Converged in "
50 << it << " iterations."
51 << std::endl;
52 }
53 }
54
55 // sincronização
56 #pragma omp barrier
57 if (r_norm < tol)
58 break;
59 }
60 }
61
62 return x;
63
64}

Verifique!

2.5.2 Exercícios

E. 2.5.1.

Desenvolva um método OpenMP/C++ para computar a multiplicação matriz-matriz AB. O código deve fazer a distribuição das threads nas colunas da matriz B. Use a cláusula omp nowait para evitar sincronizações desnecessárias.

E. 2.5.2.

Implemente uma versão OpenMP/C++ do método de Jacobi para a solução de differenças finitas da equação de Laplace

u=f (2.4)

em um domínio retangular e com condições de contorno de Dirichlet homogêneas.

E. 2.5.3.

Implemente uma versão OpenMP/C++ do método de Euler explícito para a solução de diferenças finitas da equação do calor

ut=2u (2.5)

em um domínio retangular e com condições de contorno de Dirichlet homogêneas e uma condição inicial arbitrária.


Envie seu comentário

Aproveito para agradecer a todas/os que de forma assídua ou esporádica contribuem enviando correções, sugestões e críticas!

Opcional. Preencha seu nome para que eu possa lhe contatar.
Opcional. Preencha seu e-mail para que eu possa lhe contatar.
As informações preenchidas são enviadas por e-mail para o desenvolvedor do site e tratadas de forma privada. Consulte a política de uso de dados para mais informações.

Licença Creative Commons
Este texto é disponibilizado nos termos da Licença Creative Commons Atribuição-CompartilhaIgual 4.0 Internacional. Ícones e elementos gráficos podem estar sujeitos a condições adicionais.

Pedro H A Konzen
| | | |