fdtd package

backend module

Selects the backend for the fdtd-package.

The fdtd library allows to choose a backend. The numpy backend is the default one, but there are also several additional PyTorch backends:

  • numpy (defaults to float64 arrays)

  • torch (defaults to float64 tensors)

  • torch.float32

  • torch.float64

  • torch.cuda (defaults to float64 tensors)

  • torch.cuda.float32

  • torch.cuda.float64

For example, this is how to choose the “torch” backend:

fdtd.set_backend("torch")

In general, the numpy backend is preferred for standard CPU calculations with “float64” precision. In general, float64 precision is always preferred over float32 for FDTD simulations, however, float32 might give a significant performance boost.

The cuda backends are only available for computers with a GPU.

class fdtd.backend.Backend[source]

Bases: object

Backend Base Class

pi = 3.141592653589793
class fdtd.backend.NumpyBackend[source]

Bases: fdtd.backend.Backend

Numpy Backend

static arange([start, ]stop, [step, ]dtype=None, *, like=None)

create a range of values

static array()

create an array from an array-like sequence

static asarray(a, dtype=None, order=None, *, like=None)

create an array

static bmm(arr1, arr2)[source]

batch matrix multiply two arrays

static broadcast_arrays(*args, subok=False)

broadcast arrays

static broadcast_to(array, shape, subok=False)

broadcast array into shape

cos = <ufunc 'cos'>

cosine of all elements in array

divide = <ufunc 'true_divide'>
exp = <ufunc 'exp'>

exponential of all elements in array

static fft(a, n=None, axis=- 1, norm=None)

Compute the one-dimensional discrete Fourier Transform.

This function computes the one-dimensional n-point discrete Fourier Transform (DFT) with the efficient Fast Fourier Transform (FFT) algorithm [CT].

Parameters
  • a (array_like) – Input array, can be complex.

  • n (int, optional) – Length of the transformed axis of the output. If n is smaller than the length of the input, the input is cropped. If it is larger, the input is padded with zeros. If n is not given, the length of the input along the axis specified by axis is used.

  • axis (int, optional) – Axis over which to compute the FFT. If not given, the last axis is used.

  • norm ({"backward", "ortho", "forward"}, optional) –

    New in version 1.10.0.

    Normalization mode (see numpy.fft). Default is “backward”. Indicates which direction of the forward/backward pair of transforms is scaled and with what normalization factor.

    New in version 1.20.0: The “backward”, “forward” values were added.

Returns

out – The truncated or zero-padded input, transformed along the axis indicated by axis, or the last one if axis is not specified.

Return type

complex ndarray

Raises

IndexError – If axis is not a valid axis of a.

See also

numpy.fft

for definition of the DFT and conventions used.

ifft

The inverse of fft.

fft2

The two-dimensional FFT.

fftn

The n-dimensional FFT.

rfftn

The n-dimensional FFT of real input.

fftfreq

Frequency bins for given FFT parameters.

Notes

FFT (Fast Fourier Transform) refers to a way the discrete Fourier Transform (DFT) can be calculated efficiently, by using symmetries in the calculated terms. The symmetry is highest when n is a power of 2, and the transform is therefore most efficient for these sizes.

The DFT is defined, with the conventions used in this implementation, in the documentation for the numpy.fft module.

References

CT

Cooley, James W., and John W. Tukey, 1965, “An algorithm for the machine calculation of complex Fourier series,” Math. Comput. 19: 297-301.

Examples

>>> np.fft.fft(np.exp(2j * np.pi * np.arange(8) / 8))
array([-2.33486982e-16+1.14423775e-17j,  8.00000000e+00-1.25557246e-15j,
        2.33486982e-16+2.33486982e-16j,  0.00000000e+00+1.22464680e-16j,
       -1.14423775e-17+2.33486982e-16j,  0.00000000e+00+5.20784380e-16j,
        1.14423775e-17+1.14423775e-17j,  0.00000000e+00+1.22464680e-16j])

In this example, real input has an FFT which is Hermitian, i.e., symmetric in the real part and anti-symmetric in the imaginary part, as described in the numpy.fft documentation:

>>> import matplotlib.pyplot as plt
>>> t = np.arange(256)
>>> sp = np.fft.fft(np.sin(t))
>>> freq = np.fft.fftfreq(t.shape[-1])
>>> plt.plot(freq, sp.real, freq, sp.imag)
[<matplotlib.lines.Line2D object at 0x...>, <matplotlib.lines.Line2D object at 0x...>]
>>> plt.show()
static fftfreq(n, d=1.0)

Return the Discrete Fourier Transform sample frequencies.

The returned float array f contains the frequency bin centers in cycles per unit of the sample spacing (with zero at the start). For instance, if the sample spacing is in seconds, then the frequency unit is cycles/second.

Given a window length n and a sample spacing d:

f = [0, 1, ...,   n/2-1,     -n/2, ..., -1] / (d*n)   if n is even
f = [0, 1, ..., (n-1)/2, -(n-1)/2, ..., -1] / (d*n)   if n is odd
Parameters
  • n (int) – Window length.

  • d (scalar, optional) – Sample spacing (inverse of the sampling rate). Defaults to 1.

Returns

f – Array of length n containing the sample frequencies.

Return type

ndarray

Examples

>>> signal = np.array([-2, 8, 6, 4, 1, 0, 3, 5], dtype=float)
>>> fourier = np.fft.fft(signal)
>>> n = signal.size
>>> timestep = 0.1
>>> freq = np.fft.fftfreq(n, d=timestep)
>>> freq
array([ 0.  ,  1.25,  2.5 , ..., -3.75, -2.5 , -1.25])
float

alias of numpy.float64

int

alias of numpy.int64

static is_array(arr)[source]

check if an object is an array

static linspace(start, stop, num=50, endpoint=True, retstep=False, dtype=None, axis=0)

create a linearly spaced array between two points

static max(a, axis=None, out=None, keepdims=<no value>, initial=<no value>, where=<no value>)

max element in array

static numpy()

convert the array to numpy array

static ones(shape, dtype=None, order='C', *, like=None)

create an array filled with ones

static pad(array, pad_width, mode='constant', **kwargs)

Pad an array.

Parameters
  • array (array_like of rank N) – The array to pad.

  • pad_width ({sequence, array_like, int}) – Number of values padded to the edges of each axis. ((before_1, after_1), … (before_N, after_N)) unique pad widths for each axis. ((before, after),) yields same before and after pad for each axis. (pad,) or int is a shortcut for before = after = pad width for all axes.

  • mode (str or function, optional) –

    One of the following string values or a user supplied function.

    ’constant’ (default)

    Pads with a constant value.

    ’edge’

    Pads with the edge values of array.

    ’linear_ramp’

    Pads with the linear ramp between end_value and the array edge value.

    ’maximum’

    Pads with the maximum value of all or part of the vector along each axis.

    ’mean’

    Pads with the mean value of all or part of the vector along each axis.

    ’median’

    Pads with the median value of all or part of the vector along each axis.

    ’minimum’

    Pads with the minimum value of all or part of the vector along each axis.

    ’reflect’

    Pads with the reflection of the vector mirrored on the first and last values of the vector along each axis.

    ’symmetric’

    Pads with the reflection of the vector mirrored along the edge of the array.

    ’wrap’

    Pads with the wrap of the vector along the axis. The first values are used to pad the end and the end values are used to pad the beginning.

    ’empty’

    Pads with undefined values.

    New in version 1.17.

    <function>

    Padding function, see Notes.

  • stat_length (sequence or int, optional) –

    Used in ‘maximum’, ‘mean’, ‘median’, and ‘minimum’. Number of values at edge of each axis used to calculate the statistic value.

    ((before_1, after_1), … (before_N, after_N)) unique statistic lengths for each axis.

    ((before, after),) yields same before and after statistic lengths for each axis.

    (stat_length,) or int is a shortcut for before = after = statistic length for all axes.

    Default is None, to use the entire axis.

  • constant_values (sequence or scalar, optional) –

    Used in ‘constant’. The values to set the padded values for each axis.

    ((before_1, after_1), ... (before_N, after_N)) unique pad constants for each axis.

    ((before, after),) yields same before and after constants for each axis.

    (constant,) or constant is a shortcut for before = after = constant for all axes.

    Default is 0.

  • end_values (sequence or scalar, optional) –

    Used in ‘linear_ramp’. The values used for the ending value of the linear_ramp and that will form the edge of the padded array.

    ((before_1, after_1), ... (before_N, after_N)) unique end values for each axis.

    ((before, after),) yields same before and after end values for each axis.

    (constant,) or constant is a shortcut for before = after = constant for all axes.

    Default is 0.

  • reflect_type ({'even', 'odd'}, optional) – Used in ‘reflect’, and ‘symmetric’. The ‘even’ style is the default with an unaltered reflection around the edge value. For the ‘odd’ style, the extended part of the array is created by subtracting the reflected values from two times the edge value.

Returns

pad – Padded array of rank equal to array with shape increased according to pad_width.

Return type

ndarray

Notes

New in version 1.7.0.

For an array with rank greater than 1, some of the padding of later axes is calculated from padding of previous axes. This is easiest to think about with a rank 2 array where the corners of the padded array are calculated by using padded values from the first axis.

The padding function, if used, should modify a rank 1 array in-place. It has the following signature:

padding_func(vector, iaxis_pad_width, iaxis, kwargs)

where

vectorndarray

A rank 1 array already padded with zeros. Padded values are vector[:iaxis_pad_width[0]] and vector[-iaxis_pad_width[1]:].

iaxis_pad_widthtuple

A 2-tuple of ints, iaxis_pad_width[0] represents the number of values padded at the beginning of vector where iaxis_pad_width[1] represents the number of values padded at the end of vector.

iaxisint

The axis currently being calculated.

kwargsdict

Any keyword arguments the function requires.

Examples

>>> a = [1, 2, 3, 4, 5]
>>> np.pad(a, (2, 3), 'constant', constant_values=(4, 6))
array([4, 4, 1, ..., 6, 6, 6])
>>> np.pad(a, (2, 3), 'edge')
array([1, 1, 1, ..., 5, 5, 5])
>>> np.pad(a, (2, 3), 'linear_ramp', end_values=(5, -4))
array([ 5,  3,  1,  2,  3,  4,  5,  2, -1, -4])
>>> np.pad(a, (2,), 'maximum')
array([5, 5, 1, 2, 3, 4, 5, 5, 5])
>>> np.pad(a, (2,), 'mean')
array([3, 3, 1, 2, 3, 4, 5, 3, 3])
>>> np.pad(a, (2,), 'median')
array([3, 3, 1, 2, 3, 4, 5, 3, 3])
>>> a = [[1, 2], [3, 4]]
>>> np.pad(a, ((3, 2), (2, 3)), 'minimum')
array([[1, 1, 1, 2, 1, 1, 1],
       [1, 1, 1, 2, 1, 1, 1],
       [1, 1, 1, 2, 1, 1, 1],
       [1, 1, 1, 2, 1, 1, 1],
       [3, 3, 3, 4, 3, 3, 3],
       [1, 1, 1, 2, 1, 1, 1],
       [1, 1, 1, 2, 1, 1, 1]])
>>> a = [1, 2, 3, 4, 5]
>>> np.pad(a, (2, 3), 'reflect')
array([3, 2, 1, 2, 3, 4, 5, 4, 3, 2])
>>> np.pad(a, (2, 3), 'reflect', reflect_type='odd')
array([-1,  0,  1,  2,  3,  4,  5,  6,  7,  8])
>>> np.pad(a, (2, 3), 'symmetric')
array([2, 1, 1, 2, 3, 4, 5, 5, 4, 3])
>>> np.pad(a, (2, 3), 'symmetric', reflect_type='odd')
array([0, 1, 1, 2, 3, 4, 5, 5, 6, 7])
>>> np.pad(a, (2, 3), 'wrap')
array([4, 5, 1, 2, 3, 4, 5, 1, 2, 3])
>>> def pad_with(vector, pad_width, iaxis, kwargs):
...     pad_value = kwargs.get('padder', 10)
...     vector[:pad_width[0]] = pad_value
...     vector[-pad_width[1]:] = pad_value
>>> a = np.arange(6)
>>> a = a.reshape((2, 3))
>>> np.pad(a, 2, pad_with)
array([[10, 10, 10, 10, 10, 10, 10],
       [10, 10, 10, 10, 10, 10, 10],
       [10, 10,  0,  1,  2, 10, 10],
       [10, 10,  3,  4,  5, 10, 10],
       [10, 10, 10, 10, 10, 10, 10],
       [10, 10, 10, 10, 10, 10, 10]])
>>> np.pad(a, 2, pad_with, padder=100)
array([[100, 100, 100, 100, 100, 100, 100],
       [100, 100, 100, 100, 100, 100, 100],
       [100, 100,   0,   1,   2, 100, 100],
       [100, 100,   3,   4,   5, 100, 100],
       [100, 100, 100, 100, 100, 100, 100],
       [100, 100, 100, 100, 100, 100, 100]])
static reshape(a, newshape, order='C')

reshape array into given shape

sin = <ufunc 'sin'>

sine of all elements in array

static squeeze(a, axis=None)

remove dim-1 dimensions

static stack(arrays, axis=0, out=None)

stack multiple arrays

static sum(a, axis=None, dtype=None, out=None, keepdims=<no value>, initial=<no value>, where=<no value>)

sum elements in array

static transpose(a, axes=None)

transpose array by flipping two dimensions

static zeros(shape, dtype=float, order='C', *, like=None)

create an array filled with zeros

static zeros_like(a, dtype=None, order='K', subok=True, shape=None)

create an array filled with zeros

fdtd.backend.set_backend(name: str)[source]

Set the backend for the FDTD simulations

This function monkeypatches the backend object by changing its class. This way, all methods of the backend object will be replaced.

Parameters

name – name of the backend. Allowed backend names: - numpy (defaults to float64 arrays) - torch (defaults to float64 tensors) - torch.float32 - torch.float64 - torch.cuda (defaults to float64 tensors) - torch.cuda.float32 - torch.cuda.float64

boundaries module

Boundaries for the FDTD Grid.

Available Boundaries:

  • PeriodicBoundary

  • PML

class fdtd.boundaries.Boundary(name: Optional[str] = None)[source]

Bases: object

an FDTD Boundary [base class]

__init__(name: Optional[str] = None)[source]

Create a boundary

Parameters

name – name of the boundary

update_E()[source]

Update electric field of the grid

Note

this method is called after the grid fields are updated

update_H()[source]

Update magnetic field of the grid

Note

this method is called after the grid fields are updated

update_phi_E()[source]

Update convolution [phi_E]

Note

this method is called before the electric field is updated

update_phi_H()[source]

Update convolution [phi_H]

Note

this method is called before the magnetic field is updated

fdtd.boundaries.DomainBorderPML(grid, border_cells=5)[source]

Some problem setups require a layer of PML all the way around the problem. This is a convenience function to add such a layer to an existing grid. Caution: Alters grid in-place.

class fdtd.boundaries.PML(a: float = 1e-08, name: Optional[str] = None)[source]

Bases: fdtd.boundaries.Boundary

A perfectly matched layer (PML)

a PML is an impedence-matched area at the boundary of the grid for which all fields incident perpendicular to the area are absorbed without reflection.

Note

Registering a PML to the grid will monkeypatch the PML to become one of its subclasses: _PMLXlow, _PMLYlow or _PMLZlow, _PMLXhigh, _PMLYhigh, _PMLZhigh depending on the position in the grid.

__init__(a: float = 1e-08, name: Optional[str] = None)[source]

Perfectly Matched Layer

Parameters
  • a – stability parameter

  • name – name of the PML

update_E()[source]

Update electric field of the grid

Note

this method is called after the electric field is updated

update_H()[source]

Update magnetic field of the grid

Note

this method is called after the magnetic field is updated

update_phi_E()[source]

Update convolution [phi_E]

Note

this method is called before the electric field is updated

update_phi_H()[source]

Update convolution [phi_H]

Note

this method is called before the magnetic field is updated

class fdtd.boundaries.PeriodicBoundary(name: Optional[str] = None)[source]

Bases: fdtd.boundaries.Boundary

An FDTD Periodic Boundary

Note

Registering a periodic boundary to the grid will change the periodic boundary in one of its subclasses: _PeriodicBoundaryX, _PeriodicBoundaryY or _PeriodicBoundaryY, depending on the position in the grid.

detectors module

Detectors for the FDTD Grid.

Available Detectors:

  • LineDetector

class fdtd.detectors.BlockDetector(name=None)[source]

Bases: object

A detector along a block in the FDTD grid

__init__(name=None)[source]

Create a block detector

Parameters

name – name of the Detector

detect_E()[source]

detect the electric field at a certain location in the grid

detect_H()[source]

detect the magnetic field at a certain location in the grid

detector_values()[source]

outputs what detector detects

class fdtd.detectors.CurrentDetector(name=None)[source]

Bases: object

A current detector.

__init__(name=None)[source]

Create a block detector

Parameters

name – name of the Detector

detect_E()[source]

detect the electric field at a certain location in the grid

detect_H()[source]
detector_values()[source]

outputs what detector detects

single_point_current(px, py, pz)[source]

Only Z-polarized for now. Can probably do a cross product to get arbitrary polarizations


X—->

TODO: FIXME: IMPORTANT: material magnetic permeability? find test cases!

Implements the first correction from [Fang 1994] (two cells are spatially averaged to account for Yee cell half-step inaccuracies), but not the second one (minor loss of accuracy).

Jiayuan Fang, Danwei Xue. Precautions in the calculation of impedance in FDTD computations. Proceedings of IEEE Antennas and Propagation Society International Symposium and URSI National Radio Science Meeting, vol. 3, 1994, p. 1814–7 vol.3. https://doi.org/10.1109/APS.1994.408185.

Luebbers RJ, Langdon HS. A simple feed model that reduces time steps needed for FDTD antenna and microstrip calculations. IEEE Trans Antennas Propagat 1996;44:1000–5. https://doi.org/10.1109/8.504308.

class fdtd.detectors.LineDetector(name=None)[source]

Bases: object

A detector along a line in the FDTD grid

__init__(name=None)[source]

Create a line detector

Parameters

name – name of the Detector

detect_E()[source]

detect the electric field at a certain location in the grid

detect_H()[source]

detect the magnetic field at a certain location in the grid

detector_values()[source]

outputs what detector detects

grid module

The FDTD Grid

The grid is the core of the FDTD Library. It is where everything comes together and where the biggest part of the calculations are done.

class fdtd.grid.Grid(shape: Tuple[numbers.Number, numbers.Number, numbers.Number], grid_spacing: float = 1.55e-07, permittivity: float = 1.0, permeability: float = 1.0, courant_number: Optional[float] = None)[source]

Bases: object

The FDTD Grid

The grid is the core of the FDTD Library. It is where everything comes together and where the biggest part of the calculations are done.

__init__(shape: Tuple[numbers.Number, numbers.Number, numbers.Number], grid_spacing: float = 1.55e-07, permittivity: float = 1.0, permeability: float = 1.0, courant_number: Optional[float] = None)[source]
Parameters
  • shape – shape of the FDTD grid.

  • grid_spacing – distance between the grid cells.

  • permittivity – the relative permittivity of the background.

  • permeability – the relative permeability of the background.

  • courant_number – the courant number of the FDTD simulation. Defaults to the inverse of the square root of the number of dimensions > 1 (optimal value). The timestep of the simulation will be derived from this number using the CFL-condition.

add_boundary(name, boundary)[source]

add a boundary to the grid

add_detector(name, detector)[source]

add a detector to the grid

add_object(name, obj)[source]

add an object to the grid

add_source(name, source)[source]

add a source to the grid

generate_video(delete_frames=False)[source]

Compiles frames into a video

These framed should be saved through fdtd.Grid.visualize(save=True) while having fdtd.Grid.save_simulation() enabled.

Parameters

delete_frames (optional, bool) – delete stored frames after conversion to video.

Returns

the filename of the generated video.

Note

this function requires ffmpeg to be available in your path.

reset()[source]

reset the grid by setting all fields to zero

run(total_time: numbers.Number, progress_bar: bool = True)[source]

run an FDTD simulation.

Parameters
  • total_time – the total time for the simulation to run.

  • progress_bar – choose to show a progress bar during simulation

save_data()[source]

Saves readings from all detectors in the grid into a numpy zip file. Each detector is stored in separate arrays. Electric and magnetic field field readings of each detector are also stored separately with suffix ” (E)” and ” (H)” (Example: [‘detector0 (E)’, ‘detector0 (H)’]). Therefore, the numpy zip file contains arrays twice the number of detectors. REQUIRES ‘fdtd.Grid.save_simulation()’ to be run before this function.

Parameters: None

save_simulation(sim_name=None)[source]

Creates a folder and initializes environment to store simulation or related details. saveSimulation() needs to be run before running any function that stores data (generate_video(), save_data()).

Parameters:-

(optional) sim_name (string): Preferred name for simulation

property shape: Tuple[int, int, int]

get the shape of the FDTD grid

step()[source]

do a single FDTD step by first updating the electric field and then updating the magnetic field

property time_passed: float

get the total time passed

update_E()[source]

update the electric field by using the curl of the magnetic field

update_H()[source]

update the magnetic field by using the curl of the electric field

visualize(x=None, y=None, z=None, cmap='Blues', pbcolor='C3', pmlcolor=(0, 0, 0, 0.1), objcolor=(1, 0, 0, 0.1), srccolor='C0', detcolor='C2', norm='linear', show=False, animate=False, index=None, save=False, folder=None)

visualize a projection of the grid and the optical energy inside the grid

Parameters
  • x – the x-value to make the yz-projection (leave None if using different projection)

  • y – the y-value to make the zx-projection (leave None if using different projection)

  • z – the z-value to make the xy-projection (leave None if using different projection)

  • cmap – the colormap to visualize the energy in the grid

  • pbcolor – the color to visualize the periodic boundaries

  • pmlcolor – the color to visualize the PML

  • objcolor – the color to visualize the objects in the grid

  • srccolor – the color to visualize the sources in the grid

  • detcolor – the color to visualize the detectors in the grid

  • norm – how to normalize the grid_energy color map (‘linear’ or ‘log’).

  • show – call pyplot.show() at the end of the function

  • animate – see frame by frame state of grid during simulation

  • index – index for each frame of animation (typically a loop variable is passed)

  • save – save frames in a folder

  • folder – path to folder to save frames

property x: int

get the number of grid cells in the x-direction

property y: int

get the number of grid cells in the y-direction

property z: int

get the number of grid cells in the y-direction

fdtd.grid.curl_E(E: numpy.ndarray) numpy.ndarray[source]

Transforms an E-type field into an H-type field by performing a curl operation

Parameters

E – Electric field to take the curl of (E-type field located on the edges of the grid cell [integer gridpoints])

Returns

The curl of E (H-type field located on the faces of the grid [half-integer grid points])

fdtd.grid.curl_H(H: numpy.ndarray) numpy.ndarray[source]

Transforms an H-type field into an E-type field by performing a curl operation

Parameters

H – Magnetic field to take the curl of (H-type field located on half-integer grid points)

Returns

The curl of H (E-type field located on the edges of the grid [integer grid points])

class fdtd.grid.d_[source]

Bases: object

Just a convenience function for indexing polarizations of arrays.

X = 0
Y = 1
Z = 2

objects module

The objects to place in the grid.

Objects define all the regions in the grid with a modified update equation, such as for example regions with anisotropic permittivity etc.

Available Objects:
  • Object

  • AnisotropicObject

class fdtd.objects.AbsorbingObject(permittivity: numpy.ndarray, conductivity: numpy.ndarray, name: Optional[str] = None)[source]

Bases: fdtd.objects.Object

An absorbing object takes conductivity into account

__init__(permittivity: numpy.ndarray, conductivity: numpy.ndarray, name: Optional[str] = None)[source]
Parameters
  • permittivity – permittivity tensor

  • conductivity – conductivity tensor (will introduce the loss)

  • name – name of the object (will become available as attribute to the grid)

update_E(curl_H)[source]

custom update equations for inside the absorbing object

Parameters

curl_H – the curl of magnetic field in the grid.

update_H(curl_E)[source]

custom update equations for inside the absorbing object

Parameters

curl_E – the curl of electric field in the grid.

class fdtd.objects.AnisotropicObject(permittivity: numpy.ndarray, name: Optional[str] = None)[source]

Bases: fdtd.objects.Object

An object with anisotropic permittivity tensor

update_E(curl_H)[source]

custom update equations for inside the anisotropic object

Parameters

curl_H – the curl of magnetic field in the grid.

update_H(curl_E)[source]

custom update equations for inside the anisotropic object

Parameters

curl_E – the curl of electric field in the grid.

class fdtd.objects.Object(permittivity: numpy.ndarray, name: Optional[str] = None)[source]

Bases: object

An object to place in the grid

__init__(permittivity: numpy.ndarray, name: Optional[str] = None)[source]
Parameters
  • permittivity – permittivity tensor

  • name – name of the object (will become available as attribute to the grid)

update_E(curl_H)[source]

custom update equations for inside the object

Parameters

curl_H – the curl of magnetic field in the grid.

update_H(curl_E)[source]

custom update equations for inside the object

Parameters

curl_E – the curl of electric field in the grid.

sources module

Sources are objects that inject the fields into the grid.

Available sources:

  • PointSource

  • LineSource

class fdtd.sources.LineSource(period: numbers.Number = 15, amplitude: float = 1.0, phase_shift: float = 0.0, name: Optional[str] = None, pulse: bool = False, cycle: int = 5, hanning_dt: float = 10.0)[source]

Bases: object

A source along a line in the FDTD grid

__init__(period: numbers.Number = 15, amplitude: float = 1.0, phase_shift: float = 0.0, name: Optional[str] = None, pulse: bool = False, cycle: int = 5, hanning_dt: float = 10.0)[source]

Create a LineSource with a gaussian profile

Parameters
  • period – The period of the source. The period can be specified as integer [timesteps] or as float [seconds]

  • amplitude – The amplitude of the source in simulation units

  • phase_shift – The phase offset of the source.

  • pulse – Set True to use a Hanning window pulse instead of continuous wavefunction.

  • cycle – cycles for Hanning window pulse.

  • hanning_dt – timestep used for Hanning window pulse width (optional).

update_E()[source]

Add the source to the electric field

update_H()[source]

Add the source to the magnetic field

class fdtd.sources.PlaneSource(period: numbers.Number = 15, amplitude: float = 1.0, phase_shift: float = 0.0, name: Optional[str] = None)[source]

Bases: object

A source along a plane in the FDTD grid

__init__(period: numbers.Number = 15, amplitude: float = 1.0, phase_shift: float = 0.0, name: Optional[str] = None)[source]

Create a LineSource with a gaussian profile

Parameters
  • period – The period of the source. The period can be specified as integer [timesteps] or as float [seconds]

  • amplitude – The amplitude of the source in simulation units

  • phase_shift – The phase offset of the source.

update_E()[source]

Add the source to the electric field

update_H()[source]

Add the source to the magnetic field

class fdtd.sources.PointSource(period: numbers.Number = 15, amplitude: float = 1.0, phase_shift: float = 0.0, name: Optional[str] = None, pulse: bool = False, cycle: int = 5, hanning_dt: float = 10.0)[source]

Bases: object

A source placed at a single point (grid cell) in the grid

__init__(period: numbers.Number = 15, amplitude: float = 1.0, phase_shift: float = 0.0, name: Optional[str] = None, pulse: bool = False, cycle: int = 5, hanning_dt: float = 10.0)[source]

Create a LineSource with a gaussian profile

Parameters
  • period – The period of the source. The period can be specified as integer [timesteps] or as float [seconds]

  • amplitude – The electric field amplitude in simulation units

  • phase_shift – The phase offset of the source.

  • name – name of the source.

  • pulse – Set True to use a Hanning window pulse instead of continuous wavefunction.

  • cycle – cycles for Hanning window pulse.

  • hanning_dt – timestep used for Hanning window pulse width (optional).

update_E()[source]

Add the source to the electric field

update_H()[source]

Add the source to the magnetic field

class fdtd.sources.SoftArbitraryPointSource(waveform_array: numpy.ndarray, name: Optional[str] = None, impedance: float = 0.0)[source]

Bases: object

A source placed at a single point (grid cell) in the grid. This source is special: it’s both a source and a detector.

Unlike the other sources, the input is a voltage, not an electric field. (really? why? should we convert back and forth?)

For electrical measurements I’ve only needed a single-index source, so I don’t know how the volume/line sources above work. We want the FFT function to operate over any detector. Maybe all sources should take an arbitary waveform argument?

Each index in the waveform array represents 1 value at a timestep.

There are many different geometries of “equivalent sources”. The detector/source paradigm used in /fdtd might perhaps not correspond to this in an ideal fashion.

It’s not intuitively clear to me what a “soft” source would imply in the optical case, or what impedance even means for a laser.

/fdtd/ seems to have found primary use in optical circles, so the default Z should probably be 0.

“Whilst established for microwaves and electrical circuits, this concept has only very recently been observed in the optical domain, yet is not well defined or understood.”[1]

[1]: Optical impedance of metallic nano-structures, M. Mazilu and K. Dholakia https://doi.org/10.1364/OE.14.007709

[2]: http://www.gwoptics.org/learn/02_Plane_waves/01_Fabry_Perot_cavity/02_Impedance_matched.php

//-

__init__(waveform_array: numpy.ndarray, name: Optional[str] = None, impedance: float = 0.0)[source]

Create

Parameters

waveform_array

update_E()[source]
update_H()[source]