SLiCAPmath.py

SLiCAP module with math functions.

Imported by the module SLiCAPprotos.py.

PdBm2V(p, r)

Returns the RMS value of the voltage that generates p dBm power in a resistor with resistance r.

Parameters
  • p (sympy.Symbol, sympy.Expression, int, or float) – Power in dBm

  • r (sympy.Symbol, sympy.Expression, int, or float) – Resistance

Returns

voltage

Return type

sympy.Expression

assumePosParams(expr, params='all')

Returns the sympy expression ‘expr’ in which variables, except the Laplace variable, have been redefined as real.

Parameters
  • expr (sympy.Expr, sympy.Symbol) – Sympy expression

  • params (list, str) – List with variable names (str), or ‘all’ or a variable name (str).

Returns

Expression with redefined variables.

Return type

sympy.Expr, sympy.Symbol

assumeRealParams(expr, params='all')

Returns the sympy expression ‘expr’ in which variables, except the Laplace variable, have been redefined as real.

Parameters
  • expr (sympy.Expr, sympy.Symbol) – Sympy expression

  • params (list, str) – List with variable names (str), or ‘all’ or a variable name (str).

Returns

Expression with redefined variables.

Return type

sympy.Expr, sympy.Symbol

besselPoly(n)

Returns a normalized Bessel polynomial of the n-th order of the Laplace variable.

Parameters

n (int) – order

Returns

Bessel polynomial of the n-th order of the Laplace variable

Return type

sympy.Expression

butterworthPoly(n)

Returns a narmalized Butterworth polynomial of the n-th order of the Laplace variable.

Parameters

n (int) – order

Returns

Butterworth polynomial of the n-th order of the Laplace variable

Return type

sympy.Expression

cancelPZ(poles, zeros)

Cancels poles and zeros that coincide within the displayed accuracy.

Note

The display accuracy (number of digits) is defined by ini.disp.

Parameters
  • poles (list) – List with poles (float) of a Laplace rational function.

  • zeros (list) – List with zeros (float) of a Laplace rational function.

Returns

Tuple with a list with poles (float) and a list with zeros (float).

Return type

Tuple with two lists,

checkExpression(expr)

Returns the sympy expression of expr.

Parameters

expr (str, sympy object, int, float) – argument that may represent a number or an expression.

Returns

sympy expression

Return type

int, float

checkNumber(var)

Returns a number with its value represented by var.

Parameters

var (str, sympy object, int, float) – Variable that may represent a number.

Returns

sympy expression

Return type

int, float

clearAssumptions(expr, params='all')

Returns the sympy expression ‘expr’ in which the assumtions ‘Real’ and ‘Positive’ have been deleted.

Parameters
  • expr (sympy.Expr, sympy.Symbol) – Sympy expression

  • params (list, str) – List with variable names (str), or ‘all’ or a variable name (str).

Returns

Expression with redefined variables.

Return type

sympy.Expr, sympy.Symbol

coeffsTransfer(rational, variable=s, normalize='lowest')

Returns a nested list with the coefficients of the variable of the numerator and of the denominator of ‘rational’.

The coefficients are in ascending order.

Parameters
  • rational (sympy.Expr) – Rational function of the variable.

  • variable (sympy.Symbol) – Variable of the rational function

  • normalize (str) –

    Normalization method:

    • ”highest”: the coefficients of the highest order of <variable> of the denominator will be noramalized to unity.

    • ”lowest”: the coefficients of the lowest order of <variable> of the denominator will be noramalized to unity.

Returns

Tuple with gain and two lists: [gain, numerCoeffs, denomCoeffs]

  1. gain (sympy.Expr): ratio of the nonzero coefficient of the lowest order of the numerator and the coefficient of the nonzero coefficient of the lowest order of the denominator.

  2. numerCoeffs (list): List with all coeffcients of the numerator in ascending order.

  3. denomCoeffs (list): List with all coeffcients of the denominator in ascending order.

Return type

tuple

dBmagFunc_f(LaplaceExpr, f)

Calculates the dB magnitude at the real frequency f (Fourier) from the univariate function ‘LaplaceExpr’ of the Laplace variable.

If ini.Hz == true, the Laplace variable will be replaced with 2*sp.pi*sp.I*ini.frequency.

If ini.Hz == False, the Laplace variable will be replaced with sp.I*ini.frequency.

Parameters
  • LaplaceExpr (sympy.Expr) – Univariate function of the Laplace variable.

  • f – Frequency value (float), or a numpy array with frequency values (float).

Returns

dB Magnitude at the specified frequency, or list with dB magnitudes at the specified frequencies.

Return type

float, numpy.array

delayFunc_f(LaplaceExpr, f, delta=0.0001)

Calculates the group delay at the real frequency f (Fourier) from the univariate function ‘LaplaceExpr’ of the Laplace variable.

If ini.Hz == true, the Laplace variable will be replaced with 2*sp.pi*sp.I*ini.frequency.

If ini.Hz == False, the Laplace variable will be replaced with sp.I*ini.frequency.

Parameters
  • LaplaceExpr (sympy.Expr) – Univariate function of the Laplace variable.

  • f – Frequency value (float), or a numpy array with frequency values (float).

Returns

Group delay at the specified frequency, or list with group delays at the specified frequencies.

Return type

float, numpy.array

det(M)

Returns the determinant of a square matrix ‘M’ calculated using recursive minor expansion (Laplace expansion). For large matrices with symbolic entries, this is faster than the built-in sympy.Matrix.det() method.

Parameters

M (sympy.Matrix) – Sympy matrix

Returns

Determinant of ‘M’

Return type

sympy.Expr

Example

>>> dim = 6
>>> M = sp.zeros(dim, dim)
>>> for i in range(dim):
>>>     for j in range(dim):
>>>         M[i,j] = sp.Symbol('m_' + str(i) + str(j))
>>> t1 = time()
>>> D1 = M.det()
>>> t2=time()
>>> D2 = determinant(M)
>>> t3=time()
>>> print('M.det :', t2-t1, 's')
>>> print('det(M):', t3-t2, 's')
>>> print(sp.expand(D1-D2))
M.det : 150.64841532707214 s
det(M): 0.4860818386077881 s
0
0
diff_n(expr, x, n)

Returns the n-th order derivative d^n/dx^n of expr(x).

Parameters
  • expr (sympy.Expr) – Input expression

  • x (sympy.Symbol) – Variable for differentiation

  • n (int) – Order of differentiation

Returns

d^n/dx^n of expr(x)

Return type

sympy.Expr

doCDS(noiseResult, tau)

Returns noiseResult after multiplying it with (2*sin(pi*ini.frequency*tau))^2

Parameters
  • noiseResult (sympy.Expr, sympy.Symbol, int or float) – sympy expression of a noise density spectrum in V^2/Hz or A^2/Hz

  • tau (sympy.Expr, sympy.Symbol, int or float) – Time between two samples

Returns

noiseResult*(2*sin(pi*ini.frequency*tau))^2

Return type

sympy.Expr, sympy.Symbol, int or float

doCDSint(noiseResult, tau, f_min, f_max)

Returns the integral from ini.frequency = f_min to ini.frequency = f_max, of a noise spectrum after multiplying it with (2*sin(pi*ini.frequency*tau))^2

Parameters
  • noiseResult (sympy.Expr, sympy.Symbol, int or float) – sympy expression of a noise density spectrum in V^2/Hz or A^2/Hz

  • tau (sympy.Expr, sympy.Symbol, int or float) – Time between two samples

  • f_min (sympy.Expr, sympy.Symbol, int or float) – Lower limit of the integral

  • f_max (sympy.Expr, sympy.Symbol, int or float) – Upper limit of the integral

Returns

integral of the spectrum from f_min to f_max after corelated double sampling

Return type

sympy.Expr, sympy.Symbol, int or float

equateCoeffs(protoType, transfer, noSolve=[], numeric=True)

Returns the solutions of the equation transferFunction = protoTypeFunction.

Both transfer and prototype should be Laplace rational functions. Their numerators should be polynomials of the Laplace variable of equal order and their denominators should be polynomials of the Laplace variable of equal order.

Parameters
  • protoType (sympy.Expr) – Prototype rational expression of the Laplace variable

  • transfer

Transfer fucntion of which the parameters need to be solved. The numerator and the denominator of this rational expression should be of the same order as those of the prototype.

Parameters
  • noSolve (list) – List with variables (str, sympy.core.symbol.Symbol) that do not need to be solved. These parameters will remain symbolic in the solutions.

  • numeric (bool) – True will force Maxima to use (big) floats for numeric values.

Returns

Dictionary with key-value pairs:

  • key: name of the parameter (sympy.core.symbol.Symbol)

  • value: solution of this parameter: (sympy.Expr, int, float)

Return type

dict

findServoBandwidth(loopgainRational)

Determines the intersection points of the asymptotes of the magnitude of the loopgain with unity.

Parameters

loopgainRational – Rational function of the Laplace variable, that represents the loop gain of a circuit.

Returns

Dictionary with key-value pairs:

  • hpf: frequency of high-pass intersection

  • hpo: order at high-pass intersection

  • lpf: frequency of low-pass intersection

  • lpo: order at low-pass intersection

  • mbv: mid-band value of the loopgain (highest value at order = zero)

  • mbf: lowest freqency of mbv

Return type

dict

float2rational(expr)

Converts floats in expr into rational numbers.

Parameters

expr (sympy.Expression) – Sympy expression in which floats need to be converterd into rational numbers.

Returns

expression in which floats have been replaced with rational numbers.

Return type

sympy.Expression

fullSubs(valExpr, parDefs)

Returns ‘valExpr’ after all parameters of ‘parDefs’ have been substituted into it recursively until no changes occur or until the maximum number of substitutions is achieved.

The maximum number opf recursive substitutions is set by ini.maxRexSubst.

Parameters
  • valExpr (sympy.Expr, sympy.Symbol, int, float) – Eympy expression in which the parameters should be substituted.

  • parDefs

    Dictionary with key-value pairs:

    • key (sympy.Symbol): parameter name

    • value (sympy object, int, float): value of the parameter

Returns

Expression or value obtained from recursive substitutions of parameter definitions into ‘valExpr’.

Return type

sympy object, int, float

gainValue(numer, denom)

Returns the zero frequency (s=0) value of numer/denom.

Parameters
  • numer (sympy.Expr) – Numerator of a rational function of the Laplace variable

  • denom (sympy.Expr) – Denominator of a rational function of the Laplace variable

Returns

zero frequency (s=0) value of numer/denom.

Return type

sympy.Expr

ilt(expr, s, t)

Returns the Inverse Laplace Transform f(t) of a rational expression F(s).

Parameters
  • expr (Sympy expression or str) – Rational function of the Laplace variable F(s).

  • s (sympy.Symbol) – Laplace variable

  • t (sympy.Symbol) – time variable

return: Inverse Laplace Transform f(t) :rtype: sympy.Expr

magFunc_f(LaplaceExpr, f)

Calculates the magnitude at the real frequency f (Fourier) from the univariate function ‘LaplaceExpr’ of the Laplace variable.

If ini.Hz == true, the Laplace variable will be replaced with 2*sp.pi*sp.I*ini.frequency.

If ini.Hz == False, the Laplace variable will be replaced with sp.I*ini.frequency.

Parameters
  • LaplaceExpr (sympy.Expr) – Univariate function of the Laplace variable.

  • f – Frequency value (float), or a numpy array with frequency values (float).

Returns

Magnitude at the specified frequency, or list with magnitudes at the specified frequencies.

Return type

float, numpy.array

mag_f(LaplaceExpr)

Returns the magnitude as a function of the real frequency f (Fourier) from the Laplace function ‘LaplaceExpr’.

If ini.Hz == true, the Laplace variable will be replaced with 2*sp.pi*sp.I*ini.frequency.

If ini.Hz == False, the Laplace variable will be replaced with sp.I*ini.frequency.

Parameters
  • LaplaceExpr (sympy.Expr) – Univariate function of the Laplace variable.

  • f – Frequency value (float), or a numpy array with frequency values (float).

Returns

Sympy expression representing the magnitude of the Fourier Transform.

Return type

sympy.Expr

makeNumData(yFunc, xVar, x)

Returns a list of values y, where y[i] = yFunc(x[i]).

Parameters
  • yFunc (sympy.Expr) – Function

  • xVar (sympy.Symbol) – Variable that needs to be substituted in yFunc

  • x (list) – List with values of x

Returns

list with y values: y[i] = yFunc(x[i]).

Return type

list

nonPolyCoeffs(expr, var)

Returns a dictionary with coefficients of negative and positive powers of var.

Parameters

var (sympy.Symbol) – Variable of which the coefficients will be returned

Returns

Dict with key-value pairs:

  • key: order of the variable

  • value: coefficient of that order

normalizeRational(rational, var=s, method='lowest')

Normalizes a rational expression to:

\[F(s) = gain\,s^{\ell} \frac{1+b_1s + ... + b_ms^m}{1+a_1s + ... + a_ns^n}\]
Parameters
  • Rational (sympy.Expr) – Rational function of the variable.

  • variable (sympy.Symbol) – Variable of the rational function

  • method (str) –

    Normalization method:

    • ”highest”: the coefficients of the highest order of <variable> of the denominator will be noramalized to unity.

    • ”lowest”: the coefficients of the lowest order of <variable> of the denominator will be noramalized to unity.

Returns

Normalized rational function of the variable.

Return type

sympy.Expr

numRoots(expr, var)

Returns the roots of the polynomial ‘expr’ with indeterminate ‘var’.

This function uses numpy for calculation of numeric roots.

Note

See: https://docs.scipy.org/doc/numpy/reference/generated/numpy.polynomial.polynomial.polyroots.html

Parameters
  • expr (sympy.Expr) – Univariate function.

  • var (sympy.Symbol) – Indeterminate of ‘expr’.

phaseFunc_f(LaplaceExpr, f)

Calculates the phase angle at the real frequency f (Fourier) from the univariate function ‘LaplaceExpr’ of the Laplace variable.

If ini.Hz == true, the Laplace variable will be replaced with 2*sp.pi*sp.I*ini.frequency.

If ini.Hz == False, the Laplace variable will be replaced with sp.I*ini.frequency.

Parameters
  • LaplaceExpr (sympy.Expr) – Univariate function of the Laplace variable.

  • f – Frequency value (float), or a numpy array with frequency values (float).

Returns

Angle at the specified frequency, or list with angles at the specified frequencies.

Return type

float, numpy.array

phaseMargin(LaplaceExpr)

Calculates the phase margin assuming a loop gain definition according to the asymptotic gain model.

This function uses scipy.newton() for determination of the the unity-gain frequency. It uses the function SLiCAPmath.findServoBandwidth() for the initial guess, and ini.disp for the relative accuracy.

if ini.Hz == True, the units will be degrees and Hz, else radians and radians per seconds.

Parameters

LaplaceExpr (sympy.Expr, list) – Univariate function (sympy.Expr*) or list with univariate functions (sympy.Expr*) of the Laplace variable.

Returns

Tuple with phase margin (float) and unity-gain frequency (float), or Tuple with lists with phase margins (float) and unity-gain frequencies (float).

Return type

tuple

phase_f(LaplaceExpr)

Calculates the magnitude as a function of the real frequency f (Fourier) from the Laplace function ‘LaplaceExpr’.

If ini.Hz == true, the Laplace variable will be replaced with 2*sp.pi*sp.I*ini.frequency.

If ini.Hz == False, the Laplace variable will be replaced with sp.I*ini.frequency.

Parameters
  • LaplaceExpr (sympy.Expr) – Univariate function of the Laplace variable.

  • f – Frequency value (float), or a numpy array with frequency values (float).

Returns

Sympy expression representing the phase of the Fourier Transform.

Return type

sympy.Expr

polyCoeffs(expr, var)

Returns a list with coefficients of ‘var’ in descending order.

Parameters
  • expr (sympy.Expr) – Sympy expression

  • var (sympy.Symbol) – Indeterminate of the polynomial.

Returns

List with coefficients (sympy.Expr) in descending order.

Return type

list

rational2float(expr)

Converts rational numbers in expr into floats.

Parameters

expr (sympy.Expression) – Sympy expression in which rational numbers need to be converterd intoc floats.

Returns

expression in which rational numbers have been replaced with floats.

Return type

sympy.Expression

rmsNoise(noiseResult, noise, fmin, fmax, source=None, CDS=False, tau=None)

Calculates the RMS source-referred noise or detector-referred noise, or the contribution of a specific noise source or a collection od sources to it.

Parameters
  • noiseResult (SLiCAPprotos.allResults) – Results of the execution of an instruction with data type ‘noise’.

  • noise – ‘inoise’ or ‘onoise’ for source-referred noise or detector- referred noise, respectively.

  • fmin (str, int, float, sp.Symbol) – Lower limit of the frequency range in Hz.

  • fmax (str, int, float, sp.Symbol) – Upper limit of the frequency range in Hz.

  • source (str, list) – refDes (ID) or list with IDs of noise sources of which the contribution to the RMS noise needs to be evaluated. Only IDs of current of voltage sources with a nonzero value for ‘noise’ are accepted.

  • CDS (Bool) – True if correlated double sampling is required, defaults to False If True parameter ‘tau’ must be given a nonzero finite value (can be symbolic)

  • tau (str, int, float, sp.Symbol) – CDS delay time

Returns

RMS noise over the frequency interval.

  • An expression or value if parameter stepping of the instruction is disabled.

  • A list with expressions or values if parameter stepping of the instruction is enabled.

Return type

int, float, sympy.Expr

roundN(expr, numeric=False)

Rounds all float and rational atoms in an expression to ini.disp digits, and converts integers into floats if their number of digits exceeds ini.disp

Parameters

expr (sympy.Expr) – Input expression

Returns

modified expression

Return type

sympy.Expr

routh(charPoly, eps=epsilon)

Returns the Routh array of a polynomial of the Laplace variable (ini.Laplace).

Parameters
  • charPoly (sympy.Expr) – Expression that can be written as a polynomial of the Laplace variable (ini.Laplace).

  • eps (sympy.Symbol) – Symbolic variable used to indicate marginal stability. Use a symbol that is not present in charPoly.

Returns

Routh array

Return type

sympy.Matrix

Example

>>> # ini.Laplace = sp.Symbol('s')
>>> s, eps = sp.symbols('s, epsilon')
>>> charPoly = s**4+2*s**3+(3+k)*s**2+(1+k)*s+(1+k)
>>> M = routh(charPoly, eps)
>>> print(M.col(0)) # Number of roots in the right half plane is equal to
>>>                 # the number of sign changes in the first column of the
>>>                 # Routh array
Matrix([[1], [2], [k/2 + 5/2], [(k**2 + 2*k + 1)/(k + 5)], [k + 1]])
step2PeriodicPulse(ft, t_pulse, t_period, n_periods)

Converts a step response in a periodic pulse response. Works with symbolic and numeric time functions.

For evaluation of numeric values, use the SLiCAP function: makeNumData().

Parameters
  • ft (sympy.Expr) – Time function f(t)

  • t_pulse (int, float) – Pulse width

  • t_period (int, float) – Pulse period

  • n_periods – Number of pulses

Typen_periods

int, float

Returns

modified time function

Return type

sympy.Expr