Transfer matrix method

Example 1: Using the TMM solver to calculate the reflextion of a multilayered ARC

Example 2: Looking at the effect of substrate and the no_back_reflection option in the TMM solver

Example 3: Coherent and incoherent layers in the TMM solver

The TMM interface for the solar cell solver

This is the method actually called from the solar cell solver, serving as interface between the solar cell and the lower level TMM formalism. The Beer-Lambert calculator, the RCWA calculator and the external optics calculator (where the user simply adds the reflection and the absorption profile manually) have similar interfaces.

The TMM solver can handle different input angles and polarizations, specified in the options (specifically, options.theta and options.pol). The input angle is in degrees, while pol is ‘s’, ‘p’ or ‘u’ (computationally, ‘u’ is simply the average of ‘s’ and ‘p’ and thus requires two calculations.)

Using the default options, all layers (except very thick layers, with the layer thickness more than 10 times the maximum incident wavelength) will be considered as coherent, i.e. phase information in the wave-optical calculation is retained. However, layers can also be specified to be incoherent using a coherency_list option with one entry (either ‘c’ or ‘i’) per layer.

solcore.optics.tmm.solve_tmm(solar_cell, options)[source]

Calculates the reflection, transmission and absorption of a solar cell object using the transfer matrix method. Internally, it creates an OptiStack and then it calculates the optical properties of the whole structure. A substrate can be specified in the SolarCell object, which is treated as a semi-infinite transmission medium. Shading can also be specified (as a fraction).

Relevant options are ‘wl’ (the wavelengths, in m), the incidence angle ‘theta’ (in degrees), the polarization ‘pol’ (‘s’, ‘p’ or ‘u’), ‘position’ (locations in m at which depth-dependent absorption is calculated), ‘no_back_reflection’ and ‘BL_correction’. ‘no_back_reflection’ sets whether reflections from the back surface are suppressed (if set to True, the default), or taken into account (if set to False).

If ‘BL_correction’ is set to True, thick layers (thickness > 10*maximum wavelength) are treated incoherently using the Beer-Lambert law, to avoid the calculation of unphysical interference oscillations in the R/A/T spectra.

A coherency_list option can be provided: this should have elements equal to the total number of layers (if a Junction contains multiple Layers, each should have its own entry in the coherency list). Each element is either ‘c’ for coherent treatment of that layer or ‘i’ for incoherent treatment.

Parameters:
  • solar_cell – A SolarCell object
  • options – Options for the solver
Returns:

None

solcore.optics.tmm.absorbed(self, z)[source]
solcore.optics.tmm.calculate_absorption_tmm(tmm_out, initial=1)[source]

The OptiStack and the high level TMM optical methods

This module serves as interface between solcore structures (junctions, layers, materials…) and the transfer matrix package developed by Steven Byrnes and included in the PyPi repository.

solcore.absorption_calculator.transfer_matrix.np_cache(function)[source]

This function was taken from https://stackoverflow.com/questions/52331944/cache-decorator-for-numpy-arrays/52332109#52332109

Creates a cacheable version of a function which takes a 1D numpy array as input, by using a wrapping function which converts the array to a tuple. It returns a function which returns the same output as the input function, but can be cached, avoiding a bottleneck when optical constants in a material are looked up repeatedly.

If the input argument is hashable already - eg. a single float number - then the function is not chached and is returned “as is”.

Param:function: the function of which a cacheable version is to be created
Returns:wrapper: the cacheable version of the function
class solcore.absorption_calculator.transfer_matrix.OptiStack(structure=(), no_back_reflection=False, substrate=None, incidence=None, **kwargs)[source]

Class that contains an optical structure: a sequence of layers with a thickness and a complex refractive index.

It serves as an intermediate step between solcore layers and materials and the stack of thicknesses and and n and k values necessary to run calculations involving TMM. When creating an OptiStack object, the thicknesses of all the layers forming the Solcore structure and the optical data of the materials of the layers are extracted and arranged in such a way they can be easily and fastly read by the TMM functions.

In addition to a solcore structure with Layers, it can also take a list where each element represent a layer written as a list and contains the layer thickness and the dielectrical model, the raw n and k data as a function of wavelengths, or a whole Device structure as the type used in the PDD model.

In summary, this class acepts:

  • A solcore structure with layers
  • A list where each element is [thickness, DielectricModel]
  • A list where each element is [thickness, wavelength, n, k]
  • A list mixing the above:
    [ [thickness, DielectricModel],
    [thickness, wavelength, n, k], solcore.Layer, solcore.Layer ]

This allows for maximum flexibility when creating the optical model, allowing to construct the stack with experimental data, modelled data and known material properties from the database.

Yet anther way of defining the layers mixes experimental data with a DielectricModel within the same layer but in spectrally distinct regions. The syntaxis for the layer is:

layer = [thickness, wavelength, n, k, DielectricModel, mixing]

where mixing is a list containing three elements: [the mixing point (nm), the mixing width (nm), zero or one] depending if the mixing function should be increasing with the wavelength or decreasing. If increasing (zero), the Dielectric model will be used at long wavelengths and the experimental data at short wavelengths. If decreasing (one) the oposite is done. The mixing point and mixing width control how smooth is the transition between one and the other type of data.

Extra layers such as he semi-infinite, air-like first and last medium, and a back highly absorbing layer are included at runtime to fulfill the requirements of the TMM solver or to solve some of its limitations.

get_indices(wl)[source]

Returns the complex refractive index of the stack.

Parameters:wl – Wavelength of the light in nm.
Returns:A list with the complex refractive index of each layer, including the semi-infinite front and back

layers and, opionally, the back absorbing layer used to suppress back surface relfexion.

get_widths()[source]

Returns the widths of the layers of the stack.

Returns:A list with the widths each layer, including the semi-infinite front and back layers and, opionally,

the back absorbing layer used to suppress back surface relfexion, defined as 1 mm thick.

add_layers(layers)[source]

Generic function to add layers to the OptiStack. Internally, it calls add_solcore_layer, add_modelled_layer or add_raw_nk_layer.

Parameters:layers – A list with the layers to add (even if it is just one layer) It can be one or more and it can

mixed, Solcore-based and modelled layers. :return: None

remove_layer(idx)[source]

Removes layer with index idx from the OptiStack

Parameters:idx – Index of the layer to remove
Returns:None
swap_layers(idx1, idx2)[source]

Swaps two layers in the OptiStack.

Parameters:
  • idx1 – The index of one of the layers.
  • idx2 – The index of the other.
Returns:

None

set_widths(widths)[source]

Changes the widths of the layers in the stack.

Param:widths: a list or array of widths, length equal to the number of layers
Returns:None
solcore.absorption_calculator.transfer_matrix.calculate_rat(structure, wavelength, angle=0, pol='u', coherent=True, coherency_list=None, no_back_reflection=True, **kwargs)[source]

Calculates the reflected, absorbed and transmitted intensity of the structure for the wavelengths and angles defined.

Parameters:
  • structure – A solcore Structure object with layers and materials or a OptiStack object.
  • wavelength – Wavelengths (in nm) in which calculate the data. An array.
  • angle – Angle (in degrees) of the incident light. Default: 0 (normal incidence).
  • pol – Polarisation of the light: ‘s’, ‘p’ or ‘u’. Default: ‘u’ (unpolarised).
  • coherent – If the light is coherent or not. If not, a coherency list must be added.
  • coherency_list – A list indicating in which layers light should be treated as coeherent (‘c’) and in which

incoherent (‘i’). It needs as many elements as layers in the structure. :param no_back_reflection: If reflexion from the back must be supressed. Default=True. :return: A dictionary with the R, A and T at the specified wavelengths and angle.

solcore.absorption_calculator.transfer_matrix.calculate_ellipsometry(structure, wavelength, angle, no_back_reflection=True, **kwargs)[source]

Calculates the ellipsometric parameters psi and delta. It can only deal with coherent light and the whole stack (including back surface) is considered, so caution must be taken when comparing the simulated results with experiments where the back surface is rough or layers are thick and coherent light propagation makes no sense.

The optional argument no_back_reflection can be included to add an extra layer on the back absorbing all light that reaches that position without any reflexion, to remove the reflexion from the back surface.

Parameters:
  • structure – A solcore structure with layers and materials.
  • wavelength – Wavelengths (in nm) in which calculate the data. An array.
  • angle – A tupple or list with the angles (in degrees) in which to calculate the data.
  • no_back_reflection – If reflexion from the back must be suppressed. Default=True.
Returns:

A dictionary with psi and delta at the specified wavelengths and angles (2D arrays).

solcore.absorption_calculator.transfer_matrix.calculate_absorption_profile(structure, wavelength, z_limit=None, steps_size=2, dist=None, no_back_reflection=True, angle=0, pol='u', coherent=True, coherency_list=None, zero_threshold=1e-06, RAT_out=None, **kwargs)[source]

It calculates the absorbed energy density within the material. From the documentation:

‘In principle this has units of [power]/[volume], but we can express it as a multiple of incoming light power density on the material, which has units [power]/[area], so that absorbed energy density has units of 1/[length].’

Integrating this absorption profile in the whole stack gives the same result that the absorption obtained with calculate_rat as long as the spacial mesh (controlled by steps_thinest_layer) is fine enough. If the structure is very thick and the mesh not thin enough, the calculation might diverge at short wavelengths.

Parameters:
  • structure – A solcore structure with layers and materials.
  • wavelength – Wavelengths in which calculate the data (in nm). An array
  • z_limit – Maximum value in the z direction
  • steps_size – if the dist is not specified, the step size in nm to use in the depth-dependent calculation
  • dist – the positions (in nm) at which to calculate depth-dependent absorption
  • no_back_reflection – whether to suppress reflections from the back interface (True) or not (False)
  • angle – incidence of angle in degrees
  • pol – polarization of incident light: ‘s’, ‘p’ or ‘u’ (unpolarized)
  • coherent – True if all the layers are to be treated coherently, False otherwise
  • coherency_list – if coherent is False, a list of ‘c’ (coherent) or ‘i’ (incoherent) for each layer
  • zero_threshold – when the fraction of incident light absorbed in a layer is less than this value, the absorption

profile is completely set to zero for both coherent and incoherent calculations. This is applied on a wavelength-by-wavelength basis and is intended to prevent errors where integrating a weak absorption profile in a layer over many points leads to calculated EQE > total absorption in that layer. :param RAT_out: output from calculate_rat for the same stack & options :return: A dictionary containing the positions (in nm) and a 2D array with the absorption in the structure as a function of the position and the wavelength.

Transfer matrix engine

This is mostly based on Steven Byrnes’ tmm package, with modifications to vectorize the code over wavelengths (by Diego Alonso-Alvarez), and to include depth-dependent absorption calculations in incoherent layers using the Beer-Lambert law (by Phoebe Pearce).

All credit of the algorithm, testing, etc. goes to Steven Byrnes. For more details, visit:

The most two important functions are:

coh_tmm(…) – the transfer-matrix-method calculation in the coherent case (i.e. thin films)

inc_tmm(…) – the transfer-matrix-method calculation in the incoherent case (i.e. films tens or hundreds of wavelengths thick, or whose thickness is not very uniform.

These functions are all imported into the main package (tmm) namespace, so you can call them with tmm.coh_tmm(…) etc.

solcore.absorption_calculator.tmm_core_vec.make_2x2_array(a, b, c, d, dtype=<class 'float'>)[source]

Makes a 2x2 numpy array of [[a,b],[c,d]]

Same as “numpy.array([[a,b],[c,d]], dtype=float)”, but ten times faster

solcore.absorption_calculator.tmm_core_vec.snell(n_1, n_2, th_1)[source]

return angle theta in layer 2 with refractive index n_2, assuming it has angle th_1 in layer with refractive index n_1. Use Snell’s law. Note that “angles” may be complex!!

solcore.absorption_calculator.tmm_core_vec.list_snell(n_list, th_0)[source]

return list of angle theta in each layer based on angle th_0 in layer 0, using Snell’s law. n_list is index of refraction of each layer. Note that “angles” may be complex!!

solcore.absorption_calculator.tmm_core_vec.interface_r(polarization, n_i, n_f, th_i, th_f)[source]

reflection amplitude (from Fresnel equations)

polarization is either “s” or “p” for polarization

n_i, n_f are (complex) refractive index for incident and final

th_i, th_f are (complex) propegation angle for incident and final (in radians, where 0=normal). “th” stands for “theta”.

solcore.absorption_calculator.tmm_core_vec.interface_t(polarization, n_i, n_f, th_i, th_f)[source]

transmission amplitude (frem Fresnel equations)

polarization is either “s” or “p” for polarization

n_i, n_f are (complex) refractive index for incident and final

th_i, th_f are (complex) propegation angle for incident and final (in radians, where 0=normal). “th” stands for “theta”.

solcore.absorption_calculator.tmm_core_vec.R_from_r(r)[source]

Calculate reflected power R, starting with reflection amplitude r.

solcore.absorption_calculator.tmm_core_vec.T_from_t(pol, t, n_i, n_f, th_i, th_f)[source]

Calculate transmitted power T, starting with transmission amplitude t.

n_i,n_f are refractive indices of incident and final medium.

th_i, th_f are (complex) propegation angles through incident & final medium (in radians, where 0=normal). “th” stands for “theta”.

In the case that n_i,n_f,th_i,th_f are real, formulas simplify to T=|t|^2 * (n_f cos(th_f)) / (n_i cos(th_i)).

See manual for discussion of formulas

solcore.absorption_calculator.tmm_core_vec.power_entering_from_r(pol, r, n_i, th_i)[source]

Calculate the power entering the first interface of the stack, starting with reflection amplitude r. Normally this equals 1-R, but in the unusual case that n_i is not real, it can be a bit different than 1-R. See manual.

n_i is refractive index of incident medium.

th_i is (complex) propegation angle through incident medium (in radians, where 0=normal). “th” stands for “theta”.

solcore.absorption_calculator.tmm_core_vec.interface_R(polarization, n_i, n_f, th_i, th_f)[source]

Fraction of light intensity reflected at an interface.

solcore.absorption_calculator.tmm_core_vec.interface_T(polarization, n_i, n_f, th_i, th_f)[source]

Fraction of light intensity transmitted at an interface.

solcore.absorption_calculator.tmm_core_vec.coh_tmm(pol, n_list, d_list, th_0, lam_vac)[source]

This function is vectorized. Main “coherent transfer matrix method” calc. Given parameters of a stack, calculates everything you could ever want to know about how light propagates in it. (If performance is an issue, you can delete some of the calculations without affecting the rest.)

pol is light polarization, “s” or “p”.

n_list is the list of refractive indices, in the order that the light would pass through them. The 0’th element of the list should be the semi-infinite medium from which the light enters, the last element should be the semi- infinite medium to which the light exits (if any exits).

th_0 is the angle of incidence: 0 for normal, pi/2 for glancing. Remember, for a dissipative incoming medium (n_list[0] is not real), th_0 should be complex so that n0 sin(th0) is real (intensity is constant as a function of lateral position).

d_list is the list of layer thicknesses (front to back). Should correspond one-to-one with elements of n_list. First and last elements should be “inf”.

lam_vac is vacuum wavelength of the light.

Outputs the following as a dictionary (see manual for details)

  • r–reflection amplitude
  • t–transmission amplitude
  • R–reflected wave power (as fraction of incident)
  • T–transmitted wave power (as fraction of incident)
  • power_entering–Power entering the first layer, usually (but not always) equal to 1-R (see manual).
  • vw_list– n’th element is [v_n,w_n], the forward- and backward-traveling amplitudes, respectively, in the n’th medium just after interface with (n-1)st medium.
  • kz_list–normal component of complex angular wavenumber for forward-traveling wave in each layer.
  • th_list–(complex) propagation angle (in radians) in each layer
  • pol, n_list, d_list, th_0, lam_vac–same as input
solcore.absorption_calculator.tmm_core_vec.coh_tmm_reverse(pol, n_list, d_list, th_0, lam_vac)[source]

Reverses the order of the stack then runs coh_tmm.

solcore.absorption_calculator.tmm_core_vec.ellips(n_list, d_list, th_0, lam_vac)[source]

Calculates ellipsometric parameters, in radians.

Warning: Conventions differ. You may need to subtract pi/2 or whatever.

solcore.absorption_calculator.tmm_core_vec.unpolarized_RT(n_list, d_list, th_0, lam_vac)[source]

This function is vectorized. Calculates reflected and transmitted power for unpolarized light.

solcore.absorption_calculator.tmm_core_vec.position_resolved(layer, dist, coh_tmm_data)[source]

This function is vectorized. Starting with output of coh_tmm(), calculate the Poynting vector and absorbed energy density a distance “dist” into layer number “layer”

solcore.absorption_calculator.tmm_core_vec.find_in_structure(d_list, dist)[source]

This function is vectorized. d_list is list of thicknesses of layers, all of which are finite.

dist is the distance from the front of the whole multilayer structure (i.e., from the start of layer 0.)

Function returns [layer,z], where:

layer is what number layer you’re at. (For large enough dist, layer = len(d_list), even though d_list[layer] doesn’t exist in that case.

z is the distance into that layer.

solcore.absorption_calculator.tmm_core_vec.find_in_structure_with_inf(d_list, dist)[source]

d_list is list of thicknesses of layers [inf, blah, blah, …, blah, inf]

dist is the distance from the front of the whole multilayer structure (i.e., frcom the start of layer 1.)

Function returns [layer,z], where:

layer is what number layer you’re at,

z is the distance into that layer.

solcore.absorption_calculator.tmm_core_vec.layer_starts(d_list)[source]

Gives the location of the start of any given layer, relative to the front of the whole multilayer structure. (i.e. the start of layer 1)

d_list is list of thicknesses of layers [inf, blah, blah, …, blah, inf]

class solcore.absorption_calculator.tmm_core_vec.absorp_analytic_fn[source]

This function (specifically, ‘run’) is vectorized. Absorption in a given layer is a pretty simple analytical function: The sum of four exponentials.

a(z) = A1*exp(a1*z) + A2*exp(-a1*z)
  • A3*exp(1j*a3*z) + conj(A3)*exp(-1j*a3*z)

where a(z) is absorption at depth z, with z=0 being the start of the layer, and A1,A2,a1,a3 are real numbers, with a1>0, a3>0, and A3 is complex. The class stores these five parameters, as well as d, the layer thickness.

This gives absorption as a fraction of intensity coming towards the first layer of the stack.

fill_in(coh_tmm_data, layer)[source]

fill in the absorption analytic function starting from coh_tmm_data (the output of coh_tmm), for absorption in the layer with index “layer”.

copy()[source]

Create copy of an absorp_analytic_fn object

run(z)[source]

Calculates absorption at a given depth z, where z=0 is the start of the layer.

flip()[source]

Flip the function front-to-back, to describe a(d-z) instead of a(z), where d is layer thickness.

scale(factor)[source]

multiplies the absorption at each point by “factor”.

add(b)[source]

adds another compatible absorption analytical function

solcore.absorption_calculator.tmm_core_vec.absorp_in_each_layer(coh_tmm_data)[source]

An array listing what proportion of light is absorbed in each layer.

Assumes the final layer eventually absorbs all transmitted light.

Assumes the initial layer eventually absorbs all reflected light.

Entries of array should sum to 1.

coh_tmm_data is output of coh_tmm()

solcore.absorption_calculator.tmm_core_vec.inc_group_layers(n_list, d_list, c_list)[source]

Helper function for inc_tmm. Groups and sorts layer information.

See coh_tmm for definitions of n_list, d_list.

c_list is “coherency list”. Each entry should be ‘i’ for incoherent or ‘c’ for ‘coherent’.

A “stack” is a group of one or more consecutive coherent layers. A “stack index” labels the stacks 0,1,2,…. The “within-stack index” counts the coherent layers within the stack 1,2,3… [index 0 is the incoherent layer before the stack starts]

An “incoherent layer index” labels the incoherent layers 0,1,2,…

An “alllayer index” labels all layers (all elements of d_list) 0,1,2,…

Returns info about how the layers relate:

  • stack_d_list[i] = list of thicknesses of each coherent layer in the i’th stack, plus starting and ending with “inf”
  • stack_n_list[i] = list of refractive index of each coherent layer in the i’th stack, plus the two surrounding incoherent layers
  • all_from_inc[i] = j means that the layer with incoherent index i has alllayer index j
  • inc_from_all[i] = j means that the layer with alllayer index i has incoherent index j. If j = nan then the layer is coherent.
  • all_from_stack[i1][i2] = j means that the layer with stack index i1 and within-stack index i2 has alllayer index j
  • stack_from_all[i] = [j1 j2] means that the layer with alllayer index i is part of stack j1 with withinstack-index j2. If stack_from_all[i] = nan then the layer is incoherent
  • inc_from_stack[i] = j means that the i’th stack comes after the layer with incoherent index j, and before the layer with incoherent index j+1.
  • stack_from_inc[i] = j means that the layer with incoherent index i comes immediately after the j’th stack. If j=nan, it is not immediately following a stack.
  • num_stacks = number of stacks
  • num_inc_layers = number of incoherent layers
  • num_layers = number of layers total
solcore.absorption_calculator.tmm_core_vec.inc_tmm(pol, n_list, d_list, c_list, th_0, lam_vac)[source]

This function is vectorized. Incoherent, or partly-incoherent-partly-coherent, transfer matrix method.

See coh_tmm for definitions of pol, n_list, d_list, th_0, lam_vac.

c_list is “coherency list”. Each entry should be ‘i’ for incoherent or ‘c’ for ‘coherent’.

If an incoherent layer has real refractive index (no absorption), then its thickness doesn’t affect the calculation results.

See https://arxiv.org/abs/1603.02720 for physics background and some of the definitions.

Outputs the following as a dictionary:

  • R–reflected wave power (as fraction of incident)
  • T–transmitted wave power (as fraction of incident)
  • VW_list– n’th element is [V_n,W_n], the forward- and backward-traveling intensities, respectively, at the beginning of the n’th incoherent medium.
  • coh_tmm_data_list–n’th element is coh_tmm_data[n], the output of the coh_tmm program for the n’th “stack” (group of one or more consecutive coherent layers).
  • coh_tmm_bdata_list–n’th element is coh_tmm_bdata[n], the output of the coh_tmm program for the n’th stack, but with the layers of the stack in reverse order.
  • stackFB_list–n’th element is [F,B], where F is light traveling forward towards the n’th stack and B is light traveling backwards towards the n’th stack.
  • num_layers– total number both coherent and incoherent.
  • power_entering_list–n’th element is the normalized Poynting vector crossing the interface into the n’th incoherent layer from the previous (coherent or incoherent) layer.
  • Plus, all the outputs of inc_group_layers
solcore.absorption_calculator.tmm_core_vec.inc_absorp_in_each_layer(inc_data)[source]

A list saying what proportion of light is absorbed in each layer.

Assumes all reflected light is eventually absorbed in the 0’th medium, and all transmitted light is eventually absorbed in the final medium.

Returns a list [layer0absorp, layer1absorp, …]. Entries should sum to 1.

inc_data is output of inc_tmm()

solcore.absorption_calculator.tmm_core_vec.inc_find_absorp_analytic_fn(layer, inc_data)[source]

Outputs an absorp_analytic_fn object for a coherent layer within a partly-incoherent stack.

inc_data is output of inc_tmm()

solcore.absorption_calculator.tmm_core_vec.inc_position_resolved(layer, dist, inc_tmm_data, coherency_list, alphas, zero_threshold=1e-06)[source]

This function is vectorized. Analogous to position_resolved, but for layers (incoherent or coherent) in (partly) incoherent stacks. This is a new function, not from Steven Byrnes’ tmm package. Starting with output of inc_tmm(), calculate the Poynting vector and absorbed energy density a distance “dist” into layer number “layer”

solcore.absorption_calculator.tmm_core_vec.beer_lambert(alphas, fraction, dist)[source]

Calculates absorption profile according to the Beer-Lambert law given alphas (in m-1) and a vector of distance into the layer (in m) and the fraction of incident light reaching the front of the layer. This is used to calculate the absorption profile in incoherent layers within (partly) incoherent stacks. Vectorized over wavelengths.