| | | | |

1.8 Seleção de Aplicações

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

Em revisão

1.8.1 Sistemas de Equações

Em revisão

Consideramos o seguinte problema de equações diferenciais ordinárias com valores de contorno

u0′′+u1=f0,x(0,L) (1.190)
u1′′+u0=f1,x(0,L) (1.191)
u0(0)=u00,u0(L)=u0L, (1.192)
u1(0)=u10,u1(L)=u1L, (1.193)

onde f0, f1, u00, u0L, u10, u1L são dados.

Para construirmos uma aproximação por elementos finitos podemos tomar o seguinte problema fraco associado: encontrar u=(u0,u1)V0×V1 tal que

a(u,v)=L(v),v=(v0,v1)V×V, (1.194)

onde V0={vH1(I);v0(0)=u00,v0(L)=u0L}, V1={v1H1(I);v1(0)=u10,v1(L)=u1L}, V={vH1(I);v(0)=v(L)=0}, a forma bilinear é

a(u,v)=Iu0v0𝑑x+Iu1v1𝑑x+Iu1v0𝑑x+Iu0v1𝑑x (1.195)

e a forma linear é

L(v)=If0v0𝑑x+If1v1𝑑x. (1.196)

Então, o problema de elemento finitos associado no espaço das funções lineares por partes lê-se: encontrar uh=(uh0,uh1)Vh0×Vh1 tal que

a(uh,vh)=L(vh),vh=(vh0,vh1)Vh×Vh, (1.197)

onde Vh0={vhP1(I);vh0(0)=u00,vh0(L)=u0L}, Vh1={vh1P1(I);vh1(0)=u10,vh1(L)=u1L}, Vh={vhP1(I);vh(0)=vh(L)=0}.

Exemplo 1.8.1.

Consideramos o seguinte problema de valor de contorno

u0′′+u1=sen(x)+cos(x),x(π,π) (1.198)
u1′′+u0=cos(x)sin(x),x(π,π) (1.199)
u0(π)=0,u0(π)=0, (1.200)
u1(π)=1,u1(π)=1. (1.201)

Considerando elementos lineares por partes, temos a seguinte formulação de elementos finitos: encontrar uh=(uh0,uh1)Vh0×Vh1 tal que

a(uh,vh)=L(vh),vh=(vh0,vh1)Vh×Vh, (1.202)

onde Vh0={vhP1(I);vh0(0)=vh0(L)=0}, Vh1={vh1P1(I);vh1(0)=vh1(L)=1}, Vh={vhP1(I);vh(0)=vh(L)=0}, com as formas bilinear e linear são dadas em (1.195) e (1.196), respectivamente.

A Figura 1.8 apresenta o esboço dos gráficos das soluções analíticas u0(x)=sen(x) e u1(x)=cos(x) e de suas aproximações de elementos finitos uh0 e uh1, estas construídas no espaço dos polinômios lineares por partes sobre uma malha uniforme de 5 células.

Refer to caption
Figura 1.8: Esboço dos gráficos das soluções referentes ao Exemplo 1.8.1.

Com o FEniCS, a computação do problema de elementos finitos pode ser feita com o seguinte código:

from __future__ import print_function, division
from fenics import *
import numpy as np
import matplotlib.pyplot as plt

#tolerance
tol=1e-14

# malha
mesh = IntervalMesh(10,-pi,pi)

# espaco
P1 = FiniteElement(’P’,interval,1)
element = MixedElement([P1,P1])
V = FunctionSpace(mesh, element)

#C.C.
def boundary(x,on_boundary):
    return on_boundary

bc = [DirichletBC(V.sub(0),Constant(0.0),boundary),
      DirichletBC(V.sub(1),Constant(-1.0),boundary)]
print(bc)

#MEF problem
u = TrialFunction(V)
v = TestFunction(V)
f0 = Expression(’sin(x[0]) + cos(x[0])’,
                degree=10)
f1 = Expression(’cos(x[0]) - sin(x[0])’,
                degree=10)
a  = u[0].dx(0)*v[0].dx(0)*dx
a += u[1]*v[0]*dx
a += u[1].dx(0)*v[1].dx(0)*dx
a -= u[0]*v[1]*dx
L  = f0*v[0]*dx
L += f1*v[1]*dx

#computa a sol
u = Function(V)
solve(a == L, u, bc)

#sol analitica
u0a = Expression(’sin(x[0])’,
                 degree=10)
u1a = Expression(’cos(x[0])’,
                 degree=10)

plot(u[0],mesh=mesh,marker=’.’,label=r"$u_{h0}$")
plot(u[1],mesh=mesh,marker=’.’,label=r"$u_{h1}$")
mesh = IntervalMesh(100,-pi,pi)
plot(u0a,mesh=mesh,label=r"$u_0$")
plot(u1a,mesh=mesh,label=r"$u_1$")
plt.legend(numpoints=1)
plt.grid(’on’)
plt.show()

1.8.2 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!