| | | | |

2.6 Aplicação: Equação do Calor

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

Em construção

Vamos considerar a equação do Calor com dadas condições inicial e de contorno

ut=αuxx+f,(t,x)(0,tf)×(a,b), (2.4a)
u(0,x)=u0(x),x[a,b], (2.4b)
u(t,a)=u(t,b)=0, (2.4c)

com dado α>0 e fonte f.

2.6.1 Formulação Explícita

Assumimos uma discretização no tempo com nt passos t(k)=kht, de tamanho ht=tf/nt, k=0,1,,nt. Denotando u(k)(x)u(t(k),x), temos o esquema de Euler explícito

u(0)=u0(x), (2.5a)
u(k+1)=u(k)+ht[αuxx(k)+f(k)(x)], (2.5b)
u(k+1)(a)=u(k+1)(t,b)=0.

Agora, assumimos uma malha espacial com nx+1 nodos xi=a+ihx, de tamanho hx=(ba)/nx, i=0,1,,nx, e denotando ui(k)u(t(k),xi) e usando a fórmula de diferenças finitas central, obtemos o esquema iterativo

ui(0)=u0(xi),i, (2.6a)
ui(k+1)=ui(k)+rui1(k)2rui(k)+rui+1(k)+htfi(k), (2.6b)
u0(k+1)=unx(k+1)=0, (2.6c)

onde r:=αht/hx2.

Observação 2.6.1.

O esquema de Euler explícito (2.6.1) requer a seguinte condição de estabilidade

hthx212. (2.7)
Exemplo 2.6.1.

Consideramos o seguinte problema de calor

ut=uxx+(π21)etsen(πx),(t,x)(0,1)2, (2.8a)
u(0,x)=sen(πx),x[0,1], (2.8b)
u(t,0)=u(t,1)=0. (2.8c)
Código 3: ex_mlp_calor_explicito.cc
1    #include <stdio.h>
2    #include <stdlib.h>
3
4    #include <omp.h>
5
6    #include <math.h>
7    #include <gsl/gsl_vector.h>
8
9    // fonte
10    double f(double t, double x) {
11      return (M_PI*M_PI-1.) *
12             exp(-t) * sin(M_PI * x);
13    }
14
15    int main()
16    {
17      // parâmetros no tempo
18      double tf = 1.;
19      int nt = 1000;
20      double ht = tf / nt;
21
22      // parâmetros no espaço
23      int nx = 10;
24      double hx = 1./nx;
25
26      // malha
27      gsl_vector *x = gsl_vector_alloc(nx+1);
28      #pragma omp parallel for
29      for (int i = 0; i <= nx; i++) {
30        gsl_vector_set(x, i, i*hx);
31      }
32
33      gsl_vector *u0 = gsl_vector_alloc(nx+1);
34      gsl_vector *u = gsl_vector_calloc(nx+1);
35
36      // inicialização
37      double t = 0.;
38      #pragma omp parallel for
39      for (int i = 0; i <= nx; i++) {
40        gsl_vector_set(u0, i,
41                        sin(M_PI *
42                        gsl_vector_get(x, i)));
43      }
44
45      // armazena malha
46      FILE *fx = fopen("results/x.dat", "wb");
47      gsl_vector_fwrite(fx, x);
48      fclose(fx);
49
50      // armazena solução
51      char str[21];
52      sprintf(str, "results/u_%06d.dat", 0);
53      printf("t = %g, [%s]\n", t, str);
54
55      FILE *fu = fopen(str, "wb");
56      gsl_vector_fwrite(fu, u0);
57      fclose(fu);
58
59      #pragma omp parallel
60      {
61        // iteração no tempo
62        for (int k=0; k < nt; k++) {
63          #pragma omp single
64          t += ht;
65
66          #pragma omp for
67          for (int i=1; i < nx; i++) {
68            gsl_vector_set(u, i,
69                            gsl_vector_get(u0, i) +
70                            ht/(hx*hx) * (gsl_vector_get(u0, i+1) - 2*gsl_vector_get(u0, i) + gsl_vector_get(u0, i-1)) + ht *
71                            f(t-ht, gsl_vector_get(x, i)));
72          }
73
74          #pragma omp single
75          {
76            sprintf(str, "results/u_%06d.dat", k+1);
77            printf("t = %g, [%s]\n", t, str);
78            fu = fopen(str, "wb");
79            gsl_vector_fwrite(fu, u0);
80            fclose(fu);
81
82            // u0 = u
83            gsl_vector_memcpy(u0, u);
84          }
85
86        }
87      }
88
89      gsl_vector_free(x);
90      gsl_vector_free(u);
91      gsl_vector_free(u0);
92      return 0;
93    }

2.6.2 Formulação Implícita

A condição de estabilidade (2.7) pode ser contornada aplicando-se o esquema de Euler implícito

ui(0)=u0(xi),i, (2.9a)
rui1(k+1)+(2r+1)ui(k+1)rui+1(k+1)=ui(k)+htfi(k+1), (2.9b)
u0(k+1)=unx(k+1)=0, (2.9c)

onde r:=αht/hx2.

Exemplo 2.6.2.

Implementamos o esquema de Euler implícito (2.6.2) para o problema de calor (2.6.1) do exemplo anterior. No código abaixo, o sistema linear em cada iteração é computado com o método de Jacobi.

Código 4: ex_mlp_calor_implicito.cc
1#include <stdio.h>
2#include <stdlib.h>
3
4#include <omp.h>
5
6#include <math.h>
7#include <gsl/gsl_vector.h>
8#include <gsl/gsl_blas.h>
9
10// fonte
11double f(double t, double x) {
12  return (M_PI*M_PI-1.) *
13         exp(-t) * sin(M_PI * x);
14}
15
16int main()
17{
18  // parâmetros no tempo
19  double tf = 1.;
20  int nt = 100;
21  double ht = tf / nt;
22
23  // parâmetros no espaço
24  int nx = 10;
25  double hx = 1./nx;
26
27  // parâmetros Jacobi
28  int maxit = 1000;
29  double tol = 1.49e-8;
30
31  // auxiliar
32  double r = ht/(hx*hx);
33
34  // malha
35  gsl_vector *x = gsl_vector_alloc(nx+1);
36  #pragma omp parallel for
37  for (int i = 0; i <= nx; i++) {
38    gsl_vector_set(x, i, i*hx);
39  }
40
41  gsl_vector *u0 = gsl_vector_alloc(nx+1);
42  gsl_vector *u = gsl_vector_calloc(nx+1);
43
44  // auxiliares
45  double norm;
46  gsl_vector *u00 = gsl_vector_alloc(nx+1);
47  gsl_vector *diff = gsl_vector_calloc(nx+1);
48
49  // inicialização
50  double t = 0.;
51  #pragma omp parallel for
52  for (int i = 0; i <= nx; i++) {
53    gsl_vector_set(u0, i,
54                    sin(M_PI *
55                    gsl_vector_get(x, i)));
56  }
57
58  // u00 = u0
59  gsl_vector_memcpy(u00, u0);
60
61  // armazena malha
62  FILE *fx = fopen("results/x.dat", "wb");
63  gsl_vector_fwrite(fx, x);
64  fclose(fx);
65
66  // armazena solução
67  char str[21];
68  sprintf(str, "results/u_%06d.dat", 0);
69  printf("t = %g, [%s]\n", t, str);
70
71  FILE *fu = fopen(str, "wb");
72  gsl_vector_fwrite(fu, u0);
73  fclose(fu);
74
75  #pragma omp parallel
76  {
77    // iteração no tempo
78    for (int k=0; k < nt; k++) {
79
80      #pragma omp single
81      {
82        t += ht;
83        printf("%d: t = %g\n", k+1, t);
84      }
85
86      // iteração de Jacobi
87      for (int it=0; it < maxit; it++) {
88
89        #pragma omp for
90        for (int i=1; i < nx; i++) {
91          gsl_vector_set(u, i, 1./(2.*r + 1.) *
92                             (r * (gsl_vector_get(u00, i-1) +
93                                   gsl_vector_get(u00, i+1)) +
94                             gsl_vector_get(u0, i) +
95                             + ht * f(t, gsl_vector_get(x, i))));
96
97          gsl_vector_set(diff, i, gsl_vector_get(u, i) -
98                                  gsl_vector_get(u00, i));
99        }
100
101        #pragma omp single
102        {
103          norm = gsl_blas_dnrm2(diff);
104          printf("\t%d: norm=%.2e\n", it, norm);
105
106          gsl_vector_memcpy(u00, u);
107        }
108
109        if (norm < tol) {
110            break;
111        }
112
113      }
114
115      #pragma omp single
116      {
117        sprintf(str, "results/u_%06d.dat", k+1);
118        printf("t = %g, [%s]\n", t, str);
119        fu = fopen(str, "wb");
120        gsl_vector_fwrite(fu, u0);
121        fclose(fu);
122
123        // u0 = u
124        gsl_vector_memcpy(u0, u);
125      }
126
127    }
128
129  }
130
131  gsl_vector_free(x);
132  gsl_vector_free(u);
133  gsl_vector_free(u0);
134  gsl_vector_free(u00);
135  gsl_vector_free(diff);
136
137  return 0;
138}

2.6.3 Exercícios

Em construção


Envie seu comentário

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. Aproveito para agradecer a todas/os que de forma assídua ou esporádica contribuem enviando correções, sugestões e críticas!