-
Notifications
You must be signed in to change notification settings - Fork 2
Expand file tree
/
Copy pathprojection.py
More file actions
68 lines (57 loc) · 1.9 KB
/
projection.py
File metadata and controls
68 lines (57 loc) · 1.9 KB
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
# This demo shows how the trace of a function can be projected onto a
# function space defined over the boundary of a mesh
from dolfinx import mesh, fem, io
from mpi4py import MPI
import numpy as np
import ufl
from dolfinx.fem.petsc import LinearProblem
# Create a mesh
comm = MPI.COMM_WORLD
n = 8
msh = mesh.create_unit_square(comm, n, n)
# Create a function space for the mesh function and interpolate
V = fem.functionspace(msh, ("Lagrange", 1))
u = fem.Function(V)
u.interpolate(lambda x: np.sin(2 * np.pi * x[0]))
# Create a sub-mesh of the boundary
tdim = msh.topology.dim
fdim = tdim - 1
facets = mesh.locate_entities_boundary(
msh,
fdim,
lambda x: np.isclose(x[0], 0.0)
| np.isclose(x[0], 1.0)
| np.isclose(x[1], 0.0)
| np.isclose(x[1], 1.0),
)
submsh, sm_emap = mesh.create_submesh(msh, fdim, facets)[:2]
# We take msh to be the integration domain (we will pass this mess as the domain
# when creating the measure). We need to provide entity maps relating entities in
# `msh` to each other mesh in the form (here just `submsh`)
entity_maps = [sm_emap]
# Create function space on the boundary
Vbar = fem.functionspace(submsh, ("Lagrange", 1))
ubar, vbar = ufl.TrialFunction(Vbar), ufl.TestFunction(Vbar)
# Define forms for the projection
ds = ufl.Measure("ds", domain=msh)
a = ufl.inner(ubar, vbar) * ds
L = ufl.inner(u, vbar) * ds
petsc_opts = {"ksp_type": "preonly", "pc_type": "lu", "pc_factor_mat_solver_type": "mumps"}
problem = LinearProblem(
a,
L,
bcs=[],
petsc_options_prefix="projection_",
petsc_options=petsc_opts,
entity_maps=entity_maps,
)
ubar = problem.solve()
# Compute error and check it's zero to machine precision
e = u - ubar
e_L2 = np.sqrt(
msh.comm.allreduce(fem.assemble_scalar(fem.form(ufl.inner(e, e) * ds, entity_maps=entity_maps)))
)
assert np.isclose(e_L2, 0.0)
# Write to file
with io.VTXWriter(msh.comm, "ubar.bp", ubar) as f:
f.write(0.0)