-
Notifications
You must be signed in to change notification settings - Fork 626
/
cavity-farfield.py
134 lines (109 loc) · 3.36 KB
/
cavity-farfield.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
import matplotlib
import numpy as np
import meep as mp
matplotlib.use("agg")
import matplotlib.pyplot as plt
resolution = 20 # pixels/μm
fcen = 0.25 # pulse center frequency
df = 0.2 # pulse width (in frequency)
eps = 13 # dielectric constant of waveguide
w = 1.2 # width of waveguide
r = 0.36 # radius of holes
d = 1.4 # defect spacing (ordinary spacing = 1)
N = 3 # number of holes on either side of defect
dpad = 32 # padding between last hole and PML edge
dpml = 0.5 / (fcen - 0.5 * df) # PML thickness (> half the largest wavelength)
sx = 2 * (dpad + dpml + N) + d - 1 # size of cell in x direction
d1 = 0.2 # y-distance from waveguide edge to near2far surface
d2 = 2.0 # y-distance from near2far surface to far-field line
sy = w + 2 * (d1 + d2 + dpml) # size of cell in y direction (perpendicular to wvg.)
cell = mp.Vector3(sx, sy, 0)
geometry = [
mp.Block(
center=mp.Vector3(),
size=mp.Vector3(mp.inf, w, mp.inf),
material=mp.Medium(epsilon=eps),
)
]
for i in range(N):
geometry.append(mp.Cylinder(r, center=mp.Vector3(d / 2 + i)))
geometry.append(mp.Cylinder(r, center=mp.Vector3(d / -2 - i)))
pml_layers = [mp.PML(dpml)]
sources = [
mp.Source(
src=mp.GaussianSource(fcen, fwidth=df), component=mp.Hz, center=mp.Vector3()
)
]
symmetries = [mp.Mirror(mp.X, phase=-1), mp.Mirror(mp.Y, phase=-1)]
sim = mp.Simulation(
cell_size=cell,
geometry=geometry,
sources=sources,
symmetries=symmetries,
boundary_layers=pml_layers,
resolution=resolution,
)
nearfield = sim.add_near2far(
fcen,
0,
1,
mp.Near2FarRegion(mp.Vector3(0, 0.5 * w + d1), size=mp.Vector3(sx - 2 * dpml)),
mp.Near2FarRegion(
mp.Vector3(-0.5 * sx + dpml, 0.5 * w + 0.5 * d1),
size=mp.Vector3(0, d1),
weight=-1.0,
),
mp.Near2FarRegion(
mp.Vector3(0.5 * sx - dpml, 0.5 * w + 0.5 * d1), size=mp.Vector3(0, d1)
),
)
mon = sim.add_dft_fields(
[mp.Hz],
fcen,
0,
1,
center=mp.Vector3(0, 0.5 * w + d1 + d2),
size=mp.Vector3(sx - 2 * (dpad + dpml), 0),
)
sim.run(until_after_sources=mp.stop_when_dft_decayed())
sim.plot2D()
if mp.am_master():
plt.savefig(
f"cavity_farfield_plot2D_dpad{dpad}_{d1}_{d2}.png", bbox_inches="tight", dpi=150
)
Hz_mon = sim.get_dft_array(mon, mp.Hz, 0)
(x, y, z, w) = sim.get_array_metadata(dft_cell=mon)
ff = []
for xc in x:
ff_pt = sim.get_farfield(nearfield, mp.Vector3(xc, y[0]))
ff.append(ff_pt[5])
ff = np.array(ff)
if mp.am_master():
plt.figure()
plt.subplot(1, 3, 1)
plt.plot(np.real(Hz_mon), "bo-", label="DFT")
plt.plot(np.real(ff), "ro-", label="N2F")
plt.legend()
plt.xlabel("$x$ (μm)")
plt.ylabel("real(Hz)")
plt.subplot(1, 3, 2)
plt.plot(np.imag(Hz_mon), "bo-", label="DFT")
plt.plot(np.imag(ff), "ro-", label="N2F")
plt.legend()
plt.xlabel("$x$ (μm)")
plt.ylabel("imag(Hz)")
plt.subplot(1, 3, 3)
plt.plot(np.abs(Hz_mon), "bo-", label="DFT")
plt.plot(np.abs(ff), "ro-", label="N2F")
plt.legend()
plt.xlabel("$x$ (μm)")
plt.ylabel("|Hz|")
plt.suptitle(
f"comparison of near2far and actual DFT fields\n dpad={dpad}, d1={d1}, d2={d2}"
)
plt.subplots_adjust(wspace=0.6)
plt.savefig(
f"test_Hz_dft_vs_n2f_res{resolution}_dpad{dpad}_d1{d1}_d2{d2}.png",
bbox_inches="tight",
dpi=150,
)