tidy3d.plugins.resonance.ResonanceFinder#
- class ResonanceFinder[source]#
Bases:
Tidy3dBaseModel
Tool that extracts resonance information from a time series of the form shown below. The resonance information consists of frequency \(f\), decay rate \(\alpha\), Q factor \(Q = \pi |f|/\alpha\), amplitude \(a\), and phase \(\phi\).
- Parameters:
attrs (dict = {}) โ Dictionary storing arbitrary metadata for a Tidy3D object. This dictionary can be freely used by the user for storing data without affecting the operation of Tidy3D as it is not used internally. Note that, unlike regular Tidy3D fields,
attrs
are mutable. For example, the following is allowed for setting anattr
obj.attrs['foo'] = bar
. Also note that Tidy3D` will raise aTypeError
ifattrs
contain objects that can not be serialized. One can check ifattrs
are serializable by callingobj.json()
.freq_window (Tuple[float, float]) โ [units = Hz]. Window
[fmin, fmax]
for the initial frequencies. The resonance finder is initialized with an even grid of frequencies between fmin and fmax. The resonance finder then iteratively optimizes and prunes these frequencies. Note that frequencies outside this window may be returned. A narrow frequency window that only contains a few resonances may give enhanced accuracy compared to a broad frequency window with many resonances.init_num_freqs (PositiveInt = 200) โ Number of frequencies with which the resonance finder is initialized. The resonance finder then iteratively optimizes and prunes these frequencies. The number of frequencies returned will be less than
init_num_freqs
. Making this larger can find more resonances, while making it smaller can speed up the resonance finder.rcond (NonNegativeFloat = 0.0001) โ Cutoff for eigenvalues, relative to the largest eigenvalue. The resonance finder solves a generalized eigenvalue problem of the form \(Ax = \lambda B x\). If B has small eigenvalues, this is poorly conditioned, so we eliminate eigenvalues of B less than
rcond
times the largest eigenvalue. Making this closer to zero will typically return more resonances.
Note
\[f(t) = \sum_k a_k e^{i \phi_k} e^{-2 \pi i f_k t - \alpha_k t}\]Note
We use the algorithm described in
Vladimir A. Mandelshtam and Howard S. Taylor, โHarmonic inversion of time signals and its applications,โ J. Chem. Phys. 107, 6756 (1997).
Example
>>> import numpy as np >>> from tidy3d.plugins.resonance import ResonanceFinder >>> t = np.linspace(0, 10000, 10000) >>> f1 = 2*np.pi*0.1 - 1j*0.002 >>> f2 = 2*np.pi*0.2 - 1j*0.0005 >>> sig = 2*np.exp(-1j*f1*t) + 3*1j*np.exp(-1j*f2*t) >>> resfinder = ResonanceFinder(freq_window=(0.05, 0.25)) >>> resdata = resfinder.run_raw_signal(signal=sig, time_step=1) >>> data = resdata.to_dataframe() ... # A given dataframe
Attributes
Methods
run
(signals)Finds resonances in a
FieldTimeData
or a Tuple of such.run_raw_signal
(signal,ย time_step)Finds resonances in a time series.
run_scalar_field_time
(signal)Finds resonances in a
ScalarFieldTimeDataArray
.Inherited Common Usage
- freq_window#
- init_num_freqs#
- rcond#
- run(signals)[source]#
Finds resonances in a
FieldTimeData
or a Tuple of such. The time coordinates must be uniformly spaced, and the spacing must be the same across all supplied data. The resonance finder runs on the sum of theEx
,Ey
, andEz
fields of all the supplied data, unless no electric fields are present at all, in which case it runs on the sum of theHx
,Hy
, andHz
fields. Note that the signal should start after the sources have turned off.- Parameters:
signals (
FieldTimeData
) โ Data to search for resonances- Returns:
Dataset containing the decay rate, Q, amplitude, phase, and estimation error of the resonances as a function of frequency. Modes with low Q, small amplitude, or high estimation error are likely to be spurious.
- Return type:
xr.Dataset
- run_scalar_field_time(signal)[source]#
Finds resonances in a
ScalarFieldTimeDataArray
. The time coordinates must be uniformly spaced to use the resonance finder; the time step is calculated automatically. Note that the signal should start after the sources have turned off.- Parameters:
signal (
ScalarFieldTimeDataArray
) โ Time series to search for resonances- Returns:
Dataset containing the decay rate, Q, amplitude, phase, and estimation error of the resonances as a function of frequency. Modes with low Q, small amplitude, or high estimation error are likely to be spurious.
- Return type:
xr.Dataset
- run_raw_signal(signal, time_step)[source]#
Finds resonances in a time series. Note that the signal should start after the sources have turned off.
- Parameters:
signal (List[complex]) โ One-dimensional array holding the complex-valued time series data to search for resonances.
time_step (float) โ Time step / sampling rate of the data (in seconds).
- Returns:
Dataset containing the decay rate, Q, amplitude, phase, and estimation error of the resonances as a function of frequency. Modes with low Q, small amplitude, or high estimation error are likely to be spurious.
- Return type:
xr.Dataset
- __hash__()#
Hash method.