Large-area metalens#

Run this notebook in your browser using Binder.

Introduction#

Here we use Tidy3D to simulate a very large dielectric metalens. We base this example off of the recent paper from Khorasaninejad et al. titled Metalenses at visible wavelengths: Diffraction-limited focusing and subwavelength resolution imaging, which was recently published in Science. In this paper, a 2-dimensional array of dielectric structures is used as a lens to focus transmitted light to a single position directly above the device.

Typically, these structures are simulated by simulating each dielectric unit cell individually to compute a phase and amplitude transmittance for each cell. While this approach makes for an approximation of the overall device performance, it would be useful to be able to simulate the entire device as a whole to capture the entire physics. However, a simulation of this scale requires several hours or days to do with a conventional CPU-based FDTD. With Tidy3D, we are able to complete the entire simulation in about 1 minute!

Setup#

We first perform basic imports of the packages needed.

[1]:
# standard python imports
import numpy as np
from numpy import random
import matplotlib.pyplot as plt

import tidy3d as td
from tidy3d import web

Define Simulation Parameters#

We now set the basic parameters that define the metalens. The following image (taken from the original paper) define the variables describing the unit cell of the metalens. The angle of each cell (θ) is chosen for maximum focusing effect. Note that microns are the default spatial unit in Tidy3D.

diagram

[2]:
# 1 nanometer in units of microns (for conversion)
nm = 1e-3

# free space central wavelength
wavelength = 600 * nm

# desired numerical aperture
NA = 0.5

# shape parameters of metalens unit cell (um) (refer to image above and see paper for details)
W = 85 * nm
L = 410 * nm
H = 600 * nm
S = 430 * nm

# space between bottom PML and substrate (-z)
space_below_sub = 1 * wavelength

# thickness of substrate
thickness_sub = 100 * nm

# side length of entire metalens (um)
side_length = 10

# Number of unit cells in each x and y direction (NxN grid)
N = int(side_length / S)
print(f'for diameter of {side_length:.1f} um, have {N} cells per side')
print(f'full metalens has area of {side_length**2:.1f} um^2 and {N*N} total cells')

# Define material properties at 600 nm
n_TiO2 = 2.40
n_SiO2 = 1.46
air = td.Medium(permittivity=1.0)
SiO2 = td.Medium(permittivity=n_SiO2**2)
TiO2 = td.Medium(permittivity=n_TiO2**2)
for diameter of 10.0 um, have 23 cells per side
full metalens has area of 100.0 um^2 and 529 total cells
[3]:
# using the wavelength in microns, one can use td.C_0 (um/s) to get frequency in Hz
# wavelength_meters = wavelength * meters
f0 = td.C_0 / wavelength

# Compute the domain size in x, y (note: round down from side_length)
length_xy = N * S

# focal length given diameter and numerical aperture
f = length_xy / 2 / NA * np.sqrt(1 - NA**2)

# Function describing the theoretical best angle of each box at position (x,y).  see paper for details
def theta(x, y):
    return np.pi / wavelength * (f - np.sqrt(x ** 2 + y ** 2 + f ** 2))

# total domain size in z: (space -> substrate -> unit cell -> 1.7 focal lengths)
length_z = space_below_sub + thickness_sub + H + 1.7 * f

# construct simulation size array
sim_size = (length_xy, length_xy, length_z)

Create Metalens Geometry#

Now we can automatically generate a large metalens structure using these parameters.

We will first create the substrate as a td.Box

Then, we will loop through the x and y coordinates of the lens and create each unit cell as a td.PolySlab

[4]:
# define substrate
substrate = td.Structure(
    geometry=td.Box(
        center=[0, 0, -length_z/2 + space_below_sub + thickness_sub / 2.0],
        size=[td.inf, td.inf, thickness_sub],
    ),
    medium=SiO2
)

# define coordinates of each unit cell
centers_x = S * np.arange(N) - length_xy / 2.0 + S / 2.0
centers_y = S * np.arange(N) - length_xy / 2.0 + S / 2.0
center_z = -length_z/2 + space_below_sub + thickness_sub + H / 2.0

# x, y vertices of box of size (L, W) centered at the origin
vertices_origin = np.array([[+L/2, +W/2],
                            [-L/2, +W/2],
                            [-L/2, -W/2],
                            [+L/2, -W/2]])


xs, ys = np.meshgrid(centers_x, centers_y, indexing='ij')
xs = xs.flatten()
ys = ys.flatten()

angles = theta(xs, ys)

# 2x2 rotation matrix angle `angle` with respect to x axis
rotation_matrix = np.array([
    [+np.cos(angles), -np.sin(angles)],
    [+np.sin(angles), +np.cos(angles)]
])

# rotate the origin vertices by this angle
vertices_rotated = np.einsum('ij, jkn -> nik', vertices_origin, rotation_matrix)

# shift the rotated vertices to be centered at (xs, ys)
vertices_shifted = vertices_rotated + np.stack([xs, ys], axis=-1)[:, None, :]

metalens_geometry = []
for vertices in vertices_shifted:
    # create a tidy3D PolySlab with these rotated and shifted vertices and thickness `H`
    metalens_geometry.append(
        td.PolySlab(
            vertices=vertices.tolist(),
            slab_bounds=(center_z-H/2, center_z+H/2),
            axis=2,
        ),
    )

metalens = td.Structure(
    geometry=td.GeometryGroup(geometries=metalens_geometry),
    medium=TiO2
)

Define grid specification#

We define the grid based on the properties of the geometry. The metalens is quasi-periodic in x and y, in that we have clearly defined unit cells, but each is slightly modified from its neighbors. Such structures are best resolved with a grid that matches the periodicity, which is why we use a uniform grid in x and y. In z, we use the automatic nonuniform grid that will place a higher grid density around the metalens region, and a lower one in the air region away from the metalens. To speed up the auto meshing in the region with the pillars, we put an override box in the grid specification.

[5]:
# steps per unit cell along x and y
grids_per_unit_length = 20
# uniform mesh in x and y
grid_x = td.UniformGrid(dl=S / grids_per_unit_length)
grid_y = td.UniformGrid(dl=S / grids_per_unit_length)
# in z, use an automatic nonuniform mesh with the wavelength being the "unit length"
grid_z = td.AutoGrid(min_steps_per_wvl=grids_per_unit_length)
# we need to supply the wavelength because of the automatic mesh in z
grid_spec = td.GridSpec(wavelength=wavelength, grid_x=grid_x, grid_y=grid_y, grid_z=grid_z)
# put an override box over the pillars to avoid parsing a large amount of structures in the mesher
grid_spec = grid_spec.copy(update=dict(override_structures = [
    td.Structure(
        geometry=td.Box.from_bounds(
            rmin=(-td.inf, -td.inf, -length_z/2 + space_below_sub),
            rmax=(td.inf, td.inf, center_z + H/2)
        ),
        medium=TiO2,
    )
]))

Define Sources#

Now we define the incident fields. We simply use a normally incident plane wave with Guassian time dependence centered at our central frequency. For more details, see the plane wave source documentation and the gaussian source documentation

[6]:
# Bandwidth in Hz
fwidth = f0 / 10.0

# time dependence of source
gaussian = td.GaussianPulse(freq0=f0, fwidth=fwidth, phase=0)

source = td.PlaneWave(
    source_time=gaussian,
    size=(td.inf, td.inf, 0),
    center=(0,0,-length_z/2 + space_below_sub / 10.0),
    direction='+',
    pol_angle=0)

run_time = 50 / fwidth

Define Monitors#

Now we define the monitors that measure field output from the FDTD simulation. For simplicity, we use measure the fields at the central frequency. We’ll get the fields at the following locations:

  • The y=0 cross section

  • The x=0 cross section

  • The z=f cross section at the focal plane

  • The central axis, along x=y=0

For more details on defining monitors, see the frequency-domain field monitor documentation.

[7]:
# To decrease the amount of data stored, only store the E field in 2D monitors
fields = ["Ex", "Ey", "Ez"]

# get fields along x=y=0 axis
monitor_center = td.FieldMonitor(
    center=[0., 0., 0],
    size=[0, 0, td.inf],
    freqs=[f0],
    name='center'
)

# get the fields at a few cross-sectional planes
monitor_xz = td.FieldMonitor(
    center=[0., 0., 0.],
    size=[td.inf, 0., td.inf],
    freqs=[f0],
    name='xz',
    fields=fields,
)

monitor_yz = td.FieldMonitor(
    center=[0., 0., 0.],
    size=[0., td.inf, td.inf],
    freqs=[f0],
    name='yz',
    fields=fields,
)

monitor_xy = td.FieldMonitor(
    center=[0., 0., center_z + H/2 + f],
    size=[td.inf, td.inf, 0],
    freqs=[f0],
    name='focal_plane',
    fields=fields,
)

# put them into a single list
monitors=[monitor_center, monitor_xz, monitor_yz, monitor_xy]

Create Simulation#

Now we can put everything together and define a Simulation object to be run.

We get a number of warnings about structures being too close to the PML. In FDTD simulations, this can result in instability, as PML are absorbing for propagating fields, but can be amplifying for evanescent fields. This particular simulation runs without any issues even with PML on the sides, but it is best to heed these warnings to avoid problems. There are two ways that we can fix the simulation: one is to just put some space between the last of the metalens boxes and the PML. The other is to use adiabatic absorbers on the sides, which are always stable. The only downside of the absorbers is that they are slightly thicker than the PML, making the overall simulation size slighlty larger. This is why we only put them along x and y, while we leave the PML in z.

[8]:
sim = td.Simulation(
    size=sim_size,
    grid_spec=grid_spec,
    structures=[substrate, metalens],
    sources=[source],
    monitors=monitors,
    run_time=run_time,
    boundary_spec=td.BoundarySpec(
        x=td.Boundary.absorber(),
        y=td.Boundary.absorber(),
        z=td.Boundary.pml()
    )
)

Visualize Geometry#

Lets take a look and make sure everything is defined properly

[9]:
fig, (ax1, ax2, ax3) = plt.subplots(1, 3, figsize=(14, 6))

# Time the visualization of the 2D plane
sim.plot(x=0.1, ax=ax1);
sim.plot(y=0.1, ax=ax2);
sim.plot(z=-length_z/2 + space_below_sub + thickness_sub + H / 2, ax=ax3);
<AxesSubplot:title={'center':'cross section at z=-6.93'}, xlabel='x', ylabel='y'>
<Figure size 1008x432 with 3 Axes>
../_images/notebooks_Metalens_18_2.png

Run Simulation#

Now we can run the simulation over time and measure the results to plot

[10]:
job = web.Job(simulation=sim, task_name='metalens')
sim_data = job.run(path='data/simulation_data.hdf5')
[17:06:55] INFO     Using Tidy3D credentials from stored file                      auth.py:74
[17:06:57] INFO     Uploaded task 'metalens' with task_id                       webapi.py:120
                    'abc7da22-ad8b-4111-a16e-350d080a0d55'.                                  
[17:07:01] INFO     status = queued                                             webapi.py:253
[17:07:06] INFO     Maximum flex unit cost: 1.79                                webapi.py:244
[17:07:08] INFO     status = preprocess                                         webapi.py:265
[17:07:16] INFO     starting up solver                                          webapi.py:269
[17:07:27] INFO     running solver                                              webapi.py:275
[17:13:01] INFO     early shutoff detected, exiting.                            webapi.py:286
           INFO     status = postprocess                                        webapi.py:292
[17:13:11] INFO     status = success                                            webapi.py:298
           INFO     downloading file "output/monitor_data.hdf5" to              webapi.py:574
                    "data/simulation_data.hdf5"                                              
[17:13:27] INFO     loading SimulationData from data/simulation_data.hdf5       webapi.py:398
[11]:
print(sim_data.log)
Simulation domain Nx, Ny, Nz: [541, 541, 589]
Applied symmetries: (0, 0, 0)
Number of computational grid points: 1.7426e+08.
Using subpixel averaging: True
Number of time steps: 3.4796e+04
Automatic shutoff factor: 1.00e-05
Time step (s): 2.8760e-17


Compute source modes time (s):     0.1234
Compute monitor modes time (s):    0.0725
Rest of setup time (s):            7.7395

Running solver for 34796 time steps...
- Time step    554 / time 1.59e-14s (  1 % done), field decay: 1.00e+00
- Time step   1391 / time 4.00e-14s (  4 % done), field decay: 8.68e-01
- Time step   2783 / time 8.00e-14s (  8 % done), field decay: 1.07e-01
- Time step   4175 / time 1.20e-13s ( 12 % done), field decay: 5.77e-03
- Time step   5567 / time 1.60e-13s ( 16 % done), field decay: 1.02e-03
- Time step   6959 / time 2.00e-13s ( 20 % done), field decay: 1.97e-04
- Time step   8351 / time 2.40e-13s ( 24 % done), field decay: 5.20e-05
- Time step   9742 / time 2.80e-13s ( 28 % done), field decay: 1.79e-05
- Time step  11134 / time 3.20e-13s ( 32 % done), field decay: 6.29e-06
Field decay smaller than shutoff factor, exiting solver.

Solver time (s):                   330.2838
Post-processing time (s):          1.5360



As we can see from the logs, the total time to solve this problem (not including data transfer and pre/post-processing) was on the order of a minute!

For reference, the same problem run with FDTD on a single CPU core FDTD is projected to take 40 hours!

Visualize Fields#

Let’s see the results of the simulation as captured by our field monitors.

First, we look at the field intensity along the axis of the lens

[12]:
focal_z = center_z+H/2+f
data_center_line = sim_data.at_centers('center')
I = abs(data_center_line.Ex)**2 + abs(data_center_line.Ey)**2 + abs(data_center_line.Ez)**2
I.plot()
plt.title('intensity(z)')
plt.show()
<Figure size 432x288 with 1 Axes>
../_images/notebooks_Metalens_24_1.png

We now can inspect the field patterns on the area monitors using the Tidy3D visualization nethods.

[13]:
fig, (ax1, ax2, ax3) = plt.subplots(1, 3, figsize=(14, 4))
sim_data.plot_field('xz', 'Ex', val='real', freq=f0, y=0, ax=ax1)
sim_data.plot_field('yz', 'Ex', val='real', freq=f0, x=0, ax=ax2)
sim_data.plot_field('focal_plane', 'Ex', val='real', freq=f0, z=focal_z, ax=ax3)
ax1.set_title('x-z plane')
ax2.set_title('y-z plane')
ax3.set_title('x-y (focal) plane')
plt.show()

<Figure size 1008x288 with 6 Axes>
../_images/notebooks_Metalens_26_1.png
[14]:
# plot absolute value for good measure
fig, (ax1, ax2, ax3) = plt.subplots(1, 3, figsize=(14, 4))
sim_data.plot_field('xz', 'Ex', val='abs', freq=f0, y=0, ax=ax1)
sim_data.plot_field('yz', 'Ex', val='abs', freq=f0, x=0, ax=ax2)
sim_data.plot_field('focal_plane', 'Ex', val='abs', freq=f0, z=focal_z, ax=ax3)
ax1.set_title('x-z plane')
ax2.set_title('y-z plane')
ax3.set_title('x-y (focal) plane')
plt.show()
<Figure size 1008x288 with 6 Axes>
../_images/notebooks_Metalens_27_1.png
[15]:
# and lets plot the intensites as well
fig, (ax1, ax2, ax3) = plt.subplots(1, 3, figsize=(14, 4))
sim_data.plot_field('xz', 'int', freq=f0, y=0, ax=ax1)
sim_data.plot_field('yz', 'int', freq=f0, x=0, ax=ax2)
sim_data.plot_field('focal_plane', 'int', freq=f0, z=focal_z, ax=ax3)
ax1.set_title('x-z plane')
ax2.set_title('y-z plane')
ax3.set_title('x-y (focal) plane')
plt.show()


<Figure size 1008x288 with 6 Axes>
../_images/notebooks_Metalens_28_1.png
[ ]: