Package 'planar'

Title: Multilayer Optics
Description: Solves the electromagnetic problem of reflection and transmission at a planar multilayer interface. Also computed are the decay rates and emission profile for a dipolar emitter.
Authors: Baptiste AuguiƩ [aut, cre] (<https://orcid.org/0000-0002-2749-5715>, Some functions ported from the original Matlab SPLAC code by E.C. Le Ru and P. G. Etchegoin), Steven Johnson [cph] (C code for the cubature library)
Maintainer: Baptiste AuguiĆ© <[email protected]>
License: Mozilla Public License Version 2.0
Version: 1.7.0
Built: 2024-11-23 04:20:00 UTC
Source: https://github.com/nano-optics/planar

Help Index


planar

Description

Multilayer optics

Details

R/c++ implementation of the dipole emission near a planar multilayer stack

Author(s)

baptiste Auguie [email protected]

References

Etchegoin, P. Le Ru, E., Principles of Surface-Enhanced Raman Spectroscopy, Elsevier, Amsterdam (2009).

L. Novotny, E. Hecht, Principles of Nano-optics Cambridge University Press, 2006

H. Raether. Surface Plasmons on Smooth and Rough Surfaces and on Gratings. Springer, 1988.


classify

Description

relabel factors

Usage

classify(d, id = NULL, vars = NULL, ...)

Arguments

d

data.frame

id

column id

vars

variables

...

passed on to melt

Details

Wide to long format data.frame with new factor variable(s) describing the original columns

Value

data.frame

Author(s)

Baptiste Auguie

See Also

Other helping_functions: internal_field(), invert_stack(), lfief(), modify_levels()


collection_ml

Description

Light intensity from the transmission of a bunch of plane waves at a planar interface

Usage

collection_ml(
  xyz,
  wavelength = 632.8,
  omega = c(40, 50) * pi/180,
  psi = 0,
  epsilon = c(1.5^2, epsAg(wavelength)$epsilon, 1^2, 1^2),
  thickness = c(0, 50, 10, 0),
  maxEval = 3000,
  reqAbsError = 0,
  tol = 1e-04,
  progress = FALSE
)

Arguments

xyz

position matrix

wavelength

wavelength

omega

collection angle

psi

polarisation angle

epsilon

vector of permittivities

thickness

thickness corresponding to each medium

maxEval

passed to cubature

reqAbsError

passed to cubature

tol

passed to cubature

progress

logical display progress bar

Details

Integration is performed over the solid angle defined by omega

Value

data.frame intensity at the x, y, z position

Author(s)

Baptiste Auguie


combine_layer

Description

combine layer

Usage

combine_layer(r1, r2, kd)

Arguments

r1

reflection coefficient

r2

reflection coefficient

kd

k*d

Details

reflection coefficient for a layer

Value

combined complex reflectivity

Author(s)

baptiste Auguie


dbr_analytic

Description

semi-infinite DBR

Usage

dbr_analytic(
  wavelength,
  lambda0,
  n1,
  n2,
  nleft,
  d1 = lambda0/4/n1,
  d2 = lambda0/4/n2,
  ...
)

Arguments

wavelength

in nm

lambda0

central wavelength of the stopband

n1

real refractive index for odd layers

n2

real refractive index for even layers

nleft

real refractive index for incident medium

d1

odd layer thickness in nm

d2

even layer thickness in nm

...

ignored

Details

periodic structure of dielectric layers

Value

data.frame with complex reflectivity

Note

issue at lambda0/2 needs investigating

Author(s)

baptiste Auguie

References

Amir and Vukusic, 2013, arXiv:1209.3776v2


dbr_stack

Description

DBR stack structure

Usage

dbr_stack(
  lambda0 = 630,
  n1 = 1.28,
  n2 = 1.72,
  d1 = lambda0/4/n1,
  d2 = lambda0/4/n2,
  N = 2 * pairs,
  pairs = 4,
  ...
)

Arguments

lambda0

central wavelength of the stopband

n1

real refractive index for odd layers

n2

real refractive index for even layers

d1

odd layer thickness in nm

d2

even layer thickness in nm

N

number of layers, overwrites pairs

pairs

number of pairs

...

ignored

Details

periodic structure of dielectric layers

Value

list of class 'stack'

Author(s)

baptiste Auguie

See Also

Other stack user_level: embed_stack(), layer_stack(), tamm_stack_ir(), tamm_stack_porous(), tamm_stack()


dipole

Description

Dipole decay rates near a multilayer interface

Usage

dipole(
  d = 1,
  wavelength,
  epsilon = list(incident = 1^2),
  thickness = c(0, 0),
  qcut = NULL,
  rel.err = 0.001,
  Nquadrature1 = 1000,
  Nquadrature2 = 10000,
  Nquadrature3 = 10000,
  GL = FALSE,
  show.messages = TRUE
)

Arguments

d

distance in nm

wavelength

wavelength in nm

epsilon

list of dielectric functions

thickness

list of layer thicknesses

qcut

transition between regions 2 and 3

rel.err

relative error

Nquadrature1

maximum number of quadrature points in radiative region

Nquadrature2

maximum number of quadrature points in SPPs region

Nquadrature3

maximum number of quadrature points in dipole image region

GL

logical: use Gauss Legendre quadrature, or cubature::adaptIntegrate

show.messages

logical, display integration info

Details

dipole decay rates near a multilayer interface

Author(s)

baptiste Auguie


dipole_direct

Description

Dipole total decay rate near a multilayer interface

Usage

dipole_direct(
  d = 1,
  wavelength,
  epsilon = list(incident = 1^2),
  thickness = c(0, 0),
  Nquadrature1 = 50,
  Nquadrature2 = 200,
  Nquadrature3 = 50,
  qcut = NULL,
  qmax = Inf,
  show.messages = TRUE
)

Arguments

d

distance in nm

wavelength

wavelength in nm

epsilon

list of dielectric functions

thickness

list of layer thicknesses

Nquadrature1

quadrature points in radiative region

Nquadrature2

quadrature points in SPPs region

Nquadrature3

quadrature points in dipole image region

qcut

transition between regions 2 and 3

qmax

maximum q of region 3

show.messages

logical, display integration info

Details

direct application of the textbook formula using integrand_mtot; performs poorly compared to the transformed version in dipole

Author(s)

baptiste Auguie


embed_stack

Description

Embed stack structure

Usage

embed_stack(s, nleft = 1, nright = 1, dleft = 200, dright = 200, ...)

Arguments

s

stack (finite structure)

nleft

real refractive index on the left side

nright

real refractive index on the right side

dleft

dummy layer thickness in nm

dright

dummy layer thickness in nm

...

ignored

Details

embeds a stack in semi-infinite media

Value

list of class 'stack'

Author(s)

baptiste Auguie

See Also

Other stack user_level: dbr_stack(), layer_stack(), tamm_stack_ir(), tamm_stack_porous(), tamm_stack()


epsilon_dispersion

Description

epsilon_dispersion

Usage

epsilon_dispersion(
  epsilon,
  wavelength = seq(400, 1000),
  envir = parent.frame()
)

Arguments

epsilon

list of real or complex values

wavelength

numeric vector

envir

environment to look for functions

Details

apply a function to a range of wavelength and return dielectric function

Value

list

Author(s)

baptiste Auguie


epsilon_label

Description

epsilon_label

Usage

epsilon_label(epsilon = list(3.5, 1, 3, 1, "epsAu", 3, 3.5), names = NULL)

Arguments

epsilon

list of real or complex values

names

optional unique character names in order of appearance

Details

characterise the layers of a structure with unique labels for metals and dielectrics

Value

factor

Author(s)

baptiste Auguie


gaussian_near_field_layer

Description

Electric field from the transmission of a gaussian beam at a planar interface

Usage

gaussian_near_field_layer(
  xyz,
  wavelength = 500,
  alpha = 15 * pi/180,
  psi = 0,
  w0 = 10000,
  epsilon = c(1.5^2, epsAg(wavelength)$epsilon, 1^2),
  thickness = c(0, 50, 0),
  maxEval = 3000,
  reqAbsError = 0,
  tol = 1e-04,
  progress = FALSE,
  field = FALSE
)

Arguments

xyz

position

wavelength

wavelength

alpha

beam incident angle

psi

beam polarisation angle

w0

beam waist radius

epsilon

vector of permittivities

thickness

thickness corresponding to each medium

maxEval

passed to adaptIntegrate

reqAbsError

passed to cubature

tol

passed to adaptIntegrate

progress

logical: display progress bar

field

logical: return the electric field (complex vector), or modulus squared

Details

Integration is performed over a spectrum of incident plane waves

Value

data.frame electric field at the x, y, z position

Author(s)

Baptiste Auguie

See Also

Other gaussian_beam: gaussian_near_field_ml()


gaussian_near_field_ml

Description

Electric field of a gaussian beam close to a planar interface

Usage

gaussian_near_field_ml(
  xyz,
  wavelength = 632.8,
  alpha = 15 * pi/180,
  psi = 0,
  w0 = 10000,
  epsilon = c(1.5^2, epsAg(wavelength)$epsilon, 1^2, 1^2),
  thickness = c(0, 50, 10, 0),
  maxEval = 3000,
  reqAbsError = 0,
  tol = 1e-04,
  progress = FALSE,
  field = FALSE
)

Arguments

xyz

position matrix

wavelength

wavelength

alpha

beam incident angle

psi

beam polarisation angle

w0

beam waist radius

epsilon

vector of permittivities

thickness

thickness corresponding to each medium

maxEval

passed to cubature

reqAbsError

passed to cubature

tol

passed to cubature

progress

logical display progress bar

field

logical: return the electric field (complex vector), or modulus squared

Details

Integration is performed over a spectrum of incident plane waves using integrand_gb2

Value

data.frame electric field at the x, y, z position

Author(s)

Baptiste Auguie

See Also

Other gaussian_beam: gaussian_near_field_layer()


integrand_mtot

Description

Total decay rate of a dipole near a multilayer interface

Usage

integrand_mtot(
  d = 10,
  q,
  wavelength,
  epsilon = list(incident = 1.5^2, 1^2),
  thickness = c(0, 0)
)

Arguments

d

distance in nm

q

normalised in-plane wavevector in [0, infty)

wavelength

wavelength in nm

epsilon

list of dielectric functions

thickness

list of layer thicknesses

Details

Integrand without transformation of variables

Author(s)

baptiste Auguie

See Also

Other integrands dipole: integrand_nr1(), integrand_nr2(), integrand_nr3(), integrand_rad()


integrand_nr1

Description

Dipole decay rates near a multilayer interface

Usage

integrand_nr1(
  d = 10,
  u,
  wavelength,
  epsilon = list(incident = 1.5^2, 1^2),
  thickness = c(0, 0),
  GL = FALSE
)

Arguments

d

distance in nm

u

transformed normalised in-plane wavevector sqrt(1-q^2)

wavelength

wavelength in nm

epsilon

list of dielectric functions

thickness

list of layer thicknesses

GL

logical: result formatted for use with Gauss Legendre quadrature

Details

Integrand of the dipole decay rates near a multilayer interface. Transformed part I1 (radiative) from u=0 to 1

Author(s)

baptiste Auguie

See Also

Other integrands dipole: integrand_mtot(), integrand_nr2(), integrand_nr3(), integrand_rad()


integrand_nr2

Description

Dipole decay rates near a multilayer interface

Usage

integrand_nr2(
  d = 10,
  u,
  wavelength,
  epsilon = list(incident = 1.5^2, 1^2),
  thickness = c(0, 0),
  GL = FALSE
)

Arguments

d

distance in nm

u

transformed normalised in-plane wavevector sqrt(q^2 - 1)

wavelength

wavelength in nm

epsilon

list of dielectric functions

thickness

list of layer thicknesses

GL

logical: result formatted for use with Gauss Legendre quadrature

Details

Integrand of the dipole decay rates near a multilayer interface. Transformed part I2 from u=0 to ucut

Author(s)

baptiste Auguie

See Also

Other integrands dipole: integrand_mtot(), integrand_nr1(), integrand_nr3(), integrand_rad()


integrand_nr3

Description

Dipole decay rates near a multilayer interface

Usage

integrand_nr3(
  d = 10,
  u,
  ucut,
  wavelength,
  epsilon = list(incident = 1.5^2, 1^2),
  thickness = c(0, 0),
  GL = FALSE
)

Arguments

d

distance in nm

u

transformed normalised in-plane wavevector sqrt(q^2 - 1)

ucut

limit of the integral

wavelength

wavelength in nm

epsilon

list of dielectric functions

thickness

list of layer thicknesses

GL

logical: result formatted for use with Gauss Legendre quadrature

Details

Integrand of the dipole decay rates near a multilayer interface. Transformed part III from u=ucut to infinity

Author(s)

baptiste Auguie

See Also

Other integrands dipole: integrand_mtot(), integrand_nr1(), integrand_nr2(), integrand_rad()


integrand_rad

Description

Dipole decay rates near a multilayer interface

Usage

integrand_rad(
  d = 10,
  angle,
  wavelength,
  epsilon = list(incident = 1.5^2, 1^2),
  thickness = c(0, 0),
  GL = FALSE
)

Arguments

d

distance in nm

angle

angle in radians

wavelength

wavelength in nm

epsilon

list of dielectric functions

thickness

list of layer thicknesses

GL

logical: result formatted for use with Gauss Legendre quadrature

Details

Integrand of the radiative dipole decay rates near a multilayer interface.

Author(s)

baptiste Auguie

See Also

Other integrands dipole: integrand_mtot(), integrand_nr1(), integrand_nr2(), integrand_nr3()


internal_field

Description

Internal field in a ML stack

Usage

internal_field(
  wavelength = 500,
  angle = 0,
  psi = 0,
  thickness = c(0, 20, 140, 20, 0),
  dmax = 200,
  res = 1000,
  epsilon = c(1^2, -12, 1.38^2, -12, 1.46^2),
  field = FALSE,
  ...
)

Arguments

wavelength

wavelength

angle

angle

psi

polarisation angle (0 for TM)

thickness

vector of layer thickness

dmax

maximum distance to interface

res

resolution of sampling points

epsilon

permittivities

field

logical, return complex field vector, or modulus squared

...

further args ignored

Details

returns the electric field as a function of distance inside and outside of the structure

Value

data.frame with position and electric field vector

Author(s)

baptiste Auguie

References

Principles of surface-enhanced Raman spectroscopy and related plasmonic effects

Eric C. Le Ru and Pablo G. Etchegoin, published by Elsevier, Amsterdam (2009).

See Also

Other helping_functions: classify(), invert_stack(), lfief(), modify_levels()


invert_stack

Description

invert the description of a multilayer to simulate the opposite direction of incidence

Usage

invert_stack(p)

Arguments

p

list

Details

inverts list of epsilon and thickness of layers

Value

list

Author(s)

Baptiste Auguie

See Also

Other helping_functions: classify(), internal_field(), lfief(), modify_levels()


layer_stack

Description

Single-layer stack structure

Usage

layer_stack(epsilon = "epsAu", thickness = 50, ...)

Arguments

epsilon

dielectric function (numeric, character, or complex)

thickness

layer thickness in nm

...

ignored

Details

returns a stack describing a single layer

Value

list of class 'stack'

Author(s)

baptiste Auguie

See Also

Other stack user_level: dbr_stack(), embed_stack(), tamm_stack_ir(), tamm_stack_porous(), tamm_stack()


lfief

Description

Local field intensity enhancement factors in a multilayer

Usage

lfief(
  wavelength = 500,
  angle = 0,
  polarisation = "p",
  thickness = c(0, 20, 140, 20, 0),
  dmax = 200,
  res = 1000,
  res2 = res/10,
  epsilon = list(1^2, -12, 1.38^2, -12, 1.46^2),
  displacement = FALSE,
  ...
)

Arguments

wavelength

wavelength

angle

angle

polarisation

polarisation

thickness

vector of layer thickness

dmax

maximum distance to interface, if > layer thickness

res

resolution of sampling points

res2

resolution of sampling points outside stack

epsilon

list of permittivities

displacement

logical, Mperp corresponds to displacement squared (D=epsilon x E)

...

further args passed to multilayer

Details

returns the LFIEFs as a function of distance inside and outside of the structure

Value

long format data.frame with positions and LFEF (para and perp)

Author(s)

baptiste Auguie

References

Principles of surface-enhanced Raman spectroscopy and related plasmonic effects

Eric C. Le Ru and Pablo G. Etchegoin, published by Elsevier, Amsterdam (2009).

See Also

Other helping_functions: classify(), internal_field(), invert_stack(), modify_levels()


modify_levels

Description

relabel factors

Usage

modify_levels(f, modify = list())

Arguments

f

factor

modify

named list

Value

factor

Author(s)

Baptiste Auguie

See Also

Other helping_functions: classify(), internal_field(), invert_stack(), lfief()


multilayer

Description

Multilayer Fresnel coefficients

Usage

multilayer(
  wavelength = 2 * pi/k0,
  k0 = 2 * pi/wavelength,
  angle = asin(q),
  q = sin(angle),
  epsilon = list(incident = 1.5^2, 1.33),
  thickness = c(0, 0),
  polarisation = c("p", "s"),
  d = 1,
  dout = d,
  ...
)

Arguments

wavelength

[vector] wavelength in nm

k0

[vector] wavevector in nm^-1

angle

[vector] incident angles in radians

q

[vector] normalised incident in-plane wavevector

epsilon

list of N+2 dielectric functions, each of length 1 or length(wavelength)

thickness

vector of N+2 layer thicknesses, first and last are dummy

polarisation

[character] switch between p- and s- polarisation

d

vector of distances where LFIEF are evaluated from each interface

dout

vector of distances where LFIEF are evaluated outside the stack

...

unused

Details

solves the EM problem of a multilayered interface

Value

fresnel coefficients and field profiles

Author(s)

baptiste Auguie

References

Principles of surface-enhanced Raman spectroscopy and related plasmonic effects. Eric C. Le Ru and Pablo G. Etchegoin, published by Elsevier, Amsterdam (2009).


multilayercpp

Description

Multilayer Fresnel coefficients

Usage

multilayercpp(
  wavelength = 2 * pi/k0,
  k0 = 2 * pi/wavelength,
  angle = asin(q),
  q = sin(angle),
  epsilon = list(incident = 1.5^2, 1.33),
  thickness = c(0, 0),
  ...
)

Arguments

wavelength

[vector] wavelength in nm

k0

[vector] wavevector in nm^-1

angle

[vector] incident angles in radians

q

[vector] normalised incident in-plane wavevector

epsilon

list of N+2 dielectric functions, each of length 1 or length(wavelength)

thickness

vector of N+2 layer thicknesses, first and last are dummy

...

unused

Details

solves the EM problem of a multilayered interface

Value

fresnel coefficients and field profiles

Author(s)

baptiste Auguie

Examples

library(planar)
demo(package="planar")

multilayerfull

Description

Multilayer Fresnel coefficients

Usage

multilayerfull(
  wavelength = 2 * pi/k0,
  k0 = 2 * pi/wavelength,
  angle = asin(q),
  q = sin(angle),
  epsilon = list(incident = 1.5^2, 1.33),
  thickness = c(0, 0),
  psi = 0,
  z = 0,
  ...
)

Arguments

wavelength

[vector] wavelength in nm

k0

[vector] wavevector in nm^-1

angle

[vector] incident angles in radians

q

[vector] normalised incident in-plane wavevector

epsilon

list of N+2 dielectric functions, each of length 1 or length(wavelength)

thickness

vector of N+2 layer thicknesses, first and last are dummy

psi

[numeric] polarisation angle

z

[vector] positions to calculate the electric field intensity

...

unused

Details

solves the EM problem of a multilayered interface

Value

fresnel coefficients and field profiles

Author(s)

baptiste Auguie


Colour palettes for multilayer stacks

Description

Custom palettes for DBR stacks

Custom palette for Tamm stacks

Alternative palette for Tamm stacks

Format

Colour palette (vectors of colours)

Colour palette (vectors of colours)

Colour palette (vectors of colours)

Source

See RColorBrewer package


sort_factor

Description

raman_shift

Usage

raman_shift(laser = c(514, 632.8), shift = c(520, 610))

Arguments

laser

vector of laser wavelengths in nm

shift

vector of Raman shifts in cm-1

Details

converts Raman shift to wavelength

Value

matrix of shifted wavelengths (all combinations)

Author(s)

Baptiste Auguie


recursive_fresnel

Description

Multilayer Fresnel coefficients

Usage

recursive_fresnel(
  wavelength = 2 * pi/k0,
  k0 = 2 * pi/wavelength,
  angle = NULL,
  q = sin(angle),
  epsilon = list(incident = 1.5^2, 1.33^2),
  thickness = c(0, 0),
  polarisation = c("p", "s")
)

Arguments

wavelength

[vector] wavelength in nm

k0

[vector] wavevector in nm^-1

angle

[vector] incident angles in radians

q

[vector] normalised incident in-plane wavevector

epsilon

list of N+2 dielectric functions, each of length 1 or length(wavelength)

thickness

vector of N+2 layer thicknesses, first and last are dummy

polarisation

[character] switch between p- and s- polarisation

Details

computes the reflection coefficient of a multilayered interface

Value

fresnel coefficients and field profiles

Author(s)

baptiste Auguie


recursive_fresnelcpp

Description

Multilayer Fresnel coefficients

Usage

recursive_fresnelcpp(
  wavelength = 2 * pi/k0,
  k0 = 2 * pi/wavelength,
  angle = NULL,
  q = sin(angle),
  epsilon = list(incident = 1.5^2, 1.33^2),
  thickness = c(0, 0),
  polarisation = c("p", "s")
)

Arguments

wavelength

[vector] wavelength in nm

k0

[vector] wavevector in nm^-1

angle

[vector] incident angles in radians

q

[vector] normalised incident in-plane wavevector

epsilon

list of N+2 dielectric functions, each of length 1 or length(wavelength)

thickness

vector of N+2 layer thicknesses, first and last are dummy

polarisation

[character] switch between p- and s- polarisation

Details

computes the reflection coefficient of a multilayered interface

Value

fresnel coefficients and field profiles

Author(s)

baptiste Auguie


rev.stack

Description

invert the description of a multilayer to simulate the opposite direction of incidence

Usage

## S3 method for class 'stack'
rev(x)

Arguments

x

stack

Details

inverts list of epsilon and thickness of layers

Value

stack

Author(s)

Baptiste Auguie

See Also

Other helping_functions user_level stack: simulate_ff(), simulate_nf()


simulate_ff

Description

simultate the far-field response of a multilayer stack

Usage

simulate_ff(
  ...,
  s = NULL,
  fun = tamm_stack,
  wavelength = seq(400, 1000),
  angle = 0,
  polarisation = c("p", "s")
)

Arguments

...

further arguments passed to fun

s

stack (optional)

fun

function returning a stack

wavelength

numeric vector

angle

incident angle in radians

polarisation

p or s

Details

wrapper around recursive_fresnelcpp for a stack structure

Value

data.frame

Author(s)

Baptiste Auguie

See Also

Other helping_functions user_level stack: rev.stack(), simulate_nf()


simulate_nf

Description

simultate the internal field of a multilayer stack

Usage

simulate_nf(
  ...,
  s = NULL,
  fun = tamm_stack,
  wavelength = 630,
  angle = 0,
  polarisation = c("p", "s"),
  dmax = 0,
  res = 10000,
  field = FALSE
)

Arguments

...

further arguments passed to fun

s

stack (optional)

fun

function returning a stack

wavelength

numeric vector

angle

incident angle in radians

polarisation

p or s

dmax

maximum distance from stack boundary

res

number of points

field

logical, return the real electric field

Details

wrapper around multilayer_field for a stack structure

Value

data.frame

Author(s)

Baptiste Auguie

See Also

Other helping_functions user_level stack: rev.stack(), simulate_ff()


tamm_stack

Description

DBR-metal stack structure

Usage

tamm_stack(
  lambda0 = 630,
  n1 = 1.28,
  n2 = 1.72,
  d1 = lambda0/4/n1,
  d2 = lambda0/4/n2,
  N = 2 * pairs,
  pairs = 4,
  dx1 = 0,
  dx2 = 0,
  dm = 50,
  metal = "epsAu",
  position = c("after", "before"),
  incidence = c("left", "right"),
  nleft = 1.5,
  nright = 1,
  dleft = 200,
  dright = 200,
  ...
)

Arguments

lambda0

central wavelength of the stopband

n1

real refractive index for odd layers

n2

real refractive index for even layers

d1

odd layer thickness in nm

d2

even layer thickness in nm

N

number of layers, overwrites pairs

pairs

number of pairs

dx1

variation of last odd layer thickness in nm

dx2

variation of last even layer thickness in nm

dm

thickness of metal layer

metal

character name of dielectric function

position

metal position relative to DBR

incidence

direction of incidence

nleft

refractive index of entering medium

nright

refractive index of outer medium

dleft

distance from the left side for visualisation

dright

distance from the right side for visualisation

...

ignored

Details

periodic structure of dielectric layers against metal film

Value

list of class 'stack'

Author(s)

baptiste Auguie

See Also

Other stack user_level: dbr_stack(), embed_stack(), layer_stack(), tamm_stack_ir(), tamm_stack_porous()


tamm_stack_ir

Description

DBR-metal stack structure

Usage

tamm_stack_ir(
  lambda0 = 950,
  n1 = 3,
  n2 = 3.7,
  d1 = lambda0/4/n1,
  d2 = lambda0/4/n2,
  N = 2 * pairs,
  pairs = 4,
  dx1 = 0,
  dx2 = 0,
  dm = 50,
  metal = "epsAu",
  position = "after",
  incidence = "left",
  nleft = n2,
  nright = 1,
  ...
)

Arguments

lambda0

central wavelength of the stopband

n1

real refractive index for odd layers

n2

real refractive index for even layers

d1

odd layer thickness in nm

d2

even layer thickness in nm

N

number of layers, overwrites pairs

pairs

number of pairs

dx1

variation of last odd layer thickness in nm

dx2

variation of last even layer thickness in nm

dm

thickness of metal layer

metal

character name of dielectric function

position

metal position relative to DBR

incidence

direction of incidence

nleft

refractive index of entering medium

nright

refractive index of outer medium

...

ignored

Details

periodic structure of dielectric layers against metal film

Value

list of class 'stack'

Author(s)

baptiste Auguie

See Also

Other stack user_level: dbr_stack(), embed_stack(), layer_stack(), tamm_stack_porous(), tamm_stack()


tamm_stack_porous

Description

DBR-metal stack structure

Usage

tamm_stack_porous(
  lambda0 = 600,
  n1 = 1.72,
  n2 = 1.28,
  d1 = lambda0/4/n1,
  d2 = lambda0/4/n2,
  N = 2 * pairs,
  pairs = 4,
  dx1 = 0,
  dx2 = 0,
  dm = 20,
  metal = "epsAu",
  position = "before",
  incidence = "right",
  nleft = 1.5,
  nright = 1,
  ...
)

Arguments

lambda0

central wavelength of the stopband

n1

real refractive index for odd layers

n2

real refractive index for even layers

d1

odd layer thickness in nm

d2

even layer thickness in nm

N

number of layers, overwrites pairs

pairs

number of pairs

dx1

variation of last odd layer thickness in nm

dx2

variation of last even layer thickness in nm

dm

thickness of metal layer

metal

character name of dielectric function

position

metal position relative to DBR

incidence

direction of incidence

nleft

refractive index of entering medium

nright

refractive index of outer medium

...

ignored

Details

periodic structure of dielectric layers against metal film

Value

list of class 'stack'

Author(s)

baptiste Auguie

See Also

Other stack user_level: dbr_stack(), embed_stack(), layer_stack(), tamm_stack_ir(), tamm_stack()


transmission

Description

transmission loss through a prism

Usage

transmission(n, external, polarisation = "p")

Arguments

n

prism refractive index

external

external incident angle in radians

polarisation

polarisation

Details

transmission loss through a prism

Value

transmission

Author(s)

baptiste Auguie