Data Generation for Synthetic Data

Created on Tue Mar 1 11:01:26 2022

@author: Tobias Westmeier CR/AMP4

The described scripts are used for the synthetic generation of time data for testing models. The class: ‘get_academic_data’ is used for the numerical integration of the described differential equations.

softsensor.data_gen.get_academic_data(time, Model, params, F=None, x0=None, rtol=0.001)[source]

Define Dataframe with the response of academic examples equation to given force

Parameters:
  • time (array) – A sequence of time points for which to solve for y. The initial value point should be the first element of this sequence. This sequence must be monotonically increasing or monotonically decreasing; repeated values are allowed.

  • Model (str) – either ‘Duffing’ ‘Duffing_fp’ or ‘vd_pol’ to define academic model.

  • params (dict) – Define Parameters for the models. if Model = ‘Duffing’ dict must contain c_nlin and D. if Model = ‘vd_Pol’ dict must contain epsilon.

  • F (Force, optional) – Force to apply to the ems. Must be callable by F.comp and F.comp_dt and must include self.F and self.dF with the same time steps as time The default is None.

  • x0 (array, optional) – initial conditions. The default is [0, 0].

  • rtol (float, optional) – relative tolerance for the numerical integration

Returns:

df – Dataframe of response with columns [‘F(t)’, ‘x’, ‘v’].

Return type:

pd.DataFrame

Examples

>>> import softsensor.data_gen as dg
>>> import numpy as np
>>> time = np.linspace(0, 100, 1000)
>>> params = {'D': 0.05, 'c_nlin': 0.1}
>>> F = dg.Zero(time)
>>> df = dg.get_academic_data(time, Model='Duffing', F=F, params=params)
>>> print(df.shape)
(1000, 3)
softsensor.data_gen.duffing(t, y, D, F, c_nlin=0)[source]

Duffing oscillator with external forcing as described in https://en.wikipedia.org/wiki/Duffing_equation

\[dz1 = z2\]
\[dz2 = - 2 \cdot D \cdot z2 - z1 - c_{nlin} \cdot z1^{3} + F.comp(t)\]
Parameters:
  • y (array of size 2) – initial conditions. System state (z1, z2)

  • t (float()) – time step at which the force F is computed.

  • D (float()) – damping of the system.

  • F (Force class) – Abstract class for the force excitation at the mass.

  • c_nlin (float(), optional) – indicating the nonlinearity of the system. The default is 0.

Returns:

(dz1, dz2)

Return type:

Derivatives of the duffing equation.

Examples

>>> import softsensor.data_gen as dg
>>> import numpy as np
>>> time = np.linspace(0, 100, 1000)
>>> F = dg.Zero(time)
>>> (dz1, dz2) = dg.duffing(t=0, y=(0, 0), D=0.05, F=F, c_nlin=.2)
>>> print((dz1, dz2))
(0, 0.0)
softsensor.data_gen.pendulum(t, y, D, F)[source]

Pendulum equation with external forcing and damping https://en.wikipedia.org/wiki/Pendulum_(mechanics)

\[dz1 = z2\]
\[dz2 = - 2 \cdot D \cdot z2 - sine(z1) + F.comp(t)\]
Parameters:
  • y (array of size 2) – initial conditions. System state (z1, z2)

  • t (float()) – time step at which the force F is computed.

  • D (float()) – damping of the system.

  • F (Force class) – Abstract class for the force excitation at the mass.

Returns:

(dz1, dz2)

Return type:

Derivatives of the pendulum equation.

Examples

>>> import softsensor.data_gen as dg
>>> import numpy as np
>>> time = np.linspace(0, 100, 1000)
>>> F = dg.Zero(time)
>>> (dz1, dz2) = dg.pendulum(t=0, y=(0, 0), D=0.05, F=F, c_nlin=.2)
>>> print((dz1, dz2))
(0, 0.0)
softsensor.data_gen.duffing_fp(t, y, D, c_nlin, z)[source]

Duffing oscillator with excitation at the base point as described in https://en.wikipedia.org/wiki/Duffing_equation

\[dz1 = z2\]
\[dz2 = 2 \cdot D \cdot (z.comp_dt(t) - z2) - z1 + c_{nlin} \cdot (z.comp(t) - z1)^3\]
Parameters:
  • y (array of size 2) – initial conditions. System state (z1, z2)

  • t (float()) – time step at which the Force is computed.

  • D (float()) – damping of the system.

  • c_nlin (float()) – indicating the nonlinearity of the system.

  • F (Force class) – Abstract class for the force excitation at the mass.

Returns:

(dz1, dz2)

Return type:

Derivatives of the duffing equation.

Examples

>>> import softsensor.data_gen as dg
>>> import numpy as np
>>> time = np.linspace(0, 10, 100)
>>> F = dg.Zero(time)
>>> (dz1, dz2) = dg.duffing_fp(t=0, y=(0, 0), D=0.05, c_nlin=.2, z=F)
>>> print((dz1, dz2))
(0, 0.0)
softsensor.data_gen.two_mass_oscillator(t, y, D, mue, kappa, delta, F)[source]

Two Mass oscillator with damping and Forcing on the first mass Coverning Equation:

\[\begin{split}\begin{bmatrix} 1 & 0 \\ 0 & \mu \end{bmatrix} \begin{bmatrix} q_1^{''}\\ q_2^{''} \end{bmatrix} + 2D \begin{bmatrix} 1 + \delta & -\delta \\ -\delta & \delta \end{bmatrix} \begin{bmatrix} q_1^{'}\\ q_2^{'} \end{bmatrix} + \begin{bmatrix} 1 + \kappa & -\kappa \\ -\kappa & \kappa \end{bmatrix} \begin{bmatrix} q_1\\ q_2 \end{bmatrix} = \begin{bmatrix} F.comp(t)\\ 0 \end{bmatrix}\end{split}\]

with

\[ \begin{align}\begin{aligned}z_{11} = q_{1}, z_{12} = q_{1}^{'}\\z_{21} = q_{2}, z_{22} = q_{2}^{'}\end{aligned}\end{align} \]
Parameters:
  • y (matrix of size 4) – initial conditions. System state [z11, z12, z21, z22]

  • t (float()) – time step at which the Force is computed.

  • D (float()) – damping of the system.

  • mue (float()) – relative mass of the system m2/m1

  • kappa (float()) – relative restoring force of the system c2/c1

  • delta (float()) – relative damping of the system d2/d1

  • F (Force class) – Abstract class for the force excitation at the mass.

Returns:

([dz11, dz12], [dz21, dz22])

Return type:

Derivatives of the duffing equation.

Examples

>>> import softsensor.data_gen as dg
>>> import numpy as np
>>> time = np.linspace(0, 10, 100)
>>> F = dg.Zero(time)
>>> (dz1, dz2) = dg.duffing_fp(t=0, y=(0, 0), D=0.05, c_nlin=.2, z=F)
>>> print((dz1, dz2))
(0, 0.0)
softsensor.data_gen.vd_Pol(t, y, epsilon, offset, F)[source]

van der Pol oscillation as described in https://en.wikipedia.org/wiki/Van_der_Pol_oscillator

\[dz1 = z2\]
\[dz2 = epsilon \cdot (1 - z1^2) \cdot z2 - z1 + F.comp(t)\]
Parameters:
  • y (array of size 2) – initial conditions.

  • t (float()) – time step at which the Force is computed.

  • epsilon (float()) – indicating the nonlinearity and the strength of the damping.

  • F (Force class) – forcing term of van der Pol equation.

Returns:

(dz1, dz2)

Return type:

Derivatives of the duffing equation.

Examples

>>> import softsensor.data_gen as dg
>>> import numpy as np
>>> time = np.linspace(0, 10, 100)
>>> F = dg.Zero(time)
>>> (dz1, dz2) = dg.vd_Pol(t=0, y=(0, 0), epsilon=1, F=F)
>>> print((dz1, dz2))
(0, 0)
class softsensor.data_gen.sine(time, gamma=1, w0=1)[source]

Generate Sine Wave

\[gamma \cdot sin(w_0 \cdot time)\]
Parameters:
  • time (array) – time array at which the sines are computed.

  • gamma (float(), optional) – amplification sine wave. The default is 1.

  • w0 (float(), optional) – frequency of sine wave. The default is 1.

Return type:

None.

comp(t)[source]

Compute sine at specific timestep

\[gamma \cdot sin(w_0 \cdot t)\]
Parameters:

t (float()) – time step to compute sine at.

Returns:

w – sine at the specific time step.

Return type:

float()

comp_dt(t)[source]

Compute derivative of the sine wave at specific timestep

\[gamma \cdot w_0 \cdot cos(w_0 \cdot t)\]
Parameters:

t (float()) – time step to compute derivative at.

Returns:

w – derivative at the time step.

Return type:

float()

class softsensor.data_gen.sine_2(time, gamma1=1, w01=1, gamma2=0.5, w02=26.4)[source]

Generate the addition of two sine oscillation

\[gamma1 \cdot sin(w_{01} \cdot time) + gamma2 \cdot sin(w_{02} \cdot time)\]
Parameters:
  • time (array) – time array at which the sines are computed.

  • gamma1 (float(), optional) – amplification of first sine wave. The default is 1.

  • w01 (float(), optional) – frequency of first sine wave. The default is 1.

  • gamma2 (float(), optional) – amplification of second sine wave. The default is 0.5.

  • w02 (float(), optional) – frequency of second sine wave. The default is 26.4.

Return type:

None.

comp(t)[source]

Compute sines at specific timestep

\[gamma1 \cdot sin(w_{01} \cdot t) + gamma2 \cdot sin(w_{02} \cdot t)\]
Parameters:

t (float()) – time step to compute sines at.

Returns:

w – sines at the time step.

Return type:

float()

comp_dt(t)[source]

Compute derivative of sine waves at specific timestep

\[gamma1 \cdot w_{01} \cdot cos(w_{01} \cdot t) + gamma2 \cdot w_{02} \cdot cos(w_{02} \cdot t)\]
Parameters:

t (float()) – time step to compute derivative at.

Returns:

w – derivative at the time step.

Return type:

float()

class softsensor.data_gen.sweep(time, f0, f1, t1, method='linear', direction='up', f=1)[source]

Define sweep signal and derivative based on the scipy implementation: https://docs.scipy.org/doc/scipy/reference/generated/scipy.signal.chirp.html

Parameters:
  • time (array) – time array at which the sweep is computed.

  • f0 (float()) – start frequency.

  • f1 (float()) – stop frequency.

  • t1 (float()) – time at which end frequency is met.

  • method (str, optional) – Kind of frequency sweep. The default is ‘linear’. options are: [‘linear’, ‘quadratic’, ‘logarithmic’, ‘hyperbolic’]

  • direction (str, optional) – either ‘up’ or ‘down’. The default is ‘up’.

  • f (float(), optional) – Amplification factor of the sweep. The default is 1.

Return type:

None.

comp(t)[source]

Compute sweep at specific timestep

Parameters:

t (float()) – time step to compute sweep at.

Returns:

w – sweep at the time step.

Return type:

float()

comp_dt(t)[source]

Compute derivative of the sweep at specific timestep

Parameters:

t (float()) – time step to compute derivative at.

Returns:

w – derivative at the time step.

Return type:

float()

class softsensor.data_gen.white_noise(time, f=1, seed=None)[source]

Generate White Noise signal for a defined time with numpy.random.randn method: https://numpy.org/doc/stable/reference/random/generated/numpy.random.randn.html

Parameters:
  • time (array) – time array at which the white noise is computed.

  • f (float(), optional) – Amplification factor of the white noise. The default is 1.

  • seed (int, optional) – Seed for the random number generator to ensure reproducibility. The default is None, i.e. no seed is provided.

Return type:

None.

comp(t)[source]

Compute white noise at specific timestep

Parameters:

t (float()) – time step to compute white noise at.

Returns:

w – white noise at the time step.

Return type:

float()

comp_dt(t)[source]

Compute derivative of white noise at specific timestep. Derivation is realised via interpolatable splines

Parameters:

t (float()) – time step to compute derivative at.

Returns:

w – derivative at the time step.

Return type:

float()

class softsensor.data_gen.Zero(time)[source]

Generate an array of Zeros as input Force

Parameters:

time (array) – an array of the same length as the array of zeros

Return type:

None.

comp(t)[source]

Give back zero

Parameters:

t (float()) – time step.

Returns:

zero.

Return type:

int

comp_dt(t)[source]

Give back zero

Parameters:

t (float()) – time step.

Returns:

zero.

Return type:

int