Core

core

core contains the main simulation object and the scatter routine to calculate the self-consistent fields inside the nano-structure.

simulation description class

class pyGDM2.core.simulation(struct, efield, dyads=None, dtype='f', verbose=True)

Bases: object

Main GDM simulation container class - pure-python API

Defines a linear GDM simulation. Contains information on:
  • structstructures.struct:
    • the geometry of the nano structure

    • its dielectric function

  • efieldfields.efield
    • the incident field and the wavelenghts

    • possibly further incident field configuration

  • dyadspropagators.DyadsBaseClass
    • class derived from propagators.DyadsBaseClass

    • contains set of Green’s tensors, self-terms and related functions

    • contains environment definition and polarizability calculation

Parameters:
  • struct (structures.struct) – structure object

  • efield (fields.efield) – fundamental field

  • dyads (propagators.DyadsBaseClass) – set of Green’s tensors, polarizabilities, environment info

  • dtype ((optional) str, default: ‘f’) – precision of simulation

__init__(struct, efield, dyads=None, dtype='f', verbose=True)

Initialization

scatter(**kwargs)

wrapper to scatter()

simulate electromagnetic fields inside particles

pyGDM2.core.scatter(sim, method='lu', calc_H=False, verbose=True, callback=None, **kwargs)

Perform a linear scattering GDM simulation

Calculate the electric field distribution inside a nano-structure. Optionally calculate also the internal magnetic field.

Parameters:
  • sim (simulation) – simulation description

  • method (string, default: "lu") –

    inversion method. One of [“lu”, “numpyinv”, “scipyinv”, “cupy”, “cuda”]
    • ”lu” LU-decomposition (scipy.linalg.lu_factor)

    • ”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”)

  • calc_H (bool, default: False) – if True, calculate also the internal magnetic field H. This has an impact on performance and memory requirement, so it is deactivated by default.

  • verbose (bool, default False) – Print timing info to stdout

  • callback (func, default: None) – optional callback function, which is called after each wavelength. Passes a dict to callback containing current wavelength, simulation and timing info. Available dict_keys: [‘i_wl’, ‘wavelength’, ‘sim’, ‘t_inverse’, ‘t_repropa’]. Callback needs to return True for simulation to continue. If callback returns False the simulation will be canceled, and scatter returns the simulation results until the moment when it was interrupted.

Returns:

  • 1, if finished succesfully

  • -1, if canceled

Notes

The complex electric fields inside the structure are also copied into the simulation instance as attribute simulation.E

For details on the concept of the generalized propagator, see e.g.: Martin, O. J. F. & Girard, C. & Dereux, A. Generalized Field Propagator for Electromagnetic Scattering and Light Confinement. Phys. Rev. Lett. 74, 526–529 (1995).

The scipy solvers (like ‘lu’, ‘ilu’ and ‘cg’) run parallelized if BLAS is compiled with multithreading support. See: http://scipy-cookbook.readthedocs.io/items/ParallelProgramming.html#use-parallel-primitives (last called 05/2020)

To limit the number of threads for the multi-threaded parallelized parts, you might do something like explained in: https://stackoverflow.com/questions/29559338/set-max-number-of-threads-at-runtime-on-numpy-openblas (website called 05/2020). Or use the threadpoolctl package (example for 4 threads:):

>>> from threadpoolctl import threadpool_limits
>>> threadpool_limits(limits=4, user_api='blas')
pyGDM2.core.scatter_mpi(sim, calc_H=False, verbose=False, scatter_verbose=False, **kwargs)

MPI wrapper to scatter() for embarrassingly parallel calculation of spectra

requires: mpi4py

run with “mpirun -n X python scriptname.py”, where X is the number of parallel processes (ideally an integer divisor of the number of wavelengths)

Parameters:
  • sim (simulation) – simulation description

  • calc_H (bool, default: False) – if True, calculate also the internal magnetic field H.

  • verbose (bool, default: False) – turns on some mpi-routine info printing, respectively controls verbose setting for scatter()

  • scatter_verbose (bool, default: False) – turns on some mpi-routine info printing, respectively controls verbose setting for scatter()

  • **kwargs – all kwargs are passed to scatter()

Notes

  • On single machines it is usually easier to install scipy compiled with parallel BLAS (parallel LU / CG routines). Usually the parallel BLAS is installed by default. Using both parallelisation techniques simultaneously requires proper configuration of threads / process numbers. Overloading the CPUs will usually result in **decreased* calculation speed!*

  • see scatter() for main documentation

decay rate of dipole transitions / LDOS

pyGDM2.core.decay_rate(sim, field_index=None, wavelength=None, r_probe=None, r_emitter=None, component='E', return_value='decay_rates', method='scipyinv', verbose=True)

local decay rate modification of a electric or magnetic dipole transition

Parameters:
  • sim (simulation) – simulation description

  • field_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.

  • r_probe (tuple (x,y,z) or list of 3-lists/-tuples) – (list of) coordinate(s) to evaluate the decay rate 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())

  • r_emitter (tuple (x,y,z), default: None) – optional fixed coordinate of the emitter, use to evaluate the cross-DOS. If given, the Cross-DOS due to an emitter at r_emitter is evaluated at all positions r_probe, hence the emitter will be fixed. See: Cazé et al. PRL 110, 063903 (2013)

  • component (str, default: 'E') – “E” or “H”: which LDOS component to calculate (electric or magnetic).

  • return_value (str, default: 'decay_rates') –

    Values to be returned. one of:
    • ’decay_rates’: relative decayrates for X,Y and Z oriented dipole

    • ’decay_rate_x’: relative, partial decayrate for X-oriented dipole

    • ’decay_rate_y’: relative, partial decayrate for Y-oriented dipole

    • ’decay_rate_z’: relative, partial decayrate for Z-oriented dipole

    • ’decay_rate_total’: relative, orientation-averaged decayrate (–> propto Im(Tr(Sp)) )

    • ’field_susceptibility’: complex field-susceptibility tensor Sp at each r_probe position

  • 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!

  • verbose (bool default=True) – print runtime info

Returns:

see parameter return_value

Notes

For details about the underlying formalism, see: Wiecha, P. R., Girard, C., Cuche, A., Paillard, V. & Arbouet, A. Decay Rate of Magnetic Dipoles near Non-magnetic Nanostructures. Phys. Rev. B 97(8), 085411 (2018)

For details on the concept of the Cross-DOS, see: A. Cazé, R. Pierrat & R. Carminati Spatial Coherence in Complex Photonic and Plasmonic Systems Phys. Rev. Lett. 110, 063903 (2013)

Other functions

pyGDM2.core.get_general_propagator(sim=None, wavelength=None, method='lu', fallbackmethod_gpu='lu', calc_H=False, verbose=False)

invert dipole-coupling matrix

Parameters:
  • sim (simulation) – simulation description

  • wavelength (float) – Wavelength at which to calculate susceptibility matrix (in nm)

  • method (string, default: "lu") –

    inversion method. One of [“lu”, “numpyinv”, “scipyinv”, “cupy”, “cuda”]
    • ”lu” LU-decomposition (scipy.linalg.lu_factor)

    • ”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”)

  • fallbackmethod_gpu (string, default: "lu") – method to fall back to if using cuda solver and GPU runs out of RAM

  • calc_H (bool, default: False) – if True, calculate also the internal magnetic field H. This has an impact on performance and memory requirement, so it is deactivated by default.

  • verbose (bool, default False) – Print timing info to stdout

Returns:

- K

Return type:

Generalized Propagator

Notes

For details on the concept of the generalized propagator, see e.g.: Martin, O. J. F. & Girard, C. & Dereux, A. Generalized Field Propagator for Electromagnetic Scattering and Light Confinement. Phys. Rev. Lett. 74, 526–529 (1995).

For the Electric-magnetic mixed field calculation, see e.g.: Schröter, U. Modelling of magnetic effects in near-field optics. Eur. Phys. J. B 33, 297–310 (2003).

pyGDM2.core.get_SBS_EE(sim, wavelength, invertible=True)
pyGDM2.core.get_SBS_HE(sim, wavelength)