예제 #1
0
    def __call__(
        self, simspace: SimulationSpace, wlen: float, solver, **kwargs
    ) -> Union[fdfd_tools.VecField, Tuple[fdfd_tools.VecField,
                                          fdfd_tools.Vec3d]]:
        """Creates the plane wave source.

        Args:
            simspace: Simulation space to use for the source.
            wlen: Wavelength of source.

        Returns:
            If `overwrite_bloch_vector` is `True`, a tuple containing the source
            field and the Bloch vector corresponding to the plane wave source.
            Otherwise, only the source field is returned.
        """
        space_inst = simspace(wlen)

        # Calculate the border in gridpoints and igore the border if it's larger then the simulation.
        dx = simspace.dx
        border = [int(b // dx) for b in self._params.border]
        # The plane wave is assumed to be in the z direction so the border is 0 for z.
        border.append(0)

        source, kvector = fdfd_tools.free_space_sources.build_plane_wave_source(
            omega=2 * np.pi / wlen,
            eps_grid=space_inst.eps_bg,
            mu=None,
            axis=gridlock.axisvec2axis(self._params.normal),
            polarity=gridlock.axisvec2polarity(self._params.normal),
            slices=grid_utils.create_region_slices(simspace.edge_coords,
                                                   self._params.center,
                                                   self._params.extents),
            theta=self._params.theta,
            psi=self._params.psi,
            polarization_angle=self._params.polarization_angle,
            border=border,
            power=self._params.power)

        if self._params.normalize_by_sim:
            source = fdfd_tools.free_space_sources.normalize_source_by_sim(
                omega=2 * np.pi / wlen,
                source=source,
                eps=space_inst.eps_bg.grids,
                dxes=simspace.dxes,
                pml_layers=simspace.pml_layers,
                solver=solver,
                power=self._params.power)

        if self._params.overwrite_bloch_vector:
            # TODO(logansu): Figure out what's wrong with mixing PMLs and
            # Bloch vector. It seems to create problems.
            # For now, we manually set Bloch to zero is PML is nonzero.
            if simspace.pml_layers[0] != 0 or simspace.pml_layers[1] != 0:
                kvector[0] = 0
            if simspace.pml_layers[2] != 0 or simspace.pml_layers[3] != 0:
                kvector[1] = 0
            if simspace.pml_layers[4] != 0 or simspace.pml_layers[5] != 0:
                kvector[2] = 0
            return source, kvector
        return source
예제 #2
0
    def before_sim(self, sim: FdfdSimProp) -> None:
        beam_center = self._params.beam_center
        if beam_center is None:
            beam_center = self._params.center

        eps_grid = copy.deepcopy(sim.grid)
        eps_grid.grids = sim.eps

        source, _ = fdfd_tools.free_space_sources.build_gaussian_source(
            omega=2 * np.pi / sim.wlen,
            eps_grid=eps_grid,
            mu=None,
            axis=gridlock.axisvec2axis(self._params.normal),
            polarity=gridlock.axisvec2polarity(self._params.normal),
            slices=simspace.create_region_slices(sim.grid.exyz,
                                                 self._params.center,
                                                 self._params.extents),
            theta=self._params.theta,
            psi=self._params.psi,
            polarization_angle=self._params.polarization_angle,
            w0=self._params.w0,
            center=beam_center,
            power=self._params.power)

        if self._params.normalize_by_sim:
            source = fdfd_tools.free_space_sources.normalize_source_by_sim(
                omega=2 * np.pi / sim.wlen,
                source=source,
                eps=sim.eps,
                dxes=sim.dxes,
                pml_layers=sim.pml_layers,
                solver=sim.solver,
                power=self._params.power)

        sim.source += source
예제 #3
0
    def __call__(self, simspace: SimulationSpace, wlen: float,
                 **kwargs) -> fdfd_tools.VecField:
        """Creates the source vector.

        Args:
            simspace: Simulation space object to use.
            wlen: Wavelength to operate source.

        Returns:
            The vector field corresponding to the source.
        """
        space_inst = simspace(wlen)

        return fdfd_solvers.waveguide_mode.build_waveguide_source(
            omega=2 * np.pi / wlen,
            dxes=simspace.dxes,
            eps=space_inst.eps_bg.grids,
            mu=None,
            mode_num=self._params.mode_num,
            waveguide_slice=grid_utils.create_region_slices(
                simspace.edge_coords, self._params.center,
                self._params.extents),
            axis=gridlock.axisvec2axis(self._params.normal),
            polarity=gridlock.axisvec2polarity(self._params.normal),
            power=self._params.power)
예제 #4
0
파일: poynting.py 프로젝트: yckwong/optpoly
def create_power_transmission_function(
        params: PowerTransmission,
        context: workspace.Workspace) -> PowerTransmissionFunction:
    simspace = context.get_object(params.field.simulation_space)
    return PowerTransmissionFunction(
        field=context.get_object(params.field),
        simspace=simspace,
        wlen=params.field.wavelength,
        plane_slice=grid_utils.create_region_slices(simspace.edge_coords,
                                                    params.center,
                                                    params.extents),
        axis=gridlock.axisvec2axis(params.normal),
        polarity=gridlock.axisvec2polarity(params.normal))
예제 #5
0
    def before_sim(self, sim: FdfdSimProp) -> None:
        if self._wg_mode is None:
            self._wg_mode = fdfd_solvers.waveguide_mode.build_waveguide_source(
                omega=2 * np.pi / sim.wlen,
                dxes=sim.dxes,
                eps=sim.eps,
                mu=None,
                mode_num=self._src.mode_num,
                waveguide_slice=simspace.create_region_slices(
                    sim.grid.exyz, self._src.center, self._src.extents),
                axis=gridlock.axisvec2axis(self._src.normal),
                polarity=gridlock.axisvec2polarity(self._src.normal),
                power=self._src.power)

        sim.source += self._wg_mode
예제 #6
0
파일: creator_em.py 프로젝트: zizai/spins-b
 def __call__(self, simspace: SimulationSpace, wlen: float,
              **kwargs) -> fdfd_tools.VecField:
     space_inst = simspace(wlen)
     return fdfd_solvers.waveguide_mode.build_overlap(
         omega=2 * np.pi / wlen,
         dxes=simspace.dxes,
         eps=space_inst.eps_bg.grids,
         mu=None,
         mode_num=self._params.mode_num,
         waveguide_slice=grid_utils.create_region_slices(
             simspace.edge_coords, self._params.center,
             self._params.extents),
         axis=gridlock.axisvec2axis(self._params.normal),
         polarity=gridlock.axisvec2polarity(self._params.normal),
         power=self._params.power)
예제 #7
0
    def before_sim(self, sim: FdfdSimProp) -> None:
        # Calculate the eigenmode if we have not already.
        if self._wg_overlap is None:
            self._wg_overlap = fdfd_solvers.waveguide_mode.build_overlap(
                omega=2 * np.pi / sim.wlen,
                dxes=sim.dxes,
                eps=sim.eps,
                mu=None,
                mode_num=self._overlap.mode_num,
                waveguide_slice=simspace.create_region_slices(
                    sim.grid.exyz, self._overlap.center, self._overlap.extents),
                axis=gridlock.axisvec2axis(self._overlap.normal),
                polarity=gridlock.axisvec2polarity(self._overlap.normal),
                power=self._overlap.power)

            self._wg_overlap = np.stack(self._wg_overlap, axis=0)
예제 #8
0
    def eval(self, sim: mp.Simulation) -> goos.NumericFlow:

        # Retrieve the relevant tangential fields.
        fields = []
        for comp in self._wg_mode.field_comps:
            fields.append(
                np.reshape(sim.get_dft_array(self._mon, comp, 0),
                           self._wg_mode.xyzw[3].shape))

        norms = _get_overlap_integral(self._wg_mode.mode_fields, fields,
                                      self._wg_mode.xyzw)
        if gridlock.axisvec2polarity(self._overlap.normal) > 0:
            val = 0.5 * (norms[0] + norms[1]) / self._mode_norm
        else:
            val = 0.5 * (norms[0] - norms[1]) / self._mode_norm

        return goos.NumericFlow(val * self._amp_factor)
예제 #9
0
    def test_axisvec2axis_and_axisvec2polarity(self):

        # x-dir
        vec = np.array([1, 0, 0])
        ax = axisvec2axis(vec)
        pol = axisvec2polarity(vec)
        self.assertTrue(ax == 0)
        self.assertTrue(pol == 1)

        vec = np.array([2, 1e-8, 1e-8])
        ax = axisvec2axis(vec)
        pol = axisvec2polarity(vec)
        self.assertTrue(ax == 0)
        self.assertTrue(pol == 1)

        vec = np.array([-1, 0, 0])
        ax = axisvec2axis(vec)
        pol = axisvec2polarity(vec)
        self.assertTrue(ax == 0)
        self.assertTrue(pol == -1)

        # y-dir
        vec = np.array([0, 1, 0])
        ax = axisvec2axis(vec)
        pol = axisvec2polarity(vec)
        self.assertTrue(ax == 1)
        self.assertTrue(pol == 1)

        vec = np.array([0, -1, 0])
        ax = axisvec2axis(vec)
        pol = axisvec2polarity(vec)
        self.assertTrue(ax == 1)
        self.assertTrue(pol == -1)

        # z-dir
        vec = np.array([0, 0, 1])
        ax = axisvec2axis(vec)
        pol = axisvec2polarity(vec)
        self.assertTrue(ax == 2)
        self.assertTrue(pol == 1)

        vec = np.array([0, 0, -1])
        ax = axisvec2axis(vec)
        pol = axisvec2polarity(vec)
        self.assertTrue(ax == 2)
        self.assertTrue(pol == -1)
예제 #10
0
def test_plane_power_grad():
    space = Simspace(
        TESTDATA,
        optplan.SimulationSpace(
            pml_thickness=[0, 0, 0, 0, 0, 0],
            mesh=optplan.UniformMesh(dx=40),
            sim_region=optplan.Box3d(
                center=[0, 0, 0],
                extents=[80, 80, 80],
            ),
            eps_bg=optplan.GdsEps(
                gds="straight_waveguide.gds",
                mat_stack=optplan.GdsMaterialStack(
                    background=optplan.Material(mat_name="air"),
                    stack=[
                        optplan.GdsMaterialStackLayer(
                            gds_layer=[100, 0],
                            extents=[-80, 80],
                            foreground=optplan.Material(mat_name="Si"),
                            background=optplan.Material(mat_name="air"),
                        ),
                    ],
                ),
            ),
        ))

    wlen = 1550
    power_fun = poynting.PowerTransmissionFunction(
        field=problem.Variable(1),
        simspace=space,
        wlen=wlen,
        plane_slice=grid_utils.create_region_slices(space.edge_coords,
                                                    [0, 0, 0], [40, 80, 80]),
        axis=gridlock.axisvec2axis([1, 0, 0]),
        polarity=gridlock.axisvec2polarity([1, 0, 0]))

    field = np.arange(np.prod(space.dims) * 3).astype(np.complex128) * 1j

    grad_actual = power_fun.grad([field], 1)
    fun = lambda vec: power_fun.eval([vec])
    grad_brute = eval_grad_brute_wirt(field, fun)

    np.testing.assert_array_almost_equal(grad_actual[0], grad_brute, decimal=4)
예제 #11
0
파일: creator_em.py 프로젝트: zizai/spins-b
    def __call__(self, simspace: SimulationSpace, wlen: float,
                 solver: Callable, **kwargs) -> fdfd_tools.VecField:
        """Creates the source vector.

        Args:
            simspace: Simulation space.
            wlen: Wavelength of source.
            solver: If `normalize_by_source` is `True`, `solver` will be used
                to run an EM simulation to normalize the source power.

        Returns:
            The source.
        """
        space_inst = simspace(wlen)
        source, _ = fdfd_tools.free_space_sources.build_gaussian_source(
            omega=2 * np.pi / wlen,
            eps_grid=space_inst.eps_bg,
            mu=None,
            axis=gridlock.axisvec2axis(self._params.normal),
            polarity=gridlock.axisvec2polarity(self._params.normal),
            slices=grid_utils.create_region_slices(simspace.edge_coords,
                                                   self._params.center,
                                                   self._params.extents),
            theta=self._params.theta,
            psi=self._params.psi,
            polarization_angle=self._params.polarization_angle,
            w0=self._params.w0,
            center=self._params.beam_center,
            power=self._params.power)

        if self._params.normalize_by_sim:
            source = fdfd_tools.free_space_sources.normalize_source_by_sim(
                omega=2 * np.pi / wlen,
                source=source,
                eps=space_inst.eps_bg.grids,
                dxes=simspace.dxes,
                pml_layers=simspace.pml_layers,
                solver=solver,
                power=self._params.power)

        return source
예제 #12
0
    def __call__(
        self, simspace: SimulationSpace, wlen: float, **kwargs
    ) -> Union[fdfd_tools.VecField, Tuple[fdfd_tools.VecField,
                                          fdfd_tools.Vec3d]]:
        """Creates the plane wave source.

        Args:
            simspace: Simulation space to use for the source.
            wlen: Wavelength of source.

        Returns:
            If `overwrite_bloch_vector` is `True`, a tuple containing the source
            field and the Bloch vector corresponding to the plane wave source.
            Otherwise, only the source field is returned.
        """
        space_inst = simspace(wlen)

        # Calculate the border in gridpoints and igore the border if it's larger then the simulation.
        dx = simspace.dx
        border = [int(b // dx) for b in self._params.border]
        # The plane wave is assumed to be in the z direction so the border is 0 for z.
        border.append(0)

        source, kvector = fdfd_tools.free_space_sources.build_plane_wave_source(
            omega=2 * np.pi / wlen,
            eps_grid=space_inst.eps_bg,
            mu=None,
            axis=gridlock.axisvec2axis(self._params.normal),
            polarity=gridlock.axisvec2polarity(self._params.normal),
            slices=grid_utils.create_region_slices(simspace.edge_coords,
                                                   self._params.center,
                                                   self._params.extents),
            theta=self._params.theta,
            psi=self._params.psi,
            polarization_angle=self._params.polarization_angle,
            border=border,
            power=self._params.power)

        if self._params.overwrite_bloch_vector:
            return source, kvector
        return source
예제 #13
0
def test_stored_energy_grad():
    space = Simspace(
        TESTDATA,
        optplan.SimulationSpace(
            pml_thickness=[0, 0, 0, 0, 0, 0],
            mesh=optplan.UniformMesh(dx=40),
            sim_region=optplan.Box3d(
                center=[0, 0, 0],
                extents=[80, 80, 80],
            ),
            eps_bg=optplan.GdsEps(
                gds="straight_waveguide.gds",
                mat_stack=optplan.GdsMaterialStack(
                    background=optplan.Material(mat_name="air"),
                    stack=[
                        optplan.GdsMaterialStackLayer(
                            gds_layer=[100, 0],
                            extents=[-80, 80],
                            foreground=optplan.Material(mat_name="Si"),
                            background=optplan.Material(mat_name="air"),
                        ),
                    ],
                ),
            ),
        ))

    wlen = 1550
    energy_fun = stored_energy.StoredEnergyFunction(
        input_function=problem.Variable(1),
        simspace=space,
        center=[0,0,0],
        extents=[0,0,0],
        epsilon=space._eps_bg)

        plane_slice=grid_utils.create_region_slices(space.edge_coords,
                                                    [0, 0, 0], [40, 80, 80]),
        axis=gridlock.axisvec2axis([1, 0, 0]),
        polarity=gridlock.axisvec2polarity([1, 0, 0]))
예제 #14
0
def _create_waveguide_mode_source(
        src: Union[WaveguideModeSource, WaveguideModeOverlap],
        wg_mode: MeepWaveguideMode,
        src_time: mp.SourceTime,
        amp_factor: complex = 1,
        adjoint: bool = False,
) -> List[mp.Source]:
    """Creates a TFSF waveguide mode source.

    Roughly speaking, we seek to implement the following
    ```python
    mp.EigenModeSource(
        src=mp.GaussianSource(fcen, fwidth=0.2 * fcen),
        center=src.center,
        size=src.extents,
        eig_match_freq=True,
        eig_band=src.mode_num + 1,
    )
    ```
    The issue is that MEEP seems not to have a way to control the direction
    in which the TFSF source is pointing. To handle this, we instead "manually"
    construct the effective current source ourselves through the following
    steps:
    1. Compute the eigenmode of the waveguide using MEEP.
    2. Extract the field profiles of the mode for the relevant field components
       (i.e. the ones necessary to compute the TFSF source).
    3. Add a `mp.Source` for each component of the TFSF source.

    Args:
        src: Waveguide source parameters.
        wg_mode: Waveguide mode to inject.
        src_time: Specifies the temporal source profile of the source.
        amp_factor: The source amplitude is increased by factor of `amp_factor`.
        adjoint: If `True`, conjugate the modal source.

    Returns:
        A list of MEEP sources corresponding to the TFSF source.
        If `power_wlens` is not `None`, then a list of powers emitted by the
        source is returned as a second return value.
    """
    fcen = 1 / wg_mode.wlen

    # Determine the sign of amplitude of each relevant field component.
    polarity = gridlock.axisvec2polarity(src.normal)
    if adjoint:
        polarity = -polarity
    amp_signs = [1, -1, -polarity, polarity]

    meep_sources = []
    for mode_field, mode_comp, src_comp, amp_sign in zip(
            wg_mode.mode_fields, wg_mode.field_comps,
            reversed(wg_mode.field_comps), amp_signs):
        power_factor = np.sqrt(src.power) / src_time.fourier_transform(fcen)
        if adjoint:
            mode_field = np.conj(mode_field)

        meep_sources.append(
            mp.Source(src=src_time,
                      component=src_comp,
                      center=src.center,
                      size=src.extents,
                      amplitude=amp_sign * amp_factor * power_factor,
                      amp_data=mode_field))
    return meep_sources
예제 #15
0
def test_straight_waveguide_power_poynting():
    """Tests that total through straight waveguide is unity."""
    space = Simspace(
        TESTDATA,
        optplan.SimulationSpace(
            pml_thickness=[10, 10, 10, 10, 0, 0],
            mesh=optplan.UniformMesh(dx=40),
            sim_region=optplan.Box3d(
                center=[0, 0, 0],
                extents=[5000, 5000, 40],
            ),
            eps_bg=optplan.GdsEps(
                gds="straight_waveguide.gds",
                mat_stack=optplan.GdsMaterialStack(
                    background=optplan.Material(mat_name="air"),
                    stack=[
                        optplan.GdsMaterialStackLayer(
                            gds_layer=[100, 0],
                            extents=[-80, 80],
                            foreground=optplan.Material(mat_name="Si"),
                            background=optplan.Material(mat_name="air"),
                        ),
                    ],
                ),
            ),
        ))

    source = creator_em.WaveguideModeSource(
        optplan.WaveguideModeSource(
            power=1.0,
            extents=[40, 1500, 600],
            normal=[1.0, 0.0, 0.0],
            center=[-1770, 0, 0],
            mode_num=0,
        ))

    wlen = 1550
    eps_grid = space(wlen).eps_bg.grids
    source_grid = source(space, wlen)

    eps = problem.Constant(fdfd_tools.vec(eps_grid))
    sim = creator_em.FdfdSimulation(
        eps=eps,
        solver=local_matrix_solvers.DirectSolver(),
        wlen=wlen,
        source=fdfd_tools.vec(source_grid),
        simspace=space,
    )
    power_fun = poynting.PowerTransmissionFunction(
        field=sim,
        simspace=space,
        wlen=wlen,
        plane_slice=grid_utils.create_region_slices(
            space.edge_coords, [1770, 0, 0], [40, 1500, 600]),
        axis=gridlock.axisvec2axis([1, 0, 0]),
        polarity=gridlock.axisvec2polarity([1, 0, 0]))
    # Same as `power_fun` but with opposite normal vector. Should give same
    # answer but with a negative sign.
    power_fun_back = poynting.PowerTransmissionFunction(
        field=sim,
        simspace=space,
        wlen=wlen,
        plane_slice=grid_utils.create_region_slices(
            space.edge_coords, [1770, 0, 0], [40, 1500, 600]),
        axis=gridlock.axisvec2axis([-1, 0, 0]),
        polarity=gridlock.axisvec2polarity([-1, 0, 0]))

    efield_grid = fdfd_tools.unvec(
        graph_executor.eval_fun(sim, None), eps_grid[0].shape)

    # Calculate emitted power.
    edotj = np.real(
        fdfd_tools.vec(efield_grid) * np.conj(
            fdfd_tools.vec(source_grid))) * 40**3
    power = -0.5 * np.sum(edotj)

    np.testing.assert_almost_equal(
        graph_executor.eval_fun(power_fun, None), power, decimal=4)
    np.testing.assert_almost_equal(
        graph_executor.eval_fun(power_fun_back, None), -power, decimal=4)