-
Notifications
You must be signed in to change notification settings - Fork 0
/
maxwell_eigen.py
162 lines (131 loc) · 5.05 KB
/
maxwell_eigen.py
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
# Based on [1], updated to work in FEniCSx.
# References:
# [1]: https://fenicsproject.org/olddocs/dolfin/latest/python/demos/maxwell-eigenvalues/demo_maxwell-eigenvalues.py.html
# [2]: https://bitbucket.org/fenics-project/dolfin/src/master/dolfin/la/SLEPcEigenSolver.cpp
# [3]: https://slepc.upv.es/slepc4py-current/docs/apiref/slepc4py.SLEPc-module.html
from slepc4py import SLEPc
from ufl import dx, curl, inner, TrialFunction, TestFunction
import numpy as np
from dolfinx.fem import (dirichletbc, Function, functionspace, form,
locate_dofs_topological)
from mpi4py import MPI
from dolfinx.fem.petsc import assemble_matrix
from petsc4py import PETSc
from dolfinx.mesh import (create_rectangle, locate_entities_boundary,
CellType, GhostMode, DiagonalType)
import sys
def par_print(string):
if comm.rank == 0:
print(string)
sys.stdout.flush()
def eigenvalues(n_eigs, shift, V, bcs):
# Define problem
u, v = TrialFunction(V), TestFunction(V)
a = form(inner(curl(u), curl(v)) * dx)
b = form(inner(u, v) * dx)
# Assemble matrices
A = assemble_matrix(a, bcs)
A.assemble()
# Zero rows of boundary DOFs of B. See [1]
B = assemble_matrix(b, bcs, diagonal=0.0)
B.assemble()
# Create SLEPc Eigenvalue solver
eps = SLEPc.EPS().create(comm)
eps.setOperators(A, B)
eps.setType(SLEPc.EPS.Type.KRYLOVSCHUR)
eps.setProblemType(SLEPc.EPS.ProblemType.GHEP)
eps.setWhichEigenpairs(eps.Which.TARGET_MAGNITUDE)
eps.setTarget(shift)
st = eps.getST()
st.setType(SLEPc.ST.Type.SINVERT)
st.setShift(shift)
eps.setDimensions(n_eigs, PETSc.DECIDE, PETSc.DECIDE)
eps.setFromOptions()
eps.solve()
its = eps.getIterationNumber()
eps_type = eps.getType()
n_ev, n_cv, mpd = eps.getDimensions()
tol, max_it = eps.getTolerances()
n_conv = eps.getConverged()
par_print(f"Number of iterations: {its}")
par_print(f"Solution method: {eps_type}")
par_print(f"Number of requested eigenvalues: {n_ev}")
par_print(f"Stopping condition: tol={tol}, maxit={max_it}")
par_print(f"Number of converged eigenpairs: {n_conv}")
computed_eigenvalues = []
for i in range(min(n_conv, n_eigs)):
lmbda = eps.getEigenvalue(i)
computed_eigenvalues.append(np.round(np.real(lmbda), 1))
eps.destroy()
return np.sort(computed_eigenvalues)
def boundary(x):
"Boundary marker"
return boundary_lr(x) | boundary_tb(x)
def boundary_lr(x):
"Left and right boundary marker"
return np.isclose(x[0], 0.0) | np.isclose(x[0], np.pi)
def boundary_tb(x):
"Top and bottom boundary marker"
return np.isclose(x[1], 0.0) | np.isclose(x[1], np.pi)
def print_eigenvalues(mesh):
# Nédélec
V_n = functionspace(mesh, ("N1curl", 1))
# Set boundary DOFs to 0 (u x n = 0 on \partial \Omega).
ud_n = Function(V_n)
f_dim = mesh.topology.dim - 1
boundary_facets = locate_entities_boundary(mesh, f_dim, boundary)
boundary_dofs_n = locate_dofs_topological(
V_n, f_dim, boundary_facets)
bcs_nedelec = [dirichletbc(ud_n, boundary_dofs_n)]
# Solve Maxwell eigenvalue problem
eigenvalues_nedelec = eigenvalues(n_eigs, shift, V_n, bcs_nedelec)
# Lagrange
V_l = functionspace(mesh, ("Lagrange", 1, (mesh.geometry.dim,)))
W = functionspace(mesh, ("Lagrange", 1))
ud_l = Function(W)
# Must constrain horizontal DOFs on horizontal faces and vertical DOFs
# on vertical faces
boundary_facets_tb = locate_entities_boundary(mesh, f_dim, boundary_tb)
boundary_facets_lr = locate_entities_boundary(mesh, f_dim, boundary_lr)
V_l_0, V_l_1 = V_l.sub(0), V_l.sub(1)
dofs_tb = locate_dofs_topological(
(V_l_0, W), f_dim, boundary_facets_tb)
dofs_lr = locate_dofs_topological(
(V_l_1, W), f_dim, boundary_facets_lr)
bcs_lagrange = [dirichletbc(ud_l, dofs_tb, V_l_0),
dirichletbc(ud_l, dofs_lr, V_l_1)]
# Solve Maxwell eigenvalue problem
eigenvalues_lagrange = eigenvalues(n_eigs, shift, V_l, bcs_lagrange)
# Print results
np.set_printoptions(formatter={'float': '{:5.1f}'.format})
eigenvalues_exact = np.sort(np.array([float(m**2 + n**2)
for m in range(6)
for n in range(6)]))[1:13]
par_print(f"Exact = {eigenvalues_exact}")
par_print(f"Nédélec = {eigenvalues_nedelec}")
par_print(f"Lagrange = {eigenvalues_lagrange}")
# Number of element in each direction
n = 40
# Number of eigernvalues to compute
n_eigs = 12
# Find eigenvalues near
shift = 5.5
# Domain corners
corners = ((0.0, 0.0), (np.pi, np.pi))
comm = MPI.COMM_WORLD
par_print("Right diagonal mesh:")
mesh = create_rectangle(
comm,
corners, (n, n),
CellType.triangle,
ghost_mode=GhostMode.none,
diagonal=DiagonalType.right)
print_eigenvalues(mesh)
par_print("\nCrossed diagonal mesh:")
mesh = create_rectangle(
comm,
corners, (n, n),
CellType.triangle,
ghost_mode=GhostMode.none,
diagonal=DiagonalType.crossed)
print_eigenvalues(mesh)