File: test_dft_fields.py

package info (click to toggle)
meep-openmpi 1.25.0-2
  • links: PTS, VCS
  • area: main
  • in suites: bookworm
  • size: 64,556 kB
  • sloc: cpp: 32,214; python: 27,958; lisp: 1,225; makefile: 505; sh: 249; ansic: 131; javascript: 5
file content (132 lines) | stat: -rw-r--r-- 4,737 bytes parent folder | download | duplicates (5)
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
import os
import unittest

import h5py
import numpy as np
from utils import ApproxComparisonTestCase

import meep as mp


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

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

    def init(self):
        resolution = 10
        n = 3.4
        w = 1.0
        r = 1.0
        pad = 4
        self.dpml = 2
        self.sxy = 2.0 * (r + w + pad + self.dpml)
        cell = mp.Vector3(self.sxy, self.sxy)
        pml_layers = [mp.PML(self.dpml)]

        geometry = [
            mp.Cylinder(r + w, material=mp.Medium(epsilon=n**2)),
            mp.Cylinder(r, material=mp.vacuum),
        ]

        self.fcen = 0.118
        self.df = 0.1

        src = mp.GaussianSource(self.fcen, fwidth=self.df)
        sources = [mp.Source(src=src, component=mp.Ez, center=mp.Vector3(r + 0.1))]

        return mp.Simulation(
            cell_size=cell,
            resolution=resolution,
            geometry=geometry,
            sources=sources,
            boundary_layers=pml_layers,
        )

    def test_use_centered_grid(self):
        sim = self.init()
        sim.init_sim()
        dft_fields = sim.add_dft_fields([mp.Ez], self.fcen, 0, 1, yee_grid=True)
        sim.run(until=100)

    def test_get_dft_array(self):
        sim = self.init()
        sim.init_sim()
        dft_fields = sim.add_dft_fields([mp.Ez], self.fcen, 0, 1)
        fr = mp.FluxRegion(
            mp.Vector3(), size=mp.Vector3(self.sxy, self.sxy), direction=mp.X
        )
        dft_flux = sim.add_flux(self.fcen, 0, 1, fr)

        # volumes with zero thickness in x and y directions to test collapsing
        # of empty dimensions in DFT array and HDF5 output routines
        thin_x_volume = mp.Volume(
            center=mp.Vector3(0.35 * self.sxy), size=mp.Vector3(y=0.8 * self.sxy)
        )
        thin_x_flux = sim.add_dft_fields([mp.Ez], self.fcen, 0, 1, where=thin_x_volume)
        thin_y_volume = mp.Volume(
            center=mp.Vector3(y=0.25 * self.sxy), size=mp.Vector3(x=self.sxy)
        )
        thin_y_flux = sim.add_flux(self.fcen, 0, 1, mp.FluxRegion(volume=thin_y_volume))

        sim.run(until_after_sources=100)

        # test proper collapsing of degenerate dimensions in HDF5 files and arrays
        thin_x_array = sim.get_dft_array(thin_x_flux, mp.Ez, 0)
        thin_y_array = sim.get_dft_array(thin_y_flux, mp.Ez, 0)
        np.testing.assert_equal(thin_x_array.ndim, 1)
        np.testing.assert_equal(thin_y_array.ndim, 1)

        sim.output_dft(thin_x_flux, os.path.join(self.temp_dir, "thin-x-flux"))
        sim.output_dft(thin_y_flux, os.path.join(self.temp_dir, "thin-y-flux"))

        with h5py.File(os.path.join(self.temp_dir, "thin-x-flux.h5"), "r") as thin_x:
            thin_x_h5 = mp.complexarray(thin_x["ez_0.r"][()], thin_x["ez_0.i"][()])

        with h5py.File(os.path.join(self.temp_dir, "thin-y-flux.h5"), "r") as thin_y:
            thin_y_h5 = mp.complexarray(thin_y["ez_0.r"][()], thin_y["ez_0.i"][()])

        tol = 1e-6
        self.assertClose(thin_x_array, thin_x_h5, epsilon=tol)
        self.assertClose(thin_y_array, thin_y_h5, epsilon=tol)

        # compare array data to HDF5 file content for fields and flux
        fields_arr = sim.get_dft_array(dft_fields, mp.Ez, 0)
        flux_arr = sim.get_dft_array(dft_flux, mp.Ez, 0)

        sim.output_dft(dft_fields, os.path.join(self.temp_dir, "dft-fields"))
        sim.output_dft(dft_flux, os.path.join(self.temp_dir, "dft-flux"))

        with h5py.File(
            os.path.join(self.temp_dir, "dft-fields.h5"), "r"
        ) as fields, h5py.File(os.path.join(self.temp_dir, "dft-flux.h5"), "r") as flux:
            exp_fields = mp.complexarray(fields["ez_0.r"][()], fields["ez_0.i"][()])
            exp_flux = mp.complexarray(flux["ez_0.r"][()], flux["ez_0.i"][()])

        tol = 1e-6
        self.assertClose(exp_fields, fields_arr, epsilon=tol)
        self.assertClose(exp_flux, flux_arr, epsilon=tol)

    def test_decimated_dft_fields_are_almost_equal_to_undecimated_fields(self):
        sim = self.init()
        sim.init_sim()
        undecimated_field = sim.add_dft_fields(
            [mp.Ez], self.fcen, 0, 1, decimation_factor=1
        )
        decimated_field = sim.add_dft_fields(
            [mp.Ez], self.fcen, 0, 1, decimation_factor=4
        )

        sim.run(until_after_sources=100)

        expected_dft = sim.get_dft_array(undecimated_field, mp.Ez, 0)
        actual_dft = sim.get_dft_array(decimated_field, mp.Ez, 0)
        self.assertClose(expected_dft, actual_dft, epsilon=1e-3)


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