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 1001 1002 1003 1004 1005 1006 1007 1008 1009 1010 1011 1012 1013 1014 1015 1016 1017 1018 1019 1020 1021 1022 1023 1024 1025 1026 1027 1028 1029 1030 1031 1032 1033 1034 1035 1036 1037 1038 1039 1040 1041 1042 1043 1044 1045 1046 1047 1048 1049 1050 1051 1052 1053 1054 1055 1056 1057 1058 1059 1060 1061 1062 1063 1064 1065 1066 1067 1068 1069 1070 1071 1072 1073 1074 1075 1076 1077 1078 1079 1080 1081 1082 1083 1084 1085 1086 1087 1088 1089 1090 1091 1092 1093 1094 1095 1096 1097 1098 1099 1100 1101 1102 1103 1104 1105 1106 1107 1108 1109 1110 1111 1112 1113 1114 1115 1116 1117 1118 1119 1120 1121 1122 1123 1124 1125 1126 1127 1128 1129 1130 1131 1132 1133 1134 1135 1136 1137 1138 1139 1140 1141 1142 1143 1144 1145 1146 1147 1148 1149 1150 1151 1152 1153 1154 1155 1156 1157 1158 1159 1160 1161 1162 1163 1164 1165 1166 1167 1168 1169 1170 1171 1172 1173 1174 1175 1176 1177 1178 1179 1180 1181 1182 1183 1184 1185 1186 1187 1188 1189 1190 1191 1192 1193 1194 1195 1196 1197 1198 1199 1200 1201 1202 1203 1204 1205 1206 1207 1208 1209 1210 1211 1212 1213 1214 1215 1216 1217 1218 1219 1220 1221 1222 1223 1224 1225 1226 1227 1228 1229 1230 1231 1232 1233 1234 1235 1236 1237 1238 1239 1240 1241 1242 1243 1244 1245 1246 1247 1248 1249 1250 1251 1252 1253 1254 1255 1256 1257 1258 1259 1260 1261 1262 1263 1264 1265 1266 1267 1268 1269 1270 1271 1272 1273 1274 1275 1276 1277 1278 1279 1280 1281 1282 1283 1284 1285 1286 1287 1288 1289 1290 1291 1292 1293 1294 1295 1296 1297 1298 1299 1300 1301 1302 1303 1304 1305 1306 1307 1308 1309 1310 1311 1312 1313 1314 1315 1316 1317 1318 1319 1320 1321 1322 1323 1324 1325 1326 1327 1328 1329 1330 1331 1332 1333 1334 1335 1336 1337 1338 1339 1340 1341 1342 1343 1344 1345 1346 1347 1348 1349 1350 1351 1352 1353 1354 1355 1356 1357 1358 1359 1360 1361 1362 1363 1364 1365 1366 1367 1368 1369 1370 1371 1372 1373 1374 1375 1376 1377 1378 1379 1380 1381 1382 1383 1384 1385 1386 1387 1388 1389 1390 1391 1392 1393 1394 1395 1396 1397 1398 1399 1400 1401 1402 1403 1404 1405 1406 1407 1408 1409 1410 1411 1412 1413 1414 1415 1416 1417 1418 1419 1420 1421 1422 1423 1424 1425 1426 1427 1428 1429 1430 1431 1432 1433 1434 1435 1436 1437 1438 1439 1440 1441 1442 1443 1444 1445 1446 1447 1448 1449 1450 1451 1452 1453 1454 1455 1456 1457 1458 1459 1460 1461 1462 1463 1464 1465 1466 1467 1468 1469 1470 1471 1472 1473 1474 1475 1476 1477 1478 1479 1480 1481 1482 1483 1484 1485 1486 1487 1488 1489 1490 1491 1492 1493 1494 1495 1496 1497 1498 1499 1500 1501 1502 1503 1504 1505 1506 1507 1508 1509 1510 1511 1512 1513 1514 1515 1516 1517 1518 1519 1520 1521 1522 1523 1524 1525 1526 1527 1528 1529 1530 1531 1532 1533 1534 1535 1536 1537 1538 1539 1540 1541 1542 1543 1544 1545 1546 1547 1548 1549 1550 1551 1552 1553 1554 1555 1556 1557 1558 1559 1560 1561 1562 1563 1564 1565 1566 1567 1568 1569 1570 1571 1572 1573 1574 1575 1576 1577 1578 1579 1580 1581 1582 1583 1584 1585 1586 1587 1588 1589 1590 1591 1592 1593 1594 1595 1596 1597 1598 1599 1600 1601 1602 1603 1604 1605 1606 1607 1608 1609 1610 1611 1612 1613 1614 1615 1616 1617 1618 1619 1620 1621 1622 1623 1624 1625 1626 1627 1628 1629 1630 1631 1632 1633 1634 1635 1636 1637 1638 1639 1640 1641 1642 1643 1644 1645 1646 1647 1648 1649 1650 1651 1652 1653 1654 1655 1656 1657 1658 1659 1660 1661 1662 1663 1664 1665 1666 1667 1668 1669 1670 1671 1672 1673 1674 1675 1676 1677 1678 1679 1680 1681 1682 1683 1684 1685 1686 1687 1688 1689 1690 1691 1692 1693 1694 1695 1696 1697 1698 1699 1700 1701 1702 1703 1704 1705 1706 1707 1708 1709 1710 1711 1712 1713 1714 1715 1716 1717 1718 1719 1720 1721 1722 1723 1724 1725 1726 1727 1728 1729 1730 1731 1732 1733 1734 1735 1736 1737 1738 1739 1740 1741 1742 1743 1744 1745 1746 1747 1748 1749 1750 1751 1752 1753 1754 1755 1756 1757 1758 1759 1760 1761 1762 1763 1764 1765 1766 1767 1768 1769 1770 1771 1772 1773 1774 1775 1776 1777 1778 1779 1780 1781 1782 1783 1784 1785 1786 1787 1788 1789 1790 1791 1792 1793 1794 1795 1796 1797 1798 1799 1800 1801 1802 1803 1804 1805 1806 1807 1808 1809 1810 1811 1812 1813 1814 1815 1816 1817 1818 1819 1820 1821 1822 1823 1824 1825 1826 1827 1828 1829 1830 1831 1832 1833 1834 1835 1836 1837 1838 1839 1840 1841 1842 1843 1844 1845 1846 1847 1848 1849 1850 1851 1852 1853 1854 1855 1856 1857 1858 1859 1860 1861 1862 1863 1864 1865 1866 1867 1868 1869 1870 1871 1872 1873 1874 1875 1876 1877 1878 1879 1880 1881 1882 1883 1884 1885 1886 1887 1888 1889 1890 1891 1892 1893 1894 1895 1896 1897 1898 1899 1900 1901 1902 1903 1904 1905 1906 1907 1908 1909 1910 1911 1912 1913 1914 1915 1916 1917 1918 1919 1920 1921 1922 1923 1924 1925 1926 1927 1928 1929 1930 1931 1932 1933 1934 1935 1936 1937 1938 1939 1940 1941 1942 1943 1944 1945 1946 1947 1948 1949 1950 1951 1952 1953 1954 1955 1956 1957 1958 1959 1960 1961 1962 1963 1964 1965 1966 1967 1968 1969 1970 1971 1972 1973 1974 1975 1976 1977 1978 1979 1980 1981 1982 1983 1984 1985 1986 1987 1988 1989 1990 1991 1992 1993 1994 1995 1996 1997 1998 1999 2000 2001 2002 2003 2004 2005 2006 2007 2008 2009 2010 2011 2012 2013 2014 2015 2016 2017 2018 2019 2020 2021 2022 2023 2024 2025 2026 2027 2028 2029 2030 2031 2032 2033 2034 2035 2036 2037 2038 2039 2040 2041 2042 2043 2044 2045 2046 2047 2048 2049 2050 2051 2052 2053 2054 2055 2056 2057 2058 2059 2060 2061 2062 2063 2064 2065 2066 2067 2068 2069 2070 2071 2072 2073 2074 2075 2076 2077 2078 2079 2080 2081 2082 2083 2084 2085 2086 2087 2088 2089 2090 2091 2092 2093 2094 2095 2096 2097 2098 2099 2100 2101 2102 2103 2104 2105 2106 2107 2108 2109 2110 2111 2112 2113 2114 2115 2116 2117 2118 2119 2120 2121 2122 2123 2124 2125 2126 2127 2128 2129 2130 2131 2132 2133 2134 2135 2136 2137 2138 2139 2140 2141 2142 2143 2144 2145 2146 2147 2148 2149 2150 2151 2152 2153 2154 2155 2156 2157 2158 2159 2160 2161 2162 2163 2164 2165 2166 2167 2168 2169 2170 2171 2172 2173 2174 2175 2176 2177 2178 2179 2180 2181 2182 2183 2184 2185 2186 2187 2188 2189 2190 2191 2192 2193 2194 2195 2196 2197 2198 2199 2200 2201 2202 2203 2204 2205 2206 2207 2208 2209 2210 2211 2212 2213 2214 2215 2216 2217 2218 2219 2220 2221 2222 2223 2224 2225 2226 2227 2228 2229 2230 2231 2232 2233 2234 2235 2236 2237 2238 2239 2240 2241 2242 2243 2244 2245 2246 2247 2248 2249 2250 2251 2252 2253 2254 2255 2256 2257 2258 2259 2260 2261 2262 2263 2264 2265 2266 2267 2268 2269 2270 2271 2272 2273 2274 2275 2276 2277 2278 2279 2280 2281 2282 2283 2284 2285 2286 2287 2288 2289 2290 2291 2292 2293 2294 2295 2296 2297 2298 2299 2300 2301 2302 2303 2304 2305 2306 2307 2308 2309 2310 2311 2312 2313 2314 2315 2316 2317 2318 2319 2320 2321 2322 2323 2324 2325 2326 2327 2328 2329 2330 2331 2332 2333 2334 2335 2336 2337 2338 2339 2340 2341 2342 2343 2344 2345 2346 2347 2348 2349 2350 2351 2352 2353 2354 2355 2356 2357 2358 2359 2360 2361 2362 2363 2364 2365 2366 2367 2368 2369 2370 2371 2372 2373 2374 2375 2376 2377 2378 2379 2380 2381 2382 2383 2384 2385 2386 2387 2388 2389 2390 2391 2392 2393 2394 2395 2396 2397 2398 2399 2400 2401 2402 2403 2404 2405 2406 2407 2408 2409 2410 2411 2412 2413 2414 2415 2416 2417 2418 2419 2420 2421 2422 2423 2424 2425 2426 2427 2428 2429 2430 2431 2432 2433 2434 2435 2436 2437 2438 2439 2440 2441 2442 2443 2444 2445 2446 2447 2448 2449 2450 2451 2452 2453 2454 2455 2456 2457 2458 2459 2460 2461 2462 2463 2464 2465 2466 2467 2468 2469 2470 2471 2472 2473 2474 2475 2476 2477 2478 2479 2480
|
from __future__ import division, print_function
import functools
import math
import numbers
import os
import re
import subprocess
import sys
import warnings
from collections import namedtuple
from collections import OrderedDict
from collections import Sequence
import numpy as np
import meep as mp
from meep.geom import Vector3
from meep.source import EigenModeSource, check_positive
try:
basestring
except NameError:
basestring = str
EigCoeffsResult = namedtuple('EigCoeffsResult', ['alpha', 'vgrp', 'kpoints', 'kdom'])
FluxData = namedtuple('FluxData', ['E', 'H'])
ForceData = namedtuple('ForceData', ['offdiag1', 'offdiag2', 'diag'])
NearToFarData = namedtuple('NearToFarData', ['F'])
def get_num_args(func):
if isinstance(func, Harminv):
return 2
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, v3, is_cylindrical=False):
if dims == 1:
return mp.vec(v3.z)
elif dims == 2:
if is_cylindrical:
return mp.veccyl(v3.x, v3.z)
else:
return mp.vec(v3.x, v3.y)
elif dims == 3:
return mp.vec(v3.x, v3.y, v3.z)
else:
raise ValueError("Invalid dimensions in Volume: {}".format(dims))
class PML(object):
def __init__(self, thickness,
direction=mp.ALL,
side=mp.ALL,
R_asymptotic=1e-15,
mean_stretch=1.0,
pml_profile=lambda u: u * u):
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):
if val >= 1:
self._mean_stretch = val
else:
raise ValueError("PML.mean_stretch must be >= 1. Got {}".format(val))
class Absorber(PML):
pass
class Symmetry(object):
def __init__(self, direction, phase=1):
self.direction = direction
self.phase = complex(phase)
self.swigobj = None
class Rotate2(Symmetry):
pass
class Rotate4(Symmetry):
pass
class Mirror(Symmetry):
pass
class Identity(Symmetry):
pass
class Volume(object):
def __init__(self, center, size=Vector3(), dims=2, is_cylindrical=False):
self.center = center
self.size = size
self.dims = dims
v1 = center - size.scale(0.5)
v2 = center + 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)
class FluxRegion(object):
def __init__(self, center=None, size=Vector3(), direction=mp.AUTOMATIC, weight=1.0, volume=None):
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 = center
self.size = size
self.direction = direction
self.weight = complex(weight)
ModeRegion = FluxRegion
Near2FarRegion = FluxRegion
ForceRegion = FluxRegion
class FieldsRegion(object):
def __init__(self, where=None, center=None, size=None):
if where:
self.center = where.center
self.size = where.size
else:
self.center = center
self.size = size
self.where = where
class DftObj(object):
def __init__(self, func, args):
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_min(self):
return self.swigobj_attr('freq_min')
@property
def dfreq(self):
return self.swigobj_attr('dfreq')
@property
def Nfreq(self):
return self.swigobj_attr('Nfreq')
@property
def where(self):
return self.swigobj_attr('where')
class DftFlux(DftObj):
def __init__(self, func, args):
super(DftFlux, self).__init__(func, args)
self.nfreqs = args[2]
self.regions = args[3]
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')
class DftForce(DftObj):
def __init__(self, func, args):
super(DftForce, self).__init__(func, args)
self.nfreqs = args[2]
self.regions = args[3]
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')
class DftNear2Far(DftObj):
def __init__(self, func, args):
super(DftNear2Far, self).__init__(func, args)
self.nfreqs = args[2]
self.regions = args[3]
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')
class DftFields(DftObj):
def __init__(self, func, args):
super(DftFields, self).__init__(func, args)
self.nfreqs = args[6]
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(object):
def __init__(self, band_num, freq, group_velocity, k, swigobj, kdom):
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 Harminv(object):
def __init__(self, c, pt, fcen, df, mxbands=None):
self.c = c
self.pt = pt
self.fcen = fcen
self.df = df
self.mxbands = mxbands
self.data = []
self.data_dt = 0
self.modes = []
self.spectral_density = 1.1
self.Q_thresh = 50.0
self.rel_err_thresh = mp.inf
self.err_thresh = 0.01
self.rel_amp_thresh = -1.0
self.amp_thresh = -1.0
self.step_func = self._harminv()
def __call__(self, sim, todo):
self.step_func(sim, todo)
def _collect_harminv(self):
def _collect1(c, pt):
self.t0 = 0
def _collect2(sim):
self.data_dt = sim.meep_time() - self.t0
self.t0 = sim.meep_time()
self.data.append(sim.get_field_point(c, pt))
return _collect2
return _collect1
def _check_freqs(self, sim):
source_freqs = [(s.src.frequency, 0 if s.src.width == 0 else 1 / s.src.width)
for s in sim.sources
if hasattr(s.src, 'frequency')]
harminv_max = self.fcen + 0.5 * self.df
harminv_min = self.fcen - 0.5 * self.df
for sf in source_freqs:
sf_max = sf[0] + 0.5 * sf[1]
sf_min = sf[0] - 0.5 * sf[1]
if harminv_max > sf_max:
warn_fmt = "Harminv frequency {} is outside maximum Source frequency {}"
warnings.warn(warn_fmt.format(harminv_max, sf_max), RuntimeWarning)
if harminv_min < sf_min:
warn_fmt = "Harminv frequency {} is outside minimum Source frequency {}"
warnings.warn(warn_fmt.format(harminv_min, sf_min), RuntimeWarning)
def _analyze_harminv(self, sim, maxbands):
harminv_cols = ['frequency', 'imag. freq.', 'Q', '|amp|', 'amplitude', 'error']
display_run_data(sim, 'harminv', harminv_cols)
self._check_freqs(sim)
dt = self.data_dt if self.data_dt is not None else sim.fields.dt
bands = mp.py_do_harminv(self.data, dt, self.fcen - self.df / 2, self.fcen + self.df / 2, maxbands,
self.spectral_density, self.Q_thresh, self.rel_err_thresh, self.err_thresh,
self.rel_amp_thresh, self.amp_thresh)
modes = []
for freq, amp, err in bands:
Q = freq.real / (-2 * freq.imag) if freq.imag != 0 else float('inf')
modes.append(Mode(freq.real, freq.imag, Q, amp, err))
display_run_data(sim, 'harminv', [freq.real, freq.imag, Q, abs(amp), amp, err])
return modes
def _harminv(self):
def _harm(sim):
if self.mxbands is None or self.mxbands == 0:
mb = 100
else:
mb = self.mxbands
self.modes = self._analyze_harminv(sim, mb)
f1 = self._collect_harminv()
return _combine_step_funcs(at_end(_harm), f1(self.c, self.pt))
class Simulation(object):
def __init__(self,
cell_size,
resolution,
geometry=[],
sources=[],
eps_averaging=True,
dimensions=3,
boundary_layers=[],
symmetries=[],
verbose=False,
force_complex_fields=False,
default_material=mp.Medium(),
m=0,
k_point=False,
extra_materials=[],
material_function=None,
epsilon_func=None,
epsilon_input_file='',
progress_interval=4,
subpixel_tol=1e-4,
subpixel_maxeval=100000,
ensure_periodicity=True,
num_chunks=0,
Courant=0.5,
accurate_fields_near_cylorigin=False,
filename_prefix=None,
output_volume=None,
output_single_precision=False,
load_structure='',
geometry_center=mp.Vector3()):
self.cell_size = cell_size
self.geometry = geometry
self.sources = sources
self.resolution = resolution
self.dimensions = dimensions
self.boundary_layers = boundary_layers
self.symmetries = symmetries
self.geometry_center = geometry_center
self.eps_averaging = eps_averaging
self.subpixel_tol = subpixel_tol
self.subpixel_maxeval = subpixel_maxeval
self.ensure_periodicity = ensure_periodicity
self.extra_materials = extra_materials
self.default_material = default_material
self.epsilon_input_file = epsilon_input_file
self.num_chunks = num_chunks
self.Courant = Courant
self.global_d_conductivity = 0
self.global_b_conductivity = 0
self.special_kz = False
self.k_point = k_point
self.fields = None
self.structure = None
self.accurate_fields_near_cylorigin = accurate_fields_near_cylorigin
self.m = m
self.force_complex_fields = force_complex_fields
self.verbose = verbose
self.progress_interval = progress_interval
self.init_sim_hooks = []
self.run_index = 0
self.filename_prefix = filename_prefix
self.output_append_h5 = None
self.output_single_precision = output_single_precision
self.output_volume = output_volume
self.last_eps_filename = ''
self.output_h5_hook = lambda fname: False
self.interactive = False
self.is_cylindrical = False
self.material_function = material_function
self.epsilon_func = epsilon_func
self.load_structure_file = load_structure
self.dft_objects = []
self._is_initialized = False
self._fragment_size = 10
# To prevent the user from having to specify `dims` and `is_cylindrical`
# to Volumes they create, the library will adjust them appropriately based
# on the settings in the Simulation instance. This method must be called on
# any user-defined Volume before passing it to meep via its `swigobj`.
def _fit_volume_to_simulation(self, vol):
return Volume(vol.center, vol.size, dims=self.dimensions, is_cylindrical=self.is_cylindrical)
# Every function that takes a user volume can be specified either by a volume
# (a Python Volume or a SWIG-wrapped meep::volume), or a center and a size
def _volume_from_kwargs(self, vol=None, center=None, size=None):
if vol:
if isinstance(vol, Volume):
# A pure Python Volume
return self._fit_volume_to_simulation(vol).swigobj
else:
# A SWIG-wrapped meep::volume
return vol
elif size and center:
return Volume(center, size=size, dims=self.dimensions, is_cylindrical=self.is_cylindrical).swigobj
else:
raise ValueError("Need either a Volume, or a size and center")
def _infer_dimensions(self, k):
if self.dimensions == 3:
def use_2d(self, k):
zero_z = self.cell_size.z == 0
return zero_z and (not k or self.special_kz or k.z == 0)
if use_2d(self, k):
return 2
else:
return 3
return self.dimensions
def _get_valid_material_frequencies(self):
fmin = float('-inf')
fmax = float('inf')
for mat in [go.material for go in self.geometry] + self.extra_materials:
if isinstance(mat, mp.Medium) and mat.valid_freq_range:
if mat.valid_freq_range.min > fmin:
fmin = mat.valid_freq_range.min
if mat.valid_freq_range.max < fmax:
fmax = mat.valid_freq_range.max
return fmin, fmax
def _check_material_frequencies(self):
min_freq, max_freq = self._get_valid_material_frequencies()
source_freqs = [(s.src.frequency, 0 if s.src.width == 0 else 1 / s.src.width)
for s in self.sources
if hasattr(s.src, 'frequency')]
dft_freqs = []
for dftf in self.dft_objects:
dft_freqs.append(dftf.freq_min)
dft_freqs.append(dftf.freq_min + dftf.Nfreq * dftf.dfreq)
warn_src = ('Note: your sources include frequencies outside the range of validity of the ' +
'material models. This is fine as long as you eventually only look at outputs ' +
'(fluxes, resonant modes, etc.) at valid frequencies.')
warn_dft_fmt = "DFT frequency {} is out of material's range of {}-{}"
for sf in source_freqs:
if sf[0] + 0.5 * sf[1] > max_freq or sf[0] - 0.5 * sf[1] < min_freq:
warnings.warn(warn_src, RuntimeWarning)
for dftf in dft_freqs:
if dftf > max_freq or dftf < min_freq:
warnings.warn(warn_dft_fmt.format(dftf, min_freq, max_freq), RuntimeWarning)
def _create_grid_volume(self, k):
dims = self._infer_dimensions(k)
if dims == 0 or dims == 1:
gv = mp.vol1d(self.cell_size.z, self.resolution)
elif dims == 2:
self.dimensions = 2
gv = mp.vol2d(self.cell_size.x, self.cell_size.y, self.resolution)
elif dims == 3:
gv = mp.vol3d(self.cell_size.x, self.cell_size.y, self.cell_size.z, self.resolution)
elif dims == mp.CYLINDRICAL:
gv = mp.volcyl(self.cell_size.x, self.cell_size.z, self.resolution)
self.dimensions = 2
self.is_cylindrical = True
else:
raise ValueError("Unsupported dimentionality: {}".format(dims))
gv.center_origin()
gv.shift_origin(py_v3_to_vec(self.dimensions, self.geometry_center, self.is_cylindrical))
return gv
def _create_symmetries(self, gv):
sym = mp.symmetry()
# Initialize swig objects for each symmetry and combine them into one
for s in self.symmetries:
if isinstance(s, Identity):
s.swigobj = mp.identity()
elif isinstance(s, Rotate2):
s.swigobj = mp.rotate2(s.direction, gv)
sym += s.swigobj * complex(s.phase.real, s.phase.imag)
elif isinstance(s, Rotate4):
s.swigobj = mp.rotate4(s.direction, gv)
sym += s.swigobj * complex(s.phase.real, s.phase.imag)
elif isinstance(s, Mirror):
s.swigobj = mp.mirror(s.direction, gv)
sym += s.swigobj * complex(s.phase.real, s.phase.imag)
else:
s.swigobj = mp.symmetry()
return sym
def _get_dft_volumes(self):
volumes = [self._volume_from_kwargs(vol=r.where if hasattr(r, 'where') else None,
center=r.center, size=r.size)
for dft in self.dft_objects
for r in dft.regions]
return volumes
def _boundaries_to_vols_1d(self, boundaries):
v1 = []
for bl in boundaries:
cen = mp.Vector3(z=(self.cell_size.z / 2) - (0.5 * bl.thickness))
sz = mp.Vector3(z=bl.thickness)
if bl.side == mp.High or bl.side == mp.ALL:
v1.append(self._volume_from_kwargs(center=cen, size=sz))
if bl.side == mp.Low or bl.side == mp.ALL:
v1.append(self._volume_from_kwargs(center=-1 * cen, size=sz))
return v1
def _boundaries_to_vols_2d_3d(self, boundaries, cyl=False):
side_thickness = OrderedDict()
side_thickness['top'] = 0
side_thickness['bottom'] = 0
side_thickness['left'] = 0
side_thickness['right'] = 0
side_thickness['near'] = 0
side_thickness['far'] = 0
for bl in boundaries:
d = bl.direction
s = bl.side
if d == mp.X or d == mp.ALL:
if s == mp.High or s == mp.ALL:
side_thickness['right'] = bl.thickness
if s == mp.Low or s == mp.ALL:
side_thickness['left'] = bl.thickness
if d == mp.Y or d == mp.ALL:
if s == mp.High or s == mp.ALL:
side_thickness['top'] = bl.thickness
if s == mp.Low or s == mp.ALL:
side_thickness['bottom'] = bl.thickness
if self.dimensions == 3:
if d == mp.Z or d == mp.ALL:
if s == mp.High or s == mp.ALL:
side_thickness['far'] = bl.thickness
if s == mp.Low or s == mp.ALL:
side_thickness['near'] = bl.thickness
xmax = self.cell_size.x / 2
ymax = self.cell_size.z / 2 if cyl else self.cell_size.y / 2
zmax = self.cell_size.z / 2
ytot = self.cell_size.z if cyl else self.cell_size.y
def get_overlap_0(side, d):
if side == 'top' or side == 'bottom':
ydir = 1 if side == 'top' else -1
xsz = self.cell_size.x - (side_thickness['left'] + side_thickness['right'])
ysz = d
zsz = self.cell_size.z - (side_thickness['near'] + side_thickness['far'])
xcen = xmax - side_thickness['right'] - (xsz / 2)
ycen = ydir*ymax + (-ydir*0.5*d)
zcen = zmax - side_thickness['far'] - (zsz / 2)
elif side == 'left' or side == 'right':
xdir = 1 if side == 'right' else -1
xsz = d
ysz = ytot - (side_thickness['top'] + side_thickness['bottom'])
zsz = self.cell_size.z - (side_thickness['near'] + side_thickness['far'])
xcen = xdir*xmax + (-xdir*0.5*d)
ycen = ymax - side_thickness['top'] - (ysz / 2)
zcen = zmax - side_thickness['far'] - (zsz / 2)
elif side == 'near' or side == 'far':
zdir = 1 if side == 'far' else -1
xsz = self.cell_size.x - (side_thickness['left'] + side_thickness['right'])
ysz = ytot - (side_thickness['top'] + side_thickness['bottom'])
zsz = d
xcen = xmax - side_thickness['right'] - (xsz / 2)
ycen = ymax - side_thickness['top'] - (ysz / 2)
zcen = zdir*zmax + (-zdir*0.5*d)
if cyl:
cen = mp.Vector3(xcen, 0, ycen)
sz = mp.Vector3(xsz, 0, ysz)
else:
cen = mp.Vector3(xcen, ycen, zcen)
sz = mp.Vector3(xsz, ysz, zsz)
return self._volume_from_kwargs(center=cen, size=sz)
def get_overlap_1(side1, side2, d):
if side_thickness[side2] == 0:
return []
if side1 == 'top' or side1 == 'bottom':
ydir = 1 if side1 == 'top' else -1
ysz = d
ycen = ydir*ymax + (-ydir*0.5*d)
if side2 == 'left' or side2 == 'right':
xdir = 1 if side2 == 'right' else -1
xsz = side_thickness[side2]
zsz = self.cell_size.z - (side_thickness['near'] + side_thickness['far'])
xcen = xdir*xmax + (-xdir*0.5*side_thickness[side2])
zcen = zmax - side_thickness['far'] - (zsz / 2)
elif side2 == 'near' or side2 == 'far':
zdir = 1 if side2 == 'far' else -1
xsz = self.cell_size.x - (side_thickness['left'] + side_thickness['right'])
zsz = side_thickness[side2]
xcen = xmax - side_thickness['right'] - (xsz / 2)
zcen = zdir*zmax + (-zdir*0.5*side_thickness[side2])
elif side1 == 'near' or side1 == 'far':
xdir = 1 if side2 == 'right' else -1
zdir = 1 if side1 == 'far' else -1
xsz = side_thickness[side2]
ysz = self.cell_size.y - (side_thickness['top'] + side_thickness['bottom'])
zsz = d
xcen = xdir*xmax + (-xdir*0.5*side_thickness[side2])
ycen = ymax - side_thickness['top'] - (ysz / 2)
zcen = zdir*zmax + (-zdir*0.5*d)
if cyl:
cen = mp.Vector3(xcen, 0, ycen)
sz = mp.Vector3(xsz, 0, ysz)
else:
cen = mp.Vector3(xcen, ycen, zcen)
sz = mp.Vector3(xsz, ysz, zsz)
return self._volume_from_kwargs(center=cen, size=sz)
def get_overlap_2(side1, side2, side3, d):
if side_thickness[side2] == 0 or side_thickness[side3] == 0:
return []
xdir = 1 if side2 == 'right' else -1
ydir = 1 if side1 == 'top' else -1
zdir = 1 if side3 == 'far' else -1
xsz = side_thickness[side2]
ysz = d
zsz = side_thickness[side3]
xcen = xdir*xmax + (-xdir*0.5*xsz)
ycen = ydir*ymax + (-ydir*0.5*d)
zcen = zdir*zmax + (-zdir*0.5*zsz)
cen = mp.Vector3(xcen, ycen, zcen)
sz = mp.Vector3(xsz, ysz, zsz)
return self._volume_from_kwargs(center=cen, size=sz)
v1 = []
v2 = []
v3 = []
for side, thickness in side_thickness.items():
if thickness == 0:
continue
v1.append(get_overlap_0(side, thickness))
if side == 'top' or side == 'bottom':
v2.append(get_overlap_1(side, 'left', thickness))
v2.append(get_overlap_1(side, 'right', thickness))
if self.dimensions == 3:
v2.append(get_overlap_1(side, 'near', thickness))
v2.append(get_overlap_1(side, 'far', thickness))
v3.append(get_overlap_2(side, 'left', 'near', thickness))
v3.append(get_overlap_2(side, 'right', 'near', thickness))
v3.append(get_overlap_2(side, 'left', 'far', thickness))
v3.append(get_overlap_2(side, 'right', 'far', thickness))
if side == 'near' or side == 'far':
v2.append(get_overlap_1(side, 'left', thickness))
v2.append(get_overlap_1(side, 'right', thickness))
return [v for v in v1 if v], [v for v in v2 if v], [v for v in v3 if v]
def _boundary_layers_to_vol_list(self, boundaries):
"""Returns three lists of meep::volume objects. The first represents the boundary
regions with no overlaps. The second is regions where two boundaries overlap, and
the third is regions where three boundaries overlap
"""
vols1 = []
vols2 = []
vols3 = []
if self.dimensions == 1:
vols1 = self._boundaries_to_vols_1d(boundaries)
else:
vols1, vols2, vols3 = self._boundaries_to_vols_2d_3d(boundaries, self.is_cylindrical)
return vols1, vols2, vols3
def _compute_fragment_stats(self, gv):
def convert_volumes(dft_obj):
volumes = []
for r in dft_obj.regions:
volumes.append(self._volume_from_kwargs(vol=r.where if hasattr(r, 'where') else None,
center=r.center, size=r.size))
return volumes
dft_data_list = [mp.dft_data(o.nfreqs, o.num_components, convert_volumes(o))
for o in self.dft_objects]
pmls = []
absorbers = []
for bl in self.boundary_layers:
if type(bl) is PML:
pmls.append(bl)
elif type(bl) is Absorber:
absorbers.append(bl)
pml_vols1, pml_vols2, pml_vols3 = self._boundary_layers_to_vol_list(pmls)
absorber_vols1, absorber_vols2, absorber_vols3 = self._boundary_layers_to_vol_list(absorbers)
absorber_vols = absorber_vols1 + absorber_vols2 + absorber_vols3
stats = mp.compute_fragment_stats(
self.geometry,
gv,
self.cell_size,
self.geometry_center,
self.default_material,
dft_data_list,
pml_vols1,
pml_vols2,
pml_vols3,
absorber_vols,
self.subpixel_tol,
self.subpixel_maxeval,
self.ensure_periodicity,
self._fragment_size
)
mirror_symmetries = [sym for sym in self.symmetries if isinstance(sym, Mirror)]
for sym in mirror_symmetries:
for fs in stats:
fs.num_anisotropic_eps_pixels //= 2
fs.num_anisotropic_mu_pixels //= 2
fs.num_nonlinear_pixels //= 2
fs.num_susceptibility_pixels //= 2
fs.num_nonzero_conductivity_pixels //= 2
fs.num_1d_pml_pixels //= 2
fs.num_2d_pml_pixels //= 2
fs.num_3d_pml_pixels //= 2
fs.num_pixels_in_box //= 2
return stats
def _init_structure(self, k=False):
print('-' * 11)
print('Initializing structure...')
gv = self._create_grid_volume(k)
sym = self._create_symmetries(gv)
br = _create_boundary_region_from_boundary_layers(self.boundary_layers, gv)
absorbers = [bl for bl in self.boundary_layers if type(bl) is Absorber]
if self.material_function:
self.material_function.eps = False
self.default_material = self.material_function
elif self.epsilon_func:
self.epsilon_func.eps = True
self.default_material = self.epsilon_func
elif self.epsilon_input_file:
self.default_material = self.epsilon_input_file
self.fragment_stats = self._compute_fragment_stats(gv) if isinstance(self.default_material, mp.Medium) else []
self.structure = mp.structure(gv, None, br, sym, self.num_chunks, self.Courant,
self.eps_averaging, self.subpixel_tol, self.subpixel_maxeval)
self.structure.shared_chunks = True
mp.set_materials_from_geometry(self.structure,
self.geometry,
self.geometry_center,
self.eps_averaging,
self.subpixel_tol,
self.subpixel_maxeval,
self.ensure_periodicity and not not self.k_point,
False,
self.default_material,
absorbers,
self.extra_materials)
if self.load_structure_file:
self.load_structure(self.load_structure_file)
def set_materials(self, geometry=None, default_material=None):
if self.fields:
self.fields.remove_susceptibilities()
absorbers = [bl for bl in self.boundary_layers if type(bl) is Absorber]
mp.set_materials_from_geometry(
self.structure,
geometry if geometry is not None else self.geometry,
self.geometry_center,
self.eps_averaging,
self.subpixel_tol,
self.subpixel_maxeval,
self.ensure_periodicity,
False,
default_material if default_material else self.default_material,
absorbers,
self.extra_materials
)
def load_structure(self, fname):
if self.structure is None:
raise ValueError("Fields must be initialized before calling load_structure")
self.structure.load(fname)
def dump_structure(self, fname):
if self.structure is None:
raise ValueError("Fields must be initialized before calling dump_structure")
self.structure.dump(fname)
def init_sim(self):
if self._is_initialized:
return
materials = [g.material for g in self.geometry if isinstance(g.material, mp.Medium)]
if isinstance(self.default_material, mp.Medium):
materials.append(self.default_material)
for med in materials:
if ((med.epsilon_diag.x < 1 and med.epsilon_diag.x > -mp.inf) or
(med.epsilon_diag.y < 1 and med.epsilon_diag.y > -mp.inf) or
(med.epsilon_diag.z < 1 and med.epsilon_diag.z > -mp.inf)):
eps_warning = ("Epsilon < 1 may require adjusting the Courant parameter. " +
"See the 'Numerical Stability' entry under the 'Materials' " +
"section of the documentation")
warnings.warn(eps_warning, RuntimeWarning)
if self.structure is None:
self._init_structure(self.k_point)
self.fields = mp.fields(
self.structure,
self.m if self.is_cylindrical else 0,
self.k_point.z if self.special_kz and self.k_point else 0,
not self.accurate_fields_near_cylorigin
)
if self.verbose:
self.fields.verbose()
def use_real(self):
cond1 = self.is_cylindrical and self.m != 0
cond2 = any([s.phase.imag for s in self.symmetries])
cond3 = not self.k_point
cond4 = self.special_kz and self.k_point.x == 0 and self.k_point.y == 0
cond5 = not (cond3 or cond4 or self.k_point == Vector3())
return not (self.force_complex_fields or cond1 or cond2 or cond5)
if use_real(self):
self.fields.use_real_fields()
else:
print("Meep: using complex fields.")
if self.k_point:
v = Vector3(self.k_point.x, self.k_point.y) if self.special_kz else self.k_point
self.fields.use_bloch(py_v3_to_vec(self.dimensions, v, self.is_cylindrical))
for s in self.sources:
self.add_source(s)
for hook in self.init_sim_hooks:
hook()
self._is_initialized = True
def init_fields(self):
warnings.warn('init_fields is deprecated. Please use init_sim instead', DeprecationWarning)
self.init_sim()
def require_dimensions(self):
if self.structure is None:
mp.set_dimensions(self._infer_dimensions(self.k_point))
def meep_time(self):
if self.fields is None:
self.init_sim()
return self.fields.time()
def round_time(self):
if self.fields is None:
self.init_sim()
return self.fields.round_time()
def get_field_point(self, c, pt):
v3 = py_v3_to_vec(self.dimensions, pt, self.is_cylindrical)
return self.fields.get_field_from_comp(c, v3)
def get_epsilon_point(self, pt):
v3 = py_v3_to_vec(self.dimensions, pt, self.is_cylindrical)
return self.fields.get_eps(v3)
def get_filename_prefix(self):
if isinstance(self.filename_prefix, str):
return self.filename_prefix
elif self.filename_prefix is None:
_, filename = os.path.split(sys.argv[0])
if filename == 'ipykernel_launcher.py' or filename == '__main__.py':
return ''
else:
return re.sub(r'\.py$', '', filename)
else:
raise TypeError("Expected a string for filename_prefix, or None for the default.")
def use_output_directory(self, dname=''):
if not dname:
dname = self.get_filename_prefix() + '-out'
closure = {'trashed': False}
def hook():
print("Meep: using output directory '{}'".format(dname))
self.fields.set_output_directory(dname)
if not closure['trashed']:
mp.trash_output_directory(dname)
closure['trashed'] = True
self.init_sim_hooks.append(hook)
if self.fields is not None:
hook()
self.filename_prefix = None
return dname
def _run_until(self, cond, step_funcs):
self.interactive = False
if self.fields is None:
self.init_sim()
if isinstance(cond, numbers.Number):
stop_time = cond
t0 = self.round_time()
def stop_cond(sim):
return sim.round_time() >= t0 + stop_time
cond = stop_cond
step_funcs = list(step_funcs)
step_funcs.append(display_progress(t0, t0 + stop_time, self.progress_interval))
while not cond(self):
for func in step_funcs:
_eval_step_func(self, func, 'step')
self.fields.step()
# Translating the recursive scheme version of run-until into an iterative version
# (because python isn't tail-call-optimized) means we need one extra iteration to
# be the same as scheme.
for func in step_funcs:
_eval_step_func(self, func, 'step')
for func in step_funcs:
_eval_step_func(self, func, 'finish')
print("run {} finished at t = {} ({} timesteps)".format(self.run_index, self.meep_time(), self.fields.t))
self.run_index += 1
def _run_sources_until(self, cond, step_funcs):
if self.fields is None:
self.init_sim()
ts = self.fields.last_source_time()
if isinstance(cond, numbers.Number):
new_cond = (ts - self.round_time()) + cond
else:
def f(sim):
return cond(sim) and sim.round_time() >= ts
new_cond = f
self._run_until(new_cond, step_funcs)
def _run_sources(self, step_funcs):
self._run_sources_until(self, 0, step_funcs)
def _run_k_point(self, t, k):
components = [s.component for s in self.sources]
pts = [s.center for s in self.sources]
src_freqs_min = min([s.src.frequency - 1 / s.src.width / 2 if isinstance(s.src, mp.GaussianSource) else mp.inf
for s in self.sources])
fmin = max(0, src_freqs_min)
fmax = max([s.src.frequency + 1 / s.src.width / 2 if isinstance(s.src, mp.GaussianSource) else 0
for s in self.sources])
if not components or fmin > fmax:
raise ValueError("Running with k_points requires a 'GaussianSource' source")
self.change_k_point(k)
self.restart_fields()
h = Harminv(components[0], pts[0], 0.5 * (fmin + fmax), fmax - fmin)
self.run(after_sources(h), until_after_sources=t)
return [complex(m.freq, m.decay) for m in h.modes]
def run_k_points(self, t, k_points):
k_index = 0
all_freqs = []
for k in k_points:
k_index += 1
if k_index == 1:
self.init_sim()
output_epsilon(self)
freqs = self._run_k_point(t, k)
print("freqs:, {}, {}, {}, {}, ".format(k_index, k.x, k.y, k.z), end='')
print(', '.join([str(f.real) for f in freqs]))
print("freqs-im:, {}, {}, {}, {}, ".format(k_index, k.x, k.y, k.z), end='')
print(', '.join([str(f.imag) for f in freqs]))
all_freqs.append(freqs)
return all_freqs
def set_epsilon(self, eps):
if self.fields is None:
self.init_sim()
self.structure.set_epsilon(eps, self.eps_averaging, self.subpixel_tol, self.subpixel_maxeval)
def add_source(self, src):
if self.fields is None:
self.init_sim()
where = Volume(src.center, src.size, dims=self.dimensions,
is_cylindrical=self.is_cylindrical).swigobj
if isinstance(src, EigenModeSource):
if src.direction < 0:
direction = self.fields.normal_direction(where)
else:
direction = src.direction
eig_vol = Volume(src.eig_lattice_center, src.eig_lattice_size, self.dimensions,
is_cylindrical=self.is_cylindrical).swigobj
add_eig_src_args = [
src.component,
src.src.swigobj,
direction,
where,
eig_vol,
src.eig_band,
py_v3_to_vec(self.dimensions, src.eig_kpoint, is_cylindrical=self.is_cylindrical),
src.eig_match_freq,
src.eig_parity,
src.eig_resolution,
src.eig_tolerance,
src.amplitude
]
add_eig_src = functools.partial(self.fields.add_eigenmode_source, *add_eig_src_args)
if src.amp_func is None:
add_eig_src()
else:
add_eig_src(src.amp_func)
else:
add_vol_src_args = [src.component, src.src.swigobj, where]
add_vol_src = functools.partial(self.fields.add_volume_source, *add_vol_src_args)
if src.amp_func_file:
fname_dset = src.amp_func_file.rsplit(':', 1)
if len(fname_dset) != 2:
err_msg = "Expected a string of the form 'h5filename:dataset'. Got '{}'"
raise ValueError(err_msg.format(src.amp_func_file))
fname, dset = fname_dset
if not fname.endswith('.h5'):
fname += '.h5'
add_vol_src(fname, dset, src.amplitude * 1.0,)
elif src.amp_func:
add_vol_src(src.amp_func, src.amplitude * 1.0)
elif src.amp_data is not None:
add_vol_src(src.amp_data, src.amplitude * 1.0,)
else:
add_vol_src(src.amplitude * 1.0)
def _evaluate_dft_objects(self):
for dft in self.dft_objects:
if dft.swigobj is None:
dft.swigobj = dft.func(*dft.args)
def add_dft_fields(self, components, freq_min, freq_max, nfreq, where=None, center=None, size=None):
dftf = DftFields(self._add_dft_fields, [components, where, center, size, freq_min, freq_max, nfreq])
self.dft_objects.append(dftf)
return dftf
def _add_dft_fields(self, components, where, center, size, freq_min, freq_max, nfreq):
if self.fields is None:
self.init_sim()
try:
where = self._volume_from_kwargs(where, center, size)
except ValueError:
where = self.fields.total_volume()
return self.fields.add_dft_fields(components, where, freq_min, freq_max, nfreq)
def output_dft(self, dft_fields, fname):
if self.fields is None:
self.init_sim()
if hasattr(dft_fields, 'swigobj'):
dft_fields_swigobj = dft_fields.swigobj
else:
dft_fields_swigobj = dft_fields
self.fields.output_dft(dft_fields_swigobj, fname)
def get_dft_data(self, dft_chunk):
n = mp._get_dft_data_size(dft_chunk)
arr = np.zeros(n, np.complex128)
mp._get_dft_data(dft_chunk, arr)
return arr
def add_near2far(self, fcen, df, nfreq, *near2fars):
n2f = DftNear2Far(self._add_near2far, [fcen, df, nfreq, near2fars])
self.dft_objects.append(n2f)
return n2f
def _add_near2far(self, fcen, df, nfreq, near2fars):
if self.fields is None:
self.init_sim()
return self._add_fluxish_stuff(self.fields.add_dft_near2far, fcen, df, nfreq, near2fars)
def get_farfield(self, f, v):
return mp._get_farfield(f.swigobj, py_v3_to_vec(self.dimensions, v, is_cylindrical=self.is_cylindrical))
def output_farfields(self, near2far, fname, resolution, where=None, center=None, size=None):
vol = self._volume_from_kwargs(where, center, size)
near2far.save_farfields(fname, self.get_filename_prefix(), vol, resolution)
def load_near2far(self, fname, n2f):
if self.fields is None:
self.init_sim()
n2f.load_hdf5(self.fields, fname, '', self.get_filename_prefix())
def save_near2far(self, fname, n2f):
if self.fields is None:
self.init_sim()
n2f.save_hdf5(self.fields, fname, '', self.get_filename_prefix())
def load_minus_near2far(self, fname, n2f):
self.load_near2far(fname, n2f)
n2f.scale_dfts(-1.0)
def get_near2far_data(self, n2f):
return NearToFarData(F=self.get_dft_data(n2f.F))
def load_near2far_data(self, n2f, n2fdata):
mp._load_dft_data(n2f.F, n2fdata.F)
def load_minus_near2far_data(self, n2f, n2fdata):
self.load_near2far_data(n2f, n2fdata)
n2f.scale_dfts(complex(1.0))
def add_force(self, fcen, df, nfreq, *forces):
force = DftForce(self._add_force, [fcen, df, nfreq, forces])
self.dft_objects.append(force)
return force
def _add_force(self, fcen, df, nfreq, forces):
if self.fields is None:
self.init_sim()
return self._add_fluxish_stuff(self.fields.add_dft_force, fcen, df, nfreq, forces)
def display_forces(self, *forces):
force_freqs = get_force_freqs(forces[0])
display_csv(self, 'force', zip(force_freqs, *[get_forces(f) for f in forces]))
def load_force(self, fname, force):
if self.fields is None:
self.init_sim()
force.load_hdf5(self.fields, fname, '', self.get_filename_prefix())
def save_force(self, fname, force):
if self.fields is None:
self.init_sim()
force.save_hdf5(self.fields, fname, '', self.get_filename_prefix())
def load_minus_force(self, fname, force):
self.load_force(fname, force)
force.scale_dfts(-1.0)
def get_force_data(self, force):
return ForceData(offdiag1=self.get_dft_data(force.offdiag1),
offdiag2=self.get_dft_data(force.offdiag2),
diag=self.get_dft_data(force.diag))
def load_force_data(self, force, fdata):
mp._load_dft_data(force.offdiag1, fdata.offdiag1)
mp._load_dft_data(force.offdiag2, fdata.offdiag2)
mp._load_dft_data(force.diag, fdata.diag)
def load_minus_force_data(self, force, fdata):
self.load_force_data(force, fdata)
force.scale_dfts(complex(-1.0))
def add_flux(self, fcen, df, nfreq, *fluxes):
flux = DftFlux(self._add_flux, [fcen, df, nfreq, fluxes])
self.dft_objects.append(flux)
return flux
def _add_flux(self, fcen, df, nfreq, fluxes):
if self.fields is None:
self.init_sim()
return self._add_fluxish_stuff(self.fields.add_dft_flux, fcen, df, nfreq, fluxes)
def add_mode_monitor(self, fcen, df, nfreq, *fluxes):
flux = DftFlux(self._add_mode_monitor, [fcen, df, nfreq, fluxes])
self.dft_objects.append(flux)
return flux
def _add_mode_monitor(self, fcen, df, nfreq, fluxes):
if self.fields is None:
self.init_sim()
if len(fluxes) != 1:
raise ValueError("add_mode_monitor expected just one ModeRegion. Got {}".format(len(fluxes)))
region = fluxes[0]
v = mp.Volume(region.center, region.size, dims=self.dimensions, is_cylindrical=self.is_cylindrical)
d0 = region.direction
d = self.fields.normal_direction(v.swigobj) if d0 < 0 else d0
return self.fields.add_mode_monitor(d, v.swigobj, fcen - df / 2, fcen + df / 2, nfreq)
def add_eigenmode(self, fcen, df, nfreq, *fluxes):
warnings.warn('add_eigenmode is deprecated. Please use add_mode_monitor instead.', DeprecationWarning)
return self.add_mode_monitor(fcen, df, nfreq, *fluxes)
def display_fluxes(self, *fluxes):
display_csv(self, 'flux', zip(get_flux_freqs(fluxes[0]), *[get_fluxes(f) for f in fluxes]))
def load_flux(self, fname, flux):
if self.fields is None:
self.init_sim()
flux.load_hdf5(self.fields, fname, '', self.get_filename_prefix())
load_mode = load_flux
def save_flux(self, fname, flux):
if self.fields is None:
self.init_sim()
flux.save_hdf5(self.fields, fname, '', self.get_filename_prefix())
save_mode = save_flux
def load_minus_flux(self, fname, flux):
self.load_flux(fname, flux)
flux.scale_dfts(complex(-1.0))
load_minus_mode = load_minus_flux
def get_flux_data(self, flux):
return FluxData(E=self.get_dft_data(flux.E), H=self.get_dft_data(flux.H))
get_mode_data = get_flux_data
def load_flux_data(self, flux, fdata):
mp._load_dft_data(flux.E, fdata.E)
mp._load_dft_data(flux.H, fdata.H)
load_mode_data = load_flux_data
def load_minus_flux_data(self, flux, fdata):
self.load_flux_data(flux, fdata)
flux.scale_dfts(complex(-1.0))
load_minus_mode_data = load_minus_flux_data
def flux_in_box(self, d, box=None, center=None, size=None):
if self.fields is None:
raise RuntimeError('Fields must be initialized before using flux_in_box')
box = self._volume_from_kwargs(box, center, size)
return self.fields.flux_in_box(d, box)
def electric_energy_in_box(self, box=None, center=None, size=None):
if self.fields is None:
raise RuntimeError('Fields must be initialized before using electric_energy_in_box')
box = self._volume_from_kwargs(box, center, size)
return self.fields.electric_energy_in_box(box)
def magnetic_energy_in_box(self, box=None, center=None, size=None):
if self.fields is None:
raise RuntimeError('Fields must be initialized before using magnetic_energy_in_box')
box = self._volume_from_kwargs(box, center, size)
return self.fields.magnetic_energy_in_box(box)
def field_energy_in_box(self, box=None, center=None, size=None):
if self.fields is None:
raise RuntimeError('Fields must be initialized before using field_energy_in_box')
box = self._volume_from_kwargs(box, center, size)
return self.fields.field_energy_in_box(box)
def modal_volume_in_box(self, box=None, center=None, size=None):
if self.fields is None:
raise RuntimeError('Fields must be initialized before using modal_volume_in_box')
try:
box = self._volume_from_kwargs(box, center, size)
except ValueError:
box = self.fields.total_volume()
return self.fields.modal_volume_in_box(box)
def solve_cw(self, tol=1e-8, maxiters=10000, L=2):
if self.fields is None:
raise RuntimeError('Fields must be initialized before using solve_cw')
self._evaluate_dft_objects()
return self.fields.solve_cw(tol, maxiters, L)
def _add_fluxish_stuff(self, add_dft_stuff, fcen, df, nfreq, stufflist):
vol_list = None
for s in stufflist:
v = Volume(center=s.center, size=s.size, dims=self.dimensions,
is_cylindrical=self.is_cylindrical)
d0 = s.direction
d = self.fields.normal_direction(v.swigobj) if d0 < 0 else d0
c = mp.direction_component(mp.Sx, d)
v2 = Volume(center=s.center, size=s.size, dims=self.dimensions,
is_cylindrical=self.is_cylindrical).swigobj
vol_list = mp.make_volume_list(v2, c, s.weight, vol_list)
stuff = add_dft_stuff(vol_list, fcen - df / 2, fcen + df / 2, nfreq)
vol_list.__swig_destroy__(vol_list)
return stuff
def output_component(self, c, h5file=None):
if self.fields is None:
raise RuntimeError("Fields must be initialized before calling output_component")
vol = self.fields.total_volume() if self.output_volume is None else self.output_volume
h5 = self.output_append_h5 if h5file is None else h5file
append = h5file is None and self.output_append_h5 is not None
self.fields.output_hdf5(c, vol, h5, append, self.output_single_precision, self.get_filename_prefix())
if h5file is None:
nm = self.fields.h5file_name(mp.component_name(c), self.get_filename_prefix(), True)
if c == mp.Dielectric:
self.last_eps_filename = nm
self.output_h5_hook(nm)
def output_components(self, fname, *components):
if self.fields is None:
raise RuntimeError("Fields must be initialized before calling output_component")
if self.output_append_h5 is None:
f = self.fields.open_h5file(fname, mp.h5file.WRITE, self.get_filename_prefix(), True)
else:
f = None
for c in components:
self.output_component(c, h5file=f)
if self.output_append_h5 is None:
f.prevent_deadlock()
if self.output_append_h5 is None:
self.output_h5_hook(self.fields.h5file_name(fname, self.get_filename_prefix(), True))
def h5topng(self, rm_h5, option, *step_funcs):
opts = "h5topng {}".format(option)
cmd = re.sub(r'\$EPS', self.last_eps_filename, opts)
return convert_h5(rm_h5, cmd, *step_funcs)
def get_array(self, vol=None, center=None, size=None, component=mp.Ez, cmplx=None, arr=None):
dim_sizes = np.zeros(3, dtype=np.uintp)
if vol is None and center is None and size is None:
v = self.fields.total_volume()
else:
v = self._volume_from_kwargs(vol, center, size)
self.fields.get_array_slice_dimensions(v, dim_sizes)
dims = [s for s in dim_sizes if s != 0]
if cmplx is None:
cmplx = component < mp.Dielectric and not self.fields.is_real
if arr is not None:
if cmplx and not np.iscomplexobj(arr):
raise ValueError("Requested a complex slice, but provided array of type {}.".format(arr.dtype))
for a, b in zip(arr.shape, dims):
if a != b:
fmt = "Expected dimensions {}, but got {}"
raise ValueError(fmt.format(dims, arr.shape))
arr = np.require(arr, requirements=['C', 'W'])
else:
arr = np.zeros(dims, dtype=np.complex128 if cmplx else np.float64)
if np.iscomplexobj(arr):
self.fields.get_complex_array_slice(v, component, arr)
else:
self.fields.get_array_slice(v, component, arr)
return arr
def get_dft_array(self, dft_obj, component, num_freq):
if hasattr(dft_obj, 'swigobj'):
dft_swigobj = dft_obj.swigobj
else:
dft_swigobj = dft_obj
if type(dft_swigobj) is mp.dft_fields:
return mp.get_dft_fields_array(self.fields, dft_swigobj, component, num_freq)
elif type(dft_swigobj) is mp.dft_flux:
return mp.get_dft_flux_array(self.fields, dft_swigobj, component, num_freq)
elif type(dft_swigobj) is mp.dft_force:
return mp.get_dft_force_array(self.fields, dft_swigobj, component, num_freq)
elif type(dft_swigobj) is mp.dft_near2far:
return mp.get_dft_near2far_array(self.fields, dft_swigobj, component, num_freq)
else:
raise ValueError("Invalid type of dft object: {}".format(dft_swigobj))
def get_eigenmode_coefficients(self, flux, bands, eig_parity=mp.NO_PARITY, eig_vol=None,
eig_resolution=0, eig_tolerance=1e-12, kpoint_func=None, verbose=False):
if self.fields is None:
raise ValueError("Fields must be initialized before calling get_eigenmode_coefficients")
if eig_vol is None:
eig_vol = flux.where
else:
eig_vol = self._volume_from_kwargs(vol=eig_vol)
num_bands = len(bands)
coeffs = np.zeros(2 * num_bands * flux.Nfreq, dtype=np.complex128)
vgrp = np.zeros(num_bands * flux.Nfreq)
kpoints, kdom = mp.get_eigenmode_coefficients_and_kpoints(
self.fields,
flux.swigobj,
eig_vol,
np.array(bands, dtype=np.intc),
eig_parity,
eig_resolution,
eig_tolerance,
coeffs,
vgrp,
kpoint_func,
verbose
)
return EigCoeffsResult(np.reshape(coeffs, (num_bands, flux.Nfreq, 2)), vgrp, kpoints, kdom)
def get_eigenmode(self, freq, direction, where, band_num, kpoint, eig_vol=None, match_frequency=True,
parity=mp.NO_PARITY, resolution=0, eigensolver_tol=1e-12, verbose=False):
if self.fields is None:
raise ValueError("Fields must be initialized before calling get_eigenmode")
where = self._volume_from_kwargs(vol=where)
if eig_vol is None:
eig_vol = where
else:
eig_vol = self._volume_from_kwargs(vol=eig_vol)
swig_kpoint = mp.vec(kpoint.x, kpoint.y, kpoint.z)
kdom = np.zeros(3)
emdata = mp._get_eigenmode(self.fields, freq, direction, where, eig_vol, band_num, swig_kpoint,
match_frequency, parity, resolution, eigensolver_tol, verbose, kdom)
Gk = mp._get_eigenmode_Gk(emdata)
return EigenmodeData(emdata.band_num, emdata.omega, emdata.group_velocity, Gk,
emdata, mp.Vector3(kdom[0], kdom[1], kdom[2]))
def output_field_function(self, name, cs, func, real_only=False, h5file=None):
if self.fields is None:
raise RuntimeError("Fields must be initialized before calling output_field_function")
ov = self.output_volume if self.output_volume else self.fields.total_volume()
h5 = self.output_append_h5 if h5file is None else h5file
append = h5file is None and self.output_append_h5 is not None
self.fields.output_hdf5(name, [cs, func], ov, h5, append, self.output_single_precision,
self.get_filename_prefix(), real_only)
if h5file is None:
self.output_h5_hook(self.fields.h5file_name(name, self.get_filename_prefix(), True))
def _get_field_function_volume(self, where=None, center=None, size=None):
try:
where = self._volume_from_kwargs(where, center, size)
except ValueError:
where = self.fields.total_volume()
return where
def integrate_field_function(self, cs, func, where=None, center=None, size=None):
where = self._get_field_function_volume(where, center, size)
return self.fields.integrate([cs, func], where)
def integrate2_field_function(self, fields2, cs1, cs2, func, where=None, center=None, size=None):
where = self._get_field_function_volume(where, center, size)
return self.fields.integrate2(fields2, [cs1, cs2, func], where)
def max_abs_field_function(self, cs, func, where=None, center=None, size=None):
where = self._get_field_function_volume(where, center, size)
return self.fields.max_abs([cs, func], where)
def change_k_point(self, k):
self.k_point = k
if self.fields:
needs_complex_fields = not (not self.k_point or self.k_point == mp.Vector3())
if needs_complex_fields and self.fields.is_real:
self.fields = None
self._is_initialized = False
self.init_sim()
else:
if self.k_point:
self.fields.use_bloch(py_v3_to_vec(self.dimensions, self.k_point, self.is_cylindrical))
def change_sources(self, new_sources):
self.sources = new_sources
if self.fields:
self.fields.remove_sources()
for s in self.sources:
self.add_source(s)
def reset_meep(self):
self.fields = None
self.structure = None
self.dft_objects = []
self._is_initialized = False
def restart_fields(self):
if self.fields is not None:
self.fields.t = 0
self.fields.zero_fields()
else:
self._is_initialized = False
self.init_sim()
def run(self, *step_funcs, **kwargs):
until = kwargs.pop('until', None)
until_after_sources = kwargs.pop('until_after_sources', None)
if self.fields is None:
self.init_sim()
self._evaluate_dft_objects()
self._check_material_frequencies()
if kwargs:
raise ValueError("Unrecognized keyword arguments: {}".format(kwargs.keys()))
if until_after_sources is not None:
self._run_sources_until(until_after_sources, step_funcs)
elif until is not None:
self._run_until(until, step_funcs)
else:
raise ValueError("Invalid run configuration")
def get_epsilon(self):
return self.get_array(component=mp.Dielectric)
def get_mu(self):
return self.get_array(component=mp.Permeability)
def get_hpwr(self):
return self.get_array(component=mp.H_EnergyDensity)
def get_dpwr(self):
return self.get_array(component=mp.D_EnergyDensity)
def get_tot_pwr(self):
return self.get_array(component=mp.EnergyDensity)
def get_hfield(self):
if self.is_cylindrical:
r = self.get_array(cmplx=True, component=mp.Hr)
p = self.get_array(cmplx=True, component=mp.Hp)
return np.stack([r, p], axis=-1)
else:
x = self.get_array(cmplx=True, component=mp.Hx)
y = self.get_array(cmplx=True, component=mp.Hy)
z = self.get_array(cmplx=True, component=mp.Hz)
return np.stack([x, y, z], axis=-1)
def get_hfield_x(self):
return self.get_array(cmplx=True, component=mp.Hx)
def get_hfield_y(self):
return self.get_array(cmplx=True, component=mp.Hy)
def get_hfield_z(self):
return self.get_array(cmplx=True, component=mp.Hz)
def get_hfield_r(self):
return self.get_array(cmplx=True, component=mp.Hr)
def get_hfield_p(self):
return self.get_array(cmplx=True, component=mp.Hp)
def get_bfield(self):
if self.is_cylindrical:
r = self.get_array(cmplx=True, component=mp.Br)
p = self.get_array(cmplx=True, component=mp.Bp)
return np.stack([r, p], axis=-1)
else:
x = self.get_array(cmplx=True, component=mp.Bx)
y = self.get_array(cmplx=True, component=mp.By)
z = self.get_array(cmplx=True, component=mp.Bz)
return np.stack([x, y, z], axis=-1)
def get_bfield_x(self):
return self.get_array(cmplx=True, component=mp.Bx)
def get_bfield_y(self):
return self.get_array(cmplx=True, component=mp.By)
def get_bfield_z(self):
return self.get_array(cmplx=True, component=mp.Bz)
def get_bfield_r(self):
return self.get_array(cmplx=True, component=mp.Br)
def get_bfield_p(self):
return self.get_array(cmplx=True, component=mp.Bp)
def get_efield(self):
if self.is_cylindrical:
r = self.get_array(cmplx=True, component=mp.Er)
p = self.get_array(cmplx=True, component=mp.Ep)
return np.stack([r, p], axis=-1)
else:
x = self.get_array(cmplx=True, component=mp.Ex)
y = self.get_array(cmplx=True, component=mp.Ey)
z = self.get_array(cmplx=True, component=mp.Ez)
return np.stack([x, y, z], axis=-1)
def get_efield_x(self):
return self.get_array(cmplx=True, component=mp.Ex)
def get_efield_y(self):
return self.get_array(cmplx=True, component=mp.Ey)
def get_efield_z(self):
return self.get_array(cmplx=True, component=mp.Ez)
def get_efield_r(self):
return self.get_array(cmplx=True, component=mp.Er)
def get_efield_p(self):
return self.get_array(cmplx=True, component=mp.Ep)
def get_dfield(self):
if self.is_cylindrical:
r = self.get_array(cmplx=True, component=mp.Dr)
p = self.get_array(cmplx=True, component=mp.Dp)
return np.stack([r, p], axis=-1)
else:
x = self.get_array(cmplx=True, component=mp.Dx)
y = self.get_array(cmplx=True, component=mp.Dy)
z = self.get_array(cmplx=True, component=mp.Dz)
return np.stack([x, y, z], axis=-1)
def get_dfield_x(self):
return self.get_array(cmplx=True, component=mp.Dx)
def get_dfield_y(self):
return self.get_array(cmplx=True, component=mp.Dy)
def get_dfield_z(self):
return self.get_array(cmplx=True, component=mp.Dz)
def get_dfield_r(self):
return self.get_array(cmplx=True, component=mp.Dr)
def get_dfield_p(self):
return self.get_array(cmplx=True, component=mp.Dp)
def get_sfield(self):
if self.is_cylindrical:
r = self.get_array(cmplx=True, component=mp.Sr)
p = self.get_array(cmplx=True, component=mp.Sp)
return np.stack([r, p], axis=-1)
else:
x = self.get_array(cmplx=True, component=mp.Sx)
y = self.get_array(cmplx=True, component=mp.Sy)
z = self.get_array(cmplx=True, component=mp.Sz)
return np.stack([x, y, z], axis=-1)
def get_sfield_x(self):
return self.get_array(cmplx=True, component=mp.Sx)
def get_sfield_y(self):
return self.get_array(cmplx=True, component=mp.Sy)
def get_sfield_z(self):
return self.get_array(cmplx=True, component=mp.Sz)
def get_sfield_r(self):
return self.get_array(cmplx=True, component=mp.Sr)
def get_sfield_p(self):
return self.get_array(cmplx=True, component=mp.Sp)
def _create_boundary_region_from_boundary_layers(boundary_layers, gv):
br = mp.boundary_region()
for layer in boundary_layers:
if isinstance(layer, Absorber):
continue
boundary_region_args = [
mp.boundary_region.PML,
layer.thickness,
layer.R_asymptotic,
layer.mean_stretch,
mp.py_pml_profile,
layer.pml_profile,
1 / 3,
1 / 4,
]
if layer.direction == mp.ALL:
d = mp.start_at_direction(gv.dim)
loop_stop_directi = mp.stop_at_direction(gv.dim)
while d < loop_stop_directi:
if layer.side == mp.ALL:
b = mp.High
loop_stop_bi = mp.Low
while b != loop_stop_bi:
br += mp.boundary_region(*(boundary_region_args + [d, b]))
b = (b + 1) % 2
loop_stop_bi = mp.High
else:
br += mp.boundary_region(*(boundary_region_args + [d, layer.side]))
d += 1
else:
if layer.side == mp.ALL:
b = mp.High
loop_stop_bi = mp.Low
while b != loop_stop_bi:
br += mp.boundary_region(*(boundary_region_args + [layer.direction, b]))
b = (b + 1) % 2
loop_stop_bi = mp.High
else:
br += mp.boundary_region(*(boundary_region_args + [layer.direction, layer.side]))
return br
# Private step functions
def _combine_step_funcs(*step_funcs):
def _combine(sim, todo):
for func in step_funcs:
_eval_step_func(sim, func, todo)
return _combine
def _eval_step_func(sim, func, todo):
num_args = get_num_args(func)
if num_args != 1 and num_args != 2:
raise ValueError("Step function '{}'' requires 1 or 2 arguments".format(func.__name__))
elif num_args == 1:
if todo == 'step':
func(sim)
elif num_args == 2:
func(sim, todo)
def _when_true_funcs(cond, *step_funcs):
def _true(sim, todo):
if todo == 'finish' or cond(sim):
for f in step_funcs:
_eval_step_func(sim, f, todo)
return _true
# Public step functions
def after_sources(*step_funcs):
def _after_sources(sim, todo):
time = sim.fields.last_source_time()
if sim.round_time() >= time:
for func in step_funcs:
_eval_step_func(sim, func, todo)
return _after_sources
def after_sources_and_time(t, *step_funcs):
def _after_s_and_t(sim, todo):
time = sim.fields.last_source_time() + t - sim.round_time()
if sim.round_time() >= time:
for func in step_funcs:
_eval_step_func(sim, func, todo)
return _after_s_and_t
def after_time(t, *step_funcs):
def _after_t(sim):
return sim.round_time() >= t
return _when_true_funcs(_after_t, *step_funcs)
def at_beginning(*step_funcs):
closure = {'done': False}
def _beg(sim, todo):
if not closure['done']:
for f in step_funcs:
_eval_step_func(sim, f, todo)
closure['done'] = True
return _beg
def at_end(*step_funcs):
def _end(sim, todo):
if todo == 'finish':
for func in step_funcs:
_eval_step_func(sim, func, 'step')
for func in step_funcs:
_eval_step_func(sim, func, 'finish')
return _end
def at_every(dt, *step_funcs):
closure = {'tlast': 0.0}
def _every(sim, todo):
t = sim.round_time()
if todo == 'finish' or t >= closure['tlast'] + dt + (-0.5 * sim.fields.dt):
for func in step_funcs:
_eval_step_func(sim, func, todo)
closure['tlast'] = t
return _every
def at_time(t, *step_funcs):
closure = {'done': False}
def _at_time(sim, todo):
if not closure['done'] or todo == 'finish':
for f in step_funcs:
_eval_step_func(sim, f, todo)
closure['done'] = closure['done'] or todo == 'step'
return after_time(t, _at_time)
def before_time(t, *step_funcs):
def _before_t(sim):
return sim.round_time() < t
return _when_true_funcs(_before_t, *step_funcs)
def during_sources(*step_funcs):
closure = {'finished': False}
def _during_sources(sim, todo):
time = sim.fields.last_source_time()
if sim.round_time() < time:
for func in step_funcs:
_eval_step_func(sim, func, 'step')
elif closure['finished'] is False:
for func in step_funcs:
_eval_step_func(sim, func, 'finish')
closure['finished'] = True
return _during_sources
def in_volume(v, *step_funcs):
closure = {'cur_eps': ''}
def _in_volume(sim, todo):
v_save = sim.output_volume
eps_save = sim.last_eps_filename
sim.output_volume = sim._fit_volume_to_simulation(v).swigobj
if closure['cur_eps']:
sim.last_eps_filename = closure['cur_eps']
for func in step_funcs:
_eval_step_func(sim, func, todo)
closure['cur_eps'] = sim.last_eps_filename
sim.output_volume = v_save
if eps_save:
sim.last_eps_filename = eps_save
return _in_volume
def in_point(pt, *step_funcs):
v = Volume(pt)
return in_volume(v, *step_funcs)
def to_appended(fname, *step_funcs):
closure = {'h5': None}
def _to_appended(sim, todo):
if closure['h5'] is None:
closure['h5'] = sim.fields.open_h5file(fname, mp.h5file.WRITE, sim.get_filename_prefix())
h5save = sim.output_append_h5
sim.output_append_h5 = closure['h5']
for func in step_funcs:
_eval_step_func(sim, func, todo)
if todo == 'finish':
closure['h5'] = None
sim.output_h5_hook(sim.fields.h5file_name(fname, sim.get_filename_prefix()))
sim.output_append_h5 = h5save
return _to_appended
def stop_when_fields_decayed(dt, c, pt, decay_by):
closure = {
'max_abs': 0,
'cur_max': 0,
't0': 0,
}
def _stop(sim):
fabs = abs(sim.get_field_point(c, pt)) * abs(sim.get_field_point(c, pt))
closure['cur_max'] = max(closure['cur_max'], fabs)
if sim.round_time() <= dt + closure['t0']:
return False
else:
old_cur = closure['cur_max']
closure['cur_max'] = 0
closure['t0'] = sim.round_time()
closure['max_abs'] = max(closure['max_abs'], old_cur)
if closure['max_abs'] != 0:
fmt = "field decay(t = {}): {} / {} = {}"
print(fmt.format(sim.meep_time(), old_cur, closure['max_abs'], old_cur / closure['max_abs']))
return old_cur <= closure['max_abs'] * decay_by
return _stop
def synchronized_magnetic(*step_funcs):
def _sync(sim, todo):
sim.fields.synchronize_magnetic_fields()
for f in step_funcs:
_eval_step_func(sim, f, todo)
sim.fields.restore_magnetic_fields()
return _sync
def when_true(cond, *step_funcs):
return _when_true_funcs(cond, *step_funcs)
def when_false(cond, *step_funcs):
return _when_true_funcs(lambda: not cond, *step_funcs)
def with_prefix(pre, *step_funcs):
def _with_prefix(sim, todo):
saved_pre = sim.filename_prefix
sim.filename_prefix = pre + sim.get_filename_prefix()
for f in step_funcs:
_eval_step_func(sim, f, todo)
sim.filename_prefix = saved_pre
return _with_prefix
def display_csv(sim, name, data):
for d in data:
display_run_data(sim, name, d)
def display_progress(t0, t, dt):
t_0 = mp.wall_time()
closure = {'tlast': mp.wall_time()}
def _disp(sim):
t1 = mp.wall_time()
if t1 - closure['tlast'] >= dt:
msg_fmt = "Meep progress: {}/{} = {:.1f}% done in {:.1f}s, {:.1f}s to go"
val1 = sim.meep_time() - t0
val2 = val1 / (0.01 * t)
val3 = t1 - t_0
val4 = (val3 * (t / val1) - val3) if val1 != 0 else 0
print(msg_fmt.format(val1, t, val2, val3, val4))
closure['tlast'] = t1
return _disp
def data_to_str(d):
if type(d) is complex:
sign = '+' if d.imag >= 0 else ''
return "{}{}{}i".format(d.real, sign, d.imag)
else:
return str(d)
def display_run_data(sim, data_name, data):
if isinstance(data, Sequence):
data_str = [data_to_str(f) for f in data]
else:
data_str = [data_to_str(data)]
print("{}{}:, {}".format(data_name, sim.run_index, ', '.join(data_str)))
def convert_h5(rm_h5, convert_cmd, *step_funcs):
def convert(fname):
if mp.my_rank() == 0:
cmd = convert_cmd.split()
cmd.append(fname)
ret = subprocess.call(cmd)
if ret == 0 and rm_h5:
os.remove(fname)
def _convert_h5(sim, todo):
hooksave = sim.output_h5_hook
sim.output_h5_hook = convert
for f in step_funcs:
_eval_step_func(sim, f, todo)
sim.output_h5_hook = hooksave
return _convert_h5
def output_png(compnt, options, rm_h5=True):
closure = {'maxabs': 0.0}
def _output_png(sim, todo):
if todo == 'step':
if sim.output_volume is None:
ov = sim.fields.total_volume()
else:
ov = sim.output_volume
closure['maxabs'] = max(closure['maxabs'],
sim.fields.max_abs(compnt, ov))
convert = sim.h5topng(rm_h5, "-M {} {}".format(closure['maxabs'], options),
lambda sim: sim.output_component(compnt))
convert(sim, todo)
return _output_png
def output_epsilon(sim):
sim.output_component(mp.Dielectric)
def output_mu(sim):
sim.output_component(mp.Permeability)
def output_hpwr(sim):
sim.output_component(mp.H_EnergyDensity)
def output_dpwr(sim):
sim.output_component(mp.D_EnergyDensity)
def output_tot_pwr(sim):
sim.output_component(mp.EnergyDensity)
def output_hfield(sim):
sim.output_components('h', mp.Hx, mp.Hy, mp.Hz, mp.Hr, mp.Hp)
def output_hfield_x(sim):
sim.output_component(mp.Hx)
def output_hfield_y(sim):
sim.output_component(mp.Hy)
def output_hfield_z(sim):
sim.output_component(mp.Hz)
def output_hfield_r(sim):
sim.output_component(mp.Hr)
def output_hfield_p(sim):
sim.output_component(mp.Hp)
def output_bfield(sim):
sim.output_components('b', mp.Bx, mp.By, mp.Bz, mp.Br, mp.Bp)
def output_bfield_x(sim):
sim.output_component(mp.Bx)
def output_bfield_y(sim):
sim.output_component(mp.By)
def output_bfield_z(sim):
sim.output_component(mp.Bz)
def output_bfield_r(sim):
sim.output_component(mp.Br)
def output_bfield_p(sim):
sim.output_component(mp.Bp)
def output_efield(sim):
sim.output_components('e', mp.Ex, mp.Ey, mp.Ez, mp.Er, mp.Ep)
def output_efield_x(sim):
sim.output_component(mp.Ex)
def output_efield_y(sim):
sim.output_component(mp.Ey)
def output_efield_z(sim):
sim.output_component(mp.Ez)
def output_efield_r(sim):
sim.output_component(mp.Er)
def output_efield_p(sim):
sim.output_component(mp.Ep)
def output_dfield(sim):
sim.output_components('d', mp.Dx, mp.Dy, mp.Dz, mp.Dr, mp.Dp)
def output_dfield_x(sim):
sim.output_component(mp.Dx)
def output_dfield_y(sim):
sim.output_component(mp.Dy)
def output_dfield_z(sim):
sim.output_component(mp.Dz)
def output_dfield_r(sim):
sim.output_component(mp.Dr)
def output_dfield_p(sim):
sim.output_component(mp.Dp)
# MPB compatibility
def output_poynting(sim):
sim.output_components('s', mp.Sx, mp.Sy, mp.Sz, mp.Sr, mp.Sp)
def output_poynting_x(sim):
sim.output_component(mp.Sx)
def output_poynting_y(sim):
sim.output_component(mp.Sy)
def output_poynting_z(sim):
sim.output_component(mp.Sz)
def output_poynting_r(sim):
sim.output_component(mp.Sr)
def output_poynting_p(sim):
sim.output_component(mp.Sp)
def output_sfield(sim):
sim.output_components('s', mp.Sx, mp.Sy, mp.Sz, mp.Sr, mp.Sp)
def output_sfield_x(sim):
sim.output_component(mp.Sx)
def output_sfield_y(sim):
sim.output_component(mp.Sy)
def output_sfield_z(sim):
sim.output_component(mp.Sz)
def output_sfield_r(sim):
sim.output_component(mp.Sr)
def output_sfield_p(sim):
sim.output_component(mp.Sp)
def Ldos(fcen, df, nfreq):
return mp._dft_ldos(fcen - df / 2, fcen + df / 2, nfreq)
def dft_ldos(fcen=None, df=None, nfreq=None, ldos=None):
if ldos is None:
if fcen is None or df is None or nfreq is None:
raise ValueError("Either fcen, df, and nfreq, or an Ldos is required for dft_ldos")
ldos = mp._dft_ldos(fcen - df / 2, fcen + df / 2, nfreq)
def _ldos(sim, todo):
if todo == 'step':
ldos.update(sim.fields)
else:
sim.ldos_data = mp._dft_ldos_ldos(ldos)
sim.ldos_Fdata = mp._dft_ldos_F(ldos)
sim.ldos_Jdata = mp._dft_ldos_J(ldos)
display_csv(sim, 'ldos', zip(ldos.freqs(), sim.ldos_data))
return _ldos
def scale_flux_fields(s, flux):
flux.scale_dfts(s)
def get_flux_freqs(f):
return np.linspace(f.freq_min, f.freq_min + f.dfreq * f.Nfreq, num=f.Nfreq, endpoint=False).tolist()
def get_fluxes(f):
return f.flux()
def scale_force_fields(s, force):
force.scale_dfts(s)
def get_eigenmode_freqs(f):
return np.linspace(f.freq_min, f.freq_min + f.dfreq * f.Nfreq, num=f.Nfreq, endpoint=False).tolist()
def get_force_freqs(f):
return np.linspace(f.freq_min, f.freq_min + f.dfreq * f.Nfreq, num=f.Nfreq, endpoint=False).tolist()
def get_forces(f):
return f.force()
def scale_near2far_fields(s, n2f):
n2f.scale_dfts(s)
def get_near2far_freqs(f):
return np.linspace(f.freq_min, f.freq_min + f.dfreq * f.Nfreq, num=f.Nfreq, endpoint=False).tolist()
def interpolate(n, nums):
res = []
if isinstance(nums[0], mp.Vector3):
for low, high in zip(nums, nums[1:]):
x = np.linspace(low.x, high.x, n + 1, endpoint=False).tolist()
y = np.linspace(low.y, high.y, n + 1, endpoint=False).tolist()
z = np.linspace(low.z, high.z, n + 1, endpoint=False).tolist()
for i in range(len(x)):
res.append(mp.Vector3(x[i], y[i], z[i]))
else:
for low, high in zip(nums, nums[1:]):
res.extend(np.linspace(low, high, n + 1, endpoint=False).tolist())
return res + [nums[-1]]
# extract center and size of a meep::volume
def get_center_and_size(v):
rmin = v.get_min_corner()
rmax = v.get_max_corner()
v3rmin = mp.Vector3(rmin.x(), rmin.y(), rmin.z())
v3rmax = mp.Vector3(rmax.x(), rmax.y(), rmax.z())
if v.dim == mp.D2:
v3rmin.z = 0
v3rmax.z = 0
elif v.dim == mp.D1:
v3rmin.x = 0
v3rmin.y = 0
v3rmin.y = 0
v3rmax.y = 0
center = 0.5 * (v3rmin + v3rmax)
size = v3rmax - v3rmin
return center, size
def GDSII_vol(fname, layer, zmin, zmax):
meep_vol = mp.get_GDSII_volume(fname, layer, zmin, zmax)
dims = meep_vol.dim + 1
is_cyl = False
if dims == 4:
# cylindrical
dims = 2
is_cyl = True
center, size = get_center_and_size(meep_vol)
return Volume(center, size, dims, is_cyl)
|