Response objects#
strongcoca can use response functions from several different
sources. Internally these response functions are represented via classes that
are derived from BaseResponse. The objects are usually prepared by IO
functions that provide smooth interfaces to different codes.
IO functions#
- strongcoca.response.read_gpaw_tddft(dipole_files, broadening, atoms_file=None, name='GPAWTDDFT')[source]#
Builds a
TimeDomainResponseinstance using time-dependent dipole moment data from GPAW real-time time-dependent density functional theory (RT-TDDFT) calculations.- Parameters:
dipole_files (
Union[Path,str,List[Union[Path,str]]]) – Name(s) of file(s) with time-dependent dipole moment data; each file is assumed to provide a different Cartesian direction; if there is only one file, the response is assumed to be diagonal.broadening (
Broadening) – Artificial broadening applied to the response.atoms_file (
Union[Path,str,None]) – Name of file from which to read atomic coordinates; needs to be readable by ASE.name (
str) – Name of the response.
- Returns:
Time domain response.
- Return type:
response
- Raises:
ValueError – If dipole moment files are not suitable or mutually compatible.
- strongcoca.response.read_nwchem_casida(filename, broadening=NoArtificialBroadening, name='NWChem Casida Response')[source]#
Returns a freshly initialized
CasidaResponseobject using an NWChem output file to set up the components of the Casida representation.- Parameters:
filename (
Union[Path,str]) – Name of NWChem output file; can be compressed.broadening (
Broadening) – Artificial broadening used for the continuous response; defaults to no broadening.name (
str) – Name of subsystem.
- Return type:
- strongcoca.response.read_dielec_function_file(dielectric_function_file, frequency_column=0, dielectric_function_real_column=1, dielectric_function_imag_column=2, units='eV', name='NumericDielectricFunction from file', **np_loadtxt_kwargs)[source]#
Parses a textfile with columns describing the dielectric function of a material.
- Parameters:
dielectric_function_file (
Union[Path,str]) – Path to the file containing the dielectic function.frequency_column (
int) – Index of column containing the frequency vector in units of eV.dielectric_function_real_column (
int) – Index of column containing the real part of the unitless dielectric function.dielectric_function_imag_column (
int) – Index of column containing the imaginary part of the unitless dielectric function.units (
str) – Units of the frequency column; passed tostrongcoca.response.NumericDielectricFunction.name (
str) – Passed toNumericDielectricFunction.np_loadtxt_kwargs – Parameters passed to
numpy.loadtxt(); useful parameters aredelimiterorcomments.
- Return type:
- strongcoca.response.read_dielec_function_from_refractive_index_file(dielectric_function_file, frequency_column=0, refractive_index_column=1, extinction_coefficient_column=2, units='eV', name='NumericDielectricFunction from file', **np_loadtxt_kwargs)[source]#
Parses a textfile with columns describing the refractive index and extinction coefficient of a material.
The refractive index \(n(\omega)\) and extinction coefficient \(k(\omega)\) are related to the complex dielectric function \(\varepsilon_r(\omega)\) by
\[\varepsilon_r(\omega) = n^2(\omega) - k^2(\omega) + 2jn(\omega)k(\omega)\]- Parameters:
dielectric_function_file (
Union[Path,str]) – Path to the file containing the dielectic function.frequency_column (
int) – Index of column containing the frequency vector in units of eV.refractive_index_column (
int) – Index of column containing the unitless refractive index.extinction_coefficient_column (
int) – Index of column containing the unitless extinction coefficient.units (
str) – Units of the frequency column; passed toNumericDielectricFunction.name (
str) – Passed toNumericDielectricFunction.np_loadtxt_kwargs – Parameters passed to
numpy.loadtxt(); useful parameters aredelimiterorcomments.
- Return type:
Builder functions#
- strongcoca.response.build_diagonal_casida(casida_response, name='Diagonal')[source]#
Builds a
CasidaResponseinstance in the diagonal representation from an existing non-diagonal representation.In the diagonal representation, the K matrix is zero and the eigenstates are on the diagonal.
- Parameters:
casida_response (
CasidaResponse) – Old Casida response.name (
str) – Name of the new response.
- Return type:
- strongcoca.response.build_random_casida(n_states, random_seed=42, broadening=NoArtificialBroadening, name='Random')[source]#
Builds Casida response from random numbers.
This is intended to be used for sandbox applications and testing.
- Parameters:
n_states (
int) – Number of states.random_seed (
int) – Seed for the random number generator.broadening (
Broadening) – Artificial broadening used for the continuous response; defaults to no broadening.name (
str) – Name of response function.
- Return type:
Examples
The following code snippet shows how to set up a randomized Casida matrix with 3 states.
>>> from strongcoca.response import build_random_casida >>> response = build_random_casida(n_states=3, name='Elephant', random_seed=42)
- strongcoca.response.build_time_domain_response_from_casida(casida_response, dt=None, n_steps=None, units='fs', name='PropagatedCasida')[source]#
Builds a
TimeDomainResponseinstance from an existing Casida representation via time propagation.- Parameters:
casida_response (
CasidaResponse) – Response in Casida representation.dt (
Optional[float]) – Time step; units are fs by default, optionally atomic units (seeunits).n_steps (
Optional[int]) – Number of time steps for which to propagate the equations of motion.units (
str) –fsto specify dt in fs orauto specifydtin atomic units.
- Return type:
CasidaResponse#
- class strongcoca.response.CasidaResponse(D, K, mu, broadening=NoArtificialBroadening, occ=None, units='eVA', name='Casida')[source]#
Objects of this class hold different representations of the response functions for a Casida system.
- Parameters:
D (
ndarray) –Diagonal part of decomposition of Casida matrix (initial-final energy difference of uncoupled single particle levels).
Units are eV by default, optionally atomic units (see
units).K (
ndarray) –Coupling matrix of decomposition of Casida matrix (electron-hole coupling).
Units are eV by default, optionally atomic units (see
units).mu (
ndarray) –Transition dipoles.
Units are eÅ by default, optionally atomic units (see
units).occ (
Optional[ndarray]) – Occupation number differences, defaults to one, unitless.broadening (
Broadening) – Artificial broadening used for the continuous response; defaults to no broadening.units (
str) –eVAto specify D, K in eV and mu in eÅ orauto specify inputs in atomic units.This parameter determines whether conversion should be performed during initialization and has no effect on instance methods and variables.
name (
str) – Name of response functioN.
Examples
The following snippet illustrates how a Casida response object can be constructed. Here, artificial Casida matrix components are used. In practice the Casida matrix would usually be imported from, e.g., NWChem calculations:
>>> import numpy as np >>> from strongcoca.response import CasidaResponse >>> >>> D = np.array([0.1, 0.3]) >>> K = np.zeros((2, 2)) >>> mu = np.array([[0.1, 0.0, 0.0], [0.0, 0.1, 0.0]]) >>> response = CasidaResponse(D, K, mu, name='Panda') >>> print(response) ---------------------- CasidaResponse ---------------------- name : Panda n_states : 2 D (eV) : [0.1 0.3] K (eV) : [[0. 0.] [0. 0.]] mu (eÅ) : [[0.1 0. 0. ] [0. 0.1 0. ]] broadening : NoArtificialBroadening atoms : None
- property excitations: Excitations#
Excitations as
Excitations.
- property n_states: int#
Number of states in Casida system.
MieGansResponse#
- class strongcoca.response.MieGansResponse(semiaxes, dielectric_function, name='MieGansResponse')[source]#
Objects of this class hold representations of the Mie-Gans response of a ellipsoid in vacuum.
For an ellipsoid with semi-axes \(a_x, a_y, a_z\) in the Cartesian directions, the Mie-Gans polarizability is diagonal [Sihvola, J. Nanomater. 2007, 045090]
\[\alpha_{\mu\mu}(\omega) = \varepsilon_0 V \frac{\varepsilon_r(\omega) - 1}{1 + N_\mu (\varepsilon_r(\omega) - 1)}\]where \(N_\mu\) are depolarization factors
\[N_\mu = \frac{a_x a_y a_z}{2} \int_0^\infty \frac{\mathrm{d}s}{(s+a_\mu^2)\sqrt{(s+a_x^2)(s+a_y^2)(s+a_z^2)}}.\]It is useful to keep in mind that in Hartree atomic units the vacuum permitivity
\[\varepsilon_0 = \frac{1}{4\pi}.\]Thus the expression is simplified
\[\alpha_{\mu\mu}(\omega) = \frac{a_x a_y a_z}{3} \frac{\varepsilon_r(\omega) - 1}{1 + N_\mu (\varepsilon_r(\omega) - 1)}\]- Parameters:
semiaxes (
Union[ndarray,Sequence[float]]) –Semi-axes of the ellipsoid in Cartesian coordinates in units of Å; shape {3}.
Optionally radius of sphere in units of Å; shape {1}.
dielectric_function (
DielectricFunction) – Dielectric function, which can have an analytic form or be sampled on a frequency grid.name (
str) – Name of response.
- Raises:
Examples
The following snippet illustrates how a Mie-Gans response object can be constructed. Here, the analytic Drude dielectric function of the free-electron gas is used. Instead, a numeric dielectric function can be used with the help of
read_dielec_function_file().>>> import numpy as np >>> from strongcoca.response import MieGansResponse, DrudeDielectricFunction >>> >>> semiaxes = 20 >>> omega_p = 3 >>> gamma = 0.5 >>> eps = DrudeDielectricFunction(omega_p, gamma) >>> response = MieGansResponse(semiaxes, eps, name='Leaves') >>> print(response) --------------------- MieGansResponse ---------------------- name : Leaves semiaxes (Å) : [20. 20. 20.] depol. factors : [0.33335 0.33335 0.33335] Dielectric function : DrudeDielectricFunction
Different shapes of the ellipsoid can be set by supplying three cartesian semiaxes. Here a sphere, a prolate ellipsoid and an oblate ellipsoid are created.
>>> response = MieGansResponse([30, 30, 30], eps, name='Sphere') >>> print(response) --------------------- MieGansResponse ---------------------- name : Sphere semiaxes (Å) : [30. 30. 30.] depol. factors : [0.33335 0.33335 0.33335] Dielectric function : DrudeDielectricFunction
>>> response = MieGansResponse([25, 25, 35], eps, name='Prolate ellipsoid') >>> print(response) --------------------- MieGansResponse ---------------------- name : Prolate ellipsoid semiaxes (Å) : [25. 25. 35.] depol. factors : [0.37561 0.37561 0.24881] Dielectric function : DrudeDielectricFunction
>>> response = MieGansResponse([30, 30, 20], eps, name='Oblate ellipsoid') >>> print(response) --------------------- MieGansResponse ---------------------- name : Oblate ellipsoid semiaxes (Å) : [30. 30. 20.] depol. factors : [0.27705 0.27705 0.44592] Dielectric function : DrudeDielectricFunction
- property dielectric_function: DielectricFunction#
Dielectric function of the Mie-Gans ellipsoid.
TimeDomainResponse#
- class strongcoca.response.TimeDomainResponse(times, dipole_moments, broadening, units='fseA', name='TimeDomain')[source]#
Objects of this class hold different representations of the response functions for a time domain system.
- Parameters:
times (
ndarray) – Time grid of {t} equally spaced time values units are fs by default, optionally atomic units (seeunits).dipole_moments (
ndarray) – Dipole moment tensors at the time values, shape {t}x{3}x{3} units are eÅ by default, optionally atomic units (seeunits).broadening (
Broadening) – artificial broadening applied to the responseunits (
str) –fseAto specify times in fs and dipole_moments in eÅ orauto specify inputs in atomic units.This parameter determines whether conversion should be performed during initialization and has no effect on instance methods and variables.
name (
str) – Name of response.
Examples
The following snippet illustrates how a time domain response object can be constructed. Here, artificial dipole moments are used. In practice dipole moments would usually be imported from, e.g., TDDFT calculations:
>>> import numpy as np >>> from strongcoca.response import TimeDomainResponse >>> from strongcoca.response.utilities import LorentzianBroadening >>> >>> t = np.linspace(0, 10) >>> dm = np.ones((t.shape[0], 3, 3)) >>> broadening = LorentzianBroadening(0.1) >>> response = TimeDomainResponse(t, dm, broadening, name='Eats') >>> print(response) -------------------- TimeDomainResponse -------------------- name : Eats times (fs) : [ 0. 0.20408 0.40816 ... 9.59184 9.79592 10. ] dipole_moments (eÅ) : [[[1. 1. 1.] [1. 1. 1.] [1. 1. 1.]] ...
Dielectric functions#
- class strongcoca.response.dielectric.DielectricFunction(name='DielectricFunction')[source]#
Objects of this class hold different representations of dielectric functions.
- Parameters:
name (
str) – Name of dielectric function.
- property name: str#
Name of dielectric function.
- class strongcoca.response.DrudeDielectricFunction(plasma_frequency, broadening, units='eV', name='DrudeDielectricFunction')[source]#
Objects of this class represent the Drude dielectric function, which is defined for the entire positive frequency domain. The function is implemented as
\[\varepsilon_r(\omega) = 1 - \frac{\omega_p^2}{\omega^2 + i\gamma\omega}\]where the plasma frequency \(\omega_p\) determines where the real part of the dielectric function is zero and \(\gamma\) is a broadening.
- Parameters:
plasma_frequency (
float) – The plasma frequency \(\omega_p\); units are eV by default, optionally atomic units (seeunits).broadening (
float) – Broadening of the dielectric function \(\gamma\) units are eV by default, optionally atomic units (seeunits).name (
str) – Name of the dielectric function.units (
str) –eVto specifyplasma_frequencyandbroadeningin eV, orauto specify in atomic units.This parameter determines whether conversion should be performed during initialization and has no effect on instance methods and variables.
- property broadening: float#
Broadening parameter.
- property plasma_frequency: float#
plasma frequency.
- class strongcoca.response.NumericDielectricFunction(frequencies, dielectric_function, units='eV', name='NumericDielectricFunction')[source]#
Objects of this class hold different representations of dielectric functions sampled on a discrete grid of frequencies. The dielectric function can be evaluated at any frequency within the supplied grid by linear interpolation.
- Parameters:
frequencies (
Union[ndarray,Sequence[float],float]) –List of frequencies at which the dielectric function
dielectric_functionis sampled; units are eV by default, optionally atomic units (seeunits).The list of frequencies and the corresponding dielectric function are sorted if the former is not monotonically increasing.
dielectric_function (
Union[ndarray,Sequence[complex],complex]) – Complex dielectric function sampled at the frequenciesfrequencies.units (
str) –eVto specifyfreq_win eV, orauto specify in atomic units.This parameter determines whether conversion should be performed during initialization and has no effect on instance methods and variables.
name (
str) – Name of dielectric function.
- Raises:
ValueError – If
frequenciesordielectric_functionhave zero length.ValueError – If there is a mismatch in the length of
frequenciesanddielectric_function.ValueError – If
unitsis not one of the allowed valueseVorau.
Examples
The following snippet illustrates how a NumericDielectricFunction object can be constructed. Here, the analytic Drude dielectric function of the free-electron gas is used. To construct a dielectric function object from a numeric dielectric function the function
read_dielec_function_file()can be used.>>> import numpy as np >>> from strongcoca.response.dielectric import NumericDielectricFunction >>> >>> omega_p = 3 >>> gamma = 0.5 >>> freq_w = np.linspace(0.1, 10) >>> eps_w = 1 - omega_p**2 / (freq_w**2 + 1j*gamma*freq_w) >>> dielec = NumericDielectricFunction(freq_w, eps_w, name='Drude-function') >>> print(dielec) ---------------- NumericDielectricFunction ----------------- name : Drude-function Frequency grid : 50 points between 0.10 and 10.00eV frequencies (eV) : [ 0.1 0.30204 0.50408 ... 9.59592 9.79796 10. ] dielectric_function : [-33.61538+1.73077e+02j -25.37528+4.36618e+01j -16.85366+1.77091e+01j ... 0.90253+5.07897e-03j 0.90649+4.77173e-03j 0.91022+4.48878e-03j] >>> dense_freq_w = np.linspace(3, 4) >>> print(dielec(dense_freq_w)) [0.02425453+0.16310249j 0.0367996 +0.15996402j 0.04934466+0.15682555j 0.06188972+0.15368708j 0.07443478+0.15054861j 0.08697984+0.14741014j ...
- property dielectric_function: ndarray#
Complex unitless dielectric function at frequencies
frequencies; array with shape {n}.
- property frequencies: ndarray#
Frequencies corresponding to dielectric function
dielectric_functionin units of eV; array with shape {n}.
- property maxfreq: float#
Highest frequency
frequenciesfor which there is a value of the dielectric function; in units of eV.
- property minfreq: float#
Lowest frequency
frequenciesfor which there is a value of the dielectric function; in units of eV.
- property n: int#
number of samples for the dielectric function
dielectric_function.
Utilities#
- class strongcoca.response.utilities.ArtificialBroadening(width, units='eV')[source]#
Class representing artificial broadenings.
- Parameters:
width (
float) –Broadening width; meaning interpreted by the subclasses; must be strictly positive.
Units are eV by default, atomic units can be selected via
units.units (
str) –eVto specifywidthin eV orauto specifywidthin atomic units.This parameter determines whether conversion should be performed during initialization and has no effect on instance methods and variables.
- property width: float#
Broadening width in units of eV.
- class strongcoca.response.utilities.NoArtificialBroadening[source]#
Class representing no artificial broadening.
- class strongcoca.response.utilities.LorentzianBroadening(width, units='eV')[source]#
Class for representing artificial Lorentzian broadening of spectra.
- Parameters:
width (
float) –Broadening width (\(\eta\)); must be strictly positive.
Units are eV by default, atomic units can be selected via
units.units (
str) –eVto specifywidthin eV orauto specifywidthin atomic units.This parameter determines whether conversion should be performed during initialization and has no effect on instance methods and variables.
Examples
>>> from strongcoca.response.utilities import LorentzianBroadening >>> LorentzianBroadening(0.1) LorentzianBroadening(0.1, units='eV') >>> LorentzianBroadening(0.005, units='au') LorentzianBroadening(0.136..., units='eV')
- class strongcoca.response.utilities.GaussianBroadening(width, units='eV')[source]#
Class for representing artificial Gaussian broadening of spectra.
- Parameters:
width (
float) –Broadening width (\(\sigma\)); must be strictly positive.
Units are eV by default, atomic units can be selected via
units.units (
str) –eVto specifywidthin eV orauto specifywidthin atomic units.This parameter determines whether conversion should be performed during initialization and has no effect on instance methods and variables.
Examples
>>> from strongcoca.response.utilities import GaussianBroadening >>> GaussianBroadening(0.1) GaussianBroadening(0.1, units='eV') >>> GaussianBroadening(0.005, units='au') GaussianBroadening(0.136..., units='eV')