| | | |

Método dos Elementos Finitos

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

2.3 Malha adaptativa baseada no resíduo

Em construção

Vamos estudar um método de malha adaptativa baseado no resíduo. Para isso, vamos precisar de uma estimativa a posteriori do erro.

2.3.1 Estimativa a posteriori

Teorema 2.3.1.(Estimativa a posteriori)

A solução uh do problema de elementos finitos (2.42) satisfaz

|uuh|2CK𝒦ηK2(uh), (2.110)

onde o elemento residual ηK(uh) é definido por

ηK(uh)=hKf+ΔuhL2(K)+12hK1/2[𝒏uh]L2(KΩ). (2.111)

Aqui, [𝒏uh]|K denota o salto na derivada normal de uh nos lados interiores dos elementos de 𝒦. Além disso, lembremos que Δuh=0.

Demonstração.

Denotando e:=uuh o erro entre a solução do problema forte e a solução de elementos finitos, temos

eE2=eL2(Ω)2 (2.112)
=Ωeedx (2.113)
=Ωe(eπe)dx. (2.114)

Nesta última equação, temos usado a ortogonalidade de Galerkin (Teorema 2.2.2). Daí, temos

Ωe(eπe)dx=K𝒦Ke(eπe)dx (2.115)
=K𝒦KΔe(eπe)dx
+K𝒏e(eπe)𝑑s, (2.116)
=K𝒦K(f+Δuh)(eπe)dx
+KΩ𝒏e(eπe)𝑑s, (2.117)

uma vez que Δe|K=f+Δuh|K e, ambos, e e πe se anulam em Ω.

Para computarmos o segundo termo do lado direito da ultima equação, observamos que o erro em lado E recebe contribuições dos dois elementos K± que compartilham E. Com isso, temos

K+K𝒏e(eπe)ds=E(𝒏+e+(e+πe+)
+𝒏e(eπe))ds, (2.118)

onde utilizamos a notação v±=v|K±. Lembremos que o erro e é contínuo e, portanto, (e+πe+)|E=(eπe)|E. Ainda, u é contínuo, logo (𝒏+u++𝒏u)|E=0. Entretanto, uh|E não é geralmente contínuo, sendo apenas constante por partes. Assim sendo e denotando o salto [𝒏uh]:=(𝒏+uh++𝒏uh), temos

E(𝒏+e+(eπe)+𝒏e(eπe))𝑑s
=E[𝒏uh](eπe)ds. (2.119)

Com isso, temos

K𝒦KΩ𝒏e(eπe)𝑑s=EIE[nuh](eπe)𝑑s, (2.120)

onde I é o conjunto dos lados interiores na triangularização 𝒦. Logo, retornando a (2.117), obtemos

eE2=K𝒦K(f+uh)(eπe)𝑑x
12Kω[nuh](eπe)𝑑s. (2.121)

Nos resta, agora, estimarmos estes dois termos do lado direito.

A estimativa do primeiro, segue da desigualdade de Cauchy-Schwarz seguida da estimativa padrão do erro de interpolação, i.e.

K(f+Δuh)(eπe)𝑑xf+ΔuhL2(Ω)eπeL2(Ω) (2.122)
f+ΔuhL2(Ω)ChKDeL2(Ω) (2.123)

Para estimarmos as contribuições dos lados, usamos a desigualdade do Traço [6, Subseção 7.1.6]

vL2(K)2C(hK1vL2(K)2+hKvL2(K)2). (2.124)

Com esta, a desigualdade de Cauchy-Schwarz e a estimativa padrão do erro de interpolação, temos

KΩ[nuh](eπe)𝑑s[nuh]L2(K)eπeL2(K) (2.125)
[nuh]L2(K)
C(hK1eπeL2(K)2+hKD(eπe)L2(K)2)1/2 (2.126)
[nuh]L2(K)ChK1/2DeL2(K). (2.127)

Daí, a estimativa segue das (2.123) e (2.127). ∎

2.3.2 Malha adaptativa

Em construção

1from mpi4py import MPI
2from petsc4py import PETSc
3from dolfinx import mesh
4from dolfinx import fem
5
6# malha
7domain = mesh.create_unit_square(MPI.COMM_WORLD,
8 nx = 4, ny = 4)
9
10
11# problema mef
12import ufl
13from dolfinx import default_scalar_type
14from dolfinx.fem.petsc import LinearProblem, assemble_vector, create_vector
15
16# loop over refinements
17for i in range(3):
18
19 V = fem.functionspace(domain, ('P', 1))
20
21 # condição de contorno
22 import numpy as np
23 uD = fem.Function(V)
24 uD.interpolate(lambda x: np.full(x.shape[1], 0.))
25
26 tdim = domain.topology.dim
27 fdim = tdim - 1
28 domain.topology.create_connectivity(fdim, tdim)
29 boundary_facets = mesh.exterior_facet_indices(domain.topology)
30 boundary_dofs = fem.locate_dofs_topological(V, fdim,
31 boundary_facets)
32 bc = fem.dirichletbc(uD, boundary_dofs)
33
34 u = ufl.TrialFunction(V)
35 v = ufl.TestFunction(V)
36
37 # termo de fonte
38 x = ufl.SpatialCoordinate(domain)
39 f = 2 * ufl.pi**2 * ufl.sin(ufl.pi * x[0]) * ufl.sin(ufl.pi * x[1])
40
41 # forma bilinear
42 a = ufl.dot(ufl.grad(u), ufl.grad(v)) * ufl.dx
43
44 # forma linear
45 L = f * v * ufl.dx
46
47 # problema linear
48 problem = LinearProblem(a, L, bcs=[bc],
49 petsc_options_prefix="ex_mef2d_modelo")
50 uh = problem.solve()
51
52 # gráfico
53 import pyvista
54 from dolfinx import plot
55
56 # solução numérica
57 topology, cell_types, geometry = plot.vtk_mesh(V)
58 grid_uh = pyvista.UnstructuredGrid(topology, cell_types, geometry)
59 grid_uh.point_data["uh"] = uh.x.array.real
60 grid_uh.set_active_scalars("uh")
61
62 # plot
63 p = pyvista.Plotter()
64 p.add_mesh(grid_uh, show_edges=True, scalars="uh", cmap="viridis", label="uh")
65 # p.add_mesh(contours, color="white", label="uex_interp")
66 # p.add_point_labels(np.array(pts), [f"{val:.2f}" for val in vals], font_size=10, text_color="white", shape_opacity=0.5, point_color="white", point_size=1)
67 p.add_axes()
68 p.add_legend()
69 p.view_xy()
70 p.show()
71
72 # compute the L2 norm of the residual at each cell
73 W0 = fem.functionspace(domain, ("DG", 0))
74 w0 = ufl.TestFunction(W0)
75
76 # Strong residual of -Delta(u) = f
77 r = f + ufl.div(ufl.grad(uh))
78 r2_form = fem.form((r**2) * w0 * ufl.dx)
79
80 r2_vec = create_vector(fem.extract_function_spaces(r2_form))
81 with r2_vec.localForm() as loc_r2:
82 loc_r2.set(0)
83 assemble_vector(r2_vec, r2_form)
84 r2_vec.ghostUpdate(addv=PETSc.InsertMode.ADD_VALUES,
85 mode=PETSc.ScatterMode.REVERSE)
86 r2_vec.ghostUpdate(addv=PETSc.InsertMode.INSERT_VALUES,
87 mode=PETSc.ScatterMode.FORWARD)
88
89 # cell_residual_l2 = fem.Function(W0)
90 # cell_residual_l2.x.array[:] = np.sqrt(np.maximum(r2_vec.array, 0.0))
91 # cell_residual_l2.x.scatter_forward()
92
93 # print norm of cell residuals
94 # print(f"Refinement {i}: max cell residual L2 norm = {cell_residual_l2.x.array.max()}")
95 print(f"Refinement {i}: max cell residual L2 norm = {r2_vec.array.max()}")
96
97 # plot cell residuals (DG0 is cellwise, so use mesh topology + cell_data)
98 topology, cell_types, geometry = plot.vtk_mesh(domain, domain.topology.dim)
99 grid_residual = pyvista.UnstructuredGrid(topology, cell_types, geometry)
100 # grid_residual.cell_data["cell_residual_l2"] = cell_residual_l2.x.array.real
101 grid_residual.cell_data["cell_residual_l2"] = r2_vec.array.real
102 grid_residual.set_active_scalars("cell_residual_l2")
103
104 pl = pyvista.Plotter()
105 pl.add_mesh(grid_residual, show_edges=True, scalars="cell_residual_l2", cmap="inferno", label="cell_residual_l2")
106 pl.add_axes()
107 pl.add_legend()
108 pl.view_xy()
109 pl.show()
110
111 # mark cells for refinement based on residuals
112 # DG0 dof numbering is not guaranteed to match cell numbering.
113 tdim = domain.topology.dim
114 num_local_cells = domain.topology.index_map(tdim).size_local
115 local_cells = np.arange(num_local_cells, dtype=np.int32)
116 cell_eta = np.empty(num_local_cells, dtype=np.float64)
117 for c in local_cells:
118 dof = W0.dofmap.cell_dofs(c)[0]
119 # cell_eta[c] = cell_residual_l2.x.array[dof].real
120 cell_eta[c] = r2_vec.array[dof].real
121
122 local_max = cell_eta.max() if num_local_cells > 0 else 0.0
123 global_max = domain.comm.allreduce(local_max, op=MPI.MAX)
124 threshold = 0.5 * global_max
125 cells_to_refine = local_cells[cell_eta > threshold]
126 print(f"Refinement {i}: marking {len(cells_to_refine)} cells for refinement")
127
128 # # global mesh refinement
129 # from dolfinx.mesh import refine
130 # refine_out = refine(domain)
131 # domain = refine_out[0] if isinstance(refine_out, tuple) else refine_out
132
133 # local mesh refinement
134 from dolfinx.mesh import refine
135 if len(cells_to_refine) > 0:
136 marker_entities = cells_to_refine
137 try:
138 # Some dolfinx versions expect edges, not cells, for local refine.
139 marker_entities = mesh.compute_incident_entities(domain.topology,
140 cells_to_refine,
141 tdim, 1)
142 except Exception:
143 pass
144 refine_out = refine(domain, marker_entities)
145 domain = refine_out[0] if isinstance(refine_out, tuple) else refine_out

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