skultrafast.twoD_dataset

Module Contents

Classes

SingleCLSResult

Result of a single CLS fit.

FFCFResult

Baseclass for FFCF determination methods.

CLSResult

Class holding the data of CLS-analysis. Has methods to analyze and plot them.

DiagResult

ExpFit2DResult

GaussResult

A dataclass to hold the results of the Gaussian fit.

TwoDim

Dataset for an two dimensional dataset. Requires the t- (waiting times),the

Attributes

skultrafast.twoD_dataset.PathLike[source]
class skultrafast.twoD_dataset.SingleCLSResult[source]

Result of a single CLS fit.

pump_wn: numpy.ndarray[source]

Pump wavenumbers.

max_pos: numpy.ndarray[source]

Maximum positions for each pump-slice.

max_pos_err: numpy.ndarray[source]

Standard error of the maximum positions. Nan if not available

slope: float[source]

Slope of the linear fit.

reg_result: Any[source]

Result of the linear regression using statsmodels.

recentered_pump_wn: numpy.ndarray[source]

Recentered pump wavenumbers by mean subtraction. Used in the linear regression.

linear_fit: numpy.ndarray[source]

Linear fit of the CLS data.

class skultrafast.twoD_dataset.FFCFResult[source]

Baseclass for FFCF determination methods.

For backwards compatibility, the values are always called slopes

wt: numpy.ndarray[source]

Wait times.

slopes: numpy.ndarray[source]

Values

slope_errors: numpy.ndarray | None[source]
exp_fit_result_: lmfit.model.ModelResult | None[source]
exp_fit(start_taus: List[float], use_const=True, use_weights=True)[source]
plot_cls(ax=None, model_style: Dict = {}, symlog=False, **kwargs)[source]
class skultrafast.twoD_dataset.CLSResult[source]

Bases: FFCFResult

Class holding the data of CLS-analysis. Has methods to analyze and plot them.

intercepts: numpy.ndarray[source]

Contains the intercepts, useful for plotting mostly

intercept_errors: numpy.ndarray | None[source]

Errors of the intercepts

lines: List[numpy.ndarray][source]

Contains the x and y, yerr-values used for the linear fit

exp_fit_result_: lmfit.model.ModelResult | None[source]
class skultrafast.twoD_dataset.DiagResult[source]
diag: numpy.ndarray[source]

Contains the diagonal signals of the 2d-spectra

antidiag: numpy.ndarray[source]

Contains the antidiagonal signals of the 2d-spectra

diag_coords: numpy.ndarray[source]

Contains the coordinates of the diagonal signals

antidiag_coords: numpy.ndarray[source]

Contains the coordinates of the antidiagonal signals

offset: float[source]

The offset of the diagonal signals

p: float[source]

The position crossing the diagonals

class skultrafast.twoD_dataset.ExpFit2DResult[source]
minimizer: lmfit.minimizer.MinimizerResult[source]

Lmfit minimizer result

model: TwoDim[source]

The fit data

residuals: numpy.ndarray[source]

The residuals of the fit

das: numpy.ndarray[source]

The decay amplitudes aka the DAS of the 2D-spectra

basis: numpy.ndarray[source]

The basis functions used for the fit

taus: numpy.ndarray[source]

The decay times of the fit

class skultrafast.twoD_dataset.GaussResult[source]

Bases: FFCFResult

A dataclass to hold the results of the Gaussian fit.

results: List[lmfit.model.ModelResult][source]

List of lmfit results from the gaussian fits

fit_out: TwoDim[source]

TwoDim object with the fit data, useful for plotting

class skultrafast.twoD_dataset.TwoDim[source]

Dataset for an two dimensional dataset. Requires the t- (waiting times),the probe- and pump-axes in addition to the three dimensional spec2d data.

t: numpy.ndarray[source]

Array of the waiting times

pump_wn: numpy.ndarray[source]

Array with the pump-wavenumbers

probe_wn: numpy.ndarray[source]

Array with the probe-wavenumbers

spec2d: numpy.ndarray[source]

Array with the data, shape must be (t.size, wn_probe.size, wn_pump.size)

info: Dict[source]

Meta Info

single_cls_result_: SingleCLSResult | None[source]

Contains the data from a Single CLS analysis

cls_result_: CLSResult | None[source]

Contains the data from a CLS analysis

plot: skultrafast.twoD_plotter.TwoDimPlotter[source]

Plot object offering plotting methods

interpolator_: scipy.interpolate.RegularGridInterpolator | None[source]

Contains the interpolator for the 2d-spectra

exp_fit_result_: ExpFit2DResult | None[source]

Contains the result of the exponential fit

_make_int()[source]
__attrs_post_init__()[source]
copy() TwoDim[source]

Makes a copy of the dataset.

t_idx(t: float) int[source]
t_idx(t: Iterable[float]) List[int]

Return nearest idx to nearest time value

t_d(t) numpy.ndarray[source]

Return the dataset nearest to the given time t

data_at(t: float | None = None, probe_wn: float | None = None, pump_wn: float | None = None) numpy.ndarray[source]

Extracts the data at given coordinates.

probe_idx(wn: float | Iterable[float]) int | List[int][source]

Return nearest idx to nearest probe_wn value

pump_idx(wn: float | Iterable[float]) int | List[int][source]

Return nearest idx to nearest pump_wn value

select_range(pump_range: Tuple[float, float], probe_range: Tuple[float, float], invert: bool = False) TwoDim[source]

Return a dataset containing only the selected region.

select_t_range(t_min: float = -np.inf, t_max: float = np.inf) TwoDim[source]

” Returns a dataset only containing the data within given time limits.

integrate_pump(lower: float = -np.inf, upper: float = np.inf) skultrafast.dataset.TimeResSpec[source]

Calculate and return 1D Time-resolved spectra for given range.

Uses trapezoidal integration.

Parameters:
  • lower (float) – Lower pump wl

  • upper (float) – upper pump wl

Returns:

The corresponding 1D Dataset

Return type:

TimeResSpec

single_cls(t: float, pr_range: float | Tuple[float, float] = 9.0, pu_range: float | Tuple[float, float] = 7.0, mode: Literal[neg, pos] = 'neg', method: Literal[com, quad, fit, log_quad, skew_fit] = 'com') SingleCLSResult[source]

Calculate the CLS for single 2D spectrum.

Parameters:
  • t (float) – Delay time of the spectrum to analyse

  • pr_range (float or float, optional) – How many wavenumbers away from the maximum to use for determining the exact position, by default 9, resulting a total range of 18 wavenumbers. Also accepts a tuple, which is interpreted as (lower, upper) range.

  • pu_range (float, optional) – The range around the pump-maxima used for calculating the CLS. If given a float, the range is calculated as (max - pu_range, max + pu_range). If given a tuple, it is interpreted as (lower, upper) range.

  • mode (('neg', 'pos'), optional) – negative or positive maximum, by default ‘neg’. Ignored when method is ‘nodal’.

  • method (('com', 'quad', 'fit', 'log_quad', 'skew_fit', 'nodal'), optional) – Selects the method used for determination of the maximum signal. com uses the center-of-mass, quad uses a quadratic fit and fit uses a gaussian fit. ‘log_quad’ uses a quadratic fit on the logarithm of the signal. ‘skew_fit’ uses a gaussian fit with a linear background. ‘nodal’ finds the zero-crossing of the signal.

Returns:

pump_wn max_pos max_pos_err slope reg_result recentered_pump_wn linear_fit

Return type:

Returns SingleCLSResult object with attributes

cls(joblib_kws=None, **cls_args) CLSResult[source]

Calculates the CLS for all 2d-spectra. The arguments are given to the single cls function, except joblib_kw, which is forwarded to joblib.Parallel.

Returns as CLSResult.

diag_and_antidiag(t: float, offset: float | None = None, p: float | None = None) DiagResult[source]

Extracts the diagonal and anti-diagonal.

Parameters:
  • t (float) – Waiting time of the 2d-spectra from which the data is extracted.

  • offset (float) – Offset of the diagonal, if none, it will we determined by the going through the signal minimum.

  • p (float) – The point where the anti-diagonal crosses the diagonal. If none, it also goes through the signal minimum.

Returns:

Contains the diagonals, coordinates and points.

Return type:

CLSResult

pump_slice_amp(t: float, bg_correct: bool = True) numpy.ndarray[source]

” Calculates the pump-slice-amplitude for a given delay t.

apply_filter(kind: Literal[uniform, gaussian], size, *args) TwoDim[source]

” Returns filtered dataset.

Parameters:
  • kind (str) – Which filter to use. Supported are uniform and gaussian.

  • size (tuple[float, float, float]) – Kernel of the filter

Returns:

Filtered dataset.

Return type:

TwoDim

save_txt(pname: PathLike, **kwargs)[source]

Saves 2d-spectra as a text files a directory.

Parameters:
  • pname (PathLike) – Path to the file.

  • kwargs – Additional arguments for the np.savetxt function.

save_single_txt(fname: PathLike, i: int, **kwargs)[source]

Save a single 2D spectra as a text file

Parameters:
  • fname (str) – The file name.

  • i (int) – The index of the 2D spectra.

  • kwargs – Additional arguments for the np.savetxt function.

background_correction(excluded_range: Tuple[float, float], deg: int = 3) None[source]

Fits and subtracts a background for each pump-frequency. Done for each waiting time. Does the subtraction inplace, e.g. modifies the dataset.

Parameters:
  • excluded_range (Tuple[float, float]) – The range of the pump axis which is excluded from the fit, e.g. contains the signal.

  • deg (int) – Degree of the polynomial fit.

Return type:

None

get_minmax(t: float, com: int = 3) Dict[str, float][source]

Returns the position of the minimum and maximum of the dataset at time t. If com > 0, it return the center of mass around the minimum and maximum. The com argument gives the number of points to be used for the center of mass.

integrate_reg(pump_range: Tuple[float, float], probe_range: Tuple[float, float] = None) numpy.ndarray[source]

Integrates the 2D spectra over a given range, using the trapezoidal rule.

Parameters:
  • pump_range (tuple[float, float]) – The range of the pump axis to be integrated over.

  • probe_range (tuple[float, float]) – The range of the probe axis to be integrated over. If None, the probe range is used.

Returns:

The integrated spectral signal for all waiting times.

Return type:

np.ndarray

fit_taus(taus: numpy.ndarray)[source]

Given a set of decay times, fit the data to a sum of the exponentials. Used by the fit_taus method.

fit_das(taus, fix_last_decay=False) ExpFit2DResult[source]

Fit the data to a sum of exponentials (DAS), starting from the given decay constants. The results are stored in the fit_exp_result attribute.

fit_gauss(mode='both') GaussResult[source]

Fits the 2D spectra using two gaussians peaks.