File: dispersive_eigenmode.py

package info (click to toggle)
meep-mpi-default 1.17.1-2
  • links: PTS, VCS
  • area: main
  • in suites: bullseye
  • size: 51,672 kB
  • sloc: cpp: 29,881; python: 17,210; lisp: 1,225; makefile: 477; sh: 249; ansic: 133; javascript: 5
file content (158 lines) | stat: -rw-r--r-- 5,513 bytes parent folder | download | duplicates (3)
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

# dispersive_eigenmode.py - Tests the meep eigenmode features (eigenmode source,
# eigenmode decomposition, and get_eigenmode) with dispersive materials.
# TODO: 
#  * check materials with off diagonal components
#  * check magnetic profiles
#  * once imaginary component is supported, check that

from __future__ import division

import unittest
import meep as mp
import numpy as np
from meep import mpb
import h5py
import os

class TestDispersiveEigenmode(unittest.TestCase):
    # ----------------------------------------- #
    # ----------- Helper Functions ------------ #
    # ----------------------------------------- #
    # Directly cals the C++ chi1 routine
    def call_chi1(self,material,frequency):

        sim = mp.Simulation(cell_size=mp.Vector3(1,1,1),
                    default_material=material,
                    resolution=20)

        sim.init_sim()
        v3 = mp.py_v3_to_vec(sim.dimensions, mp.Vector3(0,0,0), sim.is_cylindrical)
        chi1inv = np.zeros((3,3),dtype=np.complex128)
        for i, com in enumerate([mp.Ex,mp.Ey,mp.Ez]):
            for k, dir in enumerate([mp.X,mp.Y,mp.Z]):
                chi1inv[i,k] = sim.structure.get_chi1inv(com,dir,v3,frequency)
        n = np.real(np.sqrt(np.linalg.inv(chi1inv.astype(np.complex128))))

        n_actual = np.real(np.sqrt(material.epsilon(frequency).astype(np.complex128)))
        
        np.testing.assert_allclose(n,n_actual)

    @classmethod
    def setUpClass(cls):
        cls.temp_dir = mp.make_output_directory()

    @classmethod
    def tearDownClass(cls):
        mp.delete_directory(cls.temp_dir)

    def verify_output_and_slice(self,material,frequency):
        # Since the slice routines average the diagonals, we need to do that too:
        chi1 = material.epsilon(frequency).astype(np.complex128)
        chi1inv = np.linalg.inv(chi1)
        chi1inv = np.diag(chi1inv)
        N = chi1inv.size
        n = np.sqrt(N/np.sum(chi1inv))
        
        sim = mp.Simulation(cell_size=mp.Vector3(2,2,2),
                            default_material=material,
                            resolution=20,
                            eps_averaging=False
                            )
        sim.use_output_directory(self.temp_dir)
        sim.init_sim()
        
        # Check to make sure the get_slice routine is working with frequency
        n_slice = np.sqrt(np.max(sim.get_epsilon(frequency=frequency,snap=True)))
        self.assertAlmostEqual(n,n_slice, places=4)

        # Check to make sure h5 output is working with frequency
        filename = os.path.join(self.temp_dir, 'dispersive_eigenmode-eps-000000.00.h5')
        mp.output_epsilon(sim,frequency=frequency)
        n_h5 = 0
        mp.all_wait()
        with h5py.File(filename, 'r') as f:
            n_h5 = np.sqrt(np.max(mp.complexarray(f['eps.r'][()],f['eps.i'][()])))
        self.assertAlmostEqual(n,n_h5, places=4)
    
    # ----------------------------------------- #
    # ----------- Test Routines --------------- #
    # ----------------------------------------- #
    def test_chi1_routine(self):
        # Checks the newly implemented get_chi1inv routines within the 
        # fields and structure classes by comparing their output to the 
        # python epsilon output.

        from meep.materials import Si, Ag, LiNbO3, Au
        
        # Check Silicon
        w0 = Si.valid_freq_range.min
        w1 = Si.valid_freq_range.max
        self.call_chi1(Si,w0)
        self.call_chi1(Si,w1)

        # Check Silver
        w0 = Ag.valid_freq_range.min
        w1 = Ag.valid_freq_range.max
        self.call_chi1(Ag,w0)
        self.call_chi1(Ag,w1)

        # Check Gold
        w0 = Au.valid_freq_range.min
        w1 = Au.valid_freq_range.max
        self.call_chi1(Au,w0)
        self.call_chi1(Au,w1)   

        # Check Lithium Niobate (X,X)
        w0 = LiNbO3.valid_freq_range.min
        w1 = LiNbO3.valid_freq_range.max
        self.call_chi1(LiNbO3,w0)
        self.call_chi1(LiNbO3,w1)

        # Now let's rotate LN
        import copy
        rotLiNbO3 = copy.deepcopy(LiNbO3)
        rotLiNbO3.rotate(mp.Vector3(1,1,1),np.radians(34))
        self.call_chi1(rotLiNbO3,w0)
        self.call_chi1(rotLiNbO3,w1)
    
    def test_get_with_dispersion(self):
        # Checks the get_array_slice and output_fields method
        # with dispersive materials.

        from meep.materials import Si, Ag, LiNbO3, Au
        
        # Check Silicon
        w0 = Si.valid_freq_range.min
        w1 = Si.valid_freq_range.max
        self.verify_output_and_slice(Si,w0)
        self.verify_output_and_slice(Si,w1)

        # Check Silver
        w0 = Ag.valid_freq_range.min
        w1 = Ag.valid_freq_range.max
        self.verify_output_and_slice(Ag,w0)
        self.verify_output_and_slice(Ag,w1)

        # Check Gold
        w0 = Au.valid_freq_range.min
        w1 = Au.valid_freq_range.max
        self.verify_output_and_slice(Au,w0)
        self.verify_output_and_slice(Au,w1)       

        # Check Lithium Niobate
        w0 = LiNbO3.valid_freq_range.min
        w1 = LiNbO3.valid_freq_range.max
        self.verify_output_and_slice(LiNbO3,w0)
        self.verify_output_and_slice(LiNbO3,w1)

        # Now let's rotate LN
        import copy
        rotLiNbO3 = copy.deepcopy(LiNbO3)
        rotLiNbO3.rotate(mp.Vector3(1,1,1),np.radians(34))
        self.verify_output_and_slice(rotLiNbO3,w0)
        self.verify_output_and_slice(rotLiNbO3,w1)
        

if __name__ == '__main__':
    unittest.main()