Copyright (C) 2020 Andreas Kloeckner
from firedrake import *
import numpy as np
import numpy.linalg as la
import firedrake.mesh as fd_mesh
import matplotlib.pyplot as plt
firedrake:WARNING OMP_NUM_THREADS is not set or is set to a value greater than 1, we suggest setting OMP_NUM_THREADS=1 to improve performance
import meshpy.triangle as triangle
def round_trip_connect(start, end):
return [(i, i+1) for i in range(start, end)] + [(end, start)]
def make_mesh():
points = [(-1, -1), (1, -1), (1, 1), (-1, 1)]
facets = round_trip_connect(0, len(points)-1)
facet_markers = [1, 2, 3, 4]
circ_start = len(points)
nsegments = 30
points.extend(
(0.25 * np.cos(angle), 0.25 * np.sin(angle))
for angle in np.linspace(0, 2*np.pi, nsegments, endpoint=False))
facets.extend(round_trip_connect(circ_start, len(points)-1))
facet_markers.extend([-1] * nsegments)
def needs_refinement(vertices, area):
bary = np.sum(np.array(vertices), axis=0)/3
max_area = 0.01 + la.norm(bary, np.inf)*0.01
return bool(area > max_area)
info = triangle.MeshInfo()
info.set_points(points)
info.set_facets(facets, facet_markers=facet_markers)
built_mesh = triangle.build(info, refinement_func=needs_refinement)
plex = fd_mesh._from_cell_list(
2, np.array(built_mesh.elements), np.array(built_mesh.points), COMM_WORLD)
import firedrake.cython.dmcommon as dmcommon
v_start, v_end = plex.getDepthStratum(0) # vertices
for facet, fmarker in zip(built_mesh.facets, built_mesh.facet_markers):
vertices = [fvert + v_start for fvert in facet]
if fmarker > 0: # interior facets are negative above
join = plex.getJoin(vertices)
plex.setLabelValue(dmcommon.FACE_SETS_LABEL, join[0], fmarker)
return Mesh(plex)
mesh = make_mesh()
triplot(mesh)
plt.legend()
<matplotlib.legend.Legend at 0x7f75964a6bf0>
V = FunctionSpace(mesh, 'Lagrange', 1)
x = SpatialCoordinate(mesh)
f = conditional(le(sqrt(x[0]**2 + x[1]**2), 0.25), 25.0, 0.0)
kappa = conditional(le(sqrt(x[0]**2 + x[1]**2), 0.25), 25.0, 0.1)
g = interpolate(20.0 * (1-x[0]*x[0]), V)
bcg = DirichletBC(V, g, [1])
bcz = DirichletBC(V, 0.0, [2,3,4])
bc = [bcg, bcz]
u = TrialFunction(V)
v = TestFunction(V)
v
Argument(WithGeometry(FunctionSpace(<firedrake.mesh.MeshTopology object at 0x7f75964a6a40>, FiniteElement('Lagrange', triangle, 1), name=None), Mesh(VectorElement(FiniteElement('Lagrange', triangle, 1), dim=2), 1)), 0, None)
a = kappa*inner(grad(u), grad(v))*dx
L = f*v*dx
U = Function(V)
solve(a == L, U, bc)
tricontourf(U)
<matplotlib.tri.tricontour.TriContourSet at 0x7f758d810dc0>
fig = plt.figure(figsize=(10, 10))
axes = fig.add_subplot(111, projection='3d')
trisurf(U, axes=axes)
<mpl_toolkits.mplot3d.art3d.Poly3DCollection at 0x7f758d813940>