Materials#


The material structure in Maxwell's equations is determined by the relative permittivity and permeability .

However, is not only a function of position. In general, it also depends on frequency (material dispersion) and on the electric field itself (nonlinearity) or an applied magnetic field (gyrotropy). It may also depend on the orientation of the field (anisotropy). Material dispersion, in turn, is generally associated with absorption loss in the material, or possibly gain. All of these effects can be simulated in Meep, with certain restrictions. Similarly for the relative permeability , for which dispersion, nonlinearity, and anisotropy are all supported as well.

In this section, we describe the form of the equations and material properties that Meep can simulate. The actual user interface where these properties are specified in the simulation is described in Python Interface.

Material Dispersion#

Physically, material dispersion arises because the polarization of the material does not respond instantaneously to an applied field , and this is essentially the way that it is implemented in FDTD. In particular, is expanded to:

where , which must be positive, is the instantaneous dielectric function (the infinite-frequency response) and P is the remaining frequency-dependent polarization density in the material. P, in turn, has its own time-evolution equation, and the exact form of this equation determines the frequency-dependence .

Note: Meep's definition of uses a sign convention for the time dependence; formulas in engineering papers that use the opposite sign convention for will have a sign flip in all the imaginary terms below. If you are using parameters from the literature, you should use positive values of and as-is for loss; don't be confused by the difference in sign convention and flip the sign of the parameters.

Meep supports a Lorentzian susceptibility profile which consists of a sum of harmonic resonances plus a term for the frequency-independent electric conductivity:

where is the electric conductivity, and are user-specified constants. Actually, the numbers that one specifies are and . The is a user-specified function of position giving the strength of the -th resonance. The parameters can be anisotropic (real-symmetric) tensors, while the frequency-independent term can be an arbitrary real-symmetric tensor as well. This corresponds to evolving via the equations:

That is, we must store and evolve a set of auxiliary fields along with the electric field in order to keep track of the polarization . Essentially any could be modeled by including enough of these polarization fields — Meep allows you to specify any number of these, limited only by computer memory and time which increases with the number of polarization terms you require.

Note that the conductivity corresponds to an imaginary part of given by . This does not include the harmonic-resonance terms. When you specify frequency in Meep units, however, you are specifying without the , so the imaginary part of is .

Meep also supports polarizations of the Drude form, typically used for metals:

which corresponds to a term of the following form in 's summation:

which is equivalent to the Lorentzian model except that the term has been omitted from the denominator, and asymptotes to a conductivity as . In this case, is just a dimensional scale factor and has no interpretation as a resonance frequency.

How do I import n and k values into Meep?#

You can import into Meep any arbitrary complex permittivity profile via and values obtained via e.g. ellipsometry by fitting the wavelength- or frequency-dependent data to a sum of Lorentzian polarizability terms. In general, you have to use nonlinear optimization to perform the fit (e.g., to minimize the sum-of-squares errors or whatever error criterion you prefer). Enough Lorentzians should form a complete basis, so you should be able to fit any function given enough Lorentzians. Unfortunately, the fitting process requires some trial and error to specify the number of fitting parameters and their initial values. For a demonstration, see this script.

A wavelength-dependent, purely-real permittivity (i.e., with no loss) which can be represented using the Sellmeier equation can be directly transferred to the Lorentz model using a simple substitution of variables.

Note: Meep only does subpixel averaging of the nondispersive part of (and ).

Sellmeier Coefficients#

For a wavelength-dependent, purely-real permittivity (i.e., with no loss) which can be represented via the Sellmeier equation:

where is the vacuum wavelength, each term containing two coefficients ( and ) can be directly transferred to a Lorentzian polarization field using a simple substitution of variables: , , and . Several examples of importing Sellmeier coefficients from published fitting data including germanium (Ge) and gallium nitride (GaN) are provided in the Materials Library.

Numerical Stability#

In some cases, you may need to reduce the Courant parameter of the simulation, which relates the size of the time step to the spatial discretization: . By default, but in general you must have , where is the minimum refractive index (usually 1), so if your refractive indices are ever <1 you may need a smaller .

If a Lorentzian resonance at ω is specified at too high a frequency relative to the time discretization , the simulation becomes unstable. Essentially, the problem is that oscillates too fast compared with the time discretization for the discretization to work properly. If this happens, there are three workarounds: (1) increase the resolution which increases the resolution in both space and time, (2) decrease the Courant factor which decreases compared to , or (3) use a different model function for your dielectric response.

Roughly speaking, the equation becomes unstable for . Note that, in Meep frequency units, you specify , so this quantity should be less than .

Finally, overlapping dispersive materials with perfectly matched layer (PML) absorbing boundaries may produce instabilities. A workaround is to replace the PML with an absorber.

Loss and Gain#

If above is nonzero, then the dielectric function becomes complex, where the imaginary part is associated with absorption loss in the material if it is positive, or gain if it is negative. Alternatively, a dissipation loss or gain may be added by a positive or negative conductivity, respectively — this is often convenient if you only care about the imaginary part of in a narrow bandwidth, and is described in detail in the next section.

If you look at Maxwell's equations, then plays exactly the same role as a current . Just as is the rate of change of mechanical energy (the power expended by the electric field on moving the currents), therefore, the rate at which energy is lost to absorption is given by:

absorption rate

Meep can keep track of this energy for the Lorentzian polarizability terms but not for the conductivity terms. For gain, this gives the amount of energy expended in amplifying the field.

Conductivity and Complex ε#

Often, you only care about the absorption loss in a narrow bandwidth, where you just want to set the imaginary part of (or ) to some known experimental value, in the same way that you often just care about setting a dispersionless real that is the correct value in your bandwidth of interest.

One might hope that Meep would allow you to specify a constant, frequency-independent, imaginary part of , but this is inherently unstable: corresponds to gain (exponential growth), but every time-domain simulation will include both positive and negative frequencies ω so there is no way to choose a nonzero frequency-independent that does not yield gain at some ω. Instead, a simple approach in Meep is for you to set the conductivity (or for an imaginary part of ), chosen so that is the correct value at your frequency of interest. Note that, in Meep, you specify instead of for the frequency, however, so you need to include the factor of when computing the corresponding imaginary part of . Conductivities are absorptive at all frequencies (and are more computationally efficient than Lorentzian models).

For example, suppose you want to simulate a medium with at a frequency 0.42 (in your Meep units), and you only care about the material in a narrow bandwidth around this frequency (i.e. you don't need to simulate the full experimental frequency-dependent permittivity). Then, in Meep, you could use meep.Medium(epsilon=3.4, D_conductivity=2*math.pi*0.42*0.101/3.4) in Python or (make medium (epsilon 3.4) (D-conductivity (* 2 pi 0.42 0.101 (/ 3.4)))) in Scheme; i.e. and .

You can also use the feature to model the attenuation coefficient (units of e.g. dB/cm) obtained from experimental measurements (i.e., ellipsometry). This involves first converting into a complex refractive index (which is then converted into a complex permittivity) with imaginary part given by where is the vacuum wavelength.

Note: the "conductivity" in Meep is slightly different from the conductivity you might find in a textbook, because for computational convenience it appears as in our Maxwell equations rather than the more-conventional ; this just means that our definition is different from the usual electric conductivity by a factor of . Also, just as Meep uses the dimensionless relative permittivity for , it uses nondimensionalized units of 1/ (where is your unit of distance) for the conductivities . If you have the electric conductivity in SI units of S/m (S is siemens) and want to convert to in Meep units, you can simply use the formula: where is your unit of distance in meters, is the vacuum speed of light in m/s, is the SI vacuum permittivity, and is the real relative permittivity. The quantity in this equation is the conversion factor for frequency in SI units (s or Hz) to frequency in Meep units ().

Nonlinearity#

In general, can be changed anisotropically by the field itself, with:

where the ij is the index of the change in the 33 ε tensor and the terms are the nonlinear susceptibilities. The sum is the Pockels effect and the sum is the Kerr effect. If the above expansion is frequency-independent, then the nonlinearity is instantaneous; more generally, would depend on some average of the fields at previous times.

Meep supports instantaneous, isotropic Pockels and Kerr nonlinearities, corresponding to a frequency-independent and , respectively. Thus,

Here, "diag()" indicates the diagonal 33 matrix with the components of along the diagonal.

Normally, for nonlinear systems you will want to use real fields . This is the default. However, Meep uses complex fields if you have Bloch-periodic boundary conditions with a non-zero Bloch wavevector k, or in cylindrical coordinates with . In the C++ interface, real fields must be explicitly specified.

For complex fields in nonlinear systems, the physical interpretation of the above equations is unclear because one cannot simply obtain the physical solution by taking the real part any more. In particular, Meep simply defines the meaning of the nonlinearity for complex fields as follows: the real and imaginary parts of the fields do not interact nonlinearly. That is, the above equation should be taken to hold for the real and imaginary parts (of and ) separately: e.g., is the squared magnitude of the real part of for when computing the real part of , and conversely for the imaginary part.

For a discussion of how to relate in Meep to experimental Kerr coefficients, see Units and Nonlinearity.

Magnetic Permeability μ#

All of the above features that are supported for the electric permittivity are also supported for the magnetic permeability . That is, Meep supports with dispersion from magnetic conductivity and Lorentzian resonances, as well as magnetic and nonlinearities. The description of these is exactly the same as above, so we won't repeat it here — just take the above descriptions and replace , , , and by , , , and , respectively.

Saturable Gain and Absorption#

For some problems, simply adding gain or loss using a conductivity term does not correctly model the desired system and will lead to unphysical results. For example, attempting to model a laser by adding gain through a conductivity term will yield a diverging electric field, as the conductivity-based gain cannot saturate, and will continue to amplify arbitrarily strong electric fields within the laser cavity. Instead, such systems must be modeled using a saturable gain medium, in which the available gain is reduced for stronger electric fields, which prohibits the electric field from diverging in this manner.

Meep supports saturable gain and absorbing media through the use of a set of auxiliary equations which model both the polarization and inversion of the saturable medium. Meep's implementation of these auxiliary equations is described in arXiv:2007.09329, in which the polarization is treated using the oscillator model and the inversion is modeled as a set of rate equations for the population densities of the atomic levels of the saturable medium. The oscillator model for saturable media is a second-order differential equation for the polarization that is slightly different from the Drude-Lorentz susceptibility:

where is the inversion of the two atomic energy levels which comprise the th lasing transition, is the central frequency of the atomic transition, is the full width half-maximum of the width of the transition, and (a tensor) is the coupling strength between the electric field and the nonlinear polarization of the th transition. Note that this polarization equation is only modeling the nonlinear, saturable portion of the polarization. The atomic level population densities, , each satisfy a rate equation of the form:

where is the non-radiative decay or pumping rate from level to level . The final term in brackets is only included if level is either the upper or lower level of the th transition, where the upper/lower atomic levels are denoted by /.

In Meep, one can specify an arbitrary number of atomic levels with any number of lasing transitions between them, enabling one to realize common two- and four-level saturable media, as well as entire manifolds of levels and transitions found in realistic models of saturable media. When assigning the necessary transition frequencies, , and widths, , of the atomic transition in Meep, these are specified in units of 2π/. However, the pumping and decay rates, , are instead specified in units of . Finally, as part of initializing a saturable medium, the total atomic density, , must be specified, and Meep will ensure that .

In principle, there are a few different ways to alter the effective strength of the gain or absorption experienced by the electric field. For example, to increase the gain of a two-level medium, one could increase the density of the gain medium, (assuming ), further increase , or increase the relevant components of the coupling strength tensor, . In practice though, there are two considerations to keep in mind. First, for atomic and molecular gain media, the medium density, coupling strength, and non-radiative decay rates are all fixed by the intrinsic properties of the gain medium -- physically what can be changed is the pumping rate. However, this can occasionally be an inconvenient parameter to tune because it also effects the relaxation rate of the gain medium (see discussion about below), which can change the stability of the system, see Applied Physics Letters, Vol. 117, Num. 051102 (2020). Second, Meep will not automatically check to ensure that the non-radiative decay times, , are larger than the timestep, . If the non-radiative decay times are too short, this will lead to unphysical behavior and can result in diverging fields. As such, we would generally recommend tuning the effective gain or loss by changing .

Frequency Domain Susceptibility in Steady-State Regime#

It is possible to use saturable non-linear media to effectively implement linear gain and loss in a simulation, and treat it as a frequency-dependent linear susceptibility.

If the saturable medium only has two atomic levels (i.e. only a single non-linear polarization field) and the system is operated in the steady-state regime with only a single frequency component of the electric field, one can write the susceptibility of the non-linear saturable medium as For this two-level gain medium, , and we are assuming that the coupling matrix is proportional to the identity matrix, . Thus, if the non-linear saturation of the gain medium can be ignored, and the effectively linear susceptibility is and thus offers an independent route to realizing linear gain and loss rather than making the material conductive.

However, there are some important considerations to keep in mind when trying to use a non-linear saturable medium to approximate a linear susceptibility. First, saturable media can be a chaotic system, so one must check that the system is actually reaching a steady-state regime with only a single frequency component. Saturable non-linear media excited with more than one initial frequency can generate additional "beat" frequencies, so a sufficiently narrow spectral source must be used. Second, there can be relaxation oscillations between the electric field and the saturable medium, so even with a narrow source, one must wait for this transient to decay. Finally, keep in mind that even if your initial fields are calibrated to satisfy the above inequality to treat the medium as linear, one must still be operating below the lasing threshold, as if the system is above the lasing threshold, the fields will eventually grow to the point where the non-linear saturation is important.

Analytically, what this means is that the resonances that are present in the system must still be below the real axis in the complex plane, i.e. possess non-zero decay rates. If any added gain is sufficiently strong to move these resonances to the real axis (i.e. have a decay rate equal to zero), that is the definition of the laser threshold, and the non-linear medium is guaranteed to saturate.

Verification of the Oscillator Model Equations#

Although Meep is using an oscillator model equation for the atomic polarization and level populations, instead of the Bloch equations, Meep retains the two terms usually approximated to zero when deriving the oscillator model equations from the Bloch equations, and so these equations are exactly equivalent to the Bloch equations. For more details, see arXiv:2007.09329 and Section 6.4.1 of Nonlinear Optics (third edition) by R. W. Boyd. To verify this equivalence between the different equations for modeling the polarization, as well as confirm that saturable media have been properly implemented in Meep, we compare the results of Meep with an independent FDTD solver using the Bloch equations and the frequency domain steady-state ab initio laser theory (SALT), in a 1d, one-sided, Fabry-Perot cavity containing a two-level gain medium that exhibits steady-state, multi-mode lasing. The cavity has a length of , a background index of , an atomic transition frequency of , with width , the decay rate from level to is , the pumping rate from to is , and the total atomic density, , of the system was varied to produce different amounts of gain. These plots are given in terms of the equilibrium inversion, , which is the inversion of the saturable gain medium in the absence of an electric field, i.e. the value of when . There is agreement among the three methods close to the initial as well as past the third lasing threshold as shown in the following two figures.

Near threshold comparison

Near threshold comparison

For the two-level atomic gain model used in this example, can be calculated as:

Similar relationships can be found for systems with more than two atomic levels in Optics Express, Vol. 20, pp. 474-88 (2012).

Natural Units for Saturable Media#

A conversion chart for different choices of parameter nomenclature for the saturable medium.

There is no standard convention in the literature on lasers and saturable gain media for defining the various constants in the equations above. The following are the relationships among these constants for the three different groups of work discussed in this section:

The relationship among Meep and SALT units are:

For more details on applying SALT to atomic media with an arbitrary number of levels, see Optics Express, Vol. 23, pp. 6455-77 (2015).

An additional discussion of the natural units for saturable media is given in arXiv:2007.09329.

Gyrotropic Media#

(Experimental feature) Meep supports gyrotropic media, which break optical reciprocity and give rise to magneto-optical phenomena such as the Faraday effect. Such materials are used in devices like Faraday rotators. For a demonstration, see Tutorial/Gyrotropic Media.

In a gyrotropic medium, the polarization vector undergoes precession around a preferred direction. In the frequency domain, this corresponds to the presence of skew-symmetric off-diagonal components in the tensor (for a gyroelectric medium) or the μ tensor (for a gyromagnetic medium). Two different gyrotropy models are supported:

Gyrotropic Drude-Lorentz Model#

The first gyrotropy model is a Drude-Lorentz model with an additional precession, which is intended to describe gyroelectric materials. The polarization equation is

(Optionally, the polarization may be of Drude form, in which case the term on the left is omitted.) The third term on the left side, which breaks time-reversal symmetry, is responsible for the gyrotropy; it typically describes the deflection of electrons flowing within the material by a static external magnetic field. In the limit, the equation of motion reduces to a precession around the "bias vector" :

Hence, the magnitude of the bias vector is the angular frequency of the gyrotropic precession induced by the external field.

Gyrotropic Saturated Dipole (Linearized Landau-Lifshitz-Gilbert) Model#

The second gyrotropy model is a linearized Landau-Lifshitz-Gilbert equation, suitable for modeling gyromagnetic materials such as ferrites. Its polarization equation of motion is

Note: although the above equation is written in terms of electric susceptibilities, this model is typically used for magnetic susceptibilities. Meep places no restriction on the field type that either gyrotropy model can be applied to. As usual, electric and magnetic susceptibilities can be swapped by substituting ε with μ, with , etc.

The Landau-Lifshitz-Gilbert equation describes the precessional motion of a saturated point magnetic dipole in a magnetic field. In the above equation, the variable represents the linearized deviation of the polarization from its static equilibrium value (assumed to be much larger and aligned parallel to ). Note that this equation of motion is completely different from the Drude-Lorentz equation, though the constants σ, ω, and γ play analogous roles (σ couples the polarization to the driving field, ω is the angular frequency of precession, and γ is a damping factor).

In this model, is taken to be a unit vector (i.e., its magnitude is ignored).

Frequency Domain Susceptibility Tensors#

Suppose , and let all fields have harmonic time-dependence . Then is related to the applied field by

For the gyrotropic Lorentzian model, the components of the susceptibility tensor are

And for the gyrotropic saturated dipole (linearized Landau-Lifshitz-Gilbert) model,

Materials Library#

A materials library is available for Python and Scheme containing crystalline silicon (c-Si), amorphous silicon (a-Si) including the hydrogenated form, silicon dioxide (SiO2), indium tin oxide (ITO), alumina (Al2O3), gallium arsenide (GaAs), aluminum arsenide (AlAs), aluminum nitride (AlN), borosilicate glass (BK7), fused quartz, silicon nitride (Si3N4), germanium (Ge), indium phosphide (InP), gallium nitride (GaN), cadmium telluride (CdTe), lithium niobate (LiNbO3), barium borate (BaB2O4), calcium tungstate (CaWO4), calcium carbonate (CaCO3), yttrium oxide (Y2O3), yttrium aluminum garnet (YAG), poly(methyl methacrylate) (PMMA), polycarbonate, polystyrene, cellulose, as well as 11 elemental metals: silver (Ag), gold (Au), copper (Cu), aluminum (Al), berylium (Be), chromium (Cr), nickel (Ni), palladium (Pd), platinum (Pt), titanium (Ti), and tungsten (W).

Experimental values of the complex refractive index are fit to a Drude-Lorentz susceptibility model over various wavelength ranges. For example, the fit for crystalline silicon is based on Progress in Photovoltaics, Vol. 3, pp. 189-92 (1995) for the wavelength range of 0.4-1.0 µm as described in J. Optical Society of America A, Vol. 28, pp. 770-77 (2011). The fit for the elemental metals is over the range of 0.2-12.4 μm and is described in Applied Optics, Vol. 37, pp. 5271-83 (1998).

Fitting parameters for all materials are defined for a unit distance of 1 µm. For simulations which use a different value for the unit distance, the predefined variable um_scale (Python) or um-scale (Scheme) must be scaled by multiplying by whatever the unit distance is, in units of µm. For example, if the unit distance is 100 nm, this would require adding the line um_scale = 0.1*um_scale after the line where um_scale is defined. This change must be made directly to the materials library file.

As an example, to import aluminum from the library into a Python script requires adding the following lines:

from meep.materials import Al

Then, the material can be simply used as geometry = [meep.Cylinder(material=Al, ...].

You can inspect the frequency range over which the permittivity profile for aluminum is defined and also inspect the value of its permittivity tensor at a wavelength of 1.55 µm using the epsilon(frequency) function of the Medium class:

> Al.valid_freq_range
FreqRange(min=0.08065817067268914, max=4.839334107626791)
> Al.epsilon(1/1.55)
array([[-232.66329725+47.08102623j,  0.         +0.j        ,  0.         +0.j        ],
       [ 0.         +0.j        , -232.66329725+47.08102623j,  0.         +0.j        ],
       [ 0.         +0.j        ,  0.         +0.j       ,  -232.66329725+47.08102623j]])

In Scheme, the materials library is already included when Meep is run, so you can use it without any initialization:

(set! geometry (list (make cylinder (material Al) ...)))

As a final example, we plot the complex refractive index of SiO from the materials library over the visible wavelength range of 400 nm to 700 nm.

from meep.materials import SiO2
import numpy as np
import matplotlib.pyplot as plt

wvl_min = 0.4 # units of μm
wvl_max = 0.7 # units of μm
nwvls = 21
wvls = np.linspace(wvl_min, wvl_max, nwvls)

SiO2_epsilon = np.array([SiO2.epsilon(1/w)[0][0] for w in wvls])

plt.subplot(1,2,1)
plt.plot(wvls,np.real(SiO2_epsilon),'bo-')
plt.xlabel('wavelength (μm)')
plt.ylabel('real(ε)')

plt.subplot(1,2,2)
plt.plot(wvls,np.imag(SiO2_epsilon),'ro-')
plt.xlabel('wavelength (μm)')
plt.ylabel('imag(ε)')

plt.suptitle('SiO$_2$ from Meep materials library')
plt.subplots_adjust(wspace=0.4)
plt.show()

SiO2 from the materials library

Note: for narrowband calculations, some of the Lorentzian susceptibility terms may be unnecessary and will contribute to consuming more computational resources than are required (due to the additional storage and time stepping of the polarization fields). Computational efficiency can be improved (without significantly affecting accuracy) by removing from the material definitions those Lorentzian susceptibility terms which are far outside the spectral region of interest. In addition, when using the materials library the Courant parameter may need to be reduced to achieve meaningful results, see Numerical Stability.