Reference on kuibit.timeseries

The timeseries module provides a representation of time series and convenience functions to create TimeSeries.

TimeSeries can be evenly or unevenly sampled are rich in features. They support all the mathematical operations and operators you may expect, and have additional methods, which include ones for taking derivatives, integrals, apply windows, smooth the signal, take Fourier transform, and more. Most of these methods are available in two flavors: those that return a new TimeSeries, and those which modify the object in place. The latter have names with imperative verbs.

TimeSeries are derived from the BaseSeries, which in turn is derived from the abstract class BaseNumerical. Some of the capabilities of TimeSeries (e.g., overloading the mathematical operators) are implemented in the parent classes.

The additional functions provided in timeseries are:

  • remove_duplicated_iters() cleans the input arrays by removing duplicated times.

  • unfold_phase() takes as argument a NumPy array representing a phase and unfolds it removing all the jumps of 2 pi. This is useful in gravitational wave analysis.

  • combine_ts() takes a list of timeseries and removes all the overlapping segments.

class kuibit.timeseries.TimeSeries(t, y, guarantee_t_is_monotonic=False)[source]

This class represents real or complex valued time series.

TimeSeries are defined providing a time list or array and the corresponding values. For example,

times = np.linspace(0, 2 * np.pi, 100)
values = np.sin(times)

ts = TimeSeries(times, values)

Times cannot be empty or not monotonically increasing. Times and values must have the same length.

TimeSeries are well-behaved classed, many operations and methods are implemented. For instance, you can sum/multiply two TimeSeries.

NumPy acts on TimeSeries cleanly, eg. np.log10(TimeSeries) is a TimeSeries with log10(data).

TimeSeries have methods for smoothing, windowing, extracting phase and more.

Variables:
  • t (1D NumPy array or float) – Times.

  • y (1D NumPy array or float) – Values.

  • spline_real (tuple) – Coefficients for a spline represent of the real part of y.

  • spline_imag (tuple) – Coefficients for a spline represent of the real part of y.

Constructor.

When guarantee_t_is_monotonic is True no checks will be perform to make sure that t is monotonically increasing (increasing performance). This should is used internally whenever a new series is returned from self (since we have already checked that t is good.) or in performance critical routines.

Parameters:
  • t (1D NumPy array or list) – Sampling times, need to be strictly increasing.

  • y (1D NumPy array or list) – Data samples, can be real or complex valued.

  • guarantee_t_is_monotonic (bool) – The code will assume that t is monotonically increasing.

abs_max()

Return the maximum of the absolute value

abs_min()

Return the minimum of the absolute value

abs_nanmax()

Return the maximum of the absolute value ignoring NaNs

abs_nanmin()

Return the minimum of the absolute value ignoring NaNs

align_at_maximum()[source]

Time shift the series so that the absolute maximum is at t=0.

align_at_minimum()[source]

Time shift the series so that the absolute minimum is at t=0.

aligned_at_maximum()[source]

Return a new timeseries with absolute maximum at t=0.

Returns:

Timeseries shifted so that the maximum is a t=0.

Return type:

TimeSeries

aligned_at_minimum()[source]

Return a new timeseries with absolute minimum at t=0.

Returns:

Timeseries shifted so that the minimum is a t=0.

Return type:

TimeSeries

blackman_window()[source]

Apply Blackman window.

blackman_windowed()[source]

Return a timeseries with Blackman window applied.

clip(init=None, end=None)

Remove data outside the the interval [init, end]. If init or end are not specified or None, it does not remove anything from this side.

Parameters:
  • init (float or None) – Data with x <= init will be removed.

  • end (float or None) – Data with x >= end will be removed.

clipped(init=None, end=None)

Return a series with data removed outside the interval [init, end]. If init or end are not specified or None, it does not remove anything from this side.

Parameters:
  • init (float or None) – Data with x <= init will be removed.

  • end (float or None) – Data with x >= end will be removed.

Returns:

Series with enforced minimum and maximum

Return type:

BaseSeries or derived class

copy()

Return a deep copy.

Returns:

Deep copy of the series.

Return type:

BaseSeries or derived class

crop(init=None, end=None)

Remove data outside the the interval [init, end]. If init or end are not specified or None, it does not remove anything from this side.

Parameters:
  • init (float or None) – Data with x <= init will be removed.

  • end (float or None) – Data with x >= end will be removed.

cropped(init=None, end=None)

Return a series with data removed outside the interval [init, end]. If init or end are not specified or None, it does not remove anything from this side.

Parameters:
  • init (float or None) – Data with x <= init will be removed.

  • end (float or None) – Data with x >= end will be removed.

Returns:

Series with enforced minimum and maximum

Return type:

BaseSeries or derived class

differentiate(order=1)

Differentiate with the numerical order-differentiation.

The optional parameter order specifies the order of the derivative.

The derivative is calulated as centered differencing in the interior and one-sided derivatives at the boundaries. Higher orders are computed applying the same rule recursively.

Parameters:

order (int) – Order of derivative (e.g. 2 = second derivative).

differentiated(order=1)

Return a series that is the numerical order-differentiation of the present series.

The optional parameter order specifies the order of the derivative.

The derivative is calulated as centered differencing in the interior and one-sided derivatives at the boundaries. Higher orders are computed applying the same rule recursively.

Parameters:

order (int) – Order of derivative (e.g. 2 = second derivative).

Returns:

New series with derivative.

Return type:

BaseSeries or derived class

property dt

Return the timestep if the series is regularly sampled, otherwise raise error.

Returns:

Timestep of the series (if evenly sampled).

Return type:

float

property duration

Return the length of the covered time interval.

Returns:

Length of time covered by the timeseries (tmax - tmin).

Return type:

float

evaluate_with_spline(x, ext=2)

Evaluate the spline on the points x.

Values outside the interval are extrapolated if ext=0, set to 0 if ext=1, raise a ValueError if ext=2, or if ext=3, return the boundary value.

This method is meant to be used only if you want to use a different ext for a specific call, otherwise, just use __call__.

Parameters:
  • x (1D NumPy array of float) – Array of x where to evaluate the series or single x.

  • ext (int) – How to deal values outside the bounaries. Values outside the interval are extrapolated if ext=0, set to 0 if ext=1, raise a ValueError if ext=2, or if ext=3, return the boundary value.

Returns:

Values of the series evaluated on the input x.

Return type:

1D NumPy array or float

final_time_remove(time_end)[source]

Remove the final time_end amount of time in the timeseries.

Parameters:

time_end (float) – Amount of time to be removed from the end.

final_time_removed(time_end)[source]

Return a TimeSeries without the final time_end amount of time.

If a series goes from t=-1 to t=10 and you set time_end=2, the series will go from t=-1 to t=8.

Parameters:

time_end (float) – Amount of time to be removed from the end.

Returns:

A new TimeSeries without the final time_end.

Return type:

TimeSeries

fixed_frequency_resample(frequency)[source]

Resample the timeseries to regularly spaced times with the given frequency. The final time will change if the frequency does not lead a integer number of timesteps.

Parameters:

frequency (float) – Sampling rate.

fixed_frequency_resampled(frequency)[source]

Return a TimeSeries with same tmin and tmax but resampled at a fixed frequency. The final time will change if the frequency does not lead a integer number of timesteps.

Parameters:

frequency (float) – Sampling rate.

Returns:

Time series resampled with given frequency.

Return type:

TimeSeries

fixed_timestep_resample(timestep)[source]

Resample the timeseries to regularly spaced times with given timestep. The final time will change if the timestep does not lead a integer number of timesteps.

Parameters:

timestep (float) – New timestep.

fixed_timestep_resampled(timestep)[source]

Return a new TimeSeries with evenly spaced with the given timestep. The final time will change if the timestep does not lead a integer number of timesteps.

Parameters:

timestep (float) – New timestep.

Returns:

Time series resampled with given timestep.

Return type:

TimeSeries

hamming_window()[source]

Apply Hamming window.

hamming_windowed()[source]

Return a timeseries with Hamming window applied.

Returns:

New windowed TimeSeries.

Return type:

TimeSeries

property index

Fake pandas properties, to make Series objects plottable by matplotlib.

initial_time_remove(time_init)[source]

Remove the first time_init amount of time in the timeseries.

When tmin = 0, this is the same as cropping, otherwise the difference is that in one case the time interval is specified, whereas in the other (cropping) the new tmin is specified.

Parameters:

time_init (float) – Amount of time to be removed from the beginning.

initial_time_removed(time_init)[source]

Return a TimeSeries without the initial time_init amount of time.

When tmin = 0, this is the same as cropping, otherwise the difference is that in one case the time interval is specified, whereas in the other (cropping) the new tmin is specified.

If a series goes from t=-1 to t=10 and you set time_init=2, the series will go from t=1 to t=10.

Parameters:

time_init (float) – Amount of time to be removed from the beginning.

Returns:

A new TimeSeries without the initial time_init.

Return type:

TimeSeries

integrate(dx=None)

Integrate series with method of the rectangles.

The spacing dx can be optionally provided. If provided, it will be used (increasing performance), otherwise it will be computed internally.

integrated(dx=None)

Return a series that is the integral computed with method of the rectangles.

The spacing dx can be optionally provided. If provided, it will be used (increasing performance), otherwise it will be computed internally.

Parameters:

dx (float or None) – Delta x in the independent variable. If None it will be computed internally.

Returns:

New series with the cumulative integral.

Return type:

BaseSeries or derived class

is_complex()

Return whether the data is complex.

Returns:

True if the data is complex, false if it is not.

Return type:

bool

is_masked()

Return whether the x or y are masked.

Returns:

True if the x or y are masked, false if it is not.

Return type:

bool

is_regularly_sampled()

Return whether the series is regularly sampled.

If the series is only one point, an error is raised.

Returns:

Is the series regularly sampled?

Return type:

bool

local_maxima(*args, include_edges=True, **kwargs)

Use SciPy’s find_peaks to find the local maxima.

Unkown arguments are passed to find_peaks.

If the signal is complex, the absolute value is taken.

If include_edges is True, the edges are considered among the possible maxima.

Returns:

Coordinate and value of the peaks.

Return type:

Tuple of NumPy arrays

local_minima(*args, include_edges=True, **kwargs)

Use SciPy’s find_peaks to find the local minima.

Unkown arguments are passed to find_peaks.

If the signal is complex, the absolute value is taken.

If include_edges is True, the edges are considered among the possible minima.

Returns:

Coordinate and value of the minima.

Return type:

Tuple of NumPy arrays

property mask

Return where the data is valid (according to the mask).

Returns:

Array of True/False of the same length of the data. False where the data is valid, true where is not.

Return type:

1D array of bool

mask_applied(mask, ignore_existing=False)

Return a new series with given mask applied to the data.

If a previous mask already exists, the new mask will be added on top, unless ignore_existing is True.

Parameters:
  • mask (1D NumPy array) – Array of booleans that identify where the data is invalid. This can be obtained with the method mask().

  • ignore_existing (bool) – If True, overwrite any previously existing mask.

Returns:

New series with mask applied.

Return type:

BaseSeries

mask_apply(mask, ignore_existing=False)

Apply given mask.

If a previous mask already exists, the new mask will be added on top, unless ignore_existing is True.

Parameters:
  • mask (1D NumPy array) – Array of booleans that identify where the data is invalid. This can be obtained with the method mask().

  • ignore_existing (bool) – If True, overwrite any previously existing mask.

mask_equal(value)

Mask where data is equal to given value.

mask_greater(value)

Mask where data is greater to given value.

mask_greater_equal(value)

Mask where data is greater or equal to given value.

mask_inside(value1, value2)

Mask where data is inside the given values.

mask_invalid()

Mask where data is invalid (NaNs of infs).

mask_less(value)

Mask where data is less to given value.

mask_less_equal(value)

Mask where data is less or equal to given value.

mask_not_equal(value)

Mask where data is not equal to given value.

mask_outside(value1, value2)

Mask where data is outside the given values.

mask_remove()

Remove masked values.

mask_removed()

Remove masked value.

Return a new series with valid values only.

Returns:

A new series with only valid values.

Return type:

BaseSeries or derived class

masked_equal(value)

Return a new objected masked where data is equal to given value.

masked_greater(value)

Return a new objected masked where data is greater to given value.

masked_greater_equal(value)

Return a new objected masked where data is greater or equal to given value.

masked_inside(value1, value2)

Return a new objected masked where data is inside the given values.

masked_invalid()

Return a new objected masked where data is invalid (NaNs or infs).

masked_less(value)

Return a new objected masked where data is less to given value.

masked_less_equal(value)

Return a new objected masked where data is less or equal to given value.

masked_not_equal(value)

Return a new objected masked where data is not equal to given value.

masked_outside(value1, value2)

Return a new objected masked where data is outside the given values.

mean_remove()[source]

Remove the mean value from the data.

mean_removed()[source]

Return a TimeSeries with mean removed, so that its new total average is zero.

Returns:

A new TimeSeries with zero mean.

Return type:

TimeSeries

nans_remove()

Filter out nans/infinite values.

nans_removed()

Filter out nans/infinite values. Return a new series with finite values only.

Returns:

A new series with only finite values.

Return type:

BaseSeries or derived class

phase_angular_velocity(use_splines=True, tsmooth=None, order=3)[source]

Compute the phase angular velocity, i.e. the time derivative of the complex phase.

Optionally smooth the with a savgol filter with smoothing length tsmooth and order order. If you do so, the timeseries is resampled to regular timesteps.

Parameters:
  • use_splines (bool) – Wheter to use splines of finite differencing for the derivative.

  • tsmooth (float) – Time over which smoothing is applied.

  • order (int) – Order of the for the savgol smoothing.

Returns:

Time derivative of the complex phase.

Return type:

TimeSeries

phase_frequency(use_splines=True, tsmooth=None, order=3)[source]

Compute the phase frequency, i.e. the time derivative of the complex phase divided by 2 pi.

Optionally smooth the with a savgol filter with smoothing length tsmooth and order order. If you do so, the timeseries is resampled to regular timesteps.

Parameters:
  • use_splines (bool) – Wheter to use splines of finite differencing for the derivative.

  • tsmooth (float) – Time over which smoothing is applied.

  • order (int) – Order of the for the savgol smoothing.

Returns:

Time derivative of the complex phase divided by 2 pi

Return type:

TimeSeries

phase_shift(pshift)[source]

Shift the complex phase timeseries by pshift. If the signal is real, it is turned complex with phase of pshift.

Parameters:

pshift (float) – Amount of phase to shift.

phase_shifted(pshift)[source]

Return a new TimeSeries with complex phase shifted by pshift. If the signal is real, it is turned complex with phase of pshift.

Parameters:

pshift (float) – Amount of phase to shift.

Returns:

A new TimeSeries with phase shifted.

Return type:

TimeSeries

redshift(z)[source]

Apply redshift to the data by rescaling the time so that the frequencies are redshifted by 1 + z.

Parameters:

z (float) – Redshift factor.

redshifted(z)[source]

Return a new TimeSeries with time rescaled so that frequencies are redshifted by 1 + z.

Parameters:

z (float) – Redshift factor.

Returns:

A new redshifted TimeSeries.

Return type:

TimeSeries

regular_resample()[source]

Resample the timeseries to regularly spaced times, with the same number of points.

regular_resampled()[source]

Return a new timeseries resampled to regularly spaced times, with the same number of points.

Returns:

Regularly resampled time series.

Return type:

TimeSeries

resample(new_x, ext=2, piecewise_constant=False)

Resample the series to new independent variable new_x.

If you want to resample without using the spline, and you want a nearest neighbor resampling, pass the keyword piecewise_constant=True. This may be a good choice for data with large discontinuities, where the splines are ineffective.

Parameters:
  • new_x (1D NumPy array or list of float) – New independent variable.

  • ext (0 for extrapolation, 1 for returning zero, 2 for ValueError, 3 for extending the boundary) – How to handle points outside the interval.

  • piecewise_constant (bool) – Do not use splines, use the nearest neighbors.

resampled(new_x, ext=2, piecewise_constant=False)

Return a new series resampled from this to new_x.

You can specify the details of the spline with the method make_spline.

If you want to resample without using the spline, and you want a nearest neighbor resampling, pass the keyword piecewise_constant=True. This may be a good choice for data with large discontinuities, where the splines are ineffective.

Parameters:
  • new_x (1D NumPy array or list of float) – New independent variable.

  • ext (0 for extrapolation, 1 for returning zero, 2 for ValueError, 3 for extending the boundary) – How to handle points outside the data interval.

  • piecewise_constant (bool) – Do not use splines, use the nearest neighbors.

Returns:

Resampled series.

Return type:

BaseSeries or derived class

save(file_name, *args, **kwargs)

Saves into simple ASCII format with 2 columns (x, y) for real valued data and 3 columns (x, Re(y), Im(y)) for complex valued data.

Unknown arguments are passed to NumPy.savetxt.

Parameters:

file_name (str) – Path (with extension) of the output file.

savgol_smooth(window_size, order=3)

Smooth the series with a Savitzky-Golay filter with window of size window_size and order order.

This is just like a regular “Moving average” filter, but instead of just calculating the average, a polynomial (usually 2nd or 4th order) fit is made for every point, and only the “middle” point is chosen. Since 2nd (or 4th) order information is concerned at every point, the bias introduced in “moving average” approach at local maxima or minima, is circumvented.

Parameters:
  • window_size (int) – Number of points of the smoothing window (needs to be odd).

  • order (int) – Order of the filter.

savgol_smooth_time(tsmooth, order=3)[source]

Resample the timeseries with uniform timesteps, smooth it with savgol_smooth with a window that is tsmooth in time (as opposed to a number of points).

Parameters:
  • tsmooth (float) – Time interval over which to smooth.

  • order (int) – Order of the filter.

savgol_smoothed(window_size, order=3)

Return a smoothed series with a Savitzky-Golay filter with window of size window_size and order order.

This is just like a regular “Moving average” filter, but instead of just calculating the average, a polynomial (usually 2nd or 4th order) fit is made for every point, and only the “middle” point is chosen. Since 2nd (or 4th) order information is concerned at every point, the bias introduced in “moving average” approach at local maxima or minima, is circumvented.

Parameters:
  • window_size (int) – Number of points of the smoothing window (needs to be odd).

  • order (int) – Order of the filter.

Returns:

New smoothed series.

Return type:

BaseSeries or derived class

savgol_smoothed_time(tsmooth, order=3)[source]

Return a resampled timeseries with uniform timesteps, smoothed with savgol_smooth with a window that is tsmooth in time (as opposed to a number of points).

Parameters:
  • tsmooth (float) – Time interval over which to smooth.

  • order (int) – Order of the filter.

Returns:

New smoothed and resampled TimeSeries.

Return type:

TimeSeries

spline_differentiate(order=1)

Differentiate the series using the spline representation.

The optional parameter order specifies the order of the derivative.

Warning

The values at the boundary are typically not accurate.

Parameters:

order (int) – Order of derivative (e.g. 2 = second derivative).

spline_differentiated(order=1)

Return a series that is the derivative of the current one using the spline representation.

The optional parameter order specifies the order of the derivative.

Warning

The values at the boundary are typically not accurate.

Parameters:

order (int) – Order of derivative (e.g. 2 = second derivative).

Returns:

New series with derivative

Return type:

BaseSeries or derived class

property t

Return the time.

Returns:

Times.

Return type:

1d NumPy array.

time_at_maximum(absolute=True)[source]

Return the time at which the timeseries is maximum.

Parameters:

absolute (bool) – Whether to take the absolute value of the data.

Returns:

Time at maximum. If absolute is True, then time at absolute maximum.

Return type:

float

time_at_minimum(absolute=True)[source]

Return the time at which the timeseries is minimum.

Parameters:

absolute (bool) – Whether to take the absolute value of the data.

Returns:

Time at minimum. If absolute is True, then time at absolute minimum.

Return type:

float

property time_length

Return the length of the covered time interval.

Returns:

Length of time covered by the timeseries (tmax - tmin).

Return type:

float

time_shift(tshift)[source]

Shift the timeseries by tshift so that what was t = 0 will be tshift.

Parameters:

N (float) – Amount of time to shift.

time_shifted(tshift)[source]

Return a new timeseries with time shifted by tshift so that what was t = 0 will be tshift.

Parameters:

tshift (float) – Amount of time to shift.

Returns:

A new TimeSeries with time shifted.

Return type:

TimeSeries

time_unit_change(unit, inverse=False)[source]

Rescale time units by unit.

This amounts to sending t to t / unit. For example, if initially the units where seconds, with unit=1e-3 the new units will be milliseconds.

When inverse is True, the opposite is done and t is sent to t * unit. This is useful to convert geometrized units to physical units with unitconv. For example,

# Gravitational waves in geometrized units
gw_cu = TimeSeries(...)
# Gravitational waves in seconds, assuming a mass of 1 M_sun
CU = uc.geom_umass_msun(1)
gw_s = gw_cu.time_unit_changed(CU.time, inverse=True)
Parameters:
  • unit (float) – New time unit.

  • inverse (bool) – If True, time = 1 -> time = unit, otherwise time = unit -> 1.

time_unit_changed(unit, inverse=False)[source]

Return a new TimeSeries with time scaled by unit.

This amounts to sending t to t / unit. For example, if initially the units where seconds, with unit=1e-3 the new units will be milliseconds.

When inverse is True, the opposite is done and t is sent to t * unit. This is useful to convert geometrized units to physical units with unitconv. For example,

# Gravitational waves in geometrized units
gw_cu = TimeSeries(...)
# Gravitational waves in seconds, assuming a mass of 1 M_sun
CU = uc.geom_umass_msun(1)
gw_s = gw_cu.time_unit_changed(CU.time, inverse=True)
Parameters:
  • unit (float) – New time unit.

  • inverse (bool) – If True, time = 1 -> time = unit, otherwise time = unit -> 1.

Returns:

A TimeSeries with new time unit.

Return type:

TimeSeries

property tmax

Return the final time.

Returns:

Final time of the timeseries.

Return type:

float

property tmin

Return the starting time.

Returns:

Initial time of the timeseries.

Return type:

float

to_FrequencySeries()[source]

Return a FrequencySeries that is the Fourier transform of the timeseries.

If the signal is not complex, only positive frequencies are kept.

If the timeseries is not regularly sampled, it will be resampled before transforming.

:: warning:

To have meaningful results, you should consider removing the mean and windowing the signal before calling this method!

Returns:

Fourier Transform.

Return type:

FrequencySeries

to_numpy()

Return the data as NumPy array. Equivalent to self.y.

This function is here to enable compatibility matplotlib.

tukey_window(alpha)[source]

Apply Tukey window with parameter alpha.

Parameters:

alpha (float) – Tukey parameter.

tukey_windowed(alpha)[source]

Return a TimeSeries with Tukey window with parameter alpha applied.

Parameters:

alpha (float) – Tukey parameter.

Returns:

New windowed TimeSeries.

Return type:

TimeSeries

unfolded_phase(t_of_zero_phase=None)[source]

Compute the complex phase of a complex-valued signal such that no phase wrap-around occur, i.e. if the input is continuous, so is the output. Optionally, add a phase shift such that phase is zero at the given time.

Parameters:

t_of_zero_phase (float or None) – Time at which the phase is set to zero.

Returns:

Continuous complex phase.

Return type:

TimeSeries

property values

Fake pandas properties, to make Series objects plottable by matplotlib.

window(window_function, *args, **kwargs)[source]

Apply window_function to the data.

window_function has to be a function that takes as first argument the number of points of the signal. window_function can take additional arguments as passed by windowed. Alternatively, window_function can be a string with the name of the window function, if this is already implemented in TimeSeries (e.g., tukey).

Parameters:

window_function (callable or str) – Window function to apply to the timeseries.

windowed(window_function, *args, **kwargs)[source]

Return a TimeSeries windowed with window_function.

window_function has to be a function that takes as first argument the number of points of the signal. window_function can take additional arguments as passed by windowed. Alternatively, window_function can be a string with the name of the window function, if this is already implemented in TimeSeries (e.g., tukey).

Parameters:

window_function (callable or str) – Window function to apply to the timeseries.

Returns:

New windowed TimeSeries.

Return type:

TimeSeries

x_at_abs_maximum_y()

Return the value of x when abs(y) is maximum.

Returns:

Value of x when abs(y) is maximum.

Return type:

float

x_at_abs_minimum_y()

Return the value of x when abs(y) is minimum.

Returns:

Value of x when abs(y) is minimum.

Return type:

float

x_at_maximum_y()

Return the value of x when y is maximum.

Returns:

Value of x when y is maximum.

Return type:

float

x_at_minimum_y()

Return the value of x when y is minimum.

Returns:

Value of x when y is minimum.

Return type:

float

property xmax

Return the maximum of the independent variable x.

Rvalue:

Maximum of x

Return type:

float

property xmin

Return the minimum of the independent variable x.

Rvalue:

Minimum of x.

Return type:

float

zero_pad(N)[source]

Pad the timeseries with zeros so that it has a total of N points.

This operation will work only if the timeseries is equispaced and if N is larger than the number of points already present.

Note

N is the final number of points, not the number of points added.

Parameters:

N (int) – Total number new points with zeros at the end.

zero_padded(N)[source]

Return a TimeSeries that is zero-padded and that has in total N points.

Note

N is the final number of points, not the number of points added.

This operation will work only if the series is equispaced.

Parameters:

N (int) – Total number of points of the output.

Returns:

A new timeseries with in total N points where all the trailing ones are zero.

Return type:

TimeSeries

kuibit.timeseries.combine_ts(series, prefer_late=True)[source]

Combine several overlapping time series into one.

In intervals covered by two or more time series, which data is used depends on the parameter prefer_late. If two segments start at the same time, the longer one gets used.

Parameters:
  • series (list of TimeSeries) – The timeseries to combine.

  • prefer_late – If true, prefer data that starts later for overlapping segments, otherwise, use data from the ones that come earlier.

Returns:

The combined time series

Return type:

TimeSeries

kuibit.timeseries.remove_duplicated_iters(t, y)[source]

Remove overlapping segments from a time series in (t,y).

Only the latest of overlapping segments is kept, the rest removed.

This function is used for cleaning up simulations with multiple checkpoints.

Note, if t = [1, 2, 3, 4, 2, 3] the output will be [1, 2, 3]. The ‘4’ is discarded because it is not the last segment. The idea is that if this corresponds to a simulation restart, you may have changed the paramters, so that 4 is not anymore correct. We consider the second restart the “truth”.

Parameters:
  • t (1D NumPy array) – Times.

  • y – Values.

Returns:

Strictly monotonic time series.

Return type:

TimeSeries

kuibit.timeseries.unfold_phase(phase)[source]

Remove phase jumps to get a continuous (unfolded) phase.

Parameters:

phase (1D NumPy array) – Phase wrapped around the provided jump.

Returns:

Phase plus multiples of pi chosen to minimize jumps.

Return type:

1D NumPy array