File: mark_entities.py

package info (click to toggle)
scifem 0.16.1-1
  • links: PTS, VCS
  • area: main
  • in suites: sid
  • size: 2,404 kB
  • sloc: python: 4,293; cpp: 323; sh: 16; makefile: 3
file content (101 lines) | stat: -rw-r--r-- 2,479 bytes parent folder | download
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
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
# # Entity markers
#
# Author: Jørgen S. Dokken
#
# SPDX-License-Identifier: MIT
#
# DOLFINx gives you full control for marking entities.
# However, sometimes this can feel a bit repetative.
# In this example we will show how to use {py:func}`scifem.create_entity_markers`.

from mpi4py import MPI
import dolfinx
import numpy as np

# We start by creating a simple unit square

mesh = dolfinx.mesh.create_unit_square(MPI.COMM_WORLD, 60, 60)

# Next, we want to mark some of the cells in our domain.
# We create three marker functions below.
# Each of them takes in a ``(3, num_points)`` array, and returns a boolean array of size ``num_points``.


def left(x):
    return x[0] < 0.2


def right(x):
    return x[0] > 0.9


def inner(x):
    # We use numpy bit operator `&` for "and"
    return (x[0] > 0.3) & (x[0] < 0.7)


# We want to mark these entities with  with unique integers 1, 3 and 7.

from scifem import create_entity_markers

cell_tag = create_entity_markers(mesh, mesh.topology.dim, [(1, left), (3, right), (7, inner)])

# Next we can plot these marked entities

import pyvista

vtk_grid = dolfinx.plot.vtk_mesh(mesh, cell_tag.dim, cell_tag.indices)
grid = pyvista.UnstructuredGrid(*vtk_grid)
grid.cell_data["Marker"] = cell_tag.values

# Create plotter

plotter = pyvista.Plotter()
plotter.add_mesh(grid)
plotter.view_xy()
if not pyvista.OFF_SCREEN:
    plotter.show()


# We can also mark lower order entities, such as facets


def circle(x):
    return x[0] ** 2 + x[1] ** 2 <= 0.16**2


def top(x):
    return x[1] > 0.9


facet_tags = create_entity_markers(mesh, mesh.topology.dim - 1, [(2, top), (7, circle)])


facet_grid = dolfinx.plot.vtk_mesh(mesh, facet_tags.dim, facet_tags.indices)

fgrid = pyvista.UnstructuredGrid(*facet_grid)
fgrid.cell_data["Marker"] = facet_tags.values

fplotter = pyvista.Plotter()
fplotter.add_mesh(fgrid)
fplotter.view_xy()
if not pyvista.OFF_SCREEN:
    fplotter.show()

# We can also exclude interior facets by adding `on_boundary: True` (by default this is set to False).

boundary_facet_tags = create_entity_markers(
    mesh, mesh.topology.dim - 1, [(2, top, True), (7, circle, False)]
)


boundary_grid = dolfinx.plot.vtk_mesh(mesh, boundary_facet_tags.dim, boundary_facet_tags.indices)

bfgrid = pyvista.UnstructuredGrid(*boundary_grid)
bfgrid.cell_data["Marker"] = boundary_facet_tags.values

bfplotter = pyvista.Plotter()
bfplotter.add_mesh(bfgrid)
bfplotter.view_xy()
if not pyvista.OFF_SCREEN:
    bfplotter.show()