| | | |

Método dos Elementos Finitos

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

1.4 Aplicação: EDP não-linear

Em construção

Como um exemplo de aplicação do método de elementos finitos na solução de equações diferenciais parciais não-lineares, consideramos o problema de Troesch

u′′=μsinh(μu),x(0,1), (1.186)
u(0)=0,u(1)=1, (1.187)

onde μ>0 é um parâmetro.

A formulação fraca deste problema consiste em: encontrar uV:={vH1(I);v(0)=0,v(1)=1} tal que

a(u,v)=L(v),vV0, (1.188)

onde, V0={vH1(I);v(0)=0,v(1)=0} e as formas bilinear a(,) e linear L() são dadas por

a(u,v):=01uv𝑑x, (1.189)
L(v):=01μsinh(μu)v𝑑x. (1.190)

A formulação de elementos finitos deste problema consiste em: encontrar uhVh:={vhP1(I);vh(0)=0,vh(1)=1} tal que

a(uh,vh)=L(vh),vhVh,0, (1.191)

onde Vh,0={vhP1(I);vh(0)=0,vh(1)=0}.

Código 10: ex_mef1d_nlinear.py
1from mpi4py import MPI
2
3# malha
4from dolfinx import mesh
5domain = mesh.create_unit_interval(MPI.COMM_WORLD,
6 nx = 10)
7# espaço
8from dolfinx import fem
9V = fem.functionspace(domain, ('P', 1))
10
11# condição de contorno
12import numpy as np
13uD = fem.Function(V)
14# uD(0) = 0, uD(1) = 1
15def dirichlet_bc(x):
16 y = np.full(x.shape[1], 1.0)
17 y[np.isclose(x[0], 0.)] = 0.0
18 return y
19uD.interpolate(dirichlet_bc)
20
21tdim = domain.topology.dim
22fdim = tdim - 1
23domain.topology.create_connectivity(fdim, tdim)
24boundary_facets = mesh.exterior_facet_indices(domain.topology)
25boundary_dofs = fem.locate_dofs_topological(V, fdim,
26 boundary_facets)
27bc = fem.dirichletbc(uD, boundary_dofs)
28
29# problema mef
30import ufl
31from dolfinx import default_scalar_type
32from dolfinx.fem.petsc import NonlinearProblem
33from dolfinx.nls.petsc import NewtonSolver
34
35# unknown / test functions
36uh = fem.Function(V)
37v = ufl.TestFunction(V)
38
39# parameter
40mu = fem.Constant(domain, default_scalar_type(5))
41
42F = ufl.dot(ufl.grad(uh), ufl.grad(v)) * ufl.dx + (mu * ufl.sinh(mu * uh)) * v * ufl.dx
43
44petsc_options = {
45 "snes_type": "newtonls",
46 "snes_linesearch_type": "none",
47 "snes_atol": 1e-6,
48 "snes_rtol": 1e-6,
49 "snes_monitor": None,
50 "ksp_error_if_not_converged": True,
51 "ksp_type": "gmres",
52 "ksp_rtol": 1e-8,
53 "ksp_monitor": None,
54}
55
56problem = NonlinearProblem(
57 F,
58 uh,
59 bcs=[bc],
60 petsc_options=petsc_options,
61 petsc_options_prefix="nonlinpoisson",
62)
63
64problem.solve()
65converged = problem.solver.getConvergedReason()
66num_iter = problem.solver.getIterationNumber()
67assert converged > 0, f"Solver did not converge, got {converged}."
68print(
69 f"Solver converged after {num_iter} iterations with converged reason {converged}."
70)

1.4.1 Exercícios

Em construção


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