Processing¶
linear¶
linear contains evaluation functions for linear effects.
near-field¶
- pyGDM2.linear.nearfield(sim, field_index, r_probe, which_fields=['Es', 'Et', 'Bs', 'Bt'], val_inside_struct='field', val_outside_struct='field', N_neighbors_internal_field=4, method='')¶
Nearfield distribution in proximity of nanostructre
For a given incident field, calculate the electro-magnetic field in the proximity of the nanostructure (positions defined by MAP).
Pure python implementation.
CUDA support to run on GPU, which can be significantly faster on large problems.
- Parameters:
sim (
core.simulation
) – simulation descriptionfield_index (int) – index of evaluated self-consistent field to use for calculation. Can be obtained for specific parameter-set using
tools.get_closest_field_index()
r_probe (tuple (x,y,z) or list of 3-lists/-tuples) – (list of) coordinate(s) to evaluate nearfield on. Format: tuple (x,y,z) or list of 3 lists: [Xmap, Ymap, Zmap] (the latter can be generated e.g. using
tools.generate_NF_map()
)which_fields (list of str, default: ["Es","Et","Bs","Bt"]) –
- which fields to calculate and return. available options:
[“Es”,”Et”,”E0”, “Bs”,”Bt”,”B0”]
val_inside_struct (str, default: "field") – one of [“field”, “0”, “zero”, “NaN”, “none”, None]. value to return for positions inside structure. default “field” returns the field at the location of the closest meshpoint. Can be very coarse. Disable by setting “None”, but note that inside the structure the Green’s function will diverge in the latter case.
val_outside_struct (str, default: "field") – see val_inside_struct
N_neighbors_internal_field (int, default: 4) – Average internal field over field at N closes meshpoints. Neighbor fields are weighted by the distance of the evaluation point to the respective neighbor mesh cell.
method (str) – deprecated, has no effect
- Returns:
depending on kwarg which_fields, up to 4 lists of 6-tuples, complex –
scattered Efield (“Es”)
total Efield (“Et”, inlcuding fundamental field)
scattered Bfield (“Bs”)
total Bfield (“Bt”, inlcuding fundamental field)
the tuples are of shape (X,Y,Z, Ax,Ay,Az) with Ai the corresponding
complex field component
Notes
For details of the calculation of the scattered field outside the nano-object using the self-consistent field inside the particle, see e.g.: Girard, C. Near fields in nanostructures. Reports on Progress in Physics 68, 1883–1933 (2005).
- pyGDM2.linear.optical_chirality(sim, field_index, r_probe, which_field='t', **kwargs)¶
calculate the optical chirality of the electromagnetic field
** ——- FUNCTION STILL UNDER TESTING ——- **
Returns the normalized electromagnetic chirality C / C_LCP, as defined in [1-4]. Normalized to a left circular polarized (LCP) plane wave of amplitude |E_0|=1. Hence:
–> LCP: C = + 1 –> RCP: C = - 1
see also [3] or [4] for a short discussion about the normalization to LCP.
|C/C_LCP|>1 means that the local field is “superchiral”, hence a chiral molecule is excited with higher selectivity than it would be with circular polarized light.
kwargs are passed to
nearfield()
.- Parameters:
sim (
core.simulation
) – simulation descriptionfield_index (int) – index of evaluated self-consistent field to use for calculation. Can be obtained for specific parameter-set using
tools.get_closest_field_index()
r_probe (tuple (x,y,z) or list of 3-lists/-tuples) – (list of) coordinate(s) to evaluate nearfield on. Format: tuple (x,y,z) or list of 3 lists: [Xmap, Ymap, Zmap] (the latter can be generated e.g. using
tools.generate_NF_map()
)which_field (str, default: 't') –
- either of:
”s”, “Es” for scattered field
”t”, “Et” for total field (with illumination, i.e. Et=Es+E0; Bt=Bs+B0)
”0”, “E0” for incident field (illumination (E0, B0) only, corresponding to simulation with no nanostructure)
- Returns:
list of tuples (X,Y,Z, X). Each tuple consists of the position and the
(normalized) optical chirality *C / C_LCP.*
Notes
For details on the optical chirality, see e.g.:
[1] Tang, Y. & Cohen, A. E.: Optical Chirality and Its Interaction with Matter. Phys. Rev. Lett. 104, 163901 (2010)
[2] Meinzer, N., Hendry, E. & Barnes, W. L.: Probing the chiral nature of electromagnetic fields surrounding plasmonic nanostructures. Phys. Rev. B 88, 041407 (2013)
A discussion about the proper normalization of C can be found in:
[3] 1.Schäferling, M., Yin, X., Giessen, H. Formation of chiral fields in a symmetric environment. Opt. Express 20(24), 26326–26336 (2012)
[4] Schäferling M., Yin X., Engheta N., Giessen H. Helical Plasmonic Nanostructures as Prototypical Chiral Near-Field Sources. ACS Photonics 1(6), 530-537 (2014)
- pyGDM2.linear.internal_field_intensity(sim, field_index, which_value='average', which_field='E')¶
return internal field intensity
returns either average or peak (min or max) internal field intensity for either electric of magnetic field.
- Parameters:
sim (
core.simulation
) – simulation descriptionfield_index (int) – index of evaluated self-consistent field to use for calculation. Can be obtained for specific parameter-set using
tools.get_closest_field_index()
which_value (str, default: 'average') –
if ‘average’, ‘avg’ or ‘mean’: average field intensity
if ‘max’ or ‘maximum’: maximum field intensity
if ‘min’ or ‘minimum’: minimum field intensity
which_field (str, default: 'E') – one of ‘E’, ‘H’ or ‘both’ (either electric or magnetic field or both)
- pyGDM2.linear.poynting(sim, field_index, r_probe, which_fields='tot')¶
calculate the time average Poynting vector
calculate the time average Poynting vector of the scattered, total or fundamental (incident) field of a pyGDM simulation at one or more locations
- Parameters:
sim (
core.simulation
) – simulation descriptionfield_index (int) – index of evaluated self-consistent field to use for calculation. Can be obtained for specific parameter-set using
tools.get_closest_field_index()
r_probe (tuple (x,y,z) or list of 3-lists/-tuples) – (list of) coordinate(s) to evaluate nearfield on. Format: tuple (x,y,z) or list of 3 lists: [Xmap, Ymap, Zmap] (the latter can be generated e.g. using
tools.generate_NF_map()
)which_fields (str, default: "tot") – which type of fields to use. one of: “scat”, “tot”, “0”. scattered field, total field (=E_scat + E_0) or incident field alone.
- Returns:
P
- Return type:
list of 6-tuples (x,y,z, Px, Py, Pz)
- pyGDM2.linear.field_gradient(sim, field_index, r_probe=None, delta=None, which_fields=['Et'])¶
nearfield gradient distribution inside or in proximity of nanostructure
Calculate field-gradients (positions defined by r_probe). Using center difference for the numerical derivatives.
Implemented by C. Majorel.
- Parameters:
sim (
core.simulation
) – simulation descriptionfield_index (int) – index of evaluated self-consistent field to use for calculation. Can be obtained for specific parameter-set using
tools.get_closest_field_index()
r_probe (tuple (x,y,z) or list of 3-lists/-tuples, default: positions inside geometry) – (list of) coordinate(s) to evaluate nearfield on. Format: tuple (x,y,z) or list of 3 lists: [Xmap, Ymap, Zmap] (the latter can be generated e.g. using
tools.generate_NF_map()
) If None: calculate inside nanostructuredelta (float, default=None) – differential step for numerical center-derivative (in nanometers). If None, use stepsize, which is recommended inside a nanostructure. Finer steps are recommended for fields outside the nanostructure (1nm may be a good choice).
which_fields (list of str, default: ["Et"]) –
- for which fields to calculate the gradient. available options:
[“Es”,”Et”,”E0”, “Hs”,”Ht,”H0] (“B” can be used as equivalent to “H”)
- Returns:
depending on kwarg which_gradfields, up to 6 lists of 3 lists of 6-tuples, complex –
gradient scattered Efield (“Es”)
gradient total Efield (“Et”, includes incident field)
gradient incident Efield (“E0”)
gradient scattered Hfield (“Hs”)
gradient total Hfield (“Ht”, includes incident field)
gradient incident Hfield (“H0”)
each of the lists contains 3 lists [dAx, dAy, dAz] with dAj the differential terms
dE[0] = dE/dx = [dE0x/dx, dE0y/dx, dE0z/dx]
dE[1] = dE/dy = [dE0x/dy, dE0y/dy, dE0z/dy]
dE[2] = dE/dz = [dE0x/dz, dE0y/dz, dE0z/dz]
- pyGDM2.linear.optical_force(sim, field_index, return_value='total')¶
optical force acting on the nanostructure
** ——- FUNCTION STILL UNDER TESTING ——- **
- Parameters:
sim (
core.simulation
) – simulation descriptionfield_index (int) – index of evaluated self-consistent field to use for calculation. Can be obtained for specific parameter-set using
tools.get_closest_field_index()
return_value (str, default: 'total') –
Values to be returned. Either ‘total’ (default) or ‘structure’.
”total” : return the total force on the entire particle (3-vector [Fx, Fy, Fz])
”structure” : return spatially resolved force, in other words
the force acting on each meshpoint. –> list of tuples [x,y,z,F_ix,F_iy,F_iz])
- Returns:
F –
`return_value`==”total” : (3-vector [Fx, Fy, Fz]) total force on structure
return_value`==”structure” : (return list of tuples [xi,yi,zi,F_ix,F_iy,F_iz]) force at all mesh-cells `i of positions [xi,yi,zi].
- Return type:
tuple (Fx,Fy,Fz) or list of force vectors [x,y,z,Fx,Fy,Fz]
Notes
For details on optical force calculations see e.g.:
[1] Chaumet & Nieto-Vesperinas. Coupled dipole method determination of the electromagnetic force on a particle over a flat dielectric substrate. Phys. Rev. B 61, 14119–14127 (2000)
- pyGDM2.linear.heat(sim, field_index, power_scaling_e0=1.0, return_value='total', return_units='nw')¶
calculate the total induced heat in the nanostructure
- Parameters:
sim (
core.simulation
) – simulation descriptionfield_index (int) – index of evaluated self-consistent field to use for calculation. Can be obtained for specific parameter-set using
tools.get_closest_field_index()
power_scaling_e0 (float, default: 1) – incident laser power scaling. power_scaling_e0 = 1 corresponds to 1 mW per micron^2. See [1].
return_value (str, default: 'total') –
Values to be returned. Either ‘total’ (default) or ‘structure’.
”total” : return the total deposited heat in nW (float)
”structure” : return spatially resolved deposited heat at each meshpoint in nW (list of tuples [x,y,z,q])
return_units (str, default: "nw") – units of returned heat values, either “nw” or “uw” (nano or micro Watts)
- Returns:
Q –
`return_value`==”total” : (return float) total deposited heat in nanowatt (optionally in microwatt).
`return_value`==”structure” : (return list of tuples [x,y,z,q]) The returned quantity q is the total deposited heat at mesh-cell position [x,y,z] in nanowatt. To get the heating power-density, please divide by the mesh-cell volume.
- Return type:
float or list of tuples [x,y,z,q]
Notes
For details on heat/temperature calculations and raster-scan simulations, see:
[1] Baffou, G., Quidant, R. & Girard, C.: Heat generation in plasmonic nanostructures: Influence of morphology Applied Physics Letters 94, 153109 (2009)
[2] Teulle, A. et al.: Scanning optical microscopy modeling in nanoplasmonics Journal of the Optical Society of America B 29, 2431 (2012).
- pyGDM2.linear.temperature(sim, field_index, r_probe, kappa_env=0.6, kappa_subst=None, incident_power=1.0)¶
calculate the temperature rise at locations outside the nano-particle
Calculate the temperature increase close to a optically excited nanostructure using the approach described in [2] and [3] (Ref. [3] introduces a further correction term for particles lying on a substrate)
- Parameters:
sim (
core.simulation
) – simulation descriptionfield_index (int) – index of evaluated self-consistent field to use for calculation. Can be obtained for specific parameter-set using
tools.get_closest_field_index()
r_probe (tuple (x,y,z) or list of 3-lists/-tuples) – (list of) coordinate(s) to evaluate nearfield on. Format: tuple (x,y,z) or list of 3 lists: [Xmap, Ymap, Zmap] (the latter can be generated e.g. using
tools.generate_NF_map()
)kappa_env (float, default: 0.6) – heat conductivity of environment. default: kappa_env = 0.6 (water). (air: kappa=0.024, for more material values see e.g. [4]). In W/mK.
kappa_subst (float, default: None) – heat conductivity of substrate. default: None –> same as substrate. Using the mirror-image technique described in [3]. (glass: kappa=0.8)
incident_power (float, default: 1.0) – incident laser power density in mW per micron^2. See also [1].
- Returns:
if evaluating at a single position, D_T (float) – temperature increase in Kelvin at r_probe
if evaluating at a list of positions, list of tuples [x, y, z, D_T] – where D_T is the temperature increase in Kelvin at (x,y,z), which are the positions defined by r_probe
Notes
For details on heat/temperature calculations and raster-scan simulations, see:
[1] Baffou, G., Quidant, R. & Girard, C.: Heat generation in plasmonic nanostructures: Influence of morphology Applied Physics Letters 94, 153109 (2009)
[2] Baffou, G., Quidant, R. & Girard, C.: Thermoplasmonics modeling: A Green’s function approach Phys. Rev. B 82, 165424 (2010)
[3] Teulle, A. et al.: Scanning optical microscopy modeling in nanoplasmonics Journal of the Optical Society of America B 29, 2431 (2012).
For the thermal conductivity of common materials, see:
[4] Hugh D Young, Francis Weston Sears: University Physics, chapter 15, 8th. edition: Addison-Wesley, 1992 (see also e.g.: http://hyperphysics.phy-astr.gsu.edu/hbase/Tables/thrcn.html)
far-field¶
- pyGDM2.linear.extinct(sim, field_index, with_radiation_correction=False, normalization_E0=True)¶
Extinction, scattering and absorption cross-sections
Calculates extinction, scattering and absorption crosssections for each wavelength of the GDM simulation
Pure python implementation.
- Parameters:
sim (
core.simulation
) – simulation descriptionfield_index (int) – index of evaluated self-consistent field to use for calculation. Can be obtained for specific parameter-set using
tools.get_closest_field_index()
with_radiation_correction (bool, default: False) – Adds an optional radiative correction to the absorption section (hence it interferes also in the scattering section).See B. Draine, in: Astrophys. J. 333, 848–872 (1988), equation (3.06). Using the correction can lead to better agreement for very large discretization stepsizes, but has only a weak influence on simulations with fine meshes.
normalization_E0 (bool, default: True) – normalizes sections to peak of incident field intensity inside structure
- Returns:
ext-section (float) – extinction cross-section
scat-section (float) – scattering cross-section
abs-section (float) – apsorption cross-section
Notes
For the calculation of the cross-sections from the complex nearfield, see e.g.: Draine, B. T. & Flatau, P. J. Discrete-dipole approximation for scattering calculations. Journal of the Optical Society of America, A 11, 1491 (1994).
- pyGDM2.linear.farfield(sim, field_index, r_probe=None, r=100000.0, tetamin=0, tetamax=1.5707963267948966, Nteta=10, phimin=0, phimax=6.283185307179586, Nphi=36, polarizerangle='none', return_value='map', normalization_E0=False)¶
spatially resolved and polarization-filtered far-field scattering
For a given incident field, calculate the electro-magnetic field (E-component) in the far-field around the nanostructure (on a sphere of radius r).
Propagator for scattering into a substrate contributed by C. Majorel
Pure python implementation.
- Parameters:
sim (
core.simulation
) – simulation descriptionfield_index (int) – index of evaluated self-consistent field to use for calculation. Can be obtained for specific parameter-set using
tools.get_closest_field_index()
r_probe (tuple (x,y,z) or list of 3-lists/-tuples. optional. Default: don't use) – defaults to None, which means it is not used and a solid angle defined by a spherical coordinate range is used instead. If r_probe is given, this overrides r, tetamin, tetamax, Nteta, Nphi. (list of) coordinate(s) to evaluate farfield on. Format: tuple (x,y,z) or list of 3 lists: [Xmap, Ymap, Zmap]
r (float, default: 100000.) – radius of integration sphere (distance to coordinate origin in nm)
tetamin (float, float; defaults: 0, np.pi/2) – minimum and maximum polar angle in radians (in linear steps from tetamin to tetamax)
tetamax (float, float; defaults: 0, np.pi/2) – minimum and maximum polar angle in radians (in linear steps from tetamin to tetamax)
phimin (float, float; defaults: 0, 2*np.pi) – minimum and maximum azimuth angle in radians, excluding last position (in linear steps from phimin to phimax)
phimax (float, float; defaults: 0, 2*np.pi) – minimum and maximum azimuth angle in radians, excluding last position (in linear steps from phimin to phimax)
Nteta (int, int; defaults: 10, 36) – number of polar and azimuthal angles on sphere to calculate,
Nphi (int, int; defaults: 10, 36) – number of polar and azimuthal angles on sphere to calculate,
polarizerangle (float or 'none', default: 'none') – optional polarization filter angle **in degrees**(!). If ‘none’ (default), the total field-intensity is calculated (= no polarization filter)
return_value (str, default: 'map') –
- Values to be returned. Either ‘map’ (default) or ‘integrated’.
”map” : (default) return spatially resolved farfield intensity at each spherical coordinate (5 lists)
”efield” : return spatially resolved E-fields at each spherical coordinate (5 lists)
”int_Es” : return the integrated scattered field (as float)
”int_E0” : return the integrated fundamental field (as float)
”int_Etot” : return the integrated total field (as float)
normalization_E0 (bool, default: False) – has only effect on return_value==”int_Es”: Normalizes scattering to peak of incident field intensity inside structure
- Returns:
using r_probe for position definition –
- 3 lists of 6-tuples (x,y,z, Ex,Ey,Ez), complex :
scattered Efield or E-intensity
total Efield (inlcuding fundamental field) or total E intensity
fundamental Efield (incident field) or E0 intensity
if solid angle is defined via spherical coordinate range –
- return_value == “map”5 arrays of shape (Nteta, Nphi) :
[tetalist : teta angles - if solid angle range input]
[philist : phi angles - if solid angle range input]
I_sc : intensity of scattered field,
I_tot : intensity of total field (I_tot=|E_sc+E_0|^2),
I0 : intensity of incident field
- return_value == “efield”float
[tetalist : teta angles - if solid angle range input]
[philist : phi angles - if solid angle range input]
E_sc : complex scattered field at each pos.
E_tot : complex total field at each pos. (E_sc+E0)
E0 : complex incident field at each pos.
- return_value == “int_XX”float
integrated total intensity over specified solid angle
Notes
For details of the asymptotic (non-retarded) far-field propagators for a dipole above a substrate, see e.g.: Colas des Francs, G. & Girard, C. & Dereux, A. Theory of near-field optical imaging with a single molecule as light source. The Journal of Chemical Physics 117, 4659–4666 (2002)
multipole decomposition¶
- pyGDM2.multipole.multipole_decomposition_exact(sim, field_index, r0=None, epsilon=0.01, which_moments=['p', 'm', 'qe', 'qm'], long_wavelength_approx=False)¶
exact multipole decomposition of the nanostructure optical response
** ——- FUNCTION STILL UNDER TESTING ——- **
Multipole decomposition of electromagnetic field inside nanostructure for electric and magnetic dipole and quadrupole moments.
- Parameters:
sim (
core.simulation
) – simulation descriptionfield_index (int) – index of evaluated self-consistent field to use for calculation. Can be obtained for specific parameter-set using
tools.get_closest_field_index()
r0 (array, default: None) – [x,y,z] position of multipole decomposition development. If None, use structure’s center of gravity
epsilon (float, default: 0.01) – additional step on r0 (in nm) to avoid numerical divergence of the Bessel terms
which_moments (list of str, default: ['p', 'm', 'qe', 'qm']) –
- which multipole moments to calculate and return. supported dipole moments:
’p’: electric dipole (full)
’m’: magnetic dipole
’qe’: electric quadrupole (full)
’qm’: magnetic quadrupole
’p1’: electric dipole (first order)
’pt’: toroidal dipole
’qe1’: electric quadrupole (first order)
’qet’: toroidal quadrupole
long_wavelength_approx (bool, default: False) – if True, use long wavelength approximation
- Returns:
list of multipole moments. Default
p (3-vector) – electric dipole moment
m (3-vector) – magnetic dipole moment
Qe (3x3 tensor) – electric quadrupole moment
Qm (3x3 tensor) – magnetic quadrupole moment
Notes
For details about the method, see:
Alaee, R., Rockstuhl, C. & Fernandez-Corbaton, I. An electromagnetic multipole expansion beyond the long-wavelength approximation. Optics Communications 407, 17–21 (2018)
- pyGDM2.multipole.extinct(sim, field_index, with_toroidal=True, r0=None, eps_dd=0.001, use_generalized_polarizabilities=False, normalization_E0=True, long_wavelength_approx=False)¶
extinction cross sections from multipole decomposition
** ——- FUNCTION STILL UNDER TESTING ——- **
Returns extinction cross sections for electric and magnetic dipole and quadrupole moments of the multipole decomposition.
- Parameters:
sim (
core.simulation
) – simulation descriptionfield_index (int) – index of evaluated self-consistent field to use for calculation. Can be obtained for specific parameter-set using
tools.get_closest_field_index()
with_toroidal (bool, default: True) – whether to add toroidal moments to electric dipole and quadrupole
r0 (array, default: None) – [x,y,z] position of mulipole decomposition development. If None, use structure’s center of gravity
eps_dd (float, default: 0.1) – numerical integration step (in nm). Used for e/m quadrupole extinction.
normalization_E0 (bool, default: True) – normalizes sections to max. incident field intensity
long_wavelength_approx (bool, default: False) – if True, use long wavelength approximation
- Returns:
sigma_ext_p (float) – electric dipole extinction cross section (in nm^2)
sigma_ext_m (float) – magnetic dipole extinction cross section (in nm^2)
sigma_ext_q (float) – electric quadrupole extinction cross section (in nm^2)
sigma_ext_mq (float) – magnetic quadrupole extinction cross section (in nm^2)
Notes
For details about the extinction section of multipole moments, see:
Evlyukhin, A. B. et al. Multipole analysis of light scattering by arbitrary-shaped nanoparticles on a plane surface., JOSA B 30, 2589 (2013)
- pyGDM2.multipole.scs(sim, field_index, with_toroidal=True, use_generalized_polarizabilities=False, r0=None, normalization_E0=True, long_wavelength_approx=False)¶
total scattering cross section from multipole decomposition
** ——- FUNCTION STILL UNDER TESTING ——- **
Returns scattering cross sections for electric and magnetic dipole and quadrupole moments of the multipole decomposition.
- Parameters:
sim (
core.simulation
) – simulation descriptionfield_index (int) – index of evaluated self-consistent field to use for calculation. Can be obtained for specific parameter-set using
tools.get_closest_field_index()
with_toroidal (bool, default: True) – whether to add toroidal moments to electric dipole and quadrupole
use_generalized_polarizabilities (bool, default: False) –
if True, use generalized polarizabilities (does not require evaluation of core.scatter for new incident field
configurations, first calculation is more expensive)
r0 (array, default: None) – [x,y,z] position of mulipole decomposition development. If None, use structure’s center of gravity
normalization_E0 (bool, default: True) – normalizes sections to max. incident field intensity
long_wavelength_approx (bool, default: False) – if True, use long wavelength approximation
- Returns:
sigma_scat_p (float) – electric dipole scattering cross section (in nm^2)
sigma_scat_m (float) – magnetic dipole scattering cross section (in nm^2)
sigma_scat_q (float) – electric quadrupole scattering cross section (in nm^2)
sigma_scat_mq (float) – magnetic quadrupole scattering cross section (in nm^2)
Notes
For details about the exact multipole formalism and scs calculation, see:
Alaee, R., Rockstuhl, C. & Fernandez-Corbaton, I. An electromagnetic multipole expansion beyond the long-wavelength approximation. Optics Communications 407, 17–21 (2018)
- pyGDM2.multipole.generalized_polarizability(sim, field_index=None, wavelength=None, method='lu', epsilon=0.01, r0=None, which_moments=['p1', 'pt', 'qe', 'qt', 'm', 'qm'], long_wavelength_approx=False, verbose=1)¶
generalized electric and magnetic polarizabilities
** ——- FUNCTION STILL UNDER TESTING ——- **
Returns the generalized polarizability tensors that can be used with arbitrary, inhomogeneous illumination fields to calculate the effective electric and magnetic multipole moments, induced in the nanostructure.
- Parameters:
sim (
core.simulation
) – simulation descriptionfield_index (int, default: None) – index of evaluated self-consistent field to use for calculation. Can be obtained for specific parameter-set using
tools.get_closest_field_index()
. Either field_index or wavelength must be given.wavelength (float, default: None) – Optional wavelength (alternative to field_index) at which to calculate susceptibility matrix (in nm). Either field_index or wavelength must be given.
method (string, default: "scipyinv") –
- inversion method. One of [“lu”, “numpyinv”, “scipyinv”, “cupy”, “cuda”]
”scipyinv” scipy default inversion (scipy.linalg.inv)
”numpyinv” numpy inversion (np.linalg.inv, if numpy compiled accordingly: LAPACK’s dgesv)
”cupy” uses CUDA GPU via cupy
”cuda” (equivalent to “cupy”)
”lu” LU-decomposition (scipy.linalg.lu_factor) - inefficient for decay_rate!
epsilon (float, default: 0.01) – additional step on r0 (in nm) to avoid numerical divergence of the Bessel terms
r0 (array, default: None) – [x,y,z] position of mulipole decomposition development. If None, use structure’s center of gravity
which_moments (list of str, default: ['p1', 'pt', 'qe1', 'qt', 'm', 'qm']) –
- which generalized polarizability tensors to calculate and return. supported:
’p1’: electric dipole - only first order (rank 2)
’pt’: toroidal dipole (rank 3)
’qe1’: electric quadrupole - only first order (rank 4)
’qt’: toroidal electric quadrupole (rank 4)
’m’: magnetic dipole (rank 2)
’qm’: magnetic quadrupole (rank 3)
long_wavelength_approx (bool, default: False) – if True, use long wavelength approximation
verbose (bool default=True) – print runtime info
- Returns:
by default 6 lists of N tensors, with N the number of discretization cells (see kwarg which_moments)
K_P_E, K_T_E, K_QE_E, K_QT_E, K_M_E, K_QM_E
Notes
For details , see:
PAPER SUBMITTED
For details about the underlying exact multipole decomposition, see:
Alaee, R., Rockstuhl, C. & Fernandez-Corbaton, I. An electromagnetic multipole expansion beyond the long-wavelength approximation. Optics Communications 407, 17–21 (2018)
- pyGDM2.multipole.extract_effective_polarizability(sim, method='lu', which_moments=['p1', 'm'], long_wavelength_approx=True, illumination_mode='dipole', npoints=25, r_sphere=1000, store_simulation_object=False)¶
Extract effective electric and magnetic dipole polarizability for structure
solve inverse problem of adjusting polarizability for different illuminations via pseudoinverse
doc to be completed
- store_simulation_objectBool, default: False
optionally store full sim in output dictionnary. Caution: This contains numba compiled objects which are not necessarily reloadable on other operating systems than the one used for the simulation
nonlinear¶
nonlinear contains evaluation functions for non-linear effects.
TPL and SP-LDOS¶
- pyGDM2.nonlinear.tpl_ldos(sim, field_index, nonlin_order=2.0, beta=100000.0, verbose=False)¶
calculate the TPL for nano-object under given illumination condition
Calculate the two-photon photolouminescence (TPL). Higher order nonlinear photoluminescence can also be calculated by changing the parameter nonlin_order.
- Parameters:
sim (
core.simulation
) – simulation descriptionfield_index (int) – index of evaluated self-consistent field to use for calculation. Can be obtained for specific parameter-set using
tools.get_closest_field_index()
nonlin_order (float, default: 2.0) – order of nonlinear response. (I_npl ~ | |E|^2 | ^ nonlin_order).
beta (float, default: 1E5) – incident laser power scaling. Use 1E5 for correct LDOS values.
verbose (bool, default: False) – Enable some info printing
- Returns:
I_tpl – TPL intensity in farfield from nano-object under given illumination
- Return type:
float
Notes
For details on TPL/LDOS calculation via focused beam rasterscan simulations, see: Viarbitskaya, S. et al. Tailoring and imaging the plasmonic local density of states in crystalline nanoprisms, Nat. Mater. 12, 426–432 (2013)
fast electrons¶
electron contains evaluation functions for simulations using fast electron illumination.
electron beam incident field¶
- pyGDM2.fields.fast_electron(pos, env_dict, wavelength, electron_kinetic_energy, x0, y0, kSign=-1, avoid_div_lim_distance=10.0, returnField='E')¶
Electric field created by a fast electron moving along (OZ)
The electron beam crosses the (OXY) plane in (x0, y0)
- Parameters:
pos (np.array) – list of 3-tuple coordinates to evaluate field at: [[x1,y1,z1], [x2,y2,z2], … ]
env_dict (dict) – Must be compatible with sim.dyads.getConfigDictG typed numba dict. description of environment. Must contain either “eps_env” or [“eps1”, “eps2”].
wavelength (float) – Wavelength in nm
electron_kinetic_energy (float) – electron kinetic energy (keV)
x0 (float) – position of the electron beam (nm)
y0 (float) – position of the electron beam (nm)
kSign (int, default: -1) – sign of wavenumber. +1: electron propagation from bottom to top (towards increasing z) -1: electron propagation from top to bottom (towards smaller z, default) either kSign or k0 must be given.
avoid_div_lim_distance (float, default: 10.0) – set a min. distance (in nm) between the location where the E-field is computed and the electron trajectory to avoid divergence
returnField (str, default: 'E') – if ‘E’: returns electric field; if ‘B’: magnetic field
- Returns:
E0 (B0) – list of (complex) 3-tuples: [(Ex1, Ey1, Ez1), …]
- Return type:
Complex E-(B-)Field at each dipole position as
EELS and cathodoluminescence (CL)¶
- pyGDM2.electron.EELS(sim, field_index)¶
Electron Energy Loss probability of a fast electron moving along (OZ)
The E-field created by the fast electron is given by the function fields.fast_electron. Electron kinetic energy, beam position and propagation direction are parameters of fields.fast_electron.
- Parameters:
sim (
core.simulation
) – simulation descriptionfield_index (int) – index of evaluated self-consistent field to use for calculation. Can be obtained for specific parameter-set using
tools.get_closest_field_index()
- Returns:
eels_signal
- Return type:
electron energy loss probability per eV per electron
Notes
For details of the EELS computation: Arnaud Arbouet et al 2014 New J. Phys. 16 113012
- pyGDM2.electron.CL(sim, field_index, r=10000.0, tetamin=0, tetamax=3.141592653589793, Nteta=30, phimin=0, phimax=6.283185307179586, Nphi=36, polarizerangle='none', return_value='int_Es')¶
Cathodoluminescence emission generated by a fast electron moving along (OZ)
The fundamental E-field created by the fast electron is given by the function
pyGDM2.fields.fast_electron()
. Electron kinetic energy, beam position and propagation direction are parameters of fields.fast_electron.Using
pyGDM2.linear_py.farfield()
for calculation of farfield scattering.For further documentation on the possible kwargs, their meaning and returned values, see doc of
pyGDM2.linear_py.farfield()
- Parameters:
sim (
core.simulation
) – simulation descriptionfield_index (int) – index of evaluated self-consistent field to use for calculation. Can be obtained for specific parameter-set using
tools.get_closest_field_index()
- Returns:
CL_signal – number of emitted photons per eV per electron
- Return type:
Cathodoluminescence emission intensity
Notes
For details of the EELS computation: Arnaud Arbouet et al. Electron energy losses and cathodoluminescence from complex plasmonic nanostructures: spectra, maps and radiation patterns from a generalized field propagator, New J. Phys. 16, 113012 (2014)
tools¶
- pyGDM2.electron.electron_speed(electron_kinetic_energy)¶
Velocity of electron with kinetic energy electron_kinetic_energy
- Parameters:
electron_kinetic_energy (float) – electron kinetic energy (keV)
- Returns:
velec
- Return type:
speed of electron in units of the speed of light
- pyGDM2.electron.visu_structure_electron(struct, x0=None, y0=None, scale_e_pos=10, marker_e_pos='x', color_e_pos='r', show=True, **kwargs)¶
plot 2D structure projection together with electron beam position
all further kwargs are passed to pyGDM2.visu.structure
- Parameters:
struct (list or
core.simulation
) – either list of 3d coordinate tuples or simulation description objectx0 (float, default: None) – position of the electron beam in the (OXY) plane. If None and simulation instance is passed as structure, try to infer electron beam position(s) from incident field.
y0 (float, default: None) – position of the electron beam in the (OXY) plane. If None and simulation instance is passed as structure, try to infer electron beam position(s) from incident field.
scale_e_pos (float, default: 15) – scale of electron-position marker
marker_e_pos (str, default: 'x') – maker symbol for electron-position
color_e_pos (str, default: 'r') – color of electron-position marker
- Return type:
result returned by matplotlib’s scatter