Near to far field transformation#

To run this notebook from your browser, click this link.

This tutorial will show you how to solve for electromagnetic fields far away from your structure using field information stored on a nearby surface.

This technique is called a ‘near field to far field transformation’ and is very useful for reducing the simulation size needed for structures involving lots of empty space.

As an example, we will simulate a simple zone plate lens with a very thin domain size to get the transmitted fields measured just above the structure. Then, we’ll show how to use the Near Field to Far Field transformation feature from Tidy3D to extrapolate to the fields at the focal plane above the lens.

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

# tidy3d imports
import tidy3d as td
import tidy3d.web as web

Problem Setup#

Below is a rough sketch of the setup of a near field to far field transformation.

The transmitted near fields are measured just above the metalens on the blue line, and the near field to far field transformation is then used to project the fields to the focal plane above at the red line.

64c157a642d84e0f97de9d70fa6c0b93

Define Simulation Parameters#

As always, we first need to define our simulation parameters. As a reminder, all length units in tidy3D are specified in microns.

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

# free space central wavelength
wavelength = 1.0

# numerical aperture
NA = 0.8

# height of lens features
height_lens = 200 * nm

# space between bottom PML and substrate (-z)
# and the space between lens structure and top pml (+z)
space_below_sub = 1.5 * wavelength

# height of substrate (um)
height_sub = wavelength / 2

# side length (xy plane) of entire metalens (um)
length_xy = 20 * wavelength

# Lens and substrate refractive index
n_TiO2 = 2.40
n_SiO2 = 1.46

# define material properties
air = td.Medium(permittivity=1.0)
SiO2 = td.Medium(permittivity=n_SiO2**2)
TiO2 = td.Medium(permittivity=n_TiO2**2)

Process Geometry#

Next we perform some conversions based on these parameters to define the simulation.

[3]:
# because the wavelength is in microns, use builtin td.C_0 (um/s) to get frequency in Hz
f0 = td.C_0 / wavelength

# Define PML layers, for this application we surround the whole structure in PML to isolate the fields
boundary_spec = td.BoundarySpec.all_sides(boundary=td.PML())

# domain size in z, note, we're just simulating a thin slice: (space -> substrate -> lens height -> space)
length_z = space_below_sub + height_sub + height_lens + space_below_sub

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

Create Geometry#

Now we create the ring metalens programatically

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

# focal length
focal_length = length_xy / 2 / NA * np.sqrt(1 - NA**2)

# location from center for edge of the n-th inner ring, see https://en.wikipedia.org/wiki/Zone_plate
def edge(n):
    return np.sqrt(n * wavelength * focal_length + n**2 * wavelength**2 / 4)

# loop through the ring indeces until it's too big and add each to geometry list
n = 1
r = edge(n)
rings = []
while r < 2 * length_xy:
    # progressively wider cylinders, material alternating between air and TiO2

    cylinder = td.Structure(
        geometry=td.Cylinder(
            center=[0,0,-length_z/2  + space_below_sub + height_sub + height_lens / 2],
            axis=2,
            radius=r,
            length=height_lens),
        medium=TiO2 if n % 2 == 0 else air,
    )
    rings.append(cylinder)

    n += 1
    r = edge(n)

# reverse geometry list so that inner, smaller rings are added last and therefore override larger rings.
rings.reverse()
geometry = [substrate] + rings

Create Source#

Create a plane wave incident from below the metalens

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

# Gaussian source offset; the source peak is at time t = offset/fwidth
offset = 4.

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

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

# Simulation run time
run_time = 40 / fwidth

Create Monitor#

Create a near field monitor to measure the fields just above the metalens

[6]:
# place it halfway between top of lens and PML
pos_monitor_z = -length_z/2 + space_below_sub + height_sub + height_lens + space_below_sub / 2
monitor_near = td.FieldMonitor(
    center=[0., 0., pos_monitor_z],
    size=[td.inf, td.inf, 0],
    freqs=[f0],
    name='nearfield'
)

Create Simulation#

Put everything together and define a simulation object. A nonuniform simulation grid is generated automatically based on a given number of cells per wavelength in each material (10 by default), using the frequencies defined in the sources.

[7]:
simulation = td.Simulation(
    size=sim_size,
    grid_spec = td.GridSpec.auto(min_steps_per_wvl=20),
    structures=geometry,
    sources=[source],
    monitors=[monitor_near],
    run_time=run_time,
    boundary_spec=boundary_spec
)
[17:23:04] WARNING  Structure at structures[57] has bounds that extend      simulation.py:291
                    exactly to simulation edges. This can cause unexpected                   
                    behavior. If intending to extend the structure to                        
                    infinity along one dimension, use td.inf as a size                       
                    variable instead to make this explicit.                                  
           WARNING  Structure at structures[57] has bounds that extend      simulation.py:291
                    exactly to simulation edges. This can cause unexpected                   
                    behavior. If intending to extend the structure to                        
                    infinity along one dimension, use td.inf as a size                       
                    variable instead to make this explicit.                                  
           WARNING  Structure at structures[57] has bounds that extend      simulation.py:291
                    exactly to simulation edges. This can cause unexpected                   
                    behavior. If intending to extend the structure to                        
                    infinity along one dimension, use td.inf as a size                       
                    variable instead to make this explicit.                                  
           WARNING  Structure at structures[57] has bounds that extend      simulation.py:291
                    exactly to simulation edges. This can cause unexpected                   
                    behavior. If intending to extend the structure to                        
                    infinity along one dimension, use td.inf as a size                       
                    variable instead to make this explicit.                                  

Visualize Geometry#

Let’s take a look and make sure everything is defined properly

[8]:
fig, (ax1, ax2) = plt.subplots(1, 2, figsize=(20, 8))
simulation.plot_eps(x=0, ax=ax1);
simulation.plot_eps(z=-length_z/2  + space_below_sub + height_sub + height_lens / 2, ax=ax2);
           INFO     Auto meshing using wavelength 1.0000 defined from        grid_spec.py:472
                    sources.                                                                 
<AxesSubplot:title={'center':'cross section at z=0.25'}, xlabel='x', ylabel='y'>
<Figure size 1440x576 with 4 Axes>
../_images/notebooks_Near2Far_16_3.png

Run Simulation#

Now we can run the simulation and download the results

[9]:
import tidy3d.web as web

sim_data = web.run(simulation, task_name='near2far', path='data/simulation.hdf5')
[17:23:05] INFO     Using Tidy3D credentials from stored file                      auth.py:74
[17:23:07] INFO     Uploaded task 'near2far' with task_id                       webapi.py:120
                    '8ec6d1ad-c582-49e6-8ebd-e07486f3c22c'.                                  
[17:23:10] INFO     status = queued                                             webapi.py:253
[17:23:15] INFO     Maximum flex unit cost: 1.29                                webapi.py:244
[17:23:19] INFO     status = preprocess                                         webapi.py:265
[17:23:32] INFO     starting up solver                                          webapi.py:269
[17:23:43] INFO     running solver                                              webapi.py:275
[17:25:05] INFO     early shutoff detected, exiting.                            webapi.py:286
           INFO     status = postprocess                                        webapi.py:292
[17:25:16] INFO     status = success                                            webapi.py:298
[17:25:17] INFO     downloading file "output/monitor_data.hdf5" to              webapi.py:574
                    "data/simulation.hdf5"                                                   
[17:25:44] INFO     loading SimulationData from data/simulation.hdf5            webapi.py:398
           WARNING  Structure at structures[57] has bounds that extend      simulation.py:291
                    exactly to simulation edges. This can cause unexpected                   
                    behavior. If intending to extend the structure to                        
                    infinity along one dimension, use td.inf as a size                       
                    variable instead to make this explicit.                                  
           WARNING  Structure at structures[57] has bounds that extend      simulation.py:291
                    exactly to simulation edges. This can cause unexpected                   
                    behavior. If intending to extend the structure to                        
                    infinity along one dimension, use td.inf as a size                       
                    variable instead to make this explicit.                                  
           WARNING  Structure at structures[57] has bounds that extend      simulation.py:291
                    exactly to simulation edges. This can cause unexpected                   
                    behavior. If intending to extend the structure to                        
                    infinity along one dimension, use td.inf as a size                       
                    variable instead to make this explicit.                                  
           WARNING  Structure at structures[57] has bounds that extend      simulation.py:291
                    exactly to simulation edges. This can cause unexpected                   
                    behavior. If intending to extend the structure to                        
                    infinity along one dimension, use td.inf as a size                       
                    variable instead to make this explicit.                                  
           INFO     Auto meshing using wavelength 1.0000 defined from        grid_spec.py:472
                    sources.                                                                 

Visualization#

Let’s inspect the near field using the Tidy3D builtin field visualization methods.
For more details see the documentation of tidy3d.SimulationData.
[10]:
near_field_data = sim_data['nearfield']

fig, (ax1, ax2, ax3) = plt.subplots(1, 3, tight_layout=True, figsize=(15, 3.5))
near_field_data.Ex.interp(z=pos_monitor_z).real.plot(ax=ax1)
near_field_data.Ey.interp(z=pos_monitor_z).real.plot(ax=ax2)
near_field_data.Ez.interp(z=pos_monitor_z).real.plot(ax=ax3)
plt.show()
[17:25:45] INFO     Auto meshing using wavelength 1.0000 defined from        grid_spec.py:472
                    sources.                                                                 
<Figure size 1080x252 with 6 Axes>
../_images/notebooks_Near2Far_20_2.png

Setting Up Near 2 Far#

To set up near to far, we first need to grab the data from the nearfield monitor, for which we create a Near2FarSurface object to specify the monitor and its normal direction.

Then, we create a Near2Far object using the monitor data as follows. This object just stores near field data and provides various methods for looking at various far field quantities.

[11]:
from tidy3d.plugins import Near2Far, Near2FarSurface

n2f_surface = Near2FarSurface(monitor=monitor_near, normal_dir='+')

n2f = Near2Far(
    sim_data=sim_data,
    surfaces=[n2f_surface],
    frequency=f0
)

Getting Far Field Data#

With the Near2Far object initialized, we just need to call one of its methods to get a far field quantity.

For this example, we use Near2Far.fields_cartesian(x,y,z) to get the fields at a set of x,y,z points relative to the monitor center.

Below, we pass in an array of x and y points in a plane located at z=z0 and record the far fields.

[12]:
# points to project to
num_far = 40
xs_far = 4 * wavelength * np.linspace(-0.5, 0.5, num_far)
ys_far = 4 * wavelength * np.linspace(-0.5, 0.5, num_far)

far_fields = n2f.fields_cartesian(xs_far, ys_far, focal_length)

Plot Results#

Now we can plot the near and far fields together

[13]:
# plot everything
f, (axes_near, axes_far) =  plt.subplots(2, 3, tight_layout=True, figsize=(10, 5))

def pmesh(xs, ys, array, ax, cmap):
    im = ax.pcolormesh(xs, ys, array.T, cmap=cmap, shading='auto')
    return im

ax1, ax2, ax3 = axes_near
im = near_field_data.Ex.interp(z=pos_monitor_z).real.plot(ax=ax1)
im = near_field_data.Ey.interp(z=pos_monitor_z).real.plot(ax=ax2)
im = near_field_data.Ez.interp(z=pos_monitor_z).real.plot(ax=ax3)
# ax.set_title(f'near field E{direction}')

ax1, ax2, ax3 = axes_far
im = far_fields['Ex'].real.plot(ax=ax1)
im = far_fields['Ey'].real.plot(ax=ax2)
im = far_fields['Ez'].real.plot(ax=ax3)

plt.show()
<Figure size 720x360 with 12 Axes>
../_images/notebooks_Near2Far_26_1.png

We can also use the far field data and plot the field intensity to see the focusing effect.

[14]:
intensity_far = np.squeeze(
    np.square(np.abs(far_fields['Ex'].values)) +\
    np.square(np.abs(far_fields['Ey'].values)) +\
    np.square(np.abs(far_fields['Ez'].values))
)

_, (ax1, ax2) = plt.subplots(1, 2, figsize=(10, 5))

im1 = pmesh(xs_far, ys_far, intensity_far, ax=ax1, cmap='magma')
im2 = pmesh(xs_far, ys_far, np.sqrt(intensity_far), ax=ax2, cmap='magma')

ax1.set_title('$|E(x,y)|^2$')
ax2.set_title('$|E(x,y)|$')

plt.colorbar(im1, ax=ax1)
plt.colorbar(im2, ax=ax2)
plt.show()
<Figure size 720x360 with 4 Axes>
../_images/notebooks_Near2Far_28_1.png
[ ]: