File: 3rd-harm-1d.py

package info (click to toggle)
meep-mpi-default 1.29.0-2
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid, trixie
  • size: 79,148 kB
  • sloc: cpp: 32,541; python: 31,061; lisp: 1,225; makefile: 516; sh: 249; ansic: 131; javascript: 5
file content (146 lines) | stat: -rw-r--r-- 4,931 bytes parent folder | download | duplicates (2)
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
"""
1d simulation of a plane wave propagating through a Kerr medium
and generating the third-harmonic frequency component.

ref: https://meep.readthedocs.io/en/latest/Python_Tutorials/Third_Harmonic_Generation/
"""

from typing import Tuple, List, Union
import numpy as np
import matplotlib

matplotlib.use("agg")
import matplotlib.pyplot as plt
import meep as mp


def third_harmonic_generation(
    k: float, amp: float = 1.0, nfreq: int = 10, flux_spectrum: bool = True
) -> Union[Tuple[List[float], List[float]], Tuple[float, float]]:

    """Computes the transmission spectrum of a plane wave propagating
       through a Kerr medium.

    Args:
      k: strength of Kerr susceptibility.
      amp: amplitude of the incident planewave.
      nfreq: number of frequencies in flux spectrum.
      flux_spectrum: compute the flux spectrum over broad bandwidth (True) or
                     just the two harmonic frequencies ω and 3ω (False).

    Returns:
      The frequencies and transmitted flux over the flux spectrum or
      the transmitted flux at the harmonic frequencies ω and 3ω.
    """

    sz = 100  # size of cell in z direction
    fcen = 1 / 3.0  # center frequency of source
    df = fcen / 20.0  # frequency width of source
    dpml = 1.0  # PML thickness

    # We'll use an explicitly 1d simulation.  Setting dimensions=1 will actually
    # result in faster execution than just using two no-size dimensions.  However,
    # in this case Meep requires us to use E in the x direction (and H in y),
    # and our one no-size dimension must be z.
    dimensions = 1
    cell = mp.Vector3(0, 0, sz)
    pml_layers = [mp.PML(dpml)]
    resolution = 25

    # to put the same material in all space, we can just set the default material
    # and pass it to the Simulation constructor
    default_material = mp.Medium(index=1, chi3=k)

    sources = [
        mp.Source(
            mp.GaussianSource(fcen, fwidth=df),
            component=mp.Ex,
            center=mp.Vector3(0, 0, -0.5 * sz + dpml),
            amplitude=amp,
        )
    ]

    # frequency range for flux calculation
    fmin = fcen / 2.0
    fmax = fcen * 4

    sim = mp.Simulation(
        cell_size=cell,
        sources=sources,
        boundary_layers=pml_layers,
        default_material=default_material,
        resolution=resolution,
        dimensions=dimensions,
    )

    mon_pt = mp.Vector3(0, 0, 0.5 * sz - dpml - 0.5)

    if flux_spectrum:
        trans = sim.add_flux(
            0.5 * (fmin + fmax),
            fmax - fmin,
            nfreq,
            mp.FluxRegion(mon_pt),
        )
    else:
        trans1 = sim.add_flux(fcen, 0, 1, mp.FluxRegion(mon_pt))
        trans3 = sim.add_flux(3 * fcen, 0, 1, mp.FluxRegion(mon_pt))

    sim.run(until_after_sources=mp.stop_when_fields_decayed(50, mp.Ex, mon_pt, 1e-6))

    if flux_spectrum:
        freqs = mp.get_flux_freqs(trans)
        trans_flux = mp.get_fluxes(trans)
        return freqs, trans_flux
    else:
        print(
            f"harmonics:, {k}, {amp}, {mp.get_fluxes(trans1)[0]}, "
            f"{mp.get_fluxes(trans3)[0]}"
        )
        return mp.get_fluxes(trans1)[0], mp.get_fluxes(trans3)[0]


if __name__ == "__main__":
    # Part 1: plot transmitted power spectrum for several values of χ(3).
    nfreq = 400
    logk = range(-3, 1)
    tflux = np.zeros((nfreq, len(logk)))
    for i, lk in enumerate(logk):
        freqs, tflux[:, i] = third_harmonic_generation(
            10**lk, nfreq=nfreq, flux_spectrum=True
        )

    fig, ax = plt.subplots()
    ax.semilogy(freqs, tflux[:, 0], "bo-", label=r"$\chi^{(3)}$=0.001")
    ax.semilogy(freqs, tflux[:, 1], "ro-", label=r"$\chi^{(3)}$=0.01")
    ax.semilogy(freqs, tflux[:, 2], "go-", label=r"$\chi^{(3)}$=0.1")
    ax.semilogy(freqs, tflux[:, 3], "co-", label=r"$\chi^{(3)}$=1")
    ax.set_xlabel("frequency")
    ax.set_ylabel("transmitted power (a.u.)")
    ax.set_xlim(0.2, 1.2)
    ax.set_ylim(1e-15, 1e2)
    ax.legend()
    ax.grid(True)
    fig.savefig(
        "transmitted_power_vs_frequency_vary_chi3.png", dpi=150, bbox_inches="tight"
    )

    # Part 2: plot transmittance vs. χ(3) for frequencies ω and 3ω.
    logk = np.arange(-6.0, 0.2, 0.2)
    first_order = np.zeros(len(logk))
    third_order = np.zeros(len(logk))
    for i, lk in enumerate(logk):
        first_order[i], third_order[i] = third_harmonic_generation(
            10**lk, flux_spectrum=False
        )

    input_flux = first_order[0]
    fig, ax = plt.subplots()
    ax.loglog(10**logk, first_order / input_flux, "ro-", label=r"$\omega$")
    ax.loglog(10**logk, third_order / input_flux, "bo-", label=r"$3\omega$")
    ax.loglog(10**logk, (10**logk) ** 2, "k-", label="quadratic line")
    ax.set_xlabel(r"$\chi^{(3)}$")
    ax.set_ylabel("transmission / incident power")
    ax.legend()
    ax.grid(True, "both")
    fig.savefig("transmittance_vs_chi3.png", dpi=150, bbox_inches="tight")