File: demo_types.py

package info (click to toggle)
fenics-dolfinx 1%3A0.9.0-11
  • links: PTS, VCS
  • area: main
  • in suites: forky
  • size: 5,376 kB
  • sloc: cpp: 33,701; python: 22,338; makefile: 230; sh: 171; xml: 55
file content (236 lines) | stat: -rw-r--r-- 7,325 bytes parent folder | download | duplicates (4)
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
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
# ---
# jupyter:
#   jupytext:
#     text_representation:
#       extension: .py
#       format_name: light
#       format_version: '1.5'
#       jupytext_version: 1.13.6
# ---

# # Solving PDEs with different scalar (float) types
#
# This demo  ({download}`demo_types.py`) shows:
#
# - How to solve problems using different scalar types, .e.g. single or
#   double precision, or complex numbers
# - Interfacing with [SciPy](https://scipy.org/) sparse linear algebra
#   functionality

import sys

from mpi4py import MPI

# +
import numpy as np
import scipy.sparse
import scipy.sparse.linalg

import ufl
from dolfinx import fem, la, mesh, plot

# -

# SciPy solvers do not support MPI, so all computations will be
# performed on a single MPI rank

# +
comm = MPI.COMM_SELF
# -


# Create a function that solves the Poisson equation using different
# precision float and complex scalar types for the finite element
# solution.


def display_scalar(u, name, filter=np.real):
    """Plot the solution using pyvista"""
    try:
        import pyvista

        cells, types, x = plot.vtk_mesh(u.function_space)
        grid = pyvista.UnstructuredGrid(cells, types, x)
        grid.point_data["u"] = filter(u.x.array)
        grid.set_active_scalars("u")
        plotter = pyvista.Plotter()
        plotter.add_mesh(grid, show_edges=True)
        plotter.add_mesh(grid.warp_by_scalar())
        plotter.add_title(f"{name}: real" if filter is np.real else f"{name}: imag")
        if pyvista.OFF_SCREEN:
            pyvista.start_xvfb(wait=0.1)
            plotter.screenshot(f"u_{'real' if filter is np.real else 'imag'}.png")
        else:
            plotter.show()
    except ModuleNotFoundError:
        print("'pyvista' is required to visualise the solution")


def display_vector(u, name, filter=np.real):
    """Plot the solution using pyvista"""
    try:
        import pyvista

        V = u.function_space
        cells, types, x = plot.vtk_mesh(V)
        grid = pyvista.UnstructuredGrid(cells, types, x)
        bs = V.dofmap.index_map_bs
        grid.point_data["u"] = filter(np.insert(u.x.array.reshape(x.shape[0], bs), bs, 0, axis=1))
        plotter = pyvista.Plotter()
        plotter.add_mesh(grid.warp_by_scalar(), show_edges=True)
        plotter.add_title(f"{name}: real" if filter is np.real else f"{name}: imag")
        if pyvista.OFF_SCREEN:
            pyvista.start_xvfb(wait=0.1)
            plotter.screenshot(f"u_{'real' if filter is np.real else 'imag'}.png")
        else:
            plotter.show()
    except ModuleNotFoundError:
        print("'pyvista' is required to visualise the solution")


def poisson(dtype):
    """Poisson problem solver

    Args:
        dtype: Scalar type to use.


    """

    # Create a mesh and locate facets by a geometric condition
    msh = mesh.create_rectangle(
        comm=comm,
        points=((0.0, 0.0), (2.0, 1.0)),
        n=(32, 16),
        cell_type=mesh.CellType.triangle,
        dtype=np.real(dtype(0)).dtype,
    )
    facets = mesh.locate_entities_boundary(
        msh, dim=1, marker=lambda x: np.isclose(x[0], 0.0) | np.isclose(x[0], 2.0)
    )

    # Define a variational problem.
    V = fem.functionspace(msh, ("Lagrange", 1))
    u, v = ufl.TrialFunction(V), ufl.TestFunction(V)
    x = ufl.SpatialCoordinate(msh)
    fr = 10 * ufl.exp(-((x[0] - 0.5) ** 2 + (x[1] - 0.5) ** 2) / 0.02)
    fc = ufl.sin(2 * np.pi * x[0]) + 10 * ufl.sin(4 * np.pi * x[1]) * 1j
    gr = ufl.sin(5 * x[0])
    gc = ufl.sin(5 * x[0]) * 1j
    a = ufl.inner(ufl.grad(u), ufl.grad(v)) * ufl.dx
    L = ufl.inner(fr + fc, v) * ufl.dx + ufl.inner(gr + gc, v) * ufl.ds

    # In preparation for constructing Dirichlet boundary conditions, locate
    # facets on the constrained boundary and the corresponding
    # degrees-of-freedom.
    dofs = fem.locate_dofs_topological(V=V, entity_dim=1, entities=facets)

    # Process forms. This will compile the forms for the requested type.
    a0 = fem.form(a, dtype=dtype)
    if np.issubdtype(dtype, np.complexfloating):
        L0 = fem.form(L, dtype=dtype)
    else:
        L0 = fem.form(ufl.replace(L, {fc: 0, gc: 0}), dtype=dtype)

    # Create a Dirichlet boundary condition
    bc = fem.dirichletbc(value=dtype(0), dofs=dofs, V=V)

    # Assemble forms
    A = fem.assemble_matrix(a0, [bc])
    A.scatter_reverse()
    b = fem.assemble_vector(L0)
    fem.apply_lifting(b.array, [a0], bcs=[[bc]])
    b.scatter_reverse(la.InsertMode.add)
    bc.set(b.array)

    # Create a SciPy CSR  matrix that shares data with A and solve
    As = A.to_scipy()
    uh = fem.Function(V, dtype=dtype)
    uh.x.array[:] = scipy.sparse.linalg.spsolve(As, b.array)

    return uh


# Create a function that solves the linearised elasticity equation using
# different precision float and complex scalar types for the finite
# element solution.


def elasticity(dtype) -> fem.Function:
    """Linearised elasticity problem solver."""

    # Create a mesh and locate facets by a geometric condition
    msh = mesh.create_rectangle(
        comm=comm,
        points=((0.0, 0.0), (2.0, 1.0)),
        n=(32, 16),
        cell_type=mesh.CellType.triangle,
        dtype=np.real(dtype(0)).dtype,
    )
    facets = mesh.locate_entities_boundary(
        msh, dim=1, marker=lambda x: np.isclose(x[0], 0.0) | np.isclose(x[0], 2.0)
    )

    # Define the variational problem.
    gdim = msh.geometry.dim
    V = fem.functionspace(msh, ("Lagrange", 1, (gdim,)))
    ω, ρ = 300.0, 10.0
    x = ufl.SpatialCoordinate(msh)
    f = ufl.as_vector((ρ * ω**2 * x[0], ρ * ω**2 * x[1]))

    E, ν = 1.0e6, 0.3
    μ, λ = E / (2.0 * (1.0 + ν)), E * ν / ((1.0 + ν) * (1.0 - 2.0 * ν))

    def σ(v):
        """Return an expression for the stress σ given a displacement field"""
        return 2.0 * μ * ufl.sym(ufl.grad(v)) + λ * ufl.tr(ufl.sym(ufl.grad(v))) * ufl.Identity(
            len(v)
        )

    u, v = ufl.TrialFunction(V), ufl.TestFunction(V)
    a = ufl.inner(σ(u), ufl.grad(v)) * ufl.dx
    L = ufl.inner(f, v) * ufl.dx

    dofs = fem.locate_dofs_topological(V=V, entity_dim=1, entities=facets)

    # Process forms. This will compile the forms for the requested type.
    a0, L0 = fem.form(a, dtype=dtype), fem.form(L, dtype=dtype)

    # Create a Dirichlet boundary condition
    bc = fem.dirichletbc(np.zeros(2, dtype=dtype), dofs, V=V)

    # Assemble forms (CSR matrix)
    A = fem.assemble_matrix(a0, [bc])
    A.scatter_reverse()

    b = fem.assemble_vector(L0)
    fem.apply_lifting(b.array, [a0], bcs=[[bc]])
    b.scatter_reverse(la.InsertMode.add)
    bc.set(b.array)

    # Create a SciPy CSR  matrix that shares data with A and solve
    As = A.to_scipy()
    uh = fem.Function(V, dtype=dtype)
    uh.x.array[:] = scipy.sparse.linalg.spsolve(As, b.array)

    return uh


# Solve problems for different types


uh = poisson(dtype=np.float32)
uh = poisson(dtype=np.float64)
if not sys.platform.startswith("win32"):
    uh = poisson(dtype=np.complex64)
    uh = poisson(dtype=np.complex128)
display_scalar(uh, "poisson", np.real)
display_scalar(uh, "poisson", np.imag)


uh = elasticity(dtype=np.float32)
uh = elasticity(dtype=np.float64)
if not sys.platform.startswith("win32"):
    uh = elasticity(dtype=np.complex64)
    uh = elasticity(dtype=np.complex128)
display_vector(uh, "elasticity", np.real)