cbeam.propagator#

This page contains the documentation for the propagator submodule, which handles all things related to characterizing a Waveguide and propagating light through it.

Propagator class#

class cbeam.propagator.Propagator(wl, wvg: None | Waveguide = None, Nmax=None, save_dir=None)#

class for coupled mode propagation of tapered waveguides

__init__(wl, wvg: None | Waveguide = None, Nmax=None, save_dir=None)#

ARGS: wl: propagation wavelength wvg: waveguide to propagate through Nmax: the number of propagating modes throughout the waveguide save_dir: directory path to save output computations. default is ‘./data/’

apply_phase(u, z, zi=None)#

apply \(e^{i \beta_j z}\) phase variation to the mode amplitude of eigenmode j with propagation constant \(\beta_j\).

Parameters:
  • u – mode amplitude vector

  • z – z coord. that u is evaluated at

  • zi (opt.)

backpropagate(u0, zf=None, zi=None)#

propagate a wavefront from the back of the waveguide to the front

characterize(zi=None, zf=None, mesh=None, tag='', save=False)#

compute the modes and coupling coefficients of the currently loaded waveguide, and set things up for propagation. for reference, this function calls compute_modes() and then compute_cmat().

Parameters:
  • zi – initial z coordinate, default 0

  • zf – final z coordinate. if None, use the waveguide length

  • mesh – an initial mesh (z=zi) if you don’t want to use the auto-generated, default mesh

  • tag – a string identifier which will be attached to any saved files

  • save – set True to save function output (z values, effective indices, modes)

Returns:

a tuple containing:

  • zs: array of z coordinates

  • neffs: array of effective indices computed on zs; 0th axis is “z axis”

  • vs: eigenmodes computed on zs; 0th axis is “z axis”

  • cmats: coupling coefficients computed on zs; 0th axis is “z axis”

Return type:

(tuple)

compute_change_of_basis(newbasis, z=None, u=None)#

compute the (N x N) change of basis matrix between the current N-dimensional eigenbasis at z and a new basis

Parameters:
  • newbasis – MxN array of N eigenmodes computed over M mesh points, which we want to expand in

  • z – z coordinate at which the currently loaded eigenbasis should be evaluated

  • u – (optional) Nx1 modal vector to express in new basis

Returns:

Nmax x Nmax change of basis matrix and the vector u in the new basis, if u was provided.

Return type:

cob, _u

compute_cmats(zs=None, vs=None, mesh=None, tag='', save=False)#

compute the coupling coefficient matrices.

Parameters:
  • zs – array of z values for the mode array (next arg). if None, use self.zs

  • vs – array of eigenmodes. if None, use self.vs

  • mesh – finite element mesh of the modes. if None, use self.mesh

  • tag – a string identifier which will be attached to any saved files

  • save – set True to save function output (z values, effective indices, modes)

Returns:

coupling coefficient matrices calculated at each z value

Return type:

cmats

compute_isolated_basis(z=None)#

compute the eigenbasis corresponding to “isolated” channels of the waveguide. this only makes sense for waveguides such as PhotonicLantern, Dicoupler, and Tricoupler.

Parameters:

z – z coordinate of eigenbasis. if None, this defaults to the end of the waveguide.

Returns:

the new basis

Return type:

(array)

compute_modes(zi=None, zf=None, mesh=None, tag='', save=False)#

compute the modes through the waveguide, using an adaptive stepping scheme. this requires greater accuracy in z than get_neffs() since mode shapes may change rapidly even if eigenvalues do not. stores initial mesh to self.mesh

Parameters:
  • zi – initial z coordinate, if None zi=0

  • zf – final z coordinate. if None, zf is the waveguide length

  • mesh – an initial mesh (z=zi) if you don’t want to use the auto-generated, default mesh

  • tag – a string identifier which will be attached to any saved files

  • save – set True to save function output (z values, effective indices, modes)

  • init_modes – an initial mode basis to use for the computation

Returns:

array of z coordinates neffs: array of effective indices computed on zs vs:modes

Return type:

zs

compute_neffs(zi=0, zf=None, mesh=None, tag='', save=False)#

compute the effective refractive indices through a waveguide, using an adaptive step scheme. also saves interpolation functions to self.neffs_funcs

Parameters:
  • zi – initial z coordinate, default 0

  • zf – final z coordinate. if None, use waveguide’s length

  • mesh – a mesh object if you don’t want to use the auto-generated, default mesh

  • tag – string identifier to attach to saved files for computation results (if saving)

  • save – whether or not to save the computation results

Returns:

a tuple containing:
  • zs: array of z coordinates

  • neffs: array of effective indices computed on zs

Return type:

(tuple)

compute_transfer_matrix(channel_basis=True, zi=None, zf=None)#

compute the transfer matrix M corresponding to propagation through the waveguide. propagation of a mode vector v is equivalent to Mv.

Parameters:
  • channel_basis (bool) – set True to force the output basis of this matrix to

  • basis. (be the waveguide channel)

  • zi (float or None) – initial z value of propagation; if None, use self.zs[0]

  • zf (float or None) – ifnal z value of propagation; if None, use self.zs[-1]

Returns:

an Nmax x Nmax complex-valued transfer matrix.

Return type:

(array)

correct_degeneracy(group, v, _v, q=None)#

used least-squares minimization to transform modes _v so that they match v. mutates _v

Parameters:
  • group – a list of mode indexes setting which modes are to be considered degenerate.

  • v – the mode basis we want to match

  • _v – the mode basis we want to transform. We are allowed to apply a unitary transformation over the modes specified by group.

  • q (opt.) – if not None, q is taken as the transformation matrix. if None, q is computed through lsq min.

Returns:

a tuple containing:
  • v: the original mode basis

  • _v: the transformed mode basis

  • q: the change-of-basis matrix

Return type:

(tuple)

generate_mesh(writeto=None)#

generate a mesh for the loaded waveguide according to class attributes.

Parameters:

writeto (str or None) – set to a str (the filename) to write the mesh to file.

get_cmat(z)#

using interpolation, compute the cross-coupling matrix at z

get_int_neff(z)#

compute the antiderivative of the mode effective indices at z

get_neff(z)#

using interpolation, compute the array of mode effective indices at z

inner_product(v1, v2, B)#

compute the inner product (or overlap integral) between fields v1 & v2, which are on the same mesh. this requires mesh information, either by passing in the mesh explicitly or by passing in the B matrix. this is the matrix on the RHS of the generalized eigenvalue problem and can be calculated e.g. with wavesolve.FEsolver.construct_B()

the inner product is v1.T * B * v2 , where * is matrix mul.

load(tag='')#

load the z values, effective indices, mode profiles, coupling coefficients, and mesh points saved to files specified by <tag>.

load_init_conds(init_prop: Propagator, z=None)#

load an eigenmode basis from init_prop as the initial basis for this propagator’s calculations.

Parameters:
  • init_prop (Propagator) – the Propagator object to load modes from

  • z (float,None) – the z value of the modes to load from init_prop. if None,

  • . (use the last z value in init_prop.zs)

make_field(mode_amps, z=None, plot=False, apply_phase=True)#

construct the finite element field corresponding to the modal vector u and the eigenbasis at z.

Parameters:
  • mode_amps – the modal vector expressing the field with \(e^{ieta_j z}\) phase oscillation factored out

  • z – the z coordinate corresponding to mode_amps

  • plot (opt.) – set True to auto-plot the norm of the field.

  • apply_phase (bool) – whether the \(e^{i\beta_j z}\) phase factor needs to be applied to mode_amps; default True.

Returns:

the field evaluated on the mesh nodes.

Return type:

(array)

make_interp_funcs(zs=None, make_neff=True, make_cmat=True, make_v=True)#
construct interpolation functions for coupling matrices and mode effective indices,

loaded into self.cmats and self.neffs, which were computed on an array of z values self.zs.

Parameters:
  • zs (None or array) – the z values to use for interpolation. if none, use self.zs

  • make_neff (bool) – set True to make effective interpolation function

  • make_cmat (bool) – set True to make coupling matrix interpolation function

  • make_v (bool) – set True to make eigenmode interpolation function

make_mesh_at_z(z)#

make the mesh corresponding to the waveguide cross-section at z.

Parameters:

z – z coordinate for the mesh.

Returns:

the corresponding mesh at z.

Return type:

(meshio mesh)

make_mode_vector(field, z=None, mesh=None)#

from a field, make a complex mode amplitude vector by decomposing the field into the currently loaded basis. kind of like the opposite of make_field().

Parameters:
  • field – the finite element field to decompose into mode amplitudes

  • z – the z coordinate of field; default is the first value of self.zs, or 0 if self.zs is not found

make_sign_consistent(v, _v)#

alter the eigenmodes _v, assumed to be similar to v, so that they have consistent overall sign with v.

plot_cfield(field, z=None, mesh=None, fig=None, ax=None, show_mesh=False, res=1.0, xlim=None, ylim=None)#

plot a complex-valued field evaluated a finite element mesh. this function is a little slow because it resamples <field> onto a grid.

Parameters:
  • field – the complex field to be plotted

  • z – the z value along the waveguide corresponding to field. alternatively, explicitly pass the next arg

  • mesh – the finite element mesh corresponding to field

  • fig – matplotlib figure object to receive the plot; if None, one will be made

  • ax – matplotlib figure axis object; if None, one will be made

  • show_mesh – shows the edges of the FE mesh in the plot

  • res – the field will be evaluated on a grid with side length res.

  • xlim (tuple) – (xmin,xmax) values for the plot - useful if you want to zoom in on the field

  • ylim (tuple) – (ymin,ymax) values for the plot.

plot_coupling_coeffs(legend=True)#

plot coupling coefficient matrix vs z values.

plot_field(field, z=None, mesh=None, ax=None, show_mesh=False)#

plot a real-valued finite element field, evaluated on the points of the mesh corresponding to the waveguide cross-section at z.

Parameters:
  • field – the finite element field to plot, e.g. the output of make_field() or an eigenmode (a row of self.vs)

  • z – the z coordinate corresponding to field

  • mesh (opt.) – the mesh object defining the points that field is evaluated on. If None, a mesh will be auto-generated

  • ax (opt.) – a matplotlib axis where the plot should go. if not None, you will need to call matplotlib.pyplot.show() manually

  • show_mesh (opt.) – whether or not to draw the mesh in the plot

plot_mode_powers(zs, us)#

plot mode powers against \(z\).

Parameters:
  • zs – array of \(z\) values.

  • us – array of mode amplitudes, e.g. from propagate(). the first axis corresponds to z.

plot_neff_diffs(yscale='log')#

plot the difference between the effective index of each mode and that of the fundamental.

Parameters:
  • yscale (str) – either “lin” for linear scale or “log” for logarithmic scale. this affects only the

  • axis. (y)

plot_neffs()#

plot the effective indices of the eigenmodes

plot_wavefront(zs, us, zi=0, fig=None, ax=None)#

plot the complex-valued wavefront through the waveguide. there may be some graphical glitches. mesh lines will be visible (issue with matplotlib’s tripcolor).

Parameters:
  • zs – array of z values

  • us – array of mode amplitudes corresponding to the field at each z value

  • zi – initial z value for the plot

  • fig – a matplotlib figure; if None, one will be made

  • ax – a matplotlib axis; if None, one will be made

Returns:

the slider object. a reference to the slider needs to be kept to prevent garbage collection from removing it.

Return type:

(matplotlib.widgets.slider)

plot_waveguide_mode(i, zi=0, fig=None, ax=None)#

plot a real-valued eigenmode of the waveguide, from modes saved in self.vs. this plot comes with a slider which controls the z value.

Parameters:
  • i – the index of the mode to be plotted (mode 0,1,2…)

  • zi – starting z value for the plot

  • fig – matplotlib figure object; if None, one will be made

  • ax – matplotlib axis objects; if None, one will be made

Returns:

the slider object. a reference to the slider needs to be kept to prevent garbage collection from removing it.

Return type:

(matplotlib.widgets.slider)

propagate(u0, zi=None, zf=None)#

propagate a launch wavefront, expressed in the basis of initial eigenmodes, to z = zf

Parameters:
  • u0 – the launch field, expressed as mode amplitudes of the initial eigenmode basis.

  • zi – the initial z coordinate corresponding to u0. if None, use the initial z value used in characterize().

  • zf – the final z coordinate to propagate through to. if None, use the final z value used in characterize().

Returns:

a tuple containing:
  • zs: the array of z values used by the ODE solver

  • u: the mode amplitudes of the wavefront, evaluated along za, with \(e^{i\beta_j z}\) phase factored out.

  • uf: the final mode amplitudes of the wavefront, with \(e^{i\beta_j z}\) phase factored in.

Return type:

(tuple)

solve_at(z=0, mesh=None)#

solve for waveguide modes at given z value. returns effective indices and eigenmodes. also stores the finite element mesh used to self.mesh

Parameters:
  • z (float) – z coordinate along the waveguide for the solve, default 0

  • mesh – FE mesh to solve on. if None, we make a mesh corresponding to the given z

Returns:

a tuple containing:

  • neff: an array of effective indices

  • v: an array of eigenmodes, dimensions NxM for N modes and M mesh points

Return type:

(tuple)

swap_modes(w, _w, _v)#

permute _w so that it matches w as close as possible. permute _v in the same way.

to_channel_basis(uf, z=None)#

convert mode amplitude vector to basis of channel eigenmodes. this only makes sense if the waveguide has defined output channels (i.e. the waveguide must have the isolate() function).

Parameters:
  • uf – mode amplitude vector (phase factored in)

  • z – z coordinate. if None, z is set to the last value of self.zs

Returns:

the vector of mode amplitudes in output channel basis.

Return type:

(vector)

WKB = False#

whether to add the WKB-like correction to the coupled mode equations (usually negligible)

Type:

bool

allow_swaps = True#

whether the propagator is allowed to swap modes (to track them through eigenvalue crossings)

Type:

bool

cmat_correction_mode = 'from_interp'#

mode for coupling matrix calculation, “from_wvg” or “from_interp”

Type:

str

fixed_zstep = None#

set to a numeric value to use a fixed zstep

Type:

float or None

init_zstep = 10.0#

starting z step

Type:

float

max_zstep = inf#

maximum z step value when computing modes or effective indices

Type:

float

min_zstep = 0.625#

minimum z step value when computing modes.

Type:

float

min_zstep_neff = 10.0#

minimum z step value when computing effective indices

Type:

float

skipped_modes = []#

used to zero certain modes during calculations. useful if there’s a cladding mode which is behaving erratically and slowing down the adaptive stepping.

Type:

list

solver = 'RK45'#

solving method for scipy’s solve_ivp, default “RK45”

Type:

str

z_acc = 0.0#

controls adaptive stepping. higher = more steps (more accurate). this is base 10 logarithmic. if you want to reduce tolerance by 10x, set to 1, etc.

Type:

float

ChainPropagator class#

class cbeam.propagator.ChainPropagator(propagators: list[Propagator])#

a ChainPropagator is a series of Propagators connected end-to-end.

__init__(propagators: list[Propagator])#

ARGS: wl: propagation wavelength wvg: waveguide to propagate through Nmax: the number of propagating modes throughout the waveguide save_dir: directory path to save output computations. default is ‘./data/’

make_field(mode_amps, z, plot=False, apply_phase=True)#

construct the finite element field corresponding to the modal vector u and the eigenbasis at z.

Parameters:
  • mode_amps – the modal vector expressing the field with \(e^{ieta_j z}\) phase oscillation factored out

  • z – the z coordinate corresponding to mode_amps

  • plot (opt.) – set True to auto-plot the norm of the field.

  • apply_phase (bool) – whether the \(e^{i\beta_j z}\) phase factor needs to be applied to mode_amps; default True.

Returns:

the field evaluated on the mesh nodes.

Return type:

(array)

propagate(u0, zi=None, zf=None)#

propagate a launch wavefront, expressed in the basis of initial eigenmodes, to z = zf

Parameters:
  • u0 – the launch field, expressed as mode amplitudes of the initial eigenmode basis.

  • zi – the initial z coordinate corresponding to u0. if None, use the initial z value used in characterize().

  • zf – the final z coordinate to propagate through to. if None, use the final z value used in characterize().

Returns:

a tuple containing:
  • zs: the array of z values used by the ODE solver

  • u: the mode amplitudes of the wavefront, evaluated along za, with \(e^{i\beta_j z}\) phase factored out.

  • uf: the final mode amplitudes of the wavefront, with \(e^{i\beta_j z}\) phase factored in.

Return type:

(tuple)

to_channel_basis(uf, z=None)#

convert mode amplitude vector to basis of channel eigenmodes. this only makes sense if the waveguide has defined output channels (i.e. the waveguide must have the isolate() function).

Parameters:
  • uf – mode amplitude vector (phase factored in)

  • z – z coordinate. if None, z is set to the last value of self.zs

Returns:

the vector of mode amplitudes in output channel basis.

Return type:

(vector)