Calculators
Requiring full response
These calculators require the full response from a TDDFT calculation in the form of a BaseResponse
object.
- class rhodent.calculators.DensityCalculator(gpw_file, response, filter_occ=[], filter_unocc=[], *, times=None, pulses=None, frequencies=None, frequency_broadening=0)[source]
Calculate densities in the time or frequency domain.
The induced density (i.e. the density minus the ground state density) is to first order given by
plus PAW corrections, where
is the density of ground state Kohn-Sham pairIn the time domain, electrons and holes densities can be computed.
Refer to the documentation of
HotCarriersCalculator
for definitions of and .- Parameters:
gpw_file (
str
) – Filename of GPAW ground state file.response (
BaseResponse
) – Response object.filter_occ (
Sequence
[tuple
[float
,float
]]) – Filters for occupied states (holes). Provide a list of tuples (low, high) to compute the density of holes with energies within the interval low-high.filter_unocc (
Sequence
[tuple
[float
,float
]]) – Filters for unoccupied states (electrons). Provide a list of tuples (low, high) to compute the density of excited electrons with energies within the interval low-high.times (
Union
[list
[float
],ndarray
[tuple
[int
,...
],dtype
[float64
]],None
]) –Compute densities in the time domain, for these times (or as close to them as possible). In units of as.
May not be used together with
frequencies
orfrequency_broadening
.pulses (
Optional
[Collection
[Union
[Perturbation
,Laser
,dict
,None
]]]) –Compute densities in the time domain, in response to these pulses. If none, then no pulse convolution is performed.
May not be used together with
frequencies
orfrequency_broadening
.frequencies (
Union
[list
[float
],ndarray
[tuple
[int
,...
],dtype
[float64
]],None
]) –Compute densities in the frequency domain, for these frequencies. In units of eV.
May not be used together with
times
orpulses
.frequency_broadening (
float
) –Compute densities in the frequency domain, with Gaussian broadening of this width. In units of eV.
May not be used together with
times
orpulses
.
- property N_c: ndarray[tuple[int, ...], dtype[int64]]
Number of points in each Cartesian direction of the grid.
- calculate_and_write(out_fname, which='induced', write_extra={})[source]
Calculate density contributions.
Densities are saved in a numpy archive, ULM file or cube file depending on whether the file extension is
.npz
,.ulm
, or.cube
.If the file extension is
.cube
thenout_fname
is taken to be a format string that takes the following attributes.Accepts variables
{time}
- Time in units of as (time domain only).{freq}
- Frequency in units of eV (frequency domain only).{which}
- Thewhich
argument.{pulsefreq}
- Pulse frequency in units of eV (time domain only).{pulsefwhm}
- Pulse FWHM in units of fs (time domain only).
Examples
out_fname =
{which}_density_t{time:09.1f}.cube
.out_fname =
{which}_density_w{freq:05.2f}.cube
.
- Parameters:
out_fname (
str
) – File name of the resulting data file.which (
str
|list
[str
]) –String, or list of strings specifying the types of density to compute:
induced
- Induced density.holes
- Holes density (only allowed in the time domain).electrons
- Electrons density (only allowed in the time domain).
write_extra (
dict
[str
,Any
]) – Dictionary of extra key-value pairs to write to the data file.
- property gd: GridDescriptor
Real space grid.
- property gdshape: tuple[int, int, int]
Shape of the real space grid.
- get_density(rho_nn, nn_indices, fltn1=slice(None, None, None), fltn2=slice(None, None, None), u=0)[source]
Calculate a real space density from a density matrix in the Kohn-Sham basis.
- Parameters:
rho_nn (rhodent.typing.DistributedArray) – Density matrix
, , or .nn_indices (
str
) –Indices describing the density matrices
rho_nn
. One ofia
for induced density .ii
for holes density .aa
for electrons density .
flt_n1 – Filter selecting rows of the density matrix.
flt_n2 – Filter selecting columns of the density matrix.
u (
int
) – Kpoint index.
- Return type:
Distributed array with the density in real space on the root rank.
- get_result_keys(yield_total=True, yield_electrons=False, yield_holes=False)[source]
Get the keys that each result will contain, and dimensions thereof.
- Return type:
Object representing the data that will be present in the result objects.
- icalculate(yield_total=True, yield_electrons=False, yield_holes=False)[source]
Iteratively calculate results.
- Parameters:
yield_total (
bool
) – The results should include the total induced density.yield_holes (
bool
) – The results should include the holes densities, optionally decomposed byfilter_occ
.yield_electrons (
bool
) – The results should include the electrons densities, optionally decomposed byfilter_unocc
.
- Yields:
Tuple (work, result) on the root rank of the calculation communicator. Does not yield on non-root ranks of the calculation communicator. –
- work
An object representing the metadata (time, frequency or pulse) for the work done.
- result
Object containg the calculation results for this time, frequency or pulse.
- Return type:
Generator
[tuple
[WorkMetadata
,Result
],None
,None
]
- property occ_filters: list[slice]
List of energy filters for occupied states.
- property unocc_filters: list[slice]
List of energy filters for unoccupied states.
- class rhodent.calculators.DipoleCalculator(response, voronoi=None, *, energies_occ=None, energies_unocc=None, sigma=None, times=None, pulses=None, frequencies=None, frequency_broadening=0)[source]
Calculate the induced dipole moment in the time or frequency domain.
The induced dipole moment (i.e. the dipole moment minus the permanent component) is to first order given by
where
is the dipole matrix element of ground state Kohn-Sham pairand
the induced Kohn-Sham density matrix.In the frequency domain, this calculator calculates the polarizability, i.e. the Fourier transform of the dipole moment divided by the perturbation.
The absorption spectrum in units of dipole strength function is the imaginary part of the polarizability times a prefactor
This class can also compute projections of the above on Voronoi weights
.- Parameters:
response (
BaseResponse
) – Response object.voronoi (
Optional
[VoronoiWeights
]) – Voronoi weights object.energies_occ (
Union
[list
[float
],ndarray
[tuple
[int
,...
],dtype
[float64
]],None
]) – Energy grid in units of eV for occupied levels (holes).energies_unocc (
Union
[list
[float
],ndarray
[tuple
[int
,...
],dtype
[float64
]],None
]) – Energy grid in units of eV for unoccupied levels (electrons)sigma (
Optional
[float
]) – Gaussian broadening width for energy grid in units of eV.times (
Union
[list
[float
],ndarray
[tuple
[int
,...
],dtype
[float64
]],None
]) –Compute induced dipole in the time domain, for these times (or as close to them as possible). In units of as.
May not be used together with
frequencies
orfrequency_broadening
.pulses (
Optional
[Collection
[Union
[Perturbation
,Laser
,dict
,None
]]]) –Compute induced dipole in the time domain, in response to these pulses. If none, then no pulse convolution is performed.
May not be used together with
frequencies
orfrequency_broadening
.frequencies (
Union
[list
[float
],ndarray
[tuple
[int
,...
],dtype
[float64
]],None
]) –Compute polarizability in the frequency domain, for these frequencies. In units of eV.
May not be used together with
times
orpulses
.frequency_broadening (
float
) –Compute polarizability in the frequency domain, with Gaussian broadening of this width. In units of eV.
May not be used together with
times
orpulses
.
- calculate_and_write(out_fname, write_extra={}, save_matrix=False, only_one_pulse=True)[source]
Calculate induced dipole moments and transition contribution maps.
Dipole moments and contributions are saved in a numpy archive if the file extension is
.npz
or in an ULM file if the file extension is.ulm
.- Parameters:
out_fname (
str
) – File name of the resulting data file.write_extra (
dict
[str
,Any
]) – Dictionary of extra key-value pairs to write to the data file.save_matrix (
bool
) – Whether the transition energy distributions should be computed and saved.only_one_pulse (
bool
) – If False, group arrays by pulse. Only valid in time domain.
- get_result_keys(yield_total_ia=False, yield_proj_ia=False, yield_total_ou=False, yield_proj_ou=False, decompose_v=True, v=None)[source]
Get the keys that each result will contain, and dimensions thereof.
- Parameters:
yield_total_ia (
bool
) – The results should include the total dipole contributions in the electron-hole basis .yield_proj_ia (
bool
) – The results should include projections of the dipole contributions in the electron-hole basis .yield_total_ou (
bool
) – The results should include the total dipole contributions on the energy grid.yield_proj_ou (
bool
) – The results should include projections of the dipole contributions on the energy grid.decompose_v (
bool
) – The results should include the dipole moment and/or its contributions decomposed by Cartesian direction.v (
Optional
[int
]) – If not None, then the results should include the v:th Cartesian component of the dipole moment and its contributions.
- Return type:
- icalculate(yield_total_ia=False, yield_proj_ia=False, yield_total_ou=False, yield_proj_ou=False, decompose_v=True, v=None)[source]
Iteratively calculate dipole contributions.
- Parameters:
yield_total_ia (
bool
) – The results should include the total dipole contributions in the electron-hole basis .yield_proj_ia (
bool
) – The results should include projections of the dipole contributions in the electron-hole basis .yield_total_ou (
bool
) – The results should include the total dipole contributions on the energy grid.yield_proj_ou (
bool
) – The results should include projections of the dipole contributions on the energy grid.decompose_v (
bool
) – The results should include the dipole moment and/or its contributions decomposed by Cartesian direction.v (
Optional
[int
]) – If notNone
, then the results should include the v:th Cartesian component of the dipole moment and its contributions.
- Yields:
Tuple (work, result) on the root rank of the calculation communicator. Does not yield on non-root ranks of the calculation communicator. –
- work
An object representing the metadata (time, frequency or pulse) for the work done.
- result
Object containg the calculation results for this time, frequency or pulse.
- Return type:
Generator
[tuple
[WorkMetadata
,Result
],None
,None
]
- class rhodent.calculators.EnergyCalculator(response, voronoi=None, *, energies_occ=None, energies_unocc=None, sigma=None, times=None, pulses=None, frequencies=None, frequency_broadening=0)[source]
Calculate energy contributions in the time domain.
The total energy can be written
The contributions to the total energy are
the contributions to the Hartree energy are
and the rate of energy change is
The matrix element
can be computed from the dipole matrix elementprojected on the direction of the perturbation
, the occupation number difference and the pulse amplitude- Parameters:
response (
BaseResponse
) – Response object.voronoi (
Optional
[VoronoiWeights
]) – Voronoi weights object.energies_occ (
Union
[list
[float
],ndarray
[tuple
[int
,...
],dtype
[float64
]],None
]) – Energy grid in units of eV for occupied levels (holes).energies_unocc (
Union
[list
[float
],ndarray
[tuple
[int
,...
],dtype
[float64
]],None
]) – Energy grid in units of eV for unoccupied levels (electrons).sigma (
Optional
[float
]) – Gaussian broadening width for energy grid in units of eV.times (
Union
[list
[float
],ndarray
[tuple
[int
,...
],dtype
[float64
]],None
]) – Compute energies in the time domain, for these times (or as close to them as possible). In units of as.pulses (
Optional
[Collection
[Union
[Perturbation
,Laser
,dict
,None
]]]) – Compute energies in the time domain, in response to these pulses. If none, then no pulse convolution is performed.
- calculate_and_write(out_fname, write_extra={}, save_matrix=False, save_dist=False, only_one_pulse=True)[source]
Calculate energy contributions.
Energies are saved in a numpy archive if the file extension is
.npz
or in an ULM file if the file extension is.ulm
.- Parameters:
out_fname (
str
) – File name of the resulting data file.write_extra (
dict
[str
,Any
]) – Dictionary of extra key-value pairs to write to the data file.save_matrix (
bool
) – Whether the electron-hole map matrix should be computed and saved.save_dist (
bool
) – Whether the transition energy distributions should be computed and saved.only_one_pulse (
bool
) – If False, group arrays by pulse.
- get_result_keys(yield_total_E_ia=False, yield_proj_E_ia=False, yield_total_E_ou=False, yield_total_dists=False, direction=2)[source]
Get the keys that each result will contain, and dimensions thereof.
- Parameters:
yield_total_E_ia (
bool
) – The results should include the contributions in the electron-hole basis to the total energy and Coulomb energyyield_proj_E_ia (
bool
) – The results should include the contributions in the electron-hole basis to the total energy projected on the occupied and unoccupied Voronoi weights and .yield_total_E_ou (
bool
) – The results should include the contributions the total energy broadened on the occupied and unoccupied energy grids andyield_total_dists (
bool
) – The results should include the contributions the total energy and Coulomb energy broadened by electronic transition energy onto the unoccupied energies grid anddirection (
int
|Sequence
[int
]) – Direction of the polarization of the pulse. Integer 0, 1 or 2 to specify x, y or z, or the direction vector specified as a list of three values. Default: polarization along z.
- Return type:
- icalculate(yield_total_E_ia=False, yield_proj_E_ia=False, yield_total_E_ou=False, yield_total_dists=False, direction=2)[source]
Iteratively calculate energies.
- Parameters:
yield_total_E_ia (
bool
) – The results should include the contributions in the electron-hole basis to the total energy and Coulomb energyyield_proj_E_ia (
bool
) – The results should include the contributions in the electron-hole basis to the total energy projected on the occupied and unoccupied Voronoi weights and .yield_total_E_ou (
bool
) – The results should include the contributions the total energy broadened on the occupied and unoccupied energy grids andyield_total_dists (
bool
) – The results should include the contributions the total energy and Coulomb energy broadened by electronic transition energy onto the unoccupied energies grid anddirection (
int
|Sequence
[int
]) – Direction of the polarization of the pulse. Integer 0, 1 or 2 to specify x, y or z, or the direction vector specified as a list of three values. Default: polarization along z.
- Yields:
Tuple (work, result) on the root rank of the calculation communicator. Does not yield on non-root ranks of the calculation communicator. –
- work
An object representing the metadata (time or pulse) for the work done.
- result
Object containg the calculation results for this time and pulse.
- Return type:
Generator
[tuple
[WorkMetadata
,Result
],None
,None
]
- class rhodent.calculators.HotCarriersCalculator(response, voronoi=None, *, energies_occ=None, energies_unocc=None, sigma=None, times=None, pulses=None, frequencies=None, frequency_broadening=0)[source]
Calculate hot-carrier distributions in the time domain.
For weak perturbations, the response of the density matrix is to first order non-zero only in the occupied-unoccupied space, i.e. the block off-diagonals
The unoccupied-occupied, or electron-hole, part of the density matrix is thus linear in perturbation and can by transformed using Fourier transforms.
From the first-order response, the second order response, i.e. the hole-hole (
) and electron-electron ( ) parts can be obtained.The hole-hole part is
and the electron-hole part
where
where
is the occupation number difference of pair .Hot-carrier distributions are calculated by convolution of
(holes) and (electrons) by Gaussian broadening functions on the energy grid.where
are Kohn-Sham energies andProjected hot-carrier distributions are defined
The Voronoi weights are
where the operator
is 1 in the Voronoi region of the atomic projections and 0 outside.- Parameters:
response (
BaseResponse
) – Response object.voronoi (
Optional
[VoronoiWeights
]) – Voronoi weights object.energies_occ (
Union
[list
[float
],ndarray
[tuple
[int
,...
],dtype
[float64
]],None
]) – Energy grid in units of eV for occupied levels (holes).energies_unocc (
Union
[list
[float
],ndarray
[tuple
[int
,...
],dtype
[float64
]],None
]) – Energy grid in units of eV for unoccupied levels (electrons)sigma (
Optional
[float
]) – Gaussian broadening width for energy grid in units of eV.times (
Union
[list
[float
],ndarray
[tuple
[int
,...
],dtype
[float64
]],None
]) – Compute hot carriers in the time domain, for these times (or as close to them as possible). In units of as.pulses (
Optional
[Collection
[Union
[Perturbation
,Laser
,dict
,None
]]]) – Compute hot carriers in the time domain, in response to these pulses. IfNone
no pulse convolution is performed.
- calculate_distributions_and_write(out_fname, write_extra={}, average_times=True, only_one_pulse=True)[source]
Calculate broadened hot-carrier energy distributions, optionally averaged over time.
HC distributions are saved in a numpy archive if the file extension is
.npz
or in a comma separated file if the file extension is.dat
.- Parameters:
out_fname (
str
) – File name of the resulting data file.write_extra (
dict
[str
,Any
]) – Dictionary of extra key-value pairs to write to the data file.average_times (
bool
) – IfTrue
, an average over the given times will be taken. If false, then hot-carrier distributions are computed separately over the times, and each output is written separately for each time.only_one_pulse (
bool
) – There is only one pulse, don’t group by pulse.
- calculate_totals_by_pulse_and_write(out_fname, write_extra={})[source]
Calculate the number of generated hot carriers, projected on groups of atoms, for a list of pulses.
HC numbers are saved in a numpy archive if the file extension is
.npz
or in a comma separated file if the file extension is.dat
.- Parameters:
out_fname (
str
) – File name of the resulting data file.
- calculate_totals_by_time_and_write(out_fname, write_extra={})[source]
Calculate the number of generated hot carriers, projected on groups of atoms, for a list of times.
HC numbers are saved in a numpy archive if the file extension is
.npz
or in a comma separated file if the file extension is.dat
.- Parameters:
out_fname (
str
) – File name of the resulting data file.
- get_result_keys(yield_total_hcdists=False, yield_proj_hcdists=False, yield_total_P=False, yield_proj_P=False, yield_total_P_ou=False)[source]
Get the keys that each result will contain, and dimensions thereof.
- Parameters:
yield_total_hcdists (
bool
) – The results should include the total hot-carrier distributions on the energy grid.yield_proj_hcdists (
bool
) – The results should include the projections of the hot-carrier distributions on the energy grid.yield_total_P (
bool
) – The results should include the total hot-carrier distributions in the electron-hole basis.yield_proj_P (
bool
) – The results should include the projections of the hot-carrier distributions in the electron-hole basis.yield_total_P_ou (
bool
) – The results should include the transition matrix broadened on the energy grid.
- icalculate(yield_total_hcdists=False, yield_proj_hcdists=False, yield_total_P=False, yield_proj_P=False, yield_total_P_ou=False)[source]
Iteratively calculate second order density matrices and hot-carrier distributions.
- Parameters:
yield_total_hcdists (
bool
) – The results should include the total hot-carrier distributions on the energy grid.yield_proj_hcdists (
bool
) – The results should include the projections of the hot-carrier distributions on the energy grid.yield_total_P (
bool
) – The results should include the total hot-carrier distributions in the electron-hole basis.yield_proj_P (
bool
) – The results should include the projections of the hot-carrier distributions in the electron-hole basis.yield_total_P_ou (
bool
) – The results should include the transition matrix broadened on the energy grid.
- Yields:
Tuple (work, result) on the root rank of the calculation communicator. Does not yield on non-root ranks of the calculation communicator. –
- work
An object representing the metadata (time and pulse) for the work done.
- result
Object containg the calculation results for this time and pulse.
- Return type:
Generator
[tuple
[WorkMetadata
,Result
],None
,None
]
Other calculators
These calculators require other input files than the calculators for the full response, for example just the dipole moment file.
- class rhodent.dos.DOSCalculator(eigenvalues, fermilevel, voronoi, energies, sigma, zerofermi=False)[source]
Density of states (DOS) and partial DOS (PDOS) calculator.
Calculates DOS
and PDOS
where
are Kohn-Sham energies and a Gaussian broadening functionThe Voronoi weights are defined
where the operator
is 1 in the Voronoi region of the atomic projections and 0 outside.- Parameters:
eigenvalues (
list
[float
] |ndarray
[tuple
[int
,...
],dtype
[float64
]]) – List of eigenvalues (relative to Fermi level).fermilevel (
float
) – Fermi level.voronoi (
VoronoiWeights
|None
) –Voronoi weights object defining the atomic projections for PDOS.
Leave out if only DOS is desired.
energies (
list
[float
] |ndarray
[tuple
[int
,...
],dtype
[float64
]]) – Array of energies (in units of eV) for which the broadened PDOS is computed.sigma (
float
) – Gaussian broadening width in units of eV.zerofermi (
bool
) – Eigenvalues relative to Fermi level ifTrue
, else relative to vacuum
- calculate_dos_and_write(out_fname, write_extra={})[source]
Calculate the DOS and write to file.
The DOS is saved in a numpy archive if the file extension is
.npz
or in a comma separated file if the file extension is.dat
.- Parameters:
out_fname (
Path
|str
) – File name of the resulting data file.write_extra (
dict
[str
,Any
]) – Dictionary of additional data written to numpy archive (ignored for.dat
) files.
- calculate_pdos_and_write(out_fname, write_extra={}, write_extra_from_voronoi=False)[source]
Calculate the PDOS and write to file.
The PDOS is saved in a numpy archive if the file extension is
.npz
or in a comma separated file if the file extension is.dat
.- Parameters:
out_fname (
Path
|str
) – File name of the resulting data file.write_extra (
dict
[str
,Any
]) – Dictionary of additional data written to numpy archive (ignored for.dat
) files.write_extra_from_voronoi (
bool
) – If true, and voronoi is a ULM reader, extra key-value pairs are read from voronoi and written to the.npz
file (ignored for.dat
) files.
- classmethod from_calc(calc, voronoi, energies, sigma, zerofermi=False)[source]
Initialize from GPAW calculator.
- Parameters:
calc (gpaw.GPAW) – GPAW calculator.
voronoi (
VoronoiWeights
|None
) –Voronoi weights object defining the atomic projections for PDOS.
Leave out if only DOS is desired.
energies (
list
[float
] |ndarray
[tuple
[int
,...
],dtype
[float64
]]) – Array of energies (in units of eV) for which the broadened PDOS is computed.sigma (
float
) – Gaussian broadening width in units of eV.zerofermi (
bool
) – Eigenvalues relative to Fermi level ifTrue
, else relative to vacuum.
- classmethod from_gpw(gpw_file, voronoi, energies, sigma, zerofermi=False)[source]
Initialize from
.gpw
file.- Parameters:
gpw_file (
Path
|str
) – Filename of GPAW ground state file.voronoi (
VoronoiWeights
|None
) –Voronoi weights object defining the atomic projections for PDOS.
Leave out if only DOS is desired.
energies (
list
[float
] |ndarray
[tuple
[int
,...
],dtype
[float64
]]) – Array of energies (in units of eV) for which the broadened PDOS is computed.sigma (
float
) – Gaussian broadening width in units of eV.zerofermi (
bool
) – Eigenvalues relative to Fermi level ifTrue
, else relative to vacuum.
- icalculate_pdos()[source]
Calculate PDOS for each of the atomic projections in the
voronoi
object..
- property sigma: float
Gaussian broadening width in units of eV.
- property voronoi: VoronoiWeights | None
Voronoi weights object.
- class rhodent.spectrum.SpectrumCalculator(times, dipole_moments, perturbation)[source]
Spectrum calculator.
Calculates the dynamic polarizability and dipole strength function (the spectrum).
The polarizability
is related to the perturbation and the Fouier transform of the dipole moment asThe dipole strength function is
Both quantities can be broadened by supplying a non-zero value
. Then the convolutionis computed. When the perturbation is a delta-kick, this can efficiently be computed by multiplying the dipole moment by a Gaussian envelope before Fourier transformation
For other perturbations, a FFT and IFFT is first performed to obtain the delta-kick response.
- Parameters:
times (
ndarray
[tuple
[int
,...
],dtype
[float64
]]) – Array (length T) of times in units of as.dipole_moments (
ndarray
[tuple
[int
,...
],dtype
[float64
]]) – Array (shape (T, 3)) of dipole moments in units of eÅ.perturbation (
Union
[Perturbation
,Laser
,dict
,None
]) – The perturbation that was applied in the TDDFT calculation.
- calculate_spectrum_and_write(out_fname, frequencies, frequency_broadening=0, write_extra={})[source]
Calculate the dipole strength function (spectrum) and write to file.
The spectrum is saved in a numpy archive if the file extension is
.npz
or in a comma separated file if the file extension is.dat
.- Parameters:
out_fname (
Path
|str
) – File name of the resulting data file.frequencies (
list
[float
] |ndarray
[tuple
[int
,...
],dtype
[float64
]]) – Array of frequencies for which to calculate the spectrum; in units of eV.frequency_broadening (
float
) – Gaussian broadening width in units of eV. Default (0) is no broadening.write_extra (
dict
[str
,Any
]) – Dictionary of additional data written to numpy archive (ignored for.dat
) files.
- property dipole_moments: ndarray[tuple[int, ...], dtype[float64]]
Array of dipole moments, in units of eÅ.
- get_dipole_strength_function(frequencies, frequency_broadening=0)[source]
Calculate the dipole strength function (spectrum) in units of 1/eV.
- Parameters:
- Return type:
Array of dipole strength function in units of 1/eV. The array has shape (N, 3), where N is length of frequencies.
- get_dynamic_polarizability(frequencies, frequency_broadening=0)[source]
Calculate the dynamic polarizability.
- Parameters:
- Return type:
Array of dynamic polarizabilities in (eÅ)**2/eV. The array has shape (N, 3), where N is the length of
frequencies
.
- property perturbation: Perturbation
The perturbation that was applied in the TDDFT calculation.