-
Notifications
You must be signed in to change notification settings - Fork 626
/
simulation.py
6402 lines (5372 loc) · 252 KB
/
simulation.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
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
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
"""
A collection of objects and helper routines for setting up the simulation.
"""
import functools
import inspect
import math
import numbers
import os
import re
import signal
import subprocess
import sys
import warnings
from collections import OrderedDict, namedtuple
from typing import Callable, List, NamedTuple, Optional, Tuple, Union
from scipy.interpolate import pade
try:
from collections.abc import Sequence, Iterable
except ImportError:
from collections.abc import Sequence, Iterable
import numpy as np
from meep.geom import GeometricObject, Medium, Vector3, init_do_averaging
from meep.source import (
EigenModeSource,
GaussianBeamSource,
IndexedSource,
Source,
check_positive,
)
from meep.verbosity_mgr import Verbosity
import meep as mp
try:
basestring
except NameError:
basestring = str
try:
from IPython.display import display
from ipywidgets import FloatProgress
do_progress = True
except ImportError:
do_progress = False
from matplotlib.axes import Axes
verbosity = Verbosity(mp.cvar, "meep", 1)
mp.setup()
# Send output from Meep, ctlgeom, and MPB to Python's stdout
mp.set_meep_printf_callback(mp.py_master_printf_wrap)
mp.set_meep_printf_stderr_callback(mp.py_master_printf_stderr_wrap)
mp.set_ctl_printf_callback(mp.py_master_printf_wrap)
mp.set_mpb_printf_callback(mp.py_master_printf_wrap)
EigCoeffsResult = namedtuple(
"EigCoeffsResult", ["alpha", "vgrp", "kpoints", "kdom", "cscale"]
)
FluxData = namedtuple("FluxData", ["E", "H"])
ForceData = namedtuple("ForceData", ["offdiag1", "offdiag2", "diag"])
NearToFarData = namedtuple("NearToFarData", ["F"])
Vector3Type = Union[Vector3, Tuple[float, ...]]
def fix_dft_args(args, i):
if (
len(args) > i + 2
and isinstance(args[i], (int, float))
and isinstance(args[i + 1], (int, float))
and isinstance(args[i + 2], int)
):
fcen = args[i]
df = args[i + 1]
nfreq = args[i + 2]
freq = (
[fcen]
if nfreq == 1
else np.linspace(fcen - 0.5 * df, fcen + 0.5 * df, nfreq)
)
return args[:i] + (freq,) + args[i + 3 :]
elif not isinstance(args[i], (np.ndarray, list)):
raise TypeError(
"add_dft functions only accept fcen,df,nfreq (3 numbers) or freq (array/list)"
)
else:
return args
def get_num_args(func):
if isinstance(func, Harminv) or isinstance(func, PadeDFT):
return 2
elif inspect.ismethod(func):
return func.__code__.co_argcount - 1 # remove 'self' from count
else:
return func.__code__.co_argcount
def vec(*args):
try:
# Check for vec(x, [y, [z]])
return mp._vec(*args)
except (TypeError, NotImplementedError):
try:
# Check for vec(iterable)
if len(args) != 1:
raise TypeError
return mp._vec(*args[0])
except (TypeError, NotImplementedError):
print("Expected an iterable with three or fewer floating point values")
print(" or something of the form vec(x, [y, [z]])")
raise
def py_v3_to_vec(dims: int, iterable: Iterable, is_cylindrical: bool = False):
v3 = Vector3(*iterable)
if dims == 1:
return mp.vec(v3.z)
elif dims == 2:
if is_cylindrical:
return mp.veccyl(v3.x, v3.z)
v = mp.vec(v3.x, v3.y)
v.set_direction(mp.Z, v3.z) # for special_kz handling
return v
elif dims == 3:
return mp.vec(v3.x, v3.y, v3.z)
else:
raise ValueError(f"Invalid dimensions in Volume: {dims}")
def bands_to_diffractedplanewave(where, bands):
if bands.axis is None:
if where.in_direction(mp.X) != 0:
axis = np.array([1, 0, 0], dtype=np.float64)
elif where.in_direction(mp.Y) != 0:
axis = np.array([0, 1, 0], dtype=np.float64)
elif where.in_direction(mp.Z) != 0:
axis = np.array([0, 0, 1], dtype=np.float64)
else:
raise ValueError(
"axis parameter of DiffractedPlanewave must be a non-zero Vector3"
)
elif isinstance(bands.axis, mp.Vector3):
axis = np.array([bands.axis.x, bands.axis.y, bands.axis.z], dtype=np.float64)
else:
raise TypeError("axis parameter of DiffractedPlanewave must be a Vector3")
diffractedplanewave_args = [
np.array(bands.g, dtype=np.intc),
axis,
bands.s * 1.0,
bands.p * 1.0,
]
return mp.diffractedplanewave(*diffractedplanewave_args)
class DiffractedPlanewave:
"""
For mode decomposition or eigenmode source, specify a diffracted planewave in homogeneous media. Should be passed as the `bands` argument of `get_eigenmode_coefficients`, `band_num` of `get_eigenmode`, or `eig_band` of `EigenModeSource`.
"""
def __init__(
self,
g: List[int] = None,
axis: Vector3Type = None,
s: complex = None,
p: complex = None,
):
"""
Construct a `DiffractedPlanewave`.
+ **`g` [ list of 3 `integer`s ]** — The diffraction order $(m_x,m_y,m_z)$ corresponding to the wavevector $(k_x+2\\pi m_x/\\Lambda_x,k_y+2\\pi m_y/\\Lambda_y,k_z+2\\pi m_z/\\Lambda_z)$. $(k_x,k_y,k_z)$ is the `k_point` (wavevector specifying the Bloch-periodic boundaries) of the `Simulation` class object. The diffraction order $m_{x,y,z}$ should be non-zero only in the $d$-1 periodic directions of a $d$ dimensional cell of size $(\\Lambda_x,\\Lambda_y,\\Lambda_z)$ (e.g., a plane in 3d) in which the mode monitor or source extends the entire length of the cell.
+ **`axis` [ `Vector3` ]** — The plane of incidence for each planewave (used to define the $\\mathcal{S}$ and $\\mathcal{P}$ polarizations below) is defined to be the plane that contains the `axis` vector and the planewave's wavevector. If `None`, `axis` defaults to the first direction that lies in the plane of the monitor or source (e.g., $y$ direction for a $yz$ plane in 3d, either $x$ or $y$ in 2d).
+ **`s` [ `complex` ]** — The complex amplitude of the $\\mathcal{S}$ polarziation (i.e., electric field perpendicular to the plane of incidence).
+ **`p` [ `complex` ]** — The complex amplitude of the $\\mathcal{P}$ polarziation (i.e., electric field parallel to the plane of incidence).
"""
self._g = g
self._axis = axis
self._s = complex(s)
self._p = complex(p)
@property
def g(self):
return self._g
@property
def axis(self):
return self._axis
@property
def s(self):
return self._s
@property
def p(self):
return self._p
DefaultPMLProfile = lambda u: u * u
class PML:
"""
This class is used for specifying the PML absorbing boundary layers around the cell,
if any, via the `boundary_layers` input variable. See also [Perfectly Matched
Layers](Perfectly_Matched_Layer.md). `boundary_layers` can be zero or more `PML`
objects, with multiple objects allowing you to specify different PML layers on
different boundaries. The class represents a single PML layer specification, which
sets up one or more PML layers around the boundaries according to the following
properties.
"""
def __init__(
self,
thickness: float = None,
direction: int = mp.ALL,
side: int = mp.ALL,
R_asymptotic: float = 1e-15,
mean_stretch: float = 1.0,
pml_profile: Callable[[float], float] = DefaultPMLProfile,
):
"""
+ **`thickness` [ `number` ]** — The spatial thickness of the PML layer which
extends from the boundary towards the *inside* of the cell. The thinner it is,
the more numerical reflections become a problem. No default value.
+ **`direction` [ `direction` constant ]** — Specify the direction of the
boundaries to put the PML layers next to. e.g. if `X`, then specifies PML on the
$\\pm x$ boundaries (depending on the value of `side`, below). Default is the
special value `ALL`, which puts PML layers on the boundaries in all directions.
+ **`side` [ `side` constant ]** — Specify which side, `Low` or `High` of the
boundary or boundaries to put PML on. e.g. if side is `Low` and direction is
`meep.X`, then a PML layer is added to the $-x$ boundary. Default is the special
value `meep.ALL`, which puts PML layers on both sides.
+ **`R_asymptotic` [ `number` ]** — The asymptotic reflection in the limit of
infinite resolution or infinite PML thickness, for reflections from air (an
upper bound for other media with index > 1). For a finite resolution or
thickness, the reflection will be *much larger*, due to the discretization of
Maxwell's equation. Default value is 10<sup>−15</sup>, which should suffice for
most purposes. You want to set this to be small enough so that waves propagating
within the PML are attenuated sufficiently, but making `R_asymptotic` too small
will increase the numerical reflection due to discretization.
+ **`pml_profile` [ `function` ]** — By default, Meep turns on the PML conductivity
quadratically within the PML layer — one doesn't want to turn it on
suddenly, because that exacerbates reflections due to the discretization. More
generally, with `pml_profile` one can specify an arbitrary PML "profile"
function $f(u)$ that determines the shape of the PML absorption profile up to an
overall constant factor. *u* goes from 0 to 1 at the start and end of the PML,
and the default is $f(u) = u^2$. In some cases where a very thick PML is
required, such as in a periodic medium (where there is technically no such thing
as a true PML, only a pseudo-PML), it can be advantageous to turn on the PML
absorption more smoothly. See [Optics Express, Vol. 16, pp. 11376-92
(2008)](http://www.opticsinfobase.org/abstract.cfm?URI=oe-16-15-11376). For
example, one can use a cubic profile $f(u) = u^3$ by specifying
`pml_profile=lambda u: u*u*u`.
"""
if thickness is None:
raise ValueError("PML thickness must be specified.")
self.thickness = thickness
self.direction = direction
self.side = side
self.R_asymptotic = R_asymptotic
self.mean_stretch = mean_stretch
self.pml_profile = pml_profile
if direction == mp.ALL and side == mp.ALL:
self.swigobj = mp.pml(thickness, R_asymptotic, mean_stretch)
elif direction == mp.ALL:
self.swigobj = mp.pml(thickness, side, R_asymptotic, mean_stretch)
else:
self.swigobj = mp.pml(
thickness, direction, side, R_asymptotic, mean_stretch
)
@property
def R_asymptotic(self):
return self._R_asymptotic
@R_asymptotic.setter
def R_asymptotic(self, val):
self._R_asymptotic = check_positive("PML.R_asymptotic", val)
@property
def mean_stretch(self):
return self._mean_stretch
@mean_stretch.setter
def mean_stretch(self, val: float):
if val >= 1.0:
self._mean_stretch = val
else:
raise ValueError(f"PML.mean_stretch must be >= 1. Got {val}")
class Absorber(PML):
"""
Instead of a `PML` layer, there is an alternative class called `Absorber` which is a
**drop-in** replacement for `PML`. For example, you can do
`boundary_layers=[mp.Absorber(thickness=2)]` instead of
`boundary_layers=[mp.PML(thickness=2)]`. All the parameters are the same as for `PML`,
above. You can have a mix of `PML` on some boundaries and `Absorber` on others.
The `Absorber` class does *not* implement a perfectly matched layer (PML), however
(except in 1d). Instead, it is simply a scalar electric **and** magnetic conductivity
that turns on gradually within the layer according to the `pml_profile` (defaulting to
quadratic). Such a scalar conductivity gradient is only reflectionless in the limit as
the layer becomes sufficiently thick.
The main reason to use `Absorber` is if you have **a case in which PML fails:**
- No true PML exists for *periodic* media, and a scalar absorber is computationally
less expensive and generally just as good. See [Optics Express, Vol. 16, pp.
11376-92 (2008)](http://www.opticsinfobase.org/abstract.cfm?URI=oe-16-15-11376).
- PML can lead to *divergent* fields for certain waveguides with "backward-wave"
modes; this can readily occur in metals with surface plasmons, and a scalar
absorber is your only choice. See [Physical Review E, Vol. 79, 065601
(2009)](http://math.mit.edu/~stevenj/papers/LohOs09.pdf).
- PML can fail if you have a waveguide hitting the edge of your cell *at an angle*.
See [J. Computational Physics, Vol. 230, pp. 2369-77
(2011)](http://math.mit.edu/~stevenj/papers/OskooiJo11.pdf).
"""
class Symmetry:
"""
This class is used for the `symmetries` input variable to specify symmetries which
must preserve both the structure *and* the sources. Any number of symmetries can be
exploited simultaneously but there is no point in specifying redundant symmetries: the
cell can be reduced by at most a factor of 4 in 2d and 8 in 3d. See also [Exploiting
Symmetry](Exploiting_Symmetry.md). This is the base class of the specific symmetries
below, so normally you don't create it directly. However, it has two properties which
are shared by all symmetries:
The specific symmetry sub-classes are:
**`Mirror`** — A mirror symmetry plane. `direction` is the direction *normal* to the
mirror plane.
**`Rotate2`** — A 180° (twofold) rotational symmetry (a.k.a. $C_2$). `direction` is
the axis of the rotation.
**`Rotate4`** — A 90° (fourfold) rotational symmetry (a.k.a. $C_4$). `direction` is
the axis of the rotation.
"""
def __init__(self, direction: int = None, phase: complex = 1.0 + 0j):
"""
Construct a `Symmetry`.
+ **`direction` [ `direction` constant ]** — The direction of the symmetry (the
normal to a mirror plane or the axis for a rotational symmetry). e.g. `X`, `Y`,
or `Z` (only Cartesian/grid directions are allowed). No default value.
+ **`phase` [ `complex` ]** — An additional phase to multiply the fields by when
operating the symmetry on them. Default is +1, e.g. a phase of -1 for a mirror
plane corresponds to an *odd* mirror. Technically, you are essentially
specifying the representation of the symmetry group that your fields and sources
transform under.
"""
self.direction = direction
self.phase = complex(phase)
self.swigobj = None
class Rotate2(Symmetry):
"""
A 180° (twofold) rotational symmetry (a.k.a. $C_2$). `direction` is the axis of the
rotation.
"""
class Rotate4(Symmetry):
"""
A 90° (fourfold) rotational symmetry (a.k.a. $C_4$). `direction` is the axis of the
rotation.
"""
class Mirror(Symmetry):
"""
A mirror symmetry plane. `direction` is the direction *normal* to the mirror plane.
"""
class Identity(Symmetry):
""" """
class Volume:
"""
Many Meep functions require you to specify a volume in space, corresponding to the C++
type `meep::volume`. This class creates such a volume object, given the `center` and
`size` properties (just like e.g. a `Block` object). If the `size` is not specified,
it defaults to `(0,0,0)`, i.e. a single point. Any method that accepts such a volume
also accepts `center` and `size` keyword arguments. If these are specified instead of
the volume, the library will construct a volume for you. Alternatively, you can
specify a list of `Vector3` vertices using the `vertices` parameter. The `center` and
`size` will automatically be computed from this list.
"""
def __init__(
self,
center: Vector3Type = Vector3(),
size: Vector3Type = Vector3(),
dims: int = 2,
is_cylindrical: bool = False,
vertices: List[Vector3Type] = [],
):
"""
Construct a Volume.
"""
if len(vertices) == 0:
self.center = Vector3(*center)
self.size = Vector3(*size)
else:
vertices = np.array([np.array(i) for i in vertices])
self.center = Vector3(*np.mean(vertices, axis=0))
x_list = np.unique(vertices[:, 0])
y_list = np.unique(vertices[:, 1])
z_list = np.unique(vertices[:, 2])
x_size = 0 if x_list.size == 1 else np.abs(np.diff(x_list)[0])
y_size = 0 if y_list.size == 1 else np.abs(np.diff(y_list)[0])
z_size = 0 if z_list.size == 1 else np.abs(np.diff(z_list)[0])
self.size = Vector3(x_size, y_size, z_size)
self.dims = dims
v1 = self.center - self.size.scale(0.5)
v2 = self.center + self.size.scale(0.5)
vec1 = py_v3_to_vec(self.dims, v1, is_cylindrical)
vec2 = py_v3_to_vec(self.dims, v2, is_cylindrical)
self.swigobj = mp.volume(vec1, vec2)
def get_vertices(self):
xmin = self.center.x - self.size.x / 2
xmax = self.center.x + self.size.x / 2
ymin = self.center.y - self.size.y / 2
ymax = self.center.y + self.size.y / 2
zmin = self.center.z - self.size.z / 2
zmax = self.center.z + self.size.z / 2
# Iterate over and remove duplicates for collapsed dimensions (i.e. min=max))
return [
Vector3(x, y, z)
for x in list({xmin, xmax})
for y in list({ymin, ymax})
for z in list({zmin, zmax})
]
def get_edges(self):
vertices = self.get_vertices()
edges = []
# Useful for importing weird geometries and the sizes are slightly off
def nearly_equal(a, b, sig_fig=10):
return a == b or (abs(a - b) < 10 ** (-sig_fig))
for iter1 in range(len(vertices)):
for iter2 in range(iter1 + 1, len(vertices)):
if (
(iter1 != iter2)
and nearly_equal(
(vertices[iter1] - vertices[iter2]).norm(), self.size.x
)
or nearly_equal(
(vertices[iter1] - vertices[iter2]).norm(), self.size.y
)
or nearly_equal(
(vertices[iter1] - vertices[iter2]).norm(), self.size.z
)
):
edges.append([vertices[iter1], vertices[iter2]])
return edges
def pt_in_volume(self, pt: Vector3Type):
xmin = self.center.x - self.size.x / 2
xmax = self.center.x + self.size.x / 2
ymin = self.center.y - self.size.y / 2
ymax = self.center.y + self.size.y / 2
zmin = self.center.z - self.size.z / 2
zmax = self.center.z + self.size.z / 2
return (
pt.x >= xmin
and pt.x <= xmax
and pt.y >= ymin
and pt.y <= ymax
and pt.z >= zmin
and pt.z <= zmax
)
class FluxRegion:
"""
A `FluxRegion` object is used with [`add_flux`](#flux-spectra) to specify a region in
which Meep should accumulate the appropriate Fourier-transformed fields in order to
compute a flux spectrum. It represents a region (volume, plane, line, or point) in
which to compute the integral of the Poynting vector of the Fourier-transformed
fields. `ModeRegion` is an alias for `FluxRegion` for use with `add_mode_monitor`.
Note that the flux is always computed in the *positive* coordinate direction, although
this can effectively be flipped by using a `weight` of -1.0. This is useful, for
example, if you want to compute the outward flux through a box, so that the sides of
the box add instead of subtract.
"""
def __init__(
self,
center: Vector3Type = None,
size: Vector3Type = Vector3(),
direction: int = mp.AUTOMATIC,
weight: float = 1.0,
volume: Optional[Volume] = None,
):
"""
Construct a `FluxRegion` object.
+ **`center` [ `Vector3` ]** — The center of the flux region (no default).
+ **`size` [ `Vector3` ]** — The size of the flux region along each of the coordinate
axes. Default is `(0,0,0)`; a single point.
+ **`direction` [ `direction` constant ]** — The direction in which to compute the
flux (e.g. `mp.X`, `mp.Y`, etcetera). Default is `AUTOMATIC`, in which the
direction is determined by taking the normal direction if the flux region is a
plane (or a line, in 2d). If the normal direction is ambiguous (e.g. for a point
or volume), then you *must* specify the `direction` explicitly (not doing so
will lead to an error).
+ **`weight` [ `complex` ]** — A weight factor to multiply the flux by when it is
computed. Default is 1.0.
+ **`volume` [ `Volume` ]** — A `meep.Volume` can be used to specify the flux region
instead of a `center` and a `size`.
"""
if center is None and volume is None:
raise ValueError("Either center or volume required")
if volume:
self.center = volume.center
self.size = volume.size
else:
self.center = Vector3(*center)
self.size = Vector3(*size)
self.direction = direction
self.weight = complex(weight)
ModeRegion = FluxRegion
Near2FarRegion = FluxRegion
class ForceRegion(FluxRegion):
"""
A region (volume, plane, line, or point) in which to compute the integral of the
stress tensor of the Fourier-transformed fields. Its properties are:
+ **`center` [ `Vector3` ]** — The center of the force region (no default).
+ **`size` [ `Vector3` ]** — The size of the force region along each of the coordinate
axes. Default is `(0,0,0)` (a single point).
+ **`direction` [ `direction constant` ]** — The direction of the force that you wish
to compute (e.g. `X`, `Y`, etcetera). Unlike `FluxRegion`, you must specify this
explicitly, because there is not generally any relationship between the direction of
the force and the orientation of the force region.
+ **`weight` [ `complex` ]** — A weight factor to multiply the force by when it is
computed. Default is 1.0.
+ **`volume` [ `Volume` ]** — A `meep.Volume` can be used to specify the force region
instead of a `center` and a `size`.
In most circumstances, you should define a set of `ForceRegion`s whose union is a
closed surface lying in vacuum and enclosing the object that is experiencing the
force.
"""
class EnergyRegion(FluxRegion):
"""
A region (volume, plane, line, or point) in which to compute the integral of the
energy density of the Fourier-transformed fields. Its properties are:
+ **`center` [ `Vector3` ]** — The center of the energy region (no default).
+ **`size` [ `Vector3` ]** — The size of the energy region along each of the coordinate
axes. Default is (0,0,0): a single point.
+ **`weight` [ `complex` ]** — A weight factor to multiply the energy density by when it
is computed. Default is 1.0.
"""
class FieldsRegion:
def __init__(
self, where: Volume = None, center: Vector3Type = None, size: Vector3Type = None
):
if where:
self.center = where.center
self.size = where.size
else:
self.center = Vector3(*center) if center is not None else None
self.size = Vector3(*size) if size is not None else None
self.where = where
class DftObj:
"""Wrapper around DFT objects that allows delayed initialization of the structure.
When splitting the structure into chunks for parallel simulations, we want to know all
of the details of the simulation in order to ensure that each processor gets a similar
amount of work. The problem with DFTs is that the `add_flux` style methods immediately
initialize the structure and fields. So, if the user adds multiple DFT objects to the
simulation, the load balancing code only knows about the first one and can't split the
work up nicely. To circumvent this, we delay the execution of the `add_flux` methods
as late as possible. When `add_flux` (or `add_near2far`, etc.) is called, we:
1. Create an instance of the appropriate subclass of `DftObj` (`DftForce`, `DftFlux`, etc.).
Set its args property to the list of arguments passed to `add_flux`, and set its func
property to the 'real' `add_flux`, which is prefixed by an underscore.
2. Add this `DftObj` to the list Simulation.dft_objects. When we actually run the
simulation, we call `Simulation._evaluate_dft_objects`, which calls `dft.func(*args)`
for each dft in the list.
If the user tries to access a property or call a function on the `DftObj` before
`Simulation._evaluate_dft_objects` is called, then we initialize the C++ object through
swigobj_attr and return the property they requested.
"""
def __init__(self, func, args):
"""Construct a `DftObj`."""
self.func = func
self.args = args
self.swigobj = None
def swigobj_attr(self, attr):
if self.swigobj is None:
self.swigobj = self.func(*self.args)
return getattr(self.swigobj, attr)
@property
def save_hdf5(self):
return self.swigobj_attr("save_hdf5")
@property
def load_hdf5(self):
return self.swigobj_attr("load_hdf5")
@property
def scale_dfts(self):
return self.swigobj_attr("scale_dfts")
@property
def remove(self):
return self.swigobj_attr("remove")
@property
def freq(self):
return self.swigobj_attr("freq")
@property
def where(self):
return self.swigobj_attr("where")
class DftFlux(DftObj):
""" """
def __init__(self, func, args):
"""Construct a `DftFlux`."""
super().__init__(func, args)
self.nfreqs = len(args[0])
self.regions = args[1]
self.num_components = 4
@property
def flux(self):
return self.swigobj_attr("flux")
@property
def E(self):
return self.swigobj_attr("E")
@property
def H(self):
return self.swigobj_attr("H")
@property
def cE(self):
return self.swigobj_attr("cE")
@property
def cH(self):
return self.swigobj_attr("cH")
@property
def normal_direction(self):
return self.swigobj_attr("normal_direction")
@property
def freq(self):
return self.swigobj_attr("freq")
class DftForce(DftObj):
""" """
def __init__(self, func, args):
"""Construct a `DftForce`."""
super().__init__(func, args)
self.nfreqs = len(args[0])
self.regions = args[1]
self.num_components = 6
@property
def force(self):
return self.swigobj_attr("force")
@property
def offdiag1(self):
return self.swigobj_attr("offdiag1")
@property
def offdiag2(self):
return self.swigobj_attr("offdiag2")
@property
def diag(self):
return self.swigobj_attr("diag")
@property
def freq(self):
return self.swigobj_attr("freq")
class DftNear2Far(DftObj):
""" """
def __init__(self, func, args):
"""Construct a `DftNear2Far`."""
super().__init__(func, args)
self.nfreqs = len(args[0])
self.nperiods = args[1]
self.regions = args[2]
self.num_components = 4
@property
def farfield(self):
return self.swigobj_attr("farfield")
@property
def save_farfields(self):
return self.swigobj_attr("save_farfields")
@property
def F(self):
return self.swigobj_attr("F")
@property
def eps(self):
return self.swigobj_attr("eps")
@property
def mu(self):
return self.swigobj_attr("mu")
def flux(
self, direction: int = None, where: Volume = None, resolution: float = None
):
"""
Given a `Volume` `where` (may be 0d, 1d, 2d, or 3d) and a `resolution` (in grid
points / distance unit), compute the far fields in `where` (which may lie
*outside* the cell) in a grid with the given resolution (which may differ from the
FDTD solution) and return its Poynting flux in `direction` as a list. The dataset
is a 1d array of `nfreq` dimensions.
"""
return self.swigobj_attr("flux")(direction, where.swigobj, resolution)
@property
def freq(self):
return self.swigobj_attr("freq")
class DftEnergy(DftObj):
""" """
def __init__(self, func, args):
"""Construct a `DftEnergy`."""
super().__init__(func, args)
self.nfreqs = len(args[0])
self.regions = args[1]
self.num_components = 12
@property
def electric(self):
return self.swigobj_attr("electric")
@property
def magnetic(self):
return self.swigobj_attr("magnetic")
@property
def total(self):
return self.swigobj_attr("total")
@property
def freq(self):
return self.swigobj_attr("freq")
class DftFields(DftObj):
""" """
def __init__(self, func, args):
"""Construct a `DftFields`."""
super().__init__(func, args)
self.nfreqs = len(args[4])
self.regions = [FieldsRegion(where=args[1], center=args[2], size=args[3])]
self.num_components = len(args[0])
@property
def chunks(self):
return self.swigobj_attr("chunks")
Mode = namedtuple("Mode", ["freq", "decay", "Q", "amp", "err"])
class EigenmodeData:
def __init__(
self,
band_num,
freq: float,
group_velocity: float,
k: Vector3Type,
swigobj,
kdom: Vector3Type,
):
"""Construct an `EigenmodeData`."""
self.band_num = band_num
self.freq = freq
self.group_velocity = group_velocity
self.k = k
self.swigobj = swigobj
self.kdom = kdom
def amplitude(self, point, component):
swig_point = mp.vec(point.x, point.y, point.z)
return mp.eigenmode_amplitude(self.swigobj, swig_point, component)
class PadeDFT:
"""
Padé approximant based spectral extrapolation is implemented as a class with a [`__call__`](#PadeDFT.__call__) method,
which allows it to be used as a step function that collects field data from a given
point and runs [Padé](https://docs.scipy.org/doc/scipy/reference/generated/scipy.interpolate.pade.html)
on that data to extract an analytic rational function which approximates the frequency response.
For more information about the Padé approximant, see the [Wikipedia article](https://en.wikipedia.org/wiki/Padé_approximant).
See [`__init__`](#PadeDFT.__init__) for details about constructing a `PadeDFT`.
In particular, `PadeDFT` stores the discrete time series $\\hat{f}[n]$ corresponding to the given field
component as a function of time and expresses it as:
$$\\hat{f}(\\omega) = \\sum_n \\hat{f}[n] e^{i\\omega n \\Delta t}$$
The above is a "Taylor-like" polynomial in $n$ with a Fourier basis and
coefficients which are the sampled field data. We then compute the Padé approximant
to be the analytic form of this function as:
$$R(f) = R(\\omega / 2\\pi) = \\frac{P(f)}{Q(f)}$$
Where $P$ and $Q$ are polynomials of degree $m$ and $n$, and $m + n + 1$ is the
degree of agreement of the Padé approximant to the analytic function $f(\\omega / 2\\pi)$. This
function $R$ is stored in the callable method `pade_instance.dft`. Note that the computed polynomials
$P$ and $Q$ for each spatial point are stored as well in the instance variable `pade_instance.polys`,
as a spatial array of dicts: `[{"P": P(t), "Q": Q(t)}]` with no spectral extrapolation performed.
Be sure to save a reference to the `Pade` instance if you wish
to use the results after the simulation:
```py
sim = mp.Simulation(...)
p = mp.PadeDFT(...)
sim.run(p, until=time)
# do something with p.dft
```
"""
def __init__(
self,
c: int = None,
vol: Volume = None,
center: Vector3Type = None,
size: Vector3Type = None,
m: Optional[int] = None,
n: Optional[int] = None,
m_frac: float = 0.5,
n_frac: Optional[float] = None,
sampling_interval: int = 1,
start_time: int = 0,
stop_time: Optional[int] = None,
):
"""
Construct a Padé DFT object.
A `PadeDFT` is a step function that collects data from the field component `c`
(e.g. `meep.Ex`, etc.) at the given point `pt` (a `Vector3`). Then, at the end
of the run, it uses the scipy Padé algorithm to approximate the analytic
frequency response at the specified point.
+ **`c` [ `component` constant ]** — The field component to use for extrapolation.
No default.
+ **`vol` [ `Volume` ]** — The volume over which to accumulate the fields
(may be 0d, 1d, 2d, or 3d). No default.
+ **`center` [ `Vector3` class ]** — Alternative method for specifying volume, using a center point
+ **`size` [ `Vector3` class ]** — Alternative method for specifying volume, using a size vector
+ **`m` [ `int` ]** — The order of the numerator $P$. If not specified,
defaults to the length of aggregated field data times `m_frac`.
+ **`n` [ `int` ]** — The order of the denominator $Q$. Defaults
to length of field data - `m` - 1.
+ **`m_frac` [ `float` ]** — Method for specifying `m` as a fraction of
field samples to use as the order for numerator. Default is 0.5.
+ **`n_frac` [ `float` ]** — Fraction of field samples to use as order for
denominator. No default.
+ **`sampling_interval` [ `int` ]** — The interval at which to sample the field data.
Defaults to 1.
+ **`start_time` [ `int` ]** — The time (in increments of $\\Delta t$) at which
to start sampling the field data. Default is 0 (beginning of simulation).
+ **`stop_time` [ `int` ]** — The time (in increments of $\\Delta t$) at which
to stop sampling the field data. Default is `None` (end of simulation).
"""
self.c = c
self.vol = vol
self.center = center
self.size = size
self.m = m
self.n = n
self.m_frac = m_frac
self.n_frac = n_frac
self.sampling_interval = sampling_interval
self.start_time = start_time
self.stop_time = stop_time
self.data = []
self.data_dt = 0
self.dft: Callable = None
self.step_func = self._pade()
def __call__(self, sim, todo):
"""
Allows a PadeDFT instance to be used as a step function.
"""
self.step_func(sim, todo)
def _collect_pade(self, c, vol, center, size):
self.t0 = 0
def _collect(sim):
self.data_dt = sim.meep_time() - self.t0
self.t0 = sim.meep_time()
self.data.append(sim.get_array(c, vol, center, size))
return _collect
def _analyze_pade(self, sim):
# Ensure that the field data has dimension (# time steps x (volume dims))
self.data = np.atleast_2d(self.data)
# If the supplied volume is actually a point, np.atleast_2d will have
# shape (1, T), we want (T, 1)
if np.shape(self.data)[0] == 1:
self.data = self.data.T
# Sample the collected field data and possibly truncate the beginning or
# end to avoid transients
# TODO(Arjun-Khurana): Create a custom run function that collects field
# only at required points in time
samples = np.array(
self.data[self.start_time : self.stop_time : self.sampling_interval]
)