File: cavity-farfield.py

package info (click to toggle)
meep-openmpi 1.7.0-3
  • links: PTS, VCS
  • area: main
  • in suites: buster
  • size: 25,828 kB
  • sloc: cpp: 27,370; python: 10,574; lisp: 1,213; makefile: 437; sh: 28
file content (57 lines) | stat: -rw-r--r-- 2,237 bytes parent folder | download | duplicates (8)
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
import meep as mp


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

sy = 6                             # size of cell in y direction (perpendicular to wvg.)
pad = 2                            # padding between last hole and PML edge
dpml = 1                           # PML thickness
sx = 2 * (pad + dpml + N) + d - 1  # size of cell in x direction

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)
resolution = 20

fcen = 0.25                        # pulse center frequency
df = 0.2                           # pulse width (in frequency)

sources = mp.Source(src=mp.GaussianSource(fcen, fwidth=df), component=mp.Hz, center=mp.Vector3())

symmetries = [mp.Mirror(mp.Y, phase=-1), mp.Mirror(mp.X, phase=-1)]

d1 = 0.2

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(2 * dpml - sx)),
    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))
)

sim.run(until_after_sources=mp.stop_when_fields_decayed(50, mp.Hz, mp.Vector3(0.12, -0.37), 1e-8))

d2 = 20
h = 4

sim.output_farfields(nearfield, "spectra-{}-{}-{}".format(d1, d2, h), resolution,
                     where=mp.Volume(mp.Vector3(0, (0.5 * w) + d2 + (0.5 * h)),
                                     size=mp.Vector3(sx - 2 * dpml, h)))